phaseSystem.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) 2015-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::phaseSystem
26 
27 Description
28  Class to represent a system of phases.
29 
30 SourceFiles
31  phaseSystem.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef phaseSystem_H
36 #define phaseSystem_H
37 
38 #include "IOdictionary.H"
39 
40 #include "phaseModel.H"
41 #include "phaseInterface.H"
42 #include "phaseInterfaceKey.H"
43 #include "HashPtrTable.H"
44 #include "PtrListDictionary.H"
45 #include "hashedWordList.H"
46 
47 #include "pimpleNoLoopControl.H"
48 
49 #include "IOMRFZoneList.H"
50 #include "fvModels.H"
51 #include "fvConstraints.H"
52 
53 #include "volFields.H"
54 #include "surfaceFields.H"
55 #include "fvMatricesFwd.H"
56 #include "MULES.H"
57 
58 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59 
60 namespace Foam
61 {
62 
63 class surfaceTensionCoefficientModel;
64 class momentumTransferSystem;
65 class pressureReference;
66 class nonOrthogonalSolutionControl;
67 
68 /*---------------------------------------------------------------------------*\
69  Class phaseSystem Declaration
70 \*---------------------------------------------------------------------------*/
71 
72 class phaseSystem
73 :
74  public IOdictionary
75 {
76 public:
77 
78  // Public Types
79 
80  //- alpha solution control structure
81  struct alphaControl
82  {
83  //- Function to calculate the number of explicit MULES sub-cycles
84  // from the alpha Courant number
86 
87  //- Number of alpha equation sub-cycles
89 
90  //- Number of alpha correctors
91  // Usually used to improve the accuracy at high Courant numbers
92  // with semi-implicit MULES, MULESCorr = true
94 
95  //- Semi-implicit MULES switch
96  bool MULESCorr;
97 
98  //- Explicit correction relaxation factor (defaults to 0.5)
99  scalar MULESCorrRelax;
100 
101  //- Compressibility source stabilisation tolerance
102  scalar vDotResidualAlpha;
103 
104  //- MULES controls
106 
107  //- Optionally clip the phase-fractions
108  // between their physical limits and to sum to 1.
109  // Defaults to false
110  Switch clip;
111 
112  //- Read the alpha and MULES controls from dict
113  void read(const dictionary& dict);
114 
115  //- Correct nAlphaSubCycles for the current Courant number
116  void correct(const scalar CoNum);
117  };
118 
119 
120  // Public Type Definitions
121 
122  //- List of phase models
124 
125  //- Partial list of phase models
127 
128 
129 protected:
130 
131  // Protected Type Definitions
132 
133  //- Table of interface compression coefficients
134  typedef
135  HashTable
136  <
137  scalar,
140  >
141  cAlphaTable;
142 
143  //- Table of surface tension models
144  typedef
146  <
150  >
152 
153 
154  // Protected Data
155 
156  //- Reference to the mesh
157  const fvMesh& mesh_;
158 
159  //- Reference to pimpleNoLoopControl
161 
162  //- Optional MRF zones
164 
165  //- Name of optional reference phase which is not solved for
166  // but obtained from the sum of the other phases
168 
169  //- Phase models
171 
172  //- Moving phase models
174 
175  //- Stationary phase models
177 
178  //- Thermal phase models
180 
181  //- Multi-component phase models
183 
184  //- Total volumetric flux
186 
187  //- Rate of change of pressure
189 
190  //- Interface compression coefficients
191  const cAlphaTable cAlphas_;
192 
193  //- Stabilisation for normalisation of the interface normal
195 
196  //- Surface tension models
199 
200 
201  // Protected Member Functions
202 
203  //- Return the sum of the phase fractions of the moving phases
205 
206  //- Re-normalise the velocity of the phases
207  // around the specified mixture mean
208  void setMixtureU(const volVectorField& Um);
209 
210  //- Re-normalise the flux of the phases
211  // around the specified mixture mean
212  void setMixturePhi
213  (
214  const PtrList<surfaceScalarField>& alphafs,
215  const surfaceScalarField& phim
216  );
217 
218 
219  // Functions required for interface compression
220 
221  //- Normal to interface between two phases
222  // Used for interface compression
224  (
225  const volScalarField& alpha1,
226  const volScalarField& alpha2
227  ) const;
228 
229  //- Normal to interface between two phases dotted with face areas
230  // Used for interface compression
232  (
233  const volScalarField& alpha1,
234  const volScalarField& alpha2
235  ) const;
236 
237  //- Curvature of interface between two phases
238  // Used for interface compression
240  (
241  const phaseModel& alpha1,
242  const phaseModel& alpha2
243  ) const;
244 
245 
246 public:
247 
248  //- Runtime type information
249  TypeName("phaseSystem");
250 
251 
252  //- Default name of the phase properties dictionary
253  static const word propertiesName;
254 
255 
256  // Constructors
257 
258  //- Construct from fvMesh
259  phaseSystem(const fvMesh& mesh);
260 
261 
262  //- Destructor
263  virtual ~phaseSystem();
264 
265 
266  // Member Functions
267 
268  // Access
269 
270  //- Return the mesh
271  inline const fvMesh& mesh() const;
272 
273  //- Return pimpleNoLoopControl
274  inline const pimpleNoLoopControl& pimple() const;
275 
276  //- Return the phase models
277  inline const phaseModelList& phases() const;
278 
279  //- Access the phase models
280  inline phaseModelList& phases();
281 
282  //- Return the models for phases that are moving
283  inline const phaseModelPartialList& movingPhases() const;
284 
285  //- Access the models for phases that are moving
287 
288  //- Return the models for phases that are stationary
289  inline const phaseModelPartialList& stationaryPhases() const;
290 
291  //- Access the models for phases that are stationary
293 
294  //- Return the models for phases that have variable temperature
295  inline const phaseModelPartialList& thermalPhases() const;
296 
297  //- Access the models for phases that have variable temperature
299 
300  //- Return the models for phases that have multiple species
301  inline const phaseModelPartialList& multicomponentPhases() const;
302 
303  //- Access the models for phases that have multiple species
305 
306  //- Return the phase not given as an argument in a two-phase system
307  // An error is generated if the system is not two-phase
308  inline const phaseModel& otherPhase(const phaseModel& phase) const;
309 
310  //- Return the mixture flux
311  inline const surfaceScalarField& phi() const;
312 
313  //- Access the mixture flux
314  inline surfaceScalarField& phi();
315 
316  //- Return the rate of change of the pressure
317  inline const volScalarField& dpdt() const;
318 
319  //- Access the rate of change of the pressure
320  inline volScalarField& dpdt();
321 
322  //- Return MRF zones
323  inline const IOMRFZoneList& MRF() const;
324 
325  //- Access the fvModels
327 
328  //- Access the fvModels
329  inline const Foam::fvModels& fvModels() const;
330 
331  //- Access the fvConstraints
333 
334  //- Access the fvConstraints
335  inline const Foam::fvConstraints& fvConstraints() const;
336 
337 
338  // Sub-model lookup
339 
340  //- Check availability of a sub model for a given interface
341  template<class ModelType>
343  (
344  const phaseInterface& interface
345  ) const;
346 
347  //- Return a sub model for an interface
348  template<class ModelType>
349  const ModelType& lookupInterfacialModel
350  (
351  const phaseInterface& interface
352  ) const;
353 
354 
355  // Properties
356 
357  //- Return the mixture density
358  tmp<volScalarField> rho() const;
359 
360  //- Return the mixture velocity
361  tmp<volVectorField> U() const;
362 
363  //- Return the surface tension coefficient for an interface
365  (
366  const phaseInterfaceKey& key
367  ) const;
368 
369  //- Return the surface tension coefficient for an interface on a
370  // patch
372  (
373  const phaseInterfaceKey& key,
374  const label patchi
375  ) const;
376 
377  //- Indicator of the proximity of the interface
378  // Field values are 1 near and 0 away for the interface.
380 
381  //- Stabilisation for normalisation of the interface normal
382  inline const dimensionedScalar& deltaN() const;
383 
384  //- Return the surface tension force
386 
387  //- Return incompressibility
388  bool incompressible() const;
389 
390 
391  // Evolution
392 
393  //- Returns true if the phase pressure is treated implicitly
394  // in the phase fraction equations
395  bool implicitPhasePressure() const;
396 
397  //- Solve for the phase fractions
398  void solve
399  (
400  const alphaControl& alphaControls,
401  const PtrList<volScalarField>& rAs,
402  const momentumTransferSystem& mts
403  );
404 
405  //- Correct the fluid properties other than those listed below
406  void correct();
407 
408  //- Correct the continuity errors
410  (
412  );
413 
414  //- Correct the kinematics
415  void correctKinematics();
416 
417  //- Correct the thermodynamics
418  void correctThermo();
419 
420  //- Correct the reactions
421  void correctReactions();
422 
423  //- Correct the species mass fractions
424  void correctSpecies();
425 
426  //- Predict the momentumTransport
428 
429  //- Predict the energy transport e.g. alphat
431 
432  //- Correct the momentumTransport
434 
435  //- Correct the energy transport e.g. alphat
437 
438  //- Update the fluid properties for mesh changes
439  void meshUpdate();
440 
441  //- Correct fixed-flux BCs to be consistent with the velocity BCs
442  void correctBoundaryFlux();
443 
444  //- ...
445  void correctPhi
446  (
447  const volScalarField& p_rgh,
448  const autoPtr<volScalarField>& divU,
451  );
452 
453 
454  // IO
455 
456  //- Read base phaseProperties dictionary
457  virtual bool read();
458 };
459 
460 
461 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
462 
465 
466 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
467 
468 } // End namespace Foam
469 
470 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
471 
472 #include "phaseSystemI.H"
473 
474 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
475 
476 #ifdef NoRepository
477  #include "phaseSystemTemplates.C"
478 #endif
479 
480 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
481 
482 #endif
483 
484 // ************************************************************************* //
MULES: Multidimensional universal limiter for explicit solution.
Generic GeometricField class.
A HashTable specialisation for hashing pointers.
Definition: HashPtrTable.H:67
An STL-conforming hash table.
Definition: HashTable.H:127
List of MRF zones with IO functionality. MRF zones are specified by a list of dictionary entries,...
Definition: IOMRFZoneList.H:68
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Finite volume constraints.
Definition: fvConstraints.H:67
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
Finite volume models.
Definition: fvModels.H:65
Class which provides interfacial momentum transfer between a number of phases. Drag,...
Non-orthogonal solution control class. Provides non-orthogonal-loop control methods.
Word-pair based class used for keying interface models in hash tables.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
Class to represent a system of phases.
Definition: phaseSystem.H:74
phaseModelPartialList stationaryPhaseModels_
Stationary phase models.
Definition: phaseSystem.H:175
void correctKinematics()
Correct the kinematics.
Definition: phaseSystem.C:666
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
Definition: phaseSystem.C:543
void correctPhi(const volScalarField &p_rgh, const autoPtr< volScalarField > &divU, const pressureReference &pressureReference, nonOrthogonalSolutionControl &pimple)
...
Definition: phaseSystem.C:794
bool foundInterfacialModel(const phaseInterface &interface) const
Check availability of a sub model for a given interface.
HashTable< scalar, phaseInterfaceKey, phaseInterfaceKey::hash > cAlphaTable
Table of interface compression coefficients.
Definition: phaseSystem.H:140
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:252
const cAlphaTable cAlphas_
Interface compression coefficients.
Definition: phaseSystem.H:190
const IOMRFZoneList & MRF() const
Return MRF zones.
Definition: phaseSystemI.H:220
const phaseModel & otherPhase(const phaseModel &phase) const
Return the phase not given as an argument in a two-phase system.
Definition: phaseSystemI.H:174
void correctReactions()
Correct the reactions.
Definition: phaseSystem.C:693
const fvMesh & mesh_
Reference to the mesh.
Definition: phaseSystem.H:156
phaseModelPartialList multicomponentPhaseModels_
Multi-component phase models.
Definition: phaseSystem.H:181
phaseModelList phaseModels_
Phase models.
Definition: phaseSystem.H:169
phaseModelPartialList movingPhaseModels_
Moving phase models.
Definition: phaseSystem.H:172
void setMixtureU(const volVectorField &Um)
Re-normalise the velocity of the phases.
Definition: phaseSystem.C:83
const dimensionedScalar & deltaN() const
Stabilisation for normalisation of the interface normal.
Definition: phaseSystemI.H:250
void correctSpecies()
Correct the species mass fractions.
Definition: phaseSystem.C:702
const Foam::fvModels & fvModels() const
Access the fvModels.
Definition: phaseSystemI.H:232
tmp< surfaceScalarField > surfaceTension(const phaseModel &) const
Return the surface tension force.
Definition: phaseSystem.C:569
const surfaceScalarField & phi() const
Return the mixture flux.
Definition: phaseSystemI.H:196
void correctThermo()
Correct the thermodynamics.
Definition: phaseSystem.C:684
const pimpleNoLoopControl & pimple_
Reference to pimpleNoLoopControl.
Definition: phaseSystem.H:159
const volScalarField & dpdt() const
Return the rate of change of the pressure.
Definition: phaseSystemI.H:208
void correct()
Correct the fluid properties other than those listed below.
Definition: phaseSystem.C:621
tmp< volScalarField > sigma(const phaseInterfaceKey &key) const
Return the surface tension coefficient for an interface.
Definition: phaseSystem.C:507
const phaseModelPartialList & multicomponentPhases() const
Return the models for phases that have multiple species.
Definition: phaseSystemI.H:160
void correctContinuityError(const PtrList< volScalarField::Internal > &dmdts)
Correct the continuity errors.
Definition: phaseSystem.C:631
const phaseModelPartialList & movingPhases() const
Return the models for phases that are moving.
Definition: phaseSystemI.H:118
void predictThermophysicalTransport()
Predict the energy transport e.g. alphat.
Definition: phaseSystem.C:720
tmp< volScalarField > K(const phaseModel &alpha1, const phaseModel &alpha2) const
Curvature of interface between two phases.
Definition: phaseSystem.C:165
tmp< volVectorField > U() const
Return the mixture velocity.
Definition: phaseSystem.C:479
HashPtrTable< surfaceTensionCoefficientModel, phaseInterfaceKey, phaseInterfaceKey::hash > surfaceTensionCoefficientModelTable
Table of surface tension models.
Definition: phaseSystem.H:150
void correctMomentumTransport()
Correct the momentumTransport.
Definition: phaseSystem.C:729
void correctBoundaryFlux()
Correct fixed-flux BCs to be consistent with the velocity BCs.
Definition: phaseSystem.C:761
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:104
PtrListDictionary< phaseModel > phaseModelList
List of phase models.
Definition: phaseSystem.H:122
const surfaceTensionCoefficientModelTable surfaceTensionCoefficientModels_
Surface tension models.
Definition: phaseSystem.H:197
const dimensionedScalar deltaN_
Stabilisation for normalisation of the interface normal.
Definition: phaseSystem.H:193
const pimpleNoLoopControl & pimple() const
Return pimpleNoLoopControl.
Definition: phaseSystemI.H:97
void correctThermophysicalTransport()
Correct the energy transport e.g. alphat.
Definition: phaseSystem.C:738
phaseModelPartialList thermalPhaseModels_
Thermal phase models.
Definition: phaseSystem.H:178
void meshUpdate()
Update the fluid properties for mesh changes.
Definition: phaseSystem.C:747
const phaseModelPartialList & stationaryPhases() const
Return the models for phases that are stationary.
Definition: phaseSystemI.H:132
surfaceScalarField phi_
Total volumetric flux.
Definition: phaseSystem.H:184
bool implicitPhasePressure() const
Returns true if the phase pressure is treated implicitly.
tmp< volScalarField > sumAlphaMoving() const
Return the sum of the phase fractions of the moving phases.
Definition: phaseSystem.C:57
phaseSystem(const fvMesh &mesh)
Construct from fvMesh.
Definition: phaseSystem.C:189
void setMixturePhi(const PtrList< surfaceScalarField > &alphafs, const surfaceScalarField &phim)
Re-normalise the flux of the phases.
Definition: phaseSystem.C:104
const phaseModelPartialList & thermalPhases() const
Return the models for phases that have variable temperature.
Definition: phaseSystemI.H:146
virtual ~phaseSystem()
Destructor.
Definition: phaseSystem.C:409
const fvMesh & mesh() const
Return the mesh.
Definition: phaseSystemI.H:91
TypeName("phaseSystem")
Runtime type information.
UPtrList< phaseModel > phaseModelPartialList
Partial list of phase models.
Definition: phaseSystem.H:125
tmp< surfaceScalarField > nHatf(const volScalarField &alpha1, const volScalarField &alpha2) const
Normal to interface between two phases dotted with face areas.
Definition: phaseSystem.C:154
tmp< volScalarField > rho() const
Return the mixture density.
Definition: phaseSystem.C:452
IOMRFZoneList MRF_
Optional MRF zones.
Definition: phaseSystem.H:162
const Foam::fvConstraints & fvConstraints() const
Access the fvConstraints.
Definition: phaseSystemI.H:244
word referencePhaseName_
Name of optional reference phase which is not solved for.
Definition: phaseSystem.H:166
tmp< surfaceVectorField > nHatfv(const volScalarField &alpha1, const volScalarField &alpha2) const
Normal to interface between two phases.
Definition: phaseSystem.C:128
volScalarField dpdt_
Rate of change of pressure.
Definition: phaseSystem.H:187
void predictMomentumTransport()
Predict the momentumTransport.
Definition: phaseSystem.C:711
virtual bool read()
Read base phaseProperties dictionary.
Definition: phaseSystem.C:905
const ModelType & lookupInterfacialModel(const phaseInterface &interface) const
Return a sub model for an interface.
void solve(const alphaControl &alphaControls, const PtrList< volScalarField > &rAs, const momentumTransferSystem &mts)
Solve for the phase fractions.
bool incompressible() const
Return incompressibility.
Definition: phaseSystem.C:607
Pimple no-loop control class. Implements various option flags, but leaves loop controls to the deriva...
Provides controls for the pressure reference in closed-volume simulations.
Abstract base-class for interface surface-tension models which can be used when interface compression...
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
scalar CoNum
Forward declarations of fvMatrix specialisations.
label patchi
volScalarField sf(fieldObject, mesh)
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
tmp< volScalarField > byDt(const volScalarField &vf)
Definition: phaseSystem.C:927
dictionary dict
alpha solution control structure
Definition: phaseSystem.H:81
void correct(const scalar CoNum)
Correct nAlphaSubCycles for the current Courant number.
Definition: phaseSystem.C:446
Switch clip
Optionally clip the phase-fractions.
Definition: phaseSystem.H:109
void read(const dictionary &dict)
Read the alpha and MULES controls from dict.
Definition: phaseSystem.C:415
label nAlphaCorr
Number of alpha correctors.
Definition: phaseSystem.H:92
scalar vDotResidualAlpha
Compressibility source stabilisation tolerance.
Definition: phaseSystem.H:101
MULES::control MULES
MULES controls.
Definition: phaseSystem.H:104
autoPtr< Function1< scalar > > nAlphaSubCyclesPtr
Function to calculate the number of explicit MULES sub-cycles.
Definition: phaseSystem.H:84
bool MULESCorr
Semi-implicit MULES switch.
Definition: phaseSystem.H:95
scalar MULESCorrRelax
Explicit correction relaxation factor (defaults to 0.5)
Definition: phaseSystem.H:98
label nAlphaSubCycles
Number of alpha equation sub-cycles.
Definition: phaseSystem.H:87
Foam::surfaceFields.