Go to the documentation of this file.
41 #include "phaseModel.H"
42 #include "phasePair.H"
43 #include "orderedPhasePair.H"
56 class surfaceTensionModel;
66 public compressibleTransportModel
75 autoPtr<phasePair>, phasePairKey, phasePairKey::hash
170 template<
class modelType>
183 template<
class modelType>
197 template<
class modelType>
200 const word& modelName,
211 template<
class modelType>
214 const word& modelName,
226 template<
class modelType>
229 const word& modelName,
514 const word speciesName
521 virtual void solve() = 0;
556 template <
class modelType>
560 template <
class modelType>
phaseModelList phaseModels_
Phase models.
List< label > labelList
A List of labels.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
HashTable< autoPtr< phaseModel > > generatePhaseModels(const wordList &names) const
Generate the phases.
HashTable< autoPtr< phasePair >, phasePairKey, phasePairKey::hash > phasePairTable
interfacePorousModelTable interfacePorousModelTable_
Interface porous models.
A class for handling words, derived from Foam::string.
virtual tmp< volScalarField > alphaEff(const volScalarField &alphat) const
Effective thermal diffusivity of mixture [kg/m/s].
tmp< surfaceVectorField > nHatfv(const volScalarField &alpha1, const volScalarField &alpha2) const
Interface normal surface vector.
surfaceTensionModelTable surfaceTensionModels_
Surface tension models.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
A class for managing temporary objects.
virtual void solve()
Solve for the phase fractions.
HashTable< autoPtr< phasePair >, phasePairKey, phasePairKey::hash > phasePairTable
Ordered or unordered hashing of word pair.
virtual tmp< volScalarField > CpByCpv() const
Heat capacity ratio [].
const volScalarField & alpha2
tmp< surfaceScalarField > surfaceTensionForce() const
Calculate surface tension of the mixture.
virtual tmp< volScalarField > coeffs(const word &key) const
Return coefficients (1/rho)
const surfaceScalarField & rhoPhi() const
Constant access to the mixture mass flux.
phaseSystem(const fvMesh &mesh)
Construct from fvMesh.
void createSubModels(const dictTable &modelDicts, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables.
const surfaceScalarField & phi() const
Return the mixture flux.
Forward declarations of fvMatrix specializations.
tmp< volVectorField > U() const
Return the mixture velocity.
const modelType & lookupSubModel(const phasePair &key) const
Return a sub model between a phase pair.
void generatePairsTable()
Generate pair table.
const volScalarField & alpha1
virtual ~phaseSystem()
Destructor.
virtual tmp< volScalarField > mu() const
Dynamic viscosity of mixture [kg/m/s].
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable
virtual bool isochoric() const
Return true if the equation of state is isochoric for all phasses.
volScalarField & h
Planck constant.
virtual volScalarField & he()
Return access to the inernal energy field [J/Kg].
tmp< surfaceScalarField > nHatf(const volScalarField &alpha1, const volScalarField &alpha2) const
Interface normal volume vector.
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
const phasePairTable & totalPhasePairs() const
Constant access the total phase pairs.
virtual tmp< volScalarField > kappaEff(const volScalarField &kappat) const
Effective thermal diffusivity for temperature.
fvMatrix< vector > fvVectorMatrix
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
surfaceScalarField phi_
Total volumetric flux.
virtual tmp< scalarField > THE(const scalarField &h, const scalarField &p, const scalarField &T0, const labelList &cells) const
Temperature from enthalpy/internal energy for cell-set.
virtual tmp< volScalarField > Cv() const
Return Cv of the mixture.
void generatePairsAndSubModels(const word &modelName, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables.
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
virtual autoPtr< heatTransferTable > heatTransfer() const =0
Return the heat transfer matrices.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
phasePairTable totalPhasePairs_
Total ordered phase pairs in the system.
virtual void correct()
Correct the fluid properties other than those listed below.
const fvMesh & mesh_
Reference to the mesh.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol] of the mixture.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
virtual word thermoName() const
Return the name of the thermo physics.
const dimensionedScalar & Prt() const
Return Prandt number.
HashTable< autoPtr< surfaceTensionModel >, phasePairKey, phasePairKey::hash > surfaceTensionModelTable
virtual tmp< volScalarField > hc() const
Chemical enthalpy of the mixture [J/kg].
tmp< surfaceScalarField > generatePhi(const HashTable< autoPtr< phaseModel >> &phaseModels) const
Generate the mixture flux.
virtual tmp< volScalarField > dmdt(const phasePairKey &key) const
Return the mass transfer rate for a pair.
virtual tmp< volScalarField > Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
Mesh data needed to do the Finite Volume discretisation.
surfaceScalarField rhoPhi_
Mixture total mass flux.
wordList phaseNames_
Phase names.
const fvMesh & mesh() const
Return the mesh.
HashTable< autoPtr< porousModel >, phasePairKey, phasePairKey::hash > interfacePorousModelTable
static const word phasePropertiesName
Default name of the phase properties dictionary.
A HashTable similar to std::unordered_map.
tmp< volScalarField > K(const volScalarField &alpha1, const volScalarField &alpha2) const
Interface curvature.
dimensionedScalar Prt_
Turbulent Prandt number.
TypeName("phaseSystem")
Runtime type information.
virtual void correctTurbulence()
Correct the turbulence.
phaseModelTable phaseModels_
Phase models.
virtual tmp< volScalarField > Cp() const
Return Cp of the mixture.
virtual tmp< volScalarField > gamma() const
Gamma = Cp/Cv [].
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
HashTable< autoPtr< phaseModel > > phaseModelTable
virtual tmp< volScalarField > surfaceTensionCoeff(const phasePairKey &key) const
Return the surface tension coefficient.
void addInterfacePorosity(fvVectorMatrix &UEqn)
Add interface porosity on phasePair.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
virtual void massSpeciesTransfer(const phaseModel &phase, volScalarField::Internal &Su, volScalarField::Internal &Sp, const word speciesName)=0
Calculate mass transfer.
virtual tmp< volScalarField > alphahe() const
Thermal diffusivity for energy of mixture [kg/m/s].
tmp< volScalarField > nearInterface() const
Near Interface of alpha'n.
phasePairTable phasePairs_
Phase pairs.
const phaseModelList & phases() const
Return the phase models.
static const GeometricField< scalar, fvPatchField, volMesh > & null()
Return a null geometric field.
virtual tmp< volScalarField > kappa() const
Thermal diffusivity for temperature of mixture [J/m/s/K].
tmp< volScalarField > rho() const
Return the mixture density.
virtual bool incompressible() const
Return true if the equation of state is incompressible for all.
virtual bool read()
Read base phaseProperties dictionary.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
void generatePairs(const dictTable &modelDicts)
Generate pairs.