47#ifndef MovingPhaseModel_H
48#define MovingPhaseModel_H
50#include "phaseModel.H"
51#include "phaseCompressibleTurbulenceModel.H"
62template<
class BasePhaseModel>
121 const word& phaseName,
224 using BasePhaseModel::kappaEff;
233 using BasePhaseModel::alphaEff;
Class which represents a moving fluid phase. Holds the velocity, fluxes and turbulence model....
virtual void correctKinematics()
Correct the kinematics.
virtual tmp< volScalarField > continuityErrorSources() const
Return the continuity error due to any sources.
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
virtual tmp< surfaceScalarField > DUDtf() const
Return the substantive acceleration on the faces.
virtual tmp< volScalarField > muEff() const
Return the effective dynamic viscosity.
tmp< volScalarField > divU_
virtual void correct()
Correct the phase properties other than the thermo and turbulence.
virtual void divU(tmp< volScalarField > divU)
Set the phase dilatation rate (d(alpha)/dt + div(alpha*phi)).
virtual tmp< volScalarField > k() const
Return the turbulent kinetic energy.
virtual ~MovingPhaseModel()=default
Destructor.
virtual tmp< volVectorField > U() const
tmp< surfaceScalarField > DUDtf_
virtual tmp< volScalarField > mut() const
Return the turbulent dynamic viscosity.
volScalarField continuityErrorFlow_
virtual tmp< surfaceScalarField > alphaPhi() const
Return the volumetric flux of the phase.
autoPtr< phaseCompressibleTurbulenceModel > turbulence_
virtual surfaceScalarField & alphaRhoPhiRef()
Access the mass flux of the phase.
surfaceScalarField alphaRhoPhi_
virtual void correctTurbulence()
Correct the turbulence.
virtual surfaceScalarField & phiRef()
Access the volumetric flux.
virtual tmp< scalarField > alphaEff(const label patchi) const
Return the effective thermal conductivity on a patch.
virtual tmp< volScalarField > alphaEff() const
Return the effective thermal diffusivity.
virtual tmp< surfaceScalarField > alphaPhi() const
Constant access 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'.
MovingPhaseModel(const phaseSystem &fluid, const word &phaseName, const label index)
Construct from phase system and phase name.
virtual void correct()
Correct the phase properties other than the thermo and turbulence.
virtual tmp< volScalarField > continuityErrorFlow() const
Return the continuity error due to the flow field.
virtual tmp< volScalarField > nuEff() const
Return the effective kinematic viscosity.
MovingPhaseModel(const multiphaseInterSystem &fluid, const word &phaseName)
virtual tmp< volScalarField > divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi)).
tmp< volVectorField > DUDt_
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
Constant access the volumetric flux.
virtual tmp< fvVectorMatrix > UfEqn()
Return the momentum equation for the face-based algorithm.
virtual tmp< volScalarField > nut() const
Return the turbulent kinematic viscosity.
virtual volVectorField & URef()
Access the velocity.
virtual tmp< volVectorField > DUDt() const
Return the substantive acceleration.
virtual tmp< volScalarField > kappaEff() const
Return the effective thermal conductivity.
virtual tmp< surfaceScalarField > phi() const
Return the volumetric flux.
virtual bool stationary() const
Return whether the phase is stationary.
virtual tmp< volVectorField > U() const
Return the velocity.
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.
volScalarField continuityErrorSources_
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Class to represent a system of phases and model interfacial transfers between them.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField