twoPhaseSolver.H
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-2025 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 Class
25  Foam::solvers::twoPhaseSolver
26 
27 Description
28  Solver module base-class for 2 immiscible fluids, with optional
29  mesh motion and mesh topology changes including adaptive re-meshing.
30 
31  The momentum and other fluid properties are of the "mixture" and a single
32  momentum equation is solved.
33 
34  Either mixture or two-phase transport modelling may be selected. In the
35  mixture approach a single laminar, RAS or LES model is selected to model the
36  momentum stress. In the Euler-Euler two-phase approach separate laminar,
37  RAS or LES selected models are selected for each of the phases.
38 
39  Uses the flexible PIMPLE (PISO-SIMPLE) solution for time-resolved and
40  pseudo-transient and steady simulations.
41 
42  Optional fvModels and fvConstraints are provided to enhance the simulation
43  in many ways including adding various sources, Lagrangian
44  particles, surface film etc. and constraining or limiting the solution.
45 
46 SourceFiles
47  twoPhaseSolver.C
48 
49 See also
50  Foam::solvers::fluidSolver
51 
52 \*---------------------------------------------------------------------------*/
53 
54 #ifndef twoPhaseSolver_H
55 #define twoPhaseSolver_H
56 
57 #include "VoFSolver.H"
58 #include "twoPhaseVoFMixture.H"
59 #include "MULES.H"
60 
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62 
63 namespace Foam
64 {
65 namespace solvers
66 {
67 
68 /*---------------------------------------------------------------------------*\
69  Class twoPhaseSolver Declaration
70 \*---------------------------------------------------------------------------*/
71 
72 class twoPhaseSolver
73 :
74  public VoFSolver
75 {
76 
77 protected:
78 
79  // Phase properties
80 
81  //- Reference to the twoPhaseVoFMixture
83 
84  //- Reference to the phase1-fraction
86 
87  //- Reference to the phase2-fraction
89 
90  //- Switch indicating if this is a restart
91  bool alphaRestart;
92 
93  //- Function to calculate the number of explicit MULES sub-cycles
94  // from the alpha Courant number
96 
97  //- Number of alpha correctors
98  // Usually used to improve the accuracy at high Courant numbers
99  // with semi-implicit MULES, MULESCorr = true
101 
102  //- Semi-implicit MULES switch
103  bool MULESCorr;
104 
105  //- Apply the compression correction from the previous iteration
106  // Improves efficiency for steady-simulations but can only be applied
107  // once the alpha field is reasonably steady, i.e. fully developed
108  bool alphaApplyPrevCorr;
109 
110  //- MULES controls
112 
113 
114  // Kinematic properties
115 
116  // Phase-1 volumetric flux
118 
119  // Phase-2 volumetric flux
121 
122 
123  // Cached temporary fields
124 
125  //- MULES Correction
127 
128 
129 private:
130 
131  // Private Member Functions
132 
133  //- Solve for the phase-fractions
134  void alphaSolve(const label nAlphaSubCycles);
135 
136 
137 protected:
138 
139  // Protected Member Functions
140 
141  //- Read controls
142  virtual bool read();
143 
145  (
146  const surfaceScalarField& phi,
147  const volScalarField& alpha
148  );
149 
150  //- Solve for the phase-fractions
151  void alphaPredictor();
152 
153  //- Calculate the alpha equation sources
154  virtual void alphaSuSp
155  (
158  ) = 0;
159 
160  //- Correct the interface properties following mesh-change
161  // and phase-fraction update
162  virtual void correctInterface() = 0;
163 
164  //- Return the interface surface tension force for the momentum equation
165  virtual tmp<surfaceScalarField> surfaceTensionForce() const = 0;
166 
167  //- Construct and solve the incompressible pressure equation
168  // in the PISO loop
170 
171 
172 public:
173 
174  //- Runtime type information
175  TypeName("twoPhaseSolver");
176 
177 
178  // Constructors
179 
180  //- Construct from region mesh
182 
183  //- Disallow default bitwise copy construction
184  twoPhaseSolver(const twoPhaseSolver&) = delete;
185 
186 
187  //- Destructor
188  virtual ~twoPhaseSolver();
189 
190 
191  // Member Functions
192 
193  //- Called at the start of the time-step, before the PIMPLE loop
194  virtual void preSolve();
195 
196  //- Called at the start of the PIMPLE loop
197  virtual void prePredictor();
198 
199 
200  // Member Operators
201 
202  //- Disallow default bitwise assignment
203  void operator=(const twoPhaseSolver&) = delete;
204 };
205 
206 
207 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
208 
209 } // End namespace solvers
210 } // End namespace Foam
211 
212 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
213 
214 #endif
215 
216 // ************************************************************************* //
MULES: Multidimensional universal limiter for explicit solution.
Generic GeometricField class.
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:96
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:216
Solver module base-class for 2 immiscible fluids, with optional mesh motion and mesh topology changes...
virtual void prePredictor()
Called at the start of the PIMPLE loop.
twoPhaseSolver(fvMesh &mesh, autoPtr< twoPhaseVoFMixture >)
Construct from region mesh.
bool alphaRestart
Switch indicating if this is a restart.
virtual void alphaSuSp(tmp< volScalarField::Internal > &Su, tmp< volScalarField::Internal > &Sp)=0
Calculate the alpha equation sources.
surfaceScalarField alphaPhi1
bool alphaApplyPrevCorr
Apply the compression correction from the previous iteration.
MULES::control MULEScontrols
MULES controls.
void alphaPredictor()
Solve for the phase-fractions.
label nAlphaCorr
Number of alpha correctors.
TypeName("twoPhaseSolver")
Runtime type information.
virtual tmp< surfaceScalarField > surfaceTensionForce() const =0
Return the interface surface tension force for the momentum equation.
volScalarField & alpha1
Reference to the phase1-fraction.
void operator=(const twoPhaseSolver &)=delete
Disallow default bitwise assignment.
autoPtr< Function1< scalar > > nAlphaSubCyclesPtr
Function to calculate the number of explicit MULES sub-cycles.
bool MULESCorr
Semi-implicit MULES switch.
virtual ~twoPhaseSolver()
Destructor.
surfaceScalarField alphaPhi2
tmp< surfaceScalarField > talphaPhi1Corr0
MULES Correction.
virtual void preSolve()
Called at the start of the time-step, before the PIMPLE loop.
volScalarField & alpha2
Reference to the phase2-fraction.
virtual tmp< surfaceScalarField > alphaPhi(const surfaceScalarField &phi, const volScalarField &alpha)
twoPhaseVoFMixture & mixture
Reference to the twoPhaseVoFMixture.
void incompressiblePressureCorrector(volScalarField &p)
Construct and solve the incompressible pressure equation.
virtual void correctInterface()=0
Correct the interface properties following mesh-change.
virtual bool read()
Read controls.
A class for managing temporary objects.
Definition: tmp.H:55
Class to represent a VoF mixture.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
Definition: fvcSup.C:67
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
volScalarField & p