MovingPhaseModel< BasePhaseModel > Class Template Reference

Class which represents a moving fluid phase. Holds the velocity, fluxes and turbulence model and can generate the momentum equation. The interface is quite restrictive as it also has to support an equivalent stationary model, which does not store motion fields or a turbulence model. More...

Inheritance diagram for MovingPhaseModel< BasePhaseModel >:
Collaboration diagram for MovingPhaseModel< BasePhaseModel >:

Public Member Functions

 MovingPhaseModel (const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
 
virtual ~MovingPhaseModel ()
 Destructor. More...
 
virtual void correct ()
 Correct the phase properties other than the thermo and turbulence. More...
 
virtual void correctContinuityError (const volScalarField &source)
 Correct the continuity error. More...
 
virtual void correctKinematics ()
 Correct the kinematics. More...
 
virtual void correctTurbulence ()
 Correct the turbulence. More...
 
virtual void correctEnergyTransport ()
 Correct the energy transport e.g. alphat. More...
 
virtual bool stationary () const
 Return whether the phase is stationary. More...
 
virtual tmp< fvVectorMatrixUEqn ()
 Return the momentum equation. More...
 
virtual tmp< fvVectorMatrixUfEqn ()
 Return the momentum equation for the face-based algorithm. More...
 
virtual tmp< volVectorFieldU () const
 Return the velocity. More...
 
virtual volVectorFieldURef ()
 Access the velocity. More...
 
virtual tmp< surfaceScalarFieldphi () const
 Return the volumetric flux. More...
 
virtual surfaceScalarFieldphiRef ()
 Access the volumetric flux. More...
 
virtual tmp< surfaceScalarFieldalphaPhi () const
 Return the volumetric flux of the phase. More...
 
virtual surfaceScalarFieldalphaPhiRef ()
 Access the volumetric flux of the phase. More...
 
virtual tmp< surfaceScalarFieldalphaRhoPhi () const
 Return the mass flux of the phase. More...
 
virtual surfaceScalarFieldalphaRhoPhiRef ()
 Access the mass flux of the phase. More...
 
virtual tmp< volVectorFieldDUDt () const
 Return the substantive acceleration. More...
 
virtual tmp< surfaceScalarFieldDUDtf () const
 Return the substantive acceleration on the faces. More...
 
virtual tmp< volScalarFieldcontinuityError () const
 Return the continuity error. More...
 
virtual tmp< volScalarFieldK () const
 Return the phase kinetic energy. More...
 
virtual tmp< volScalarFielddivU () const
 Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi)) More...
 
virtual void divU (tmp< volScalarField > divU)
 Set the phase dilatation rate (d(alpha)/dt + div(alpha*phi)) More...
 
virtual tmp< scalarFieldkappaEff (const label patchi) const
 Return the effective thermal conductivity on a patch. More...
 
virtual tmp< volScalarFieldk () const
 Return the turbulent kinetic energy. More...
 
virtual tmp< volScalarFieldpPrime () const
 Return the phase-pressure'. More...
 
virtual tmp< fvScalarMatrixdivq (volScalarField &he) const
 Return the source term for the energy equation. More...
 
virtual tmp< fvScalarMatrixdivj (volScalarField &Yi) const
 Return the source term for the given specie mass-fraction. More...
 

Protected Attributes

volVectorField U_
 Velocity field. More...
 
surfaceScalarField phi_
 Flux. More...
 
surfaceScalarField alphaPhi_
 Volumetric flux. More...
 
surfaceScalarField alphaRhoPhi_
 Mass flux. More...
 
tmp< volVectorFieldDUDt_
 Lagrangian acceleration field (needed for virtual-mass) More...
 
tmp< surfaceScalarFieldDUDtf_
 Lagrangian acceleration field on the faces (needed for virtual-mass) More...
 
tmp< volScalarFielddivU_
 Dilatation rate. More...
 
autoPtr< phaseCompressibleMomentumTransportModelturbulence_
 Turbulence model. More...
 
autoPtr< PhaseThermophysicalTransportModel< phaseCompressibleMomentumTransportModel, typename BasePhaseModel::thermoModel > > thermophysicalTransport_
 Thermophysical transport model. More...
 
volScalarField continuityError_
 Continuity error. More...
 
tmp< volScalarFieldK_
 Kinetic Energy. More...
 

Detailed Description

template<class BasePhaseModel>
class Foam::MovingPhaseModel< BasePhaseModel >

Class which represents a moving fluid phase. Holds the velocity, fluxes and turbulence model and can generate the momentum equation. The interface is quite restrictive as it also has to support an equivalent stationary model, which does not store motion fields or a turbulence model.

Possible future extensions include separating the turbulent functionality into another layer.

See also
StationaryPhaseModel
Source files

Definition at line 60 of file MovingPhaseModel.H.

Constructor & Destructor Documentation

◆ MovingPhaseModel()

MovingPhaseModel ( const phaseSystem fluid,
const word phaseName,
const bool  referencePhase,
const label  index 
)

◆ ~MovingPhaseModel()

virtual ~MovingPhaseModel ( )
virtual

Destructor.

Member Function Documentation

◆ correct()

virtual void correct ( )
virtual

Correct the phase properties other than the thermo and turbulence.

◆ correctContinuityError()

virtual void correctContinuityError ( const volScalarField source)
virtual

Correct the continuity error.

◆ correctKinematics()

virtual void correctKinematics ( )
virtual

Correct the kinematics.

◆ correctTurbulence()

virtual void correctTurbulence ( )
virtual

Correct the turbulence.

◆ correctEnergyTransport()

virtual void correctEnergyTransport ( )
virtual

Correct the energy transport e.g. alphat.

◆ stationary()

virtual bool stationary ( ) const
virtual

Return whether the phase is stationary.

◆ UEqn()

virtual tmp<fvVectorMatrix> UEqn ( )
virtual

Return the momentum equation.

◆ UfEqn()

virtual tmp<fvVectorMatrix> UfEqn ( )
virtual

Return the momentum equation for the face-based algorithm.

◆ U()

virtual tmp<volVectorField> U ( ) const
virtual

Return the velocity.

◆ URef()

virtual volVectorField& URef ( )
virtual

Access the velocity.

◆ phi()

virtual tmp<surfaceScalarField> phi ( ) const
virtual

Return the volumetric flux.

◆ phiRef()

virtual surfaceScalarField& phiRef ( )
virtual

Access the volumetric flux.

◆ alphaPhi()

virtual tmp<surfaceScalarField> alphaPhi ( ) const
virtual

Return the volumetric flux of the phase.

◆ alphaPhiRef()

virtual surfaceScalarField& alphaPhiRef ( )
virtual

Access the volumetric flux of the phase.

◆ alphaRhoPhi()

virtual tmp<surfaceScalarField> alphaRhoPhi ( ) const
virtual

Return the mass flux of the phase.

◆ alphaRhoPhiRef()

virtual surfaceScalarField& alphaRhoPhiRef ( )
virtual

Access the mass flux of the phase.

◆ DUDt()

virtual tmp<volVectorField> DUDt ( ) const
virtual

Return the substantive acceleration.

◆ DUDtf()

virtual tmp<surfaceScalarField> DUDtf ( ) const
virtual

Return the substantive acceleration on the faces.

◆ continuityError()

virtual tmp<volScalarField> continuityError ( ) const
virtual

Return the continuity error.

◆ K()

virtual tmp<volScalarField> K ( ) const
virtual

Return the phase kinetic energy.

◆ divU() [1/2]

virtual tmp<volScalarField> divU ( ) const
virtual

Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))

◆ divU() [2/2]

virtual void divU ( tmp< volScalarField divU)
virtual

Set the phase dilatation rate (d(alpha)/dt + div(alpha*phi))

◆ kappaEff()

virtual tmp<scalarField> kappaEff ( const label  patchi) const
virtual

Return the effective thermal conductivity on a patch.

◆ k()

virtual tmp<volScalarField> k ( ) const
virtual

Return the turbulent kinetic energy.

◆ pPrime()

virtual tmp<volScalarField> pPrime ( ) const
virtual

Return the phase-pressure'.

(derivative of phase-pressure w.r.t. phase-fraction)

◆ divq()

virtual tmp<fvScalarMatrix> divq ( volScalarField he) const
virtual

Return the source term for the energy equation.

◆ divj()

virtual tmp<fvScalarMatrix> divj ( volScalarField Yi) const
virtual

Return the source term for the given specie mass-fraction.

equation

Member Data Documentation

◆ U_

volVectorField U_
protected

Velocity field.

Definition at line 69 of file MovingPhaseModel.H.

◆ phi_

surfaceScalarField phi_
protected

Flux.

Definition at line 72 of file MovingPhaseModel.H.

◆ alphaPhi_

surfaceScalarField alphaPhi_
protected

Volumetric flux.

Definition at line 75 of file MovingPhaseModel.H.

◆ alphaRhoPhi_

surfaceScalarField alphaRhoPhi_
protected

Mass flux.

Definition at line 78 of file MovingPhaseModel.H.

◆ DUDt_

tmp<volVectorField> DUDt_
mutableprotected

Lagrangian acceleration field (needed for virtual-mass)

Definition at line 81 of file MovingPhaseModel.H.

◆ DUDtf_

tmp<surfaceScalarField> DUDtf_
mutableprotected

Lagrangian acceleration field on the faces (needed for virtual-mass)

Definition at line 84 of file MovingPhaseModel.H.

◆ divU_

tmp<volScalarField> divU_
protected

Dilatation rate.

Definition at line 87 of file MovingPhaseModel.H.

◆ turbulence_

Turbulence model.

Definition at line 90 of file MovingPhaseModel.H.

◆ thermophysicalTransport_

autoPtr< PhaseThermophysicalTransportModel < phaseCompressibleMomentumTransportModel, typename BasePhaseModel::thermoModel > > thermophysicalTransport_
protected

Thermophysical transport model.

Definition at line 100 of file MovingPhaseModel.H.

◆ continuityError_

volScalarField continuityError_
protected

Continuity error.

Definition at line 103 of file MovingPhaseModel.H.

◆ K_

tmp<volScalarField> K_
mutableprotected

Kinetic Energy.

Definition at line 106 of file MovingPhaseModel.H.


The documentation for this class was generated from the following file: