41 #include "phaseCompressibleTurbulenceModelFwd.H" 63 const phaseSystem& fluid_;
79 autoPtr<diameterModel> diameterModel_;
96 const phaseSystem&
fluid,
97 const word& phaseName,
100 (fluid, phaseName, index)
108 const phaseSystem& fluid,
109 const word& phaseName,
114 autoPtr<phaseModel>
clone()
const;
119 static autoPtr<phaseModel>
New 121 const phaseSystem& fluid,
122 const word& phaseName,
130 const phaseSystem& fluid_;
131 mutable label indexCounter_;
241 virtual bool pure()
const = 0;
virtual tmp< volScalarField > mut() const =0
Return the turbulent dynamic viscosity.
virtual surfaceScalarField & phiRef()=0
Access the volumetric flux.
virtual UPtrList< volScalarField > & YActiveRef()=0
Access the active species mass fractions.
void correctInflowOutflow(surfaceScalarField &alphaPhi) const
Ensure that the flux at inflow/outflow BCs is preserved.
label index() const
Return the index of the phase.
static autoPtr< phaseModel > New(const phaseSystem &fluid, const word &phaseName, const label index)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
virtual tmp< volScalarField > continuityError() const =0
Return the continuity error.
virtual surfaceScalarField & alphaRhoPhiRef()=0
Access the mass flux of the phase.
virtual tmp< volScalarField > divU() const =0
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
virtual bool compressible() const =0
Return true if the phase is compressible otherwise false.
virtual tmp< volScalarField > alpha() const =0
Thermal diffusivity for enthalpy of mixture [kg/m/s].
virtual tmp< volScalarField > alphaEff() const =0
Effective thermal turbulent diffusivity of mixture [kg/m/s].
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
autoPtr< phaseModel > operator()(Istream &is) const
virtual void correctKinematics()
Correct the kinematics.
autoPtr< phaseModel > clone() const
Return clone.
const word & name() const
virtual tmp< volVectorField > DUDt() const =0
Return the substantive acceleration.
virtual tmp< volScalarField > muEff() const =0
Return the effective dynamic viscosity.
tmp< volScalarField > d() const
phaseModel(const word &phaseName, const volScalarField &p, const volScalarField &T)
Construct from components.
virtual bool isothermal() const =0
Return whether the phase is isothermal.
const surfaceScalarField & phi() const
virtual tmp< volScalarField > continuityErrorSources() const =0
Return the continuity error due to any sources.
iNew(const volScalarField &p, const volScalarField &T)
virtual tmp< fvScalarMatrix > YiEqn(volScalarField &Yi)=0
Return the species fraction equation.
const surfaceScalarField & alphaPhi() const
GeometricField< scalar, fvPatchField, volMesh > volScalarField
virtual tmp< fvScalarMatrix > heEqn()=0
Return the enthalpy equation.
const volVectorField & U() const
virtual ~phaseModel()
Destructor.
virtual tmp< volScalarField > mu() const =0
Return the laminar dynamic viscosity.
virtual PtrList< volScalarField > & YRef()=0
Access the species mass fractions.
virtual bool read()
Read phase properties dictionary.
Class to represent a system of phases and model interfacial transfers between them.
virtual tmp< volScalarField > k() const =0
Return the turbulent kinetic energy.
A class for handling words, derived from string.
virtual const UPtrList< volScalarField > & YActive() const =0
Return the active species mass fractions.
const rhoThermo & thermo() const
Return const-access to phase rhoThermo.
virtual bool pure() const =0
Return whether the phase is pure (i.e., not multi-component)
virtual tmp< volScalarField > continuityErrorFlow() const =0
Return the continuity error due to the flow field.
const word & keyword() const
virtual tmp< volScalarField > K() const =0
Return the phase kinetic energy.
virtual tmp< volScalarField > pPrime() const =0
Return the phase-pressure'.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
void correct()
Correct the laminar viscosity.
virtual const PtrList< volScalarField > & Y() const =0
Return the species mass fractions.
const dimensionedScalar & kappa() const
scalar alphaMax() const
Return the maximum phase-fraction (e.g. packing limit)
const dimensionedScalar & nu() const
Return the laminar viscosity.
virtual tmp< fvVectorMatrix > UEqn()=0
Return the momentum equation.
virtual tmp< surfaceScalarField > alphaRhoPhi() const =0
Return the mass flux of the phase.
virtual tmp< surfaceScalarField > DUDtf() const =0
Return the substantive acceleration on the faces.
Forward declarations of fvMatrix specializations.
virtual surfaceScalarField & alphaPhiRef()=0
Access the volumetric flux of the phase.
ClassName("phaseModel")
Runtime type information.
const dimensionedScalar & rho() const
virtual tmp< volScalarField > nuEff() const =0
Return the effective kinematic viscosity.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
const phaseSystem & fluid() const
Return the system to which this phase belongs.
virtual tmp< volScalarField > alphahe() const =0
Thermal diffusivity for energy of mixture [kg/m/s].
const autoPtr< diameterModel > & dPtr() const
Return const-reference to diameterModel of the phase.
virtual tmp< volScalarField > kappaEff() const =0
Effective thermal turbulent diffusivity for temperature.
virtual bool stationary() const =0
Return whether the phase is stationary.
Basic thermodynamic properties based on density.
virtual void correctThermo()
Correct the thermodynamics.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Macros to ease declaration of run-time selection tables.
A class for managing temporary objects.
virtual tmp< volScalarField > nut() const =0
Return the turbulent kinematic viscosity.
virtual void correctEnergyTransport()
Correct the energy transport.
virtual void correctTurbulence()
Correct the turbulence.
declareRunTimeSelectionTable(autoPtr, phaseModel, phaseSystem,(const phaseSystem &fluid, const word &phaseName, const label index),(fluid, phaseName, index))
virtual tmp< fvVectorMatrix > UfEqn()=0
Return the momentum equation for the face-based algorithm.
virtual volVectorField & URef()=0
Access the velocity.
virtual rhoThermo & thermoRef()=0
Access the thermophysical model.