35 template<
class BasePhaseModel>
39 const tmp<volScalarField>& pressureWork
44 scalar pressureWorkAlphaLimit =
45 this->thermo_->properties()
46 .lookupOrDefault(
"pressureWorkAlphaLimit", 0.0);
48 if (pressureWorkAlphaLimit > 0)
52 max(
alpha - pressureWorkAlphaLimit, scalar(0))
53 /
max(
alpha - pressureWorkAlphaLimit, pressureWorkAlphaLimit)
65 template<
class BasePhaseModel>
69 const word& phaseName,
70 const bool referencePhase,
74 BasePhaseModel(fluid, phaseName, referencePhase, index),
76 thermophysicalTransport_
82 >::
New(this->momentumTransport_, this->thermo_)
89 template<
class BasePhaseModel>
96 template<
class BasePhaseModel>
99 BasePhaseModel::correctThermo();
101 this->thermo_->correct();
105 template<
class BasePhaseModel>
112 template<
class BasePhaseModel>
116 BasePhaseModel::predictThermophysicalTransport();
117 thermophysicalTransport_->predict();
121 template<
class BasePhaseModel>
125 BasePhaseModel::correctThermophysicalTransport();
126 thermophysicalTransport_->correct();
130 template<
class BasePhaseModel>
134 return thermophysicalTransport_->kappaEff(
patchi);
138 template<
class BasePhaseModel>
142 return thermophysicalTransport_->divq(
he);
146 template<
class BasePhaseModel>
150 return thermophysicalTransport_->divj(Yi);
154 template<
class BasePhaseModel>
191 if (he.
name() == this->thermo_->phasePropertyName(
"e"))
193 tEEqn.
ref() += filterPressureWork
203 else if (this->thermo_->dpdt())
205 tEEqn.
ref() -= filterPressureWork(
alpha*this->fluid().dpdt());
Class which represents a phase for which the temperature (strictly energy) varies....
virtual ~AnisothermalPhaseModel()
Destructor.
AnisothermalPhaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
virtual void correctThermo()
Correct the thermodynamics.
virtual void predictThermophysicalTransport()
Predict the energy transport e.g. alphat.
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction.
virtual void correctThermophysicalTransport()
Correct the energy transport e.g. alphat.
virtual bool isothermal() const
Return whether the phase is isothermal.
virtual tmp< fvScalarMatrix > heEqn()
Return the enthalpy equation.
virtual tmp< scalarField > kappaEff(const label patchi) const
Return the effective thermal conductivity on a patch.
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
Generic GeometricField class.
const word & name() const
Return name.
Templated base class for multiphase thermophysical transport models.
Base-class for fluid thermodynamic properties.
Abstract base class for turbulence models (RAS, LES and laminar).
Class to represent a system of phases and model interfacial transfers between them.
A class for managing temporary objects.
T & ref() const
Return non-const reference or generate a fatal error.
A class for handling words, derived from string.
Calculate the first temporal derivative.
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for implicit and explicit sources.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
VolField< scalar > volScalarField
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)