45 #ifndef MovingPhaseModel_H
46 #define MovingPhaseModel_H
61 template<
class ThermoModel>
80 template<
class BasePhaseModel>
93 typename BasePhaseModel::thermoModel
159 const word& phaseName,
160 const bool referencePhase,
Generic GeometricField class.
Class which represents a moving fluid phase. Holds the velocity, fluxes and momentumTransport model a...
virtual void correctUf()
Correct the face velocity for moving meshes.
virtual void correctKinematics()
Correct the kinematics.
autoPtr< volScalarField > divU_
Dilatation rate.
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
virtual tmp< surfaceScalarField > DUDtf() const
Return the substantive acceleration on the faces.
MovingPhaseModelTransportThermoModel< typename BasePhaseModel::thermoModel >::type transportThermoModel
Thermo type for the thermophysical transport model.
volVectorField U_
Velocity field.
virtual void correct()
Correct the phase properties other than the thermo.
virtual tmp< volScalarField > k() const
Return the turbulent kinetic energy.
virtual void predictThermophysicalTransport()
Predict the energy transport e.g. alphat.
virtual ~MovingPhaseModel()
Destructor.
virtual tmp< volVectorField > U() const
Return the velocity.
virtual void correctMomentumTransport()
Correct the momentumTransport.
tmp< surfaceScalarField > DUDtf_
Lagrangian acceleration field on the faces (needed for virtual-mass)
autoPtr< PhaseThermophysicalTransportModel< phaseCompressible::momentumTransportModel, transportThermoModel > > thermophysicalTransport_
Thermophysical transport model.
volScalarField continuityError_
Continuity error.
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction.
virtual surfaceScalarField & alphaRhoPhiRef()
Access the mass flux of the phase.
surfaceScalarField alphaRhoPhi_
Mass flux.
virtual const autoPtr< volScalarField > & divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
virtual void correctThermophysicalTransport()
Correct the energy transport e.g. alphat.
virtual surfaceScalarField & phiRef()
Access the volumetric flux.
virtual tmp< surfaceScalarField > alphaPhi() const
Return the volumetric flux of the phase.
virtual tmp< volScalarField > K() const
Return the phase kinetic energy.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure'.
surfaceScalarField phi_
Flux.
MovingPhaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
surfaceScalarField alphaPhi_
Volumetric flux.
virtual const autoPtr< surfaceVectorField > & Uf() const
Return the face velocity.
virtual surfaceVectorField & UfRef()
Access the face velocity.
tmp< volVectorField > DUDt_
Lagrangian acceleration field (needed for virtual-mass)
virtual surfaceScalarField & alphaPhiRef()
Access the volumetric flux of the phase.
virtual tmp< surfaceScalarField > alphaRhoPhi() const
Return the mass flux of the phase.
virtual tmp< surfaceScalarField > phi() const
Return the volumetric flux.
virtual tmp< fvVectorMatrix > UfEqn()
Return the momentum equation for the face-based algorithm.
virtual volVectorField & URef()
Access the velocity.
virtual tmp< volVectorField > DUDt() const
Return the substantive acceleration.
virtual bool stationary() const
Return whether the phase is stationary.
tmp< volScalarField > K_
Kinetic Energy.
virtual void predictMomentumTransport()
Predict the momentumTransport.
virtual tmp< volScalarField > continuityError() const
Return the continuity error.
virtual tmp< scalarField > kappaEff(const label patchi) const
Return the effective thermal conductivity on a patch.
autoPtr< phaseCompressible::momentumTransportModel > momentumTransport_
Turbulence model.
virtual void correctContinuityError(const volScalarField &source)
Correct the continuity error.
autoPtr< surfaceVectorField > Uf_
Face velocity field.
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
Templated base class for multiphase thermophysical transport models.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Base-class for multi-component fluid thermodynamic properties.
Base-class for fluid thermodynamic properties.
Class to represent a system of phases and model interfacial transfers between them.
Base-class for multi-component fluid thermodynamic properties based on density.
Base-class for fluid thermodynamic properties based on density.
A class for managing temporary objects.
A class for handling words, derived from string.
phaseCompressibleMomentumTransportModel momentumTransportModel
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.