twoPhaseVoFSolver.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 "twoPhaseVoFSolver.H"
27 #include "fvcAverage.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace solvers
34 {
36 }
37 }
38 
39 
40 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41 
42 void Foam::solvers::twoPhaseVoFSolver::correctCoNum()
43 {
45 
46  const scalarField sumPhi
47  (
48  interface.nearInterface()().primitiveField()
49  *fvc::surfaceSum(mag(phi))().primitiveField()
50  );
51 
52  alphaCoNum =
53  0.5*gMax(sumPhi/mesh.V().primitiveField())*runTime.deltaTValue();
54 
55  const scalar meanAlphaCoNum =
56  0.5
57  *(gSum(sumPhi)/gSum(mesh.V().primitiveField()))
59 
60  Info<< "Interface Courant Number mean: " << meanAlphaCoNum
61  << " max: " << alphaCoNum << endl;
62 }
63 
64 
65 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
66 
68 {
69  interface.correct();
70 }
71 
72 
75 {
76  return interface.surfaceTensionForce();
77 }
78 
79 
80 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
81 
83 (
84  fvMesh& mesh,
86 )
87 :
88  twoPhaseSolver(mesh, mixturePtr),
89 
90  interface(mixture, alpha1, alpha2, U)
91 {
92  if (transient())
93  {
94  correctCoNum();
95  }
96 }
97 
98 
99 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
100 
102 {}
103 
104 
105 // ************************************************************************* //
const scalar meanAlphaCoNum
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:34
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.
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
const surfaceScalarField & phi
Reference to the mass-flux field.
Definition: VoFSolver.H:212
scalar alphaCoNum
Phase-fraction flux Courant number.
Definition: VoFSolver.H:88
virtual void correctCoNum()=0
Correct the cached Courant numbers.
Definition: VoFSolver.C:75
Solver module base-class for 2 immiscible fluids, with optional mesh motion and mesh topology changes...
Solver module base-class for 2 immiscible fluids using a VOF (volume of fluid) phase-fraction based i...
twoPhaseVoFSolver(fvMesh &mesh, autoPtr< twoPhaseVoFMixture >)
Construct from region mesh.
virtual tmp< surfaceScalarField > surfaceTensionForce() const
Return the interface surface tension force for the momentum equation.
virtual void correctInterface()
Correct the interface properties following mesh-change.
virtual ~twoPhaseVoFSolver()
Destructor.
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)
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)