42 #ifndef StationaryPhaseModel_H 43 #define StationaryPhaseModel_H 45 #include "phaseModel.H" 56 template<
class BasePhaseModel>
66 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
71 const bool cache =
false 80 const bool cache =
false 89 const bool cache =
false 100 const word& phaseName,
190 using BasePhaseModel::kappaEff;
virtual tmp< volScalarField > mut() const
Return the turbulent dynamic viscosity.
virtual tmp< surfaceScalarField > alphaPhi() const
Return the volumetric flux of the phase.
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 surfaceScalarField & phiRef()
Access the volumetric flux.
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
virtual tmp< fvVectorMatrix > UfEqn()
Return the momentum equation for the face-based algorithm.
virtual tmp< surfaceScalarField > phi() const
Return the volumetric flux.
virtual tmp< surfaceScalarField > DUDtf() const
Return the substantive acceleration on the faces.
virtual tmp< volScalarField > continuityErrorSources() const
Return the continuity error due to any sources.
virtual tmp< volScalarField > divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
virtual tmp< volScalarField > nuEff() const
Return the effective kinematic viscosity.
virtual tmp< volScalarField > kappaEff() const
Return the effective thermal conductivity.
virtual tmp< volVectorField > U() const
Return the velocity.
Dimension set for the base types.
Class which represents a stationary (and therefore probably solid) phase. Generates, but does not store, zero velocity and flux field and turbulent qauantities. Throws an error when non-const access is requested to the motion fields or when the momentum equation is requested. Usage must must protect against such calls.
Class to represent a system of phases and model interfacial transfers between them.
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
A class for handling words, derived from string.
virtual tmp< surfaceScalarField > alphaRhoPhi() const
Return the mass flux of the phase.
virtual surfaceScalarField & alphaRhoPhiRef()
Access the mass flux of the phase.
virtual tmp< volScalarField > continuityErrorFlow() const
Return the continuity error due to the flow field.
virtual tmp< volScalarField > k() const
Return the turbulent kinetic energy.
word name(const complex &)
Return a string representation of a complex.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure'.
virtual surfaceScalarField & alphaPhiRef()
Access the volumetric flux of the phase.
virtual volVectorField & URef()
Access the velocity.
virtual tmp< volScalarField > muEff() const
Return the effective dynamic viscosity.
virtual tmp< volScalarField > K() const
Return the phase kinetic energy.
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
A class for managing temporary objects.
virtual tmp< volScalarField > alphaEff() const
Return the effective thermal diffusivity.
virtual tmp< volScalarField > continuityError() const
Return the continuity error.
virtual ~StationaryPhaseModel()
Destructor.
virtual tmp< volScalarField > nut() const
Return the turbulent kinematic viscosity.
virtual bool stationary() const
Return whether the phase is stationary.
virtual tmp< volVectorField > DUDt() const
Return the substantive acceleration.
StationaryPhaseModel(const phaseSystem &fluid, const word &phaseName, const label index)