Go to the documentation of this file.
46 #ifndef MovingPhaseModel_H
47 #define MovingPhaseModel_H
49 #include "phaseModel.H"
61 template<
class BasePhaseModel>
110 const word& phaseName,
volVectorField U_
Velocity field.
tmp< volScalarField > divU_
Dilatation rate.
surfaceScalarField alphaRhoPhi_
Mass flux.
A class for handling words, derived from string.
surfaceScalarField phi_
Flux.
virtual tmp< surfaceScalarField > alphaPhi() const
Constant access the volumetric flux of the phase.
virtual tmp< volVectorField > U() const
Constant access the velocity.
A class for managing temporary objects.
volScalarField continuityError_
Continuity error.
surfaceScalarField alphaPhi_
Volumetric flux.
tmp< surfaceScalarField > DbyA_
Phase diffusivity divided by the momentum coefficient.
Class which represents a moving fluid phase. Holds the velocity, fluxes and turbulence model....
virtual ~MovingPhaseModel()
Destructor.
virtual tmp< surfaceScalarField > phi() const
Constant access the volumetric flux.
virtual const surfaceScalarField & DbyA() const
Return the phase diffusivity divided by the momentum coefficient.
virtual void correctKinematics()
Correct the kinematics.
virtual tmp< volVectorField > DUDt() const
Return the substantive acceleration.
MovingPhaseModel(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 void correctTurbulence()
Correct the turbulence.
autoPtr< phaseCompressibleTurbulenceModel > turbulence_
Turbulence model.
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
virtual const phaseCompressibleTurbulenceModel & turbulence() const
Return the turbulence model.
virtual tmp< surfaceScalarField > alphaRhoPhi() const
Constant access the mass flux of the phase.
virtual const tmp< volScalarField > & divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
virtual void correct()
Correct the phase properties other than the thermo and turbulence.
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
virtual tmp< volScalarField > continuityError() const
Constant access the continuity error.
Class to represent a system of phases and model interfacial transfers between them.
Generic GeometricField class.
volVectorField DUDt_
Lagrangian acceleration field (needed for virtual-mass)