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-2022 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 and model interfacial transfers
29  between them.
30 
31 SourceFiles
32  phaseSystem.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef phaseSystem_H
37 #define phaseSystem_H
38 
39 #include "IOdictionary.H"
40 
41 #include "phaseModel.H"
42 #include "phaseInterface.H"
43 #include "phaseInterfaceKey.H"
44 #include "HashPtrTable.H"
45 #include "PtrListDictionary.H"
46 #include "hashedWordList.H"
47 
48 #include "IOMRFZoneList.H"
49 #include "fvModels.H"
50 #include "fvConstraints.H"
51 
52 #include "volFields.H"
53 #include "surfaceFields.H"
54 #include "fvMatricesFwd.H"
55 
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58 namespace Foam
59 {
60 
61 class surfaceTensionModel;
62 class pressureReference;
63 class nonOrthogonalSolutionControl;
64 
65 /*---------------------------------------------------------------------------*\
66  Class phaseSystem Declaration
67 \*---------------------------------------------------------------------------*/
68 
69 class phaseSystem
70 :
71  public IOdictionary
72 {
73 public:
74 
75  // Public Typedefs
76 
78 
80 
82 
84 
86 
87  typedef
89  <
91  phaseInterfaceKey,
93  >
94  dmdtfTable;
95 
96  typedef
98  <
99  HashPtrTable<volScalarField>,
100  phaseInterfaceKey,
102  >
103  dmidtfTable;
104 
105 
106 
107 protected:
108 
109  // Protected typedefs
110 
111  typedef
112  HashTable
113  <
114  autoPtr<surfaceTensionModel>,
115  phaseInterfaceKey,
117  >
119 
121  cAlphaTable;
122 
123 
124  // Protected data
125 
126  //- Reference to the mesh
127  const fvMesh& mesh_;
128 
129  //- Name of optional reference phase which is not solved for
130  // but obtained from the sum of the other phases
132 
133  //- Phase models
134  phaseModelList phaseModels_;
135 
136  //- Moving phase models
137  phaseModelPartialList movingPhaseModels_;
138 
139  //- Stationary phase models
140  phaseModelPartialList stationaryPhaseModels_;
141 
142  //- Anisothermal phase models
143  phaseModelPartialList anisothermalPhaseModels_;
144 
145  //- Multi-component phase models
146  phaseModelPartialList multiComponentPhaseModels_;
147 
148  //- Total volumetric flux
150 
151  //- Rate of change of pressure
152  volScalarField dpdt_;
153 
154  //- Optional MRF zones
156 
157  //- Interface compression coefficients
158  cAlphaTable cAlphas_;
159 
160  //- Stabilisation for normalisation of the interface normal
162 
163 
164  // Sub Models
165 
166  //- Surface tension models
168 
169 
170  //- Flag to indicate that returned lists of fields are "complete"; i.e.,
171  // that an absence of force is returned as a constructed list of zeros,
172  // rather than a null pointer
173  static const bool fillFields_ = false;
174 
175 
176  // Protected member functions
177 
178  //- Calculate and return the mixture flux
180  (
181  const phaseModelList& phaseModels
182  ) const;
183 
184  //- Check that mass transfer is supported across the given interface
185  template<class ModelType>
186  void validateMassTransfer(const phaseInterface& interface) const;
187 
188  //- Return the sum of the phase fractions of the moving phases
190 
191  //- Re-normalise the velocity of the phases
192  // around the specified mixture mean
193  void setMixtureU(const volVectorField& Um);
194 
195 
196  // Functions required for interface compression
197 
198  //- Normal to interface between two phases
199  // Used for interface compression
201  (
202  const volScalarField& alpha1,
203  const volScalarField& alpha2
204  ) const;
205 
206  //- Normal to interface between two phases dotted with face areas
207  // Used for interface compression
209  (
210  const volScalarField& alpha1,
211  const volScalarField& alpha2
212  ) const;
213 
214  //- Curvature of interface between two phases
215  // Used for interface compression
217  (
218  const phaseModel& alpha1,
219  const phaseModel& alpha2
220  ) const;
221 
222 
223  // Sub-model construction
224 
225  //- Return the dictionary containing interfacial model or value
226  // settings for the given name. Performs backwards compatibility
227  // conversions.
228  template<class Type>
229  dictionary interfacialDict(const word& name) const;
230 
231 
232 public:
233 
234  //- Runtime type information
235  TypeName("phaseSystem");
236 
237  //- Default name of the phase properties dictionary
238  static const word propertiesName;
239 
240 
241  // Declare runtime construction
242 
244  (
245  autoPtr,
246  phaseSystem,
247  dictionary,
248  (
249  const fvMesh& mesh
250  ),
251  (mesh)
252  );
253 
254 
255  // Constructors
256 
257  //- Construct from fvMesh
258  phaseSystem(const fvMesh& mesh);
259 
260 
261  // Selectors
262 
264  (
265  const fvMesh& mesh
266  );
267 
268 
269  //- Destructor
270  virtual ~phaseSystem();
271 
272 
273  // Member Functions
274 
275  // Access
276 
277  //- Return the mesh
278  inline const fvMesh& mesh() const;
279 
280  //- Return the phase models
281  inline const phaseModelList& phases() const;
282 
283  //- Access the phase models
284  inline phaseModelList& phases();
285 
286  //- Return the models for phases that are moving
287  inline const phaseModelPartialList& movingPhases() const;
288 
289  //- Access the models for phases that are moving
290  inline phaseModelPartialList& movingPhases();
291 
292  //- Return the models for phases that are stationary
293  inline const phaseModelPartialList& stationaryPhases() const;
294 
295  //- Access the models for phases that are stationary
296  inline phaseModelPartialList& stationaryPhases();
297 
298  //- Return the models for phases that have variable temperature
299  inline const phaseModelPartialList& anisothermalPhases() const;
300 
301  //- Access the models for phases that have variable temperature
302  inline phaseModelPartialList& anisothermalPhases();
303 
304  //- Return the models for phases that have multiple species
305  inline const phaseModelPartialList& multiComponentPhases() const;
306 
307  //- Access the models for phases that have multiple species
308  inline phaseModelPartialList& multiComponentPhases();
309 
310  //- Return the phase not given as an argument in a two-phase system
311  // An error is generated if the system is not two-phase
312  inline const phaseModel& otherPhase(const phaseModel& phase) const;
313 
314  //- Return the mixture flux
315  inline const surfaceScalarField& phi() const;
316 
317  //- Access the mixture flux
318  inline surfaceScalarField& phi();
319 
320  //- Return the rate of change of the pressure
321  inline const volScalarField& dpdt() const;
322 
323  //- Access the rate of change of the pressure
324  inline volScalarField& dpdt();
325 
326  //- Return MRF zones
327  inline const IOMRFZoneList& MRF() const;
328 
329  //- Access the fvModels
330  inline Foam::fvModels& fvModels(fvMesh& mesh);
331 
332  //- Access the fvModels
333  inline const Foam::fvModels& fvModels() const;
334 
335  //- Access the fvConstraints
337 
338  //- Access the fvConstraints
339  inline const Foam::fvConstraints& fvConstraints() const;
340 
341 
342  // Sub-model construction
343 
344  //- Return the model name. This is the same as the model's typename
345  // but without "Model" on the end.
346  template<class ModelType>
347  word modelName() const;
348 
349  //- Generate interfacial-model lists
350  template<class ModelType, class ... InterfaceTypes>
352  (
353  const dictionary& dict,
354  const phaseInterface& interface,
355  PtrList<phaseInterface>& interfaces,
356  PtrList<ModelType>& models
357  ) const;
358 
359  //- Generate interfacial-model tables
360  template<class ModelType>
362  (
363  const dictionary& dict,
364  HashTable
365  <
367  phaseInterfaceKey,
369  >& models
370  ) const;
371 
372  //- Generate interfacial-model tables
373  template<class ModelType>
375  (
376  HashTable
377  <
379  phaseInterfaceKey,
381  >& models
382  ) const;
383 
384  //- Generate interfacial-model tables
385  template<class ValueType>
387  (
388  const dictionary& dict,
389  HashTable
390  <
391  ValueType,
392  phaseInterfaceKey,
394  >& values
395  ) const;
396 
397  //- Generate interfacial-model tables
398  template<class ValueType>
400  (
401  const word& valueName,
402  HashTable
403  <
404  ValueType,
405  phaseInterfaceKey,
407  >& values
408  ) const;
409 
410  //- Return the dictionary from which to construct a low-level
411  // sub-model. Checks that there is just one sub-dictionary then
412  // returns it.
413  template<class ModelType>
414  static const dictionary& modelSubDict(const dictionary& dict);
415 
416 
417  // Sub-model lookup
418 
419  //- Check availability of a sub model for a given interface
420  template<class ModelType>
421  bool foundInterfacialModel(const phaseInterface& interface) const;
422 
423  //- Return a sub model for an interface
424  template<class ModelType>
425  const ModelType& lookupInterfacialModel
426  (
427  const phaseInterface& interface
428  ) const;
429 
430 
431  // Field construction
432 
433  //- Fill up gaps in a phase-indexed list of fields with zeros
434  template
435  <
436  class Type,
437  template<class> class PatchField,
438  class GeoMesh
439  >
440  void fillFields
441  (
442  const word& name,
443  const dimensionSet& dims,
445  ) const;
446 
447  //- Fill up gaps in a phase-indexed table of fields with zeros
448  template
449  <
450  class Type,
451  template<class> class PatchField,
452  class GeoMesh
453  >
454  void fillFields
455  (
456  const word& name,
457  const dimensionSet& dims,
459  fieldTable
460  ) const;
461 
462 
463  // Properties
464 
465  //- Return the mixture density
466  tmp<volScalarField> rho() const;
467 
468  //- Return the mixture velocity
469  tmp<volVectorField> U() const;
470 
471  //- Return the surface tension coefficient for an interface
472  tmp<volScalarField> sigma(const phaseInterfaceKey& key) const;
473 
474  //- Return the surface tension coefficient for an interface on a
475  // patch
477  (
478  const phaseInterfaceKey& key,
479  const label patchi
480  ) const;
481 
482  //- Indicator of the proximity of the interface
483  // Field values are 1 near and 0 away for the interface.
485 
486  //- Stabilisation for normalisation of the interface normal
487  inline const dimensionedScalar& deltaN() const;
488 
489  //- Return the mass transfer rate for an interface
490  virtual tmp<volScalarField> dmdtf
491  (
492  const phaseInterfaceKey& key
493  ) const;
494 
495  //- Return the mass transfer rates for each phase
496  virtual PtrList<volScalarField> dmdts() const;
497 
498  //- Return the mass transfer pressure implicit coefficients
499  // for each phase
500  virtual PtrList<volScalarField> d2mdtdps() const;
501 
502  //- Return incompressibility
503  bool incompressible() const;
504 
505 
506  // Transfers
507 
508  //- Return the momentum transfer matrices for the cell-based
509  // algorithm
511 
512  //- Return the momentum transfer matrices for the face-based
513  // algorithm
515 
516  //- Return the implicit force coefficients for the face-based
517  // algorithm
518  virtual PtrList<surfaceScalarField> AFfs() const = 0;
519 
520  //- Return the force fluxes for the cell-based algorithm
522  (
524  ) = 0;
525 
526  //- Return the force fluxes for the face-based algorithm
528  (
530  ) = 0;
531 
532  //- Return the force fluxes for the cell-based algorithm
534  (
536  ) const = 0;
537 
538  //- Return the force fluxes for the face-based algorithm
540  (
542  ) const = 0;
543 
544  //- Return the explicit part of the drag force
546  (
548  ) const = 0;
549 
550  //- Returns true if the phase pressure is treated implicitly
551  // in the phase fraction equation
552  virtual bool implicitPhasePressure(const phaseModel& phase) const;
553 
554  //- Returns true if the phase pressure is treated implicitly
555  // in the phase fraction equation for any phase
556  virtual bool implicitPhasePressure() const;
557 
558  //- Return the phase diffusivity
559  // divided by the momentum central coefficient
561  (
564  ) const = 0;
565 
566  //- Solve the drag system for the new velocities and fluxes
567  virtual void partialElimination
568  (
573  ) = 0;
574 
575  //- Solve the drag system for the new fluxes
576  virtual void partialEliminationf
577  (
581  ) = 0;
582 
583  //- Re-normalise the flux of the phases
584  // around the specified mixture mean
585  void setMixturePhi
586  (
588  const surfaceScalarField& phim
589  );
590 
591  //- Return the flux corrections for the cell-based algorithm
593  (
595  const bool includeVirtualMass = false
596  ) const = 0;
597 
598  //- Return the heat transfer matrices
599  virtual autoPtr<heatTransferTable> heatTransfer() const = 0;
600 
601  //- Return the specie transfer matrices
602  virtual autoPtr<specieTransferTable> specieTransfer() const = 0;
603 
604  //- Return the surface tension force
606  (
607  const phaseModel& phase
608  ) const;
609 
610 
611  // Evolution
612 
613  //- Solve for the phase fractions
614  virtual void solve
615  (
618  );
619 
620  //- Correct the fluid properties other than those listed below
621  virtual void correct();
622 
623  //- Correct the continuity errors
624  virtual void correctContinuityError();
625 
626  //- Correct the kinematics
627  virtual void correctKinematics();
628 
629  //- Correct the thermodynamics
630  virtual void correctThermo();
631 
632  //- Correct the reactions
633  virtual void correctReactions();
634 
635  //- Correct the species mass fractions
636  virtual void correctSpecies();
637 
638  //- Correct the turbulence
639  virtual void correctTurbulence();
640 
641  //- Correct the energy transport e.g. alphat
642  virtual void correctEnergyTransport();
643 
644  //- Update the fluid properties for mesh changes
645  virtual void meshUpdate();
646 
647  //- Correct fixed-flux BCs to be consistent with the velocity BCs
648  void correctBoundaryFlux();
649 
650  void correctPhi
651  (
652  const volScalarField& p_rgh,
653  const tmp<volScalarField>& divU,
656  );
657 
658 
659  // IO
660 
661  //- Read base phaseProperties dictionary
662  virtual bool read();
663 };
664 
665 
666 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
667 
668 tmp<volScalarField> byDt(const volScalarField& vf);
670 
671 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
672 
673 } // End namespace Foam
674 
675 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
676 
677 #include "phaseSystemI.H"
678 
679 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
680 
681 #ifdef NoRepository
682  #include "phaseSystemTemplates.C"
683 #endif
684 
685 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
686 
687 #endif
688 
689 // ************************************************************************* //
Foam::surfaceFields.
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:237
pressureReference & pressureReference
volScalarField & p_rgh
virtual PtrList< surfaceScalarField > AFfs() const =0
Return the implicit force coefficients for the face-based.
dictionary dict
virtual void partialEliminationf(const PtrList< surfaceScalarField > &rAUfs, const PtrList< surfaceScalarField > &alphafs, const PtrList< surfaceScalarField > &phiKdPhifs)=0
Solve the drag system for the new fluxes.
const phaseModelPartialList & movingPhases() const
Return the models for phases that are moving.
Definition: phaseSystemI.H:49
surfaceTensionModelTable surfaceTensionModels_
Surface tension models.
Definition: phaseSystem.H:166
virtual void correctTurbulence()
Correct the turbulence.
pimpleNoLoopControl & pimple
const Foam::fvModels & fvModels() const
Access the fvModels.
Definition: phaseSystemI.H:163
virtual PtrList< surfaceScalarField > phiFfs(const PtrList< surfaceScalarField > &rAUfs)=0
Return the force fluxes for the face-based algorithm.
virtual autoPtr< specieTransferTable > specieTransfer() const =0
Return the specie transfer matrices.
tmp< surfaceScalarField > nHatf(const volScalarField &alpha1, const volScalarField &alpha2) const
Normal to interface between two phases dotted with face areas.
virtual PtrList< volScalarField > d2mdtdps() const
Return the mass transfer pressure implicit coefficients.
const Foam::fvConstraints & fvConstraints() const
Access the fvConstraints.
Definition: phaseSystemI.H:175
declareRunTimeSelectionTable(autoPtr, phaseSystem, dictionary,(const fvMesh &mesh),(mesh))
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
const dimensionedScalar deltaN_
Stabilisation for normalisation of the interface normal.
Definition: phaseSystem.H:160
PtrList< surfaceScalarField > alphafs(phases.size())
HashPtrTable< fvVectorMatrix > momentumTransferTable
Definition: phaseSystem.H:76
void setMixturePhi(const PtrList< surfaceScalarField > &alphafs, const surfaceScalarField &phim)
Re-normalise the flux of the phases.
virtual autoPtr< heatTransferTable > heatTransfer() const =0
Return the heat transfer matrices.
HashPtrTable< HashPtrTable< volScalarField >, phaseInterfaceKey, phaseInterfaceKey::hash > dmidtfTable
Definition: phaseSystem.H:102
volScalarField dpdt_
Rate of change of pressure.
Definition: phaseSystem.H:151
volScalarField & alpha1(mixture.alpha1())
TypeName("phaseSystem")
Runtime type information.
virtual autoPtr< momentumTransferTable > momentumTransferf()=0
Return the momentum transfer matrices for the face-based.
tmp< volScalarField > K(const phaseModel &alpha1, const phaseModel &alpha2) const
Curvature of interface between two phases.
HashPtrTable< fvScalarMatrix > heatTransferTable
Definition: phaseSystem.H:78
phaseModelPartialList stationaryPhaseModels_
Stationary phase models.
Definition: phaseSystem.H:139
void correctPhi(const volScalarField &p_rgh, const tmp< volScalarField > &divU, const pressureReference &pressureReference, nonOrthogonalSolutionControl &pimple)
HashPtrTable< fvScalarMatrix > specieTransferTable
Definition: phaseSystem.H:80
virtual PtrList< surfaceScalarField > ddtCorrByAs(const PtrList< volScalarField > &rAUs, const bool includeVirtualMass=false) const =0
Return the flux corrections for the cell-based algorithm.
A HashTable specialisation for hashing pointers.
Definition: HashPtrTable.H:50
Provides controls for the pressure reference in closed-volume simulations.
virtual autoPtr< momentumTransferTable > momentumTransfer()=0
Return the momentum transfer matrices for the cell-based.
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
void correctBoundaryFlux()
Correct fixed-flux BCs to be consistent with the velocity BCs.
HashTable< autoPtr< surfaceTensionModel >, phaseInterfaceKey, phaseInterfaceKey::hash > surfaceTensionModelTable
Definition: phaseSystem.H:117
phaseModelList phaseModels_
Phase models.
Definition: phaseSystem.H:133
alpha2
Definition: alphaEqn.H:115
static const dictionary & modelSubDict(const dictionary &dict)
Return the dictionary from which to construct a low-level.
phaseModelPartialList multiComponentPhaseModels_
Multi-component phase models.
Definition: phaseSystem.H:145
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
const phaseModelPartialList & stationaryPhases() const
Return the models for phases that are stationary.
Definition: phaseSystemI.H:63
virtual void solve(const PtrList< volScalarField > &rAUs, const PtrList< surfaceScalarField > &rAUfs)
Solve for the phase fractions.
tmp< surfaceScalarField > surfaceTension(const phaseModel &phase) const
Return the surface tension force.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:35
virtual void correctContinuityError()
Correct the continuity errors.
virtual PtrList< volVectorField > KdUByAs(const PtrList< volScalarField > &rAUs) const =0
Return the explicit part of the drag force.
virtual void correctKinematics()
Correct the kinematics.
Dimension set for the base types.
Definition: dimensionSet.H:121
const fileName & name() const
Return the dictionary name.
Definition: dictionary.H:109
void generateInterfacialModels(const dictionary &dict, const phaseInterface &interface, PtrList< phaseInterface > &interfaces, PtrList< ModelType > &models) const
Generate interfacial-model lists.
phaseSystem(const fvMesh &mesh)
Construct from fvMesh.
Info<< "Reading field p_rgh\"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);if(fluid.found("pMin")){ IOWarningInFunction(fluid)<< "Pressure limits, pMin and pMax, are now read from "<< pimple.dict().name()<< endl;}pressureReference pressureReference(p, p_rgh, pimple.dict(), fluid.incompressible());if(fluid.incompressible()){ p=p_rgh+fluid.rho() *gh;}if(p_rgh.needReference() &&fluid.incompressible()){ p+=dimensionedScalar("p", p.dimensions(), pressureReference.refValue() - getRefCellValue(p, pressureReference.refCell()));}p_rgh=p - fluid.rho() *gh;mesh.schemes().setFluxRequired(p_rgh.name());PtrList< volScalarField > rAUs
Definition: createFields.H:67
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:68
A class for handling words, derived from string.
Definition: word.H:59
virtual tmp< volScalarField > dmdtf(const phaseInterfaceKey &key) const
Return the mass transfer rate for an interface.
tmp< volVectorField > U() const
Return the mixture velocity.
Finite volume models.
Definition: fvModels.H:60
const IOMRFZoneList & MRF() const
Return MRF zones.
Definition: phaseSystemI.H:151
const phaseModelPartialList & anisothermalPhases() const
Return the models for phases that have variable temperature.
Definition: phaseSystemI.H:77
virtual void correctThermo()
Correct the thermodynamics.
PtrListDictionary< phaseModel > phaseModelList
Definition: phaseSystem.H:82
const phaseModel & otherPhase(const phaseModel &phase) const
Return the phase not given as an argument in a two-phase system.
Definition: phaseSystemI.H:105
word referencePhaseName_
Name of optional reference phase which is not solved for.
Definition: phaseSystem.H:130
bool foundInterfacialModel(const phaseInterface &interface) const
Check availability of a sub model for a given interface.
UPtrList< phaseModel > phaseModelPartialList
Definition: phaseSystem.H:84
const fvMesh & mesh_
Reference to the mesh.
Definition: phaseSystem.H:126
static autoPtr< phaseSystem > New(const fvMesh &mesh)
virtual PtrList< surfaceScalarField > phiFs(const PtrList< volScalarField > &rAUs)=0
Return the force fluxes for the cell-based algorithm.
An STL-conforming hash table.
Definition: HashTable.H:61
phaseModelPartialList anisothermalPhaseModels_
Anisothermal phase models.
Definition: phaseSystem.H:142
cAlphaTable cAlphas_
Interface compression coefficients.
Definition: phaseSystem.H:157
tmp< surfaceVectorField > nHatfv(const volScalarField &alpha1, const volScalarField &alpha2) const
Normal to interface between two phases.
phaseModelPartialList movingPhaseModels_
Moving phase models.
Definition: phaseSystem.H:136
void generateInterfacialValues(const dictionary &dict, HashTable< ValueType, phaseInterfaceKey, phaseInterfaceKey::hash > &values) const
Generate interfacial-model tables.
virtual void correct()
Correct the fluid properties other than those listed below.
tmp< volScalarField > rho() const
Return the mixture density.
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
tmp< volScalarField > sumAlphaMoving() const
Return the sum of the phase fractions of the moving phases.
IOMRFZoneList MRF_
Optional MRF zones.
Definition: phaseSystem.H:154
PtrList< surfaceScalarField > rAUfs
Definition: createFields.H:68
HashPtrTable< volScalarField, phaseInterfaceKey, phaseInterfaceKey::hash > dmdtfTable
Definition: phaseSystem.H:93
dictionary interfacialDict(const word &name) const
Return the dictionary containing interfacial model or value.
tmp< volScalarField > divU
Definition: alphaSuSp.H:9
volScalarField sf(fieldObject, mesh)
const volScalarField & dpdt() const
Return the rate of change of the pressure.
Definition: phaseSystemI.H:139
bool incompressible() const
Return incompressibility.
tmp< volScalarField > byDt(const volScalarField &vf)
tmp< volScalarField > sigma(const phaseInterfaceKey &key) const
Return the surface tension coefficient for an interface.
virtual void correctReactions()
Correct the reactions.
virtual void partialElimination(const PtrList< volScalarField > &rAUs, const PtrList< volVectorField > &KdUByAs, const PtrList< surfaceScalarField > &alphafs, const PtrList< surfaceScalarField > &phiKdPhis)=0
Solve the drag system for the new velocities and fluxes.
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
void fillFields(const word &name, const dimensionSet &dims, PtrList< GeometricField< Type, PatchField, GeoMesh >> &fieldList) const
Fill up gaps in a phase-indexed list of fields with zeros.
Forward declarations of fvMatrix specialisations.
label patchi
const ModelType & lookupInterfacialModel(const phaseInterface &interface) const
Return a sub model for an interface.
tmp< surfaceScalarField > calcPhi(const phaseModelList &phaseModels) const
Calculate and return the mixture flux.
const surfaceScalarField & phi() const
Return the mixture flux.
Definition: phaseSystemI.H:127
virtual ~phaseSystem()
Destructor.
Finite volume constraints.
Definition: fvConstraints.H:57
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
virtual bool read()
Read base phaseProperties dictionary.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
virtual PtrList< surfaceScalarField > DByAfs(const PtrList< volScalarField > &rAUs, const PtrList< surfaceScalarField > &rAUfs) const =0
Return the phase diffusivity.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:52
void setMixtureU(const volVectorField &Um)
Re-normalise the velocity of the phases.
void validateMassTransfer(const phaseInterface &interface) const
Check that mass transfer is supported across the given interface.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
surfaceScalarField phi_
Total volumetric flux.
Definition: phaseSystem.H:148
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:46
virtual PtrList< surfaceScalarField > phiKdPhifs(const PtrList< surfaceScalarField > &rAUfs) const =0
Return the force fluxes for the face-based algorithm.
A class for managing temporary objects.
Definition: PtrList.H:53
virtual PtrList< surfaceScalarField > phiKdPhis(const PtrList< volScalarField > &rAUs) const =0
Return the force fluxes for the cell-based algorithm.
const phaseModelPartialList & multiComponentPhases() const
Return the models for phases that have multiple species.
Definition: phaseSystemI.H:91
const dimensionedScalar & deltaN() const
Stabilisation for normalisation of the interface normal.
Definition: phaseSystemI.H:181
const fvMesh & mesh() const
Return the mesh.
Definition: phaseSystemI.H:28
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
Non-orthogonal solution control class. Provides non-orthogonal-loop control methods.
virtual void meshUpdate()
Update the fluid properties for mesh changes.
List of MRF zones with IO functionality. MRF zones are specified by a list of dictionary entries...
Definition: IOMRFZoneList.H:64
virtual bool implicitPhasePressure() const
Returns true if the phase pressure is treated implicitly.
word modelName() const
Return the model name. This is the same as the model&#39;s typename.
HashTable< scalar, phaseInterfaceKey, phaseInterfaceKey::hash > cAlphaTable
Definition: phaseSystem.H:120
virtual void correctSpecies()
Correct the species mass fractions.
Namespace for OpenFOAM.