41 #include "phaseModel.H" 42 #include "phasePair.H" 43 #include "orderedPhasePair.H" 60 template<
class modelType>
class BlendedInterfacialModel;
61 class surfaceTensionModel;
62 class aspectRatioModel;
174 volScalarField
dpdt_;
197 const phaseModelList& phaseModels
203 const dictTable& modelDicts
207 template<
class modelType>
210 const dictTable& modelDicts,
220 template<
class modelType>
223 const word& modelName,
233 template<
class modelType>
236 const word& modelName,
246 template<
class modelType>
249 const word& modelName,
284 inline const phaseModelList&
phases()
const;
287 inline phaseModelList&
phases();
290 inline const phasePairTable&
phasePairs()
const;
305 inline const volScalarField&
dpdt()
const;
308 inline volScalarField&
dpdt();
323 template<
class modelType>
327 template<
class modelType>
335 virtual void solve();
static const word propertiesName
Default name of the phase properties dictionary.
Hashing function class, shared by all the derived classes.
void generatePairs(const dictTable &modelDicts)
Generate pairs.
fv::options & fvOptions() const
Optional FV-options.
fvMatrix< scalar > fvScalarMatrix
surfaceTensionModelTable surfaceTensionModels_
Surface tension models.
tmp< volScalarField > E(const phasePairKey &key) const
Return the aspect-ratio.
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > VmTable
virtual void correctTurbulence()
Correct the turbulence.
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > KdTable
HashTable< autoPtr< phasePair >, phasePairKey, phasePairKey::hash > phasePairTable
virtual void solve()
Solve for the phase fractions.
phasePairTable phasePairs_
Phase pairs.
const modelType & lookupSubModel(const phasePair &key) const
Access a sub model between a phase pair.
tmp< volScalarField > sigma(const phasePairKey &key) const
Return the surface tension coefficient.
volScalarField dpdt_
Rate of change of pressure.
TypeName("phaseSystem")
Runtime type information.
A HashTable specialization for hashing pointers.
void createSubModels(const dictTable &modelDicts, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables.
HashTable< autoPtr< surfaceTensionModel >, phasePairKey, phasePairKey::hash > surfaceTensionModelTable
phaseModelList phaseModels_
Phase models.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const phasePairTable & phasePairs() const
Constant access the phase pairs.
HashPtrTable< fvScalarMatrix, word, string::hash > heatTransferTable
const phaseModelList & phases() const
Constant access the phase models.
virtual void correctKinematics()
Correct the kinematics.
phaseSystem(const fvMesh &mesh)
Construct from fvMesh.
Class to represent a system of phases and model interfacial transfers between them.
A class for handling words, derived from string.
tmp< volVectorField > U() const
Return the mixture velocity.
HashPtrTable< fvVectorMatrix, word, string::hash > momentumTransferTable
const IOMRFZoneList & MRF() const
Return MRF zones.
virtual void correctThermo()
Correct the thermodynamics.
PtrListDictionary< phaseModel > phaseModelList
const fvMesh & mesh_
Reference to the mesh.
An STL-conforming hash table.
virtual void correct()
Correct the fluid properties other than the thermo and turbulence.
tmp< volScalarField > rho() const
Return the mixture density.
aspectRatioModelTable aspectRatioModels_
Aspect ratio models.
IOMRFZoneList MRF_
Optional MRF zones.
const volScalarField & dpdt() const
Constant access the rate of change of the pressure.
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
void generatePairsAndSubModels(const word &modelName, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables.
Forward declarations of fvMatrix specializations.
tmp< surfaceScalarField > calcPhi(const phaseModelList &phaseModels) const
Calculate and return the mixture flux.
const surfaceScalarField & phi() const
Constant access the mixture flux.
virtual ~phaseSystem()
Destructor.
virtual bool read()
Read base phaseProperties dictionary.
blendingMethodTable blendingMethods_
Blending methods.
HashTable< autoPtr< blendingMethod >, word, word::hash > blendingMethodTable
Mesh data needed to do the Finite Volume discretisation.
fvMatrix< vector > fvVectorMatrix
HashTable< autoPtr< aspectRatioModel >, phasePairKey, phasePairKey::hash > aspectRatioModelTable
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
HashPtrTable< fvScalarMatrix, word, string::hash > massTransferTable
surfaceScalarField phi_
Total volumetric flux.
A class for managing temporary objects.
const fvMesh & mesh() const
Constant access the mesh.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
List of MRF zones with IO functionality. MRF zones are specified by a list of dictionary entries...
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable