Go to the documentation of this file.
43 #include "phaseCompressibleTurbulenceModelFwd.H"
99 const word& phaseName,
111 const word& phaseName,
124 const word& phaseName,
133 mutable label indexCounter_;
243 virtual bool pure()
const = 0;
virtual tmp< surfaceScalarField > DUDtf() const =0
Return the substantive acceleration on the faces.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
virtual void correctTurbulence()
Correct the turbulence.
virtual bool compressible() const =0
Return true if the phase is compressible otherwise false.
label index() const
Return the index of the phase.
virtual tmp< fvScalarMatrix > heEqn()=0
Return the enthalpy equation.
A class for handling words, derived from Foam::string.
virtual tmp< volScalarField > k() const =0
Return the turbulent kinetic energy.
virtual tmp< fvScalarMatrix > YiEqn(volScalarField &Yi)=0
Return the species fraction equation.
virtual tmp< volScalarField > kappaEff() const =0
Effective thermal turbulent diffusivity for temperature.
virtual tmp< volScalarField > rho() const =0
Return the density field.
tmp< volScalarField > d() const
Return the Sauter-mean diameter.
scalar alphaMax() const
Return the maximum phase-fraction (e.g. packing limit)
virtual void correctEnergyTransport()
Correct the energy transport.
const word & keyword() const
Return the name of the phase for use as the keyword in PtrDictionary.
virtual tmp< volScalarField > continuityErrorFlow() const =0
Return the continuity error due to the flow field.
virtual void correct()
Correct the phase properties.
Forward declarations of fvMatrix specializations.
virtual const rhoThermo & thermo() const =0
Return the thermophysical model.
virtual tmp< fvVectorMatrix > UEqn()=0
Return the momentum equation.
virtual volVectorField & URef()=0
Access the velocity.
static autoPtr< phaseModel > New(const phaseSystem &fluid, const word &phaseName, const label index)
virtual const PtrList< volScalarField > & Y() const =0
Return the species mass fractions.
declareRunTimeSelectionTable(autoPtr, phaseModel, phaseSystem,(const phaseSystem &fluid, const word &phaseName, const label index),(fluid, phaseName, index))
virtual tmp< volScalarField > divU() const =0
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
virtual ~phaseModel()
Destructor.
virtual rhoThermo & thermoRef()=0
Access the thermophysical model.
virtual tmp< volScalarField > kappa() const =0
Thermal diffusivity for temperature of mixture [J/m/s/K].
virtual void correctThermo()
Correct the thermodynamics.
Basic thermodynamic properties based on density.
virtual tmp< volScalarField > nu() const =0
Return the laminar kinematic viscosity.
virtual tmp< volScalarField > alpha() const =0
Thermal diffusivity for enthalpy of mixture [kg/m/s].
virtual tmp< surfaceScalarField > alphaRhoPhi() const =0
Return the mass flux of the phase.
virtual tmp< volScalarField > continuityErrorSources() const =0
Return the continuity error due to any sources.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
virtual tmp< volScalarField > nuEff() const =0
Return the effective kinematic viscosity.
autoPtr< phaseModel > operator()(Istream &is) const
virtual tmp< volScalarField > alphahe() const =0
Thermal diffusivity for energy of mixture [kg/m/s].
ClassName("phaseModel")
Runtime type information.
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Return a pointer to a new phase created on freestore.
virtual bool isothermal() const =0
Return whether the phase is isothermal.
virtual tmp< volScalarField > mu() const =0
Return the laminar dynamic viscosity.
virtual PtrList< volScalarField > & YRef()=0
Access the species mass fractions.
virtual const UPtrList< volScalarField > & YActive() const =0
Return the active species mass fractions.
const autoPtr< diameterModel > & dPtr() const
Return const-reference to diameterModel of the phase.
virtual tmp< surfaceScalarField > phi() const =0
Return the volumetric flux.
virtual tmp< surfaceScalarField > alphaPhi() const =0
Return the volumetric flux of the phase.
const word & name() const
Return the name of this phase.
virtual tmp< volScalarField > mut() const =0
Return the turbulent dynamic viscosity.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
virtual surfaceScalarField & phiRef()=0
Access the volumetric flux.
virtual tmp< fvVectorMatrix > UfEqn()=0
Return the momentum equation for the face-based algorithm.
void correctInflowOutflow(surfaceScalarField &alphaPhi) const
Ensure that the flux at inflow/outflow BCs is preserved.
Macros to ease declaration of run-time selection tables.
virtual surfaceScalarField & alphaRhoPhiRef()=0
Access the mass flux of the phase.
virtual bool read()
Read phase properties dictionary.
virtual tmp< volScalarField > continuityError() const =0
Return the continuity error.
const phaseSystem & fluid() const
Return the system to which this phase belongs.
virtual tmp< volVectorField > U() const =0
Return the velocity.
iNew(const phaseSystem &fluid)
virtual tmp< volScalarField > alphaEff() const =0
Effective thermal turbulent diffusivity of mixture [kg/m/s].
virtual bool pure() const =0
Return whether the phase is pure (i.e., not multi-component)
Class to represent a system of phases and model interfacial transfers between them.
virtual surfaceScalarField & alphaPhiRef()=0
Access the volumetric flux of the phase.
virtual tmp< volScalarField > K() const =0
Return the phase kinetic energy.
virtual tmp< volScalarField > pPrime() const =0
Return the phase-pressure'.
autoPtr< phaseModel > clone() const
Return clone.
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
virtual tmp< volScalarField > muEff() const =0
Return the effective dynamic viscosity.
virtual UPtrList< volScalarField > & YActiveRef()=0
Access the active species mass fractions.
virtual bool stationary() const =0
Return whether the phase is stationary.
virtual void correctKinematics()
Correct the kinematics.
phaseModel(const phaseSystem &fluid, const word &phaseName, const label index)
virtual tmp< volVectorField > DUDt() const =0
Return the substantive acceleration.
virtual tmp< volScalarField > nut() const =0
Return the turbulent kinematic viscosity.