41 #include "phaseModel.H" 42 #include "phasePair.H" 43 #include "orderedPhasePair.H" 60 template<
class modelType>
class BlendedInterfacialModel;
61 class surfaceTensionModel;
62 class aspectRatioModel;
172 const phaseModelList& phaseModels
178 const dictTable& modelDicts
182 template<
class modelType>
185 const dictTable& modelDicts,
195 template<
class modelType>
198 const word& modelName,
208 template<
class modelType>
211 const word& modelName,
218 const bool correctFixedFluxBCs =
true 222 template<
class modelType>
225 const word& modelName,
232 const bool correctFixedFluxBCs =
true 237 template<
class GeoField>
241 const word& fieldName,
248 template<
class GeoField>
252 const word& fieldName,
253 const GeoField& field,
259 template<
class GeoField>
263 const word& fieldName,
270 template<
class GeoField>
274 const word& fieldName,
275 const GeoField& field,
307 inline const phaseModelList&
phases()
const;
310 inline phaseModelList&
phases();
313 inline const phaseModelPartialList&
movingPhases()
const;
337 inline const phasePairTable&
phasePairs()
const;
361 template<
class modelType>
365 template<
class modelType>
373 template<
class modelType>
384 template<
class>
class PatchField,
398 template<
class>
class PatchField,
491 const bool includeVirtualMass =
false 508 virtual void solve();
const BlendedInterfacialModel< modelType > & lookupBlendedSubModel(const phasePair &key) const
Return a blended sub model between a phase pair.
A simple container for copying or transferring objects of type <T>.
static const word propertiesName
Default name of the phase properties dictionary.
Hashing function class, shared by all the derived classes.
virtual Xfer< PtrList< surfaceScalarField > > ddtCorrByAs(const PtrList< volScalarField > &rAUs, const bool includeVirtualMass=false) const =0
Return the flux corrections for the cell-based algorithm.
void generatePairs(const dictTable &modelDicts)
Generate pairs.
fv::options & fvOptions() const
Access the fvOptions.
const phaseModelPartialList & movingPhases() const
Return the models for phases that are moving.
surfaceTensionModelTable surfaceTensionModels_
Surface tension models.
tmp< volScalarField > E(const phasePairKey &key) const
Return the aspect-ratio for a pair.
virtual void correctTurbulence()
Correct the turbulence.
virtual Xfer< PtrList< volVectorField > > KdUByAs(const PtrList< volScalarField > &rAUs) const =0
Return the explicit part of the drag force.
PtrList< volScalarField > rAUs
HashPtrTable< fvScalarMatrix > massTransferTable
HashTable< autoPtr< phasePair >, phasePairKey, phasePairKey::hash > phasePairTable
HashPtrTable< fvVectorMatrix > momentumTransferTable
virtual void solve()
Solve for the phase fractions.
virtual Xfer< PtrList< surfaceScalarField > > phiKdPhifs(const PtrList< surfaceScalarField > &rAUfs) const =0
Return the force fluxes for the face-based algorithm.
phasePairTable phasePairs_
Phase pairs.
const modelType & lookupSubModel(const phasePair &key) const
Return a sub model between a phase pair.
virtual void partialElimination(const PtrList< volScalarField > &rAUs)=0
Solve the drag system for the new velocities and fluxes.
virtual autoPtr< heatTransferTable > heatTransfer() const =0
Return the heat transfer matrices.
virtual Xfer< PtrList< surfaceScalarField > > AFfs() const =0
Return the implicit force coefficients for the face-based.
tmp< volScalarField > sigma(const phasePairKey &key) const
Return the surface tension coefficient for a pair.
volScalarField dpdt_
Rate of change of pressure.
TypeName("phaseSystem")
Runtime type information.
virtual autoPtr< momentumTransferTable > momentumTransferf()=0
Return the momentum transfer matrices for the face-based.
HashPtrTable< fvScalarMatrix > heatTransferTable
phaseModelPartialList stationaryPhaseModels_
Stationary phase models.
A HashTable specialization for hashing pointers.
void createSubModels(const dictTable &modelDicts, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables.
virtual tmp< volScalarField > dmdt(const phasePairKey &key) const
Return the mass transfer rate for a pair.
HashTable< autoPtr< surfaceTensionModel >, phasePairKey, phasePairKey::hash > surfaceTensionModelTable
virtual autoPtr< momentumTransferTable > momentumTransfer()=0
Return the momentum transfer matrices for the cell-based.
phaseModelList phaseModels_
Phase models.
phaseModelPartialList multiComponentPhaseModels_
Multi-component phase models.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
const phaseModelPartialList & stationaryPhases() const
Return the models for phases that are stationary.
const phasePairTable & phasePairs() const
Return the phase pairs.
const phaseModelList & phases() const
Return the phase models.
virtual void correctKinematics()
Correct the kinematics.
Dimension set for the base types.
An ordered pair of two objects of type <T> with first() and second() elements.
phaseSystem(const fvMesh &mesh)
Construct from fvMesh.
Class to represent a system of phases and model interfacial transfers between them.
virtual void partialEliminationf(const PtrList< surfaceScalarField > &rAUfs)=0
Solve the drag system for the new fluxes.
A class for handling words, derived from string.
tmp< volVectorField > U() const
Return the mixture velocity.
const IOMRFZoneList & MRF() const
Return MRF zones.
const phaseModelPartialList & anisothermalPhases() const
Return the models for phases that have variable temperature.
virtual void correctThermo()
Correct the thermodynamics.
PtrListDictionary< phaseModel > phaseModelList
UPtrList< phaseModel > phaseModelPartialList
const fvMesh & mesh_
Reference to the mesh.
const word & name() const
Name function is needed to disambiguate those inherited.
An STL-conforming hash table.
PtrList< surfaceScalarField > rAUfs
virtual Xfer< PtrList< surfaceScalarField > > phiKdPhis(const PtrList< volScalarField > &rAUs) const =0
Return the force fluxes for the cell-based algorithm.
phaseModelPartialList anisothermalPhaseModels_
Anisothermal phase models.
phaseModelPartialList movingPhaseModels_
Moving phase models.
virtual void correct()
Correct the fluid properties other than those listed below.
tmp< volScalarField > rho() const
Return the mixture density.
virtual Xfer< PtrList< surfaceScalarField > > phiFfs(const PtrList< surfaceScalarField > &rAUfs)=0
Return the force fluxes for the face-based algorithm.
aspectRatioModelTable aspectRatioModels_
Aspect ratio models.
IOMRFZoneList MRF_
Optional MRF zones.
volScalarField sf(fieldObject, mesh)
const volScalarField & dpdt() const
Return the rate of change of the pressure.
tmp< volScalarField > byDt(const volScalarField &vf)
virtual Xfer< PtrList< surfaceScalarField > > phiFs(const PtrList< volScalarField > &rAUs)=0
Return the force fluxes for the cell-based algorithm.
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.
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
Return the mixture flux.
virtual ~phaseSystem()
Destructor.
void addField(const phaseModel &phase, const word &fieldName, tmp< GeoField > field, PtrList< GeoField > &fieldList) const
Add the field to a phase-indexed list, with the given name,.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
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.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
virtual autoPtr< massTransferTable > massTransfer() const =0
Return the mass transfer matrices.
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...
surfaceScalarField phi_
Total volumetric flux.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
A class for managing temporary objects.
const phaseModelPartialList & multiComponentPhases() const
Return the models for phases that have multiple species.
const fvMesh & mesh() const
Return 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
virtual const HashPtrTable< surfaceScalarField > & DByAfs() const =0
Return the phase diffusivities divided by the momentum.
virtual Xfer< PtrList< volScalarField > > dmdts() const
Return the mass transfer rates for each phase.