multiphaseVoFSolver.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2023-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "multiphaseVoFSolver.H"
27 #include "localEulerDdtScheme.H"
28 #include "fvcAverage.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace solvers
35 {
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
42 
43 void Foam::solvers::multiphaseVoFSolver::correctCoNum()
44 {
46 
47  const scalarField sumPhi
48  (
49  mixture.nearInterface()().primitiveField()
50  *fvc::surfaceSum(mag(phi))().primitiveField()
51  );
52 
53  alphaCoNum =
54  0.5*gMax(sumPhi/mesh.V().primitiveField())*runTime.deltaTValue();
55 
56  const scalar meanAlphaCoNum =
57  0.5
58  *(gSum(sumPhi)/gSum(mesh.V().primitiveField()))
60 
61  Info<< "Interface Courant Number mean: " << meanAlphaCoNum
62  << " max: " << alphaCoNum << endl;
63 }
64 
65 
66 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
67 
69 {}
70 
71 
74 {
75  return mixture.surfaceTensionForce(U);
76 }
77 
78 
79 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
80 
82 (
83  fvMesh& mesh,
85 )
86 :
87  VoFSolver(mesh, autoPtr<VoFMixture>(mixturePtr.ptr())),
88 
90 
91  phases(mixture.phases())
92 {
93  if (transient())
94  {
95  correctCoNum();
96  }
97 }
98 
99 
100 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
101 
103 {}
104 
105 
106 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
107 
109 {
111 }
112 
113 
115 {
117 
118  alphaPredictor();
119 }
120 
121 
122 // ************************************************************************* //
const scalar meanAlphaCoNum
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:34
Two-phase VoF mixture.
Definition: VoFMixture.H:52
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Multiphase VoF mixture with support for interface properties.
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
const Time & runTime
Time.
Definition: solver.H:104
const fvMesh & mesh
Region mesh.
Definition: solver.H:101
Base solver module base-class for the solution of immiscible fluids using a VOF (volume of fluid) pha...
Definition: VoFSolver.H:73
const surfaceScalarField & phi
Reference to the mass-flux field.
Definition: VoFSolver.H:212
virtual void prePredictor()=0
Called at the start of the PIMPLE loop.
Definition: VoFSolver.C:250
scalar alphaCoNum
Phase-fraction flux Courant number.
Definition: VoFSolver.H:88
virtual void preSolve()
Called at the start of the time-step, before the PIMPLE loop.
Definition: VoFSolver.C:198
virtual void correctCoNum()=0
Correct the cached Courant numbers.
Definition: VoFSolver.C:75
Base solver module for the solution of multiple immiscible fluids using a VOF (volume of fluid) phase...
virtual tmp< surfaceScalarField > surfaceTensionForce() const
Return the interface surface tension force for the momentum equation.
virtual void prePredictor()
Called at the start of the PIMPLE loop.
multiphaseVoFSolver(fvMesh &mesh, autoPtr< multiphaseVoFMixture >)
Construct from region mesh.
multiphaseVoFMixture & mixture
Reference to the multiphaseVoFMixture.
virtual void preSolve()
Called at the start of the time-step, before the PIMPLE loop.
virtual void correctInterface()
Correct the interface properties following mesh-change.
A class for managing temporary objects.
Definition: tmp.H:55
Area-weighted average a surfaceField creating a volField.
U
Definition: pEqn.H:72
tmp< VolField< Type > > surfaceSum(const SurfaceField< Type > &ssf)
defineTypeNameAndDebug(compressibleMultiphaseVoF, 0)
Namespace for OpenFOAM.
Type gSum(const FieldField< Field, Type > &f)
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:129
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
messageStream Info
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensioned< scalar > mag(const dimensioned< Type > &)
Type gMax(const FieldField< Field, Type > &f)