MovingPhaseModel< BasePhaseModel > Class Template Reference

Class which represents a moving fluid phase. Holds the velocity, fluxes and momentumTransport 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 momentumTransport 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. More...
 
virtual void correctContinuityError (const volScalarField &source)
 Correct the continuity error. More...
 
virtual void correctKinematics ()
 Correct the kinematics. More...
 
virtual void predictMomentumTransport ()
 Predict the momentumTransport. More...
 
virtual void correctMomentumTransport ()
 Correct the momentumTransport. More...
 
virtual void correctUf ()
 Correct the face velocity for moving meshes. 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 const volVectorFieldURef () const
 Access the velocity. More...
 
virtual tmp< surfaceScalarFieldphi () const
 Return the volumetric flux. More...
 
virtual surfaceScalarFieldphiRef ()
 Access the volumetric flux. More...
 
virtual const surfaceScalarFieldphiRef () const
 Access the volumetric flux. More...
 
virtual const autoPtr< surfaceVectorField > & Uf () const
 Return the face velocity. More...
 
virtual surfaceVectorFieldUfRef ()
 Access the face velocity. More...
 
virtual const surfaceVectorFieldUfRef () const
 Access the face velocity. More...
 
virtual tmp< surfaceScalarFieldalphaPhi () const
 Return the volumetric flux of the phase. More...
 
virtual surfaceScalarFieldalphaPhiRef ()
 Access the volumetric flux of the phase. More...
 
virtual const surfaceScalarFieldalphaPhiRef () const
 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 const surfaceScalarFieldalphaRhoPhiRef () const
 Access the mass flux of the phase. More...
 
virtual tmp< fvVectorMatrixUgradU () const
 Return the velocity transport matrix. More...
 
virtual tmp< fvVectorMatrixDUDt () const
 Return the substantive acceleration matrix. More...
 
virtual tmp< volScalarFieldcontinuityError () const
 Return the continuity error. More...
 
virtual tmp< volScalarFieldK () const
 Return the phase kinetic energy. More...
 
virtual const autoPtr< volScalarField > & divU () 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< volScalarFieldk () const
 Return the turbulent kinetic energy. More...
 
virtual tmp< surfaceScalarFieldpPrimef () const
 Return the face-phase-pressure'. More...
 

Protected Attributes

volVectorField U_
 Velocity field. More...
 
surfaceScalarField phi_
 Flux. More...
 
surfaceScalarField alphaPhi_
 Volumetric flux. More...
 
surfaceScalarField alphaRhoPhi_
 Mass flux. More...
 
autoPtr< surfaceVectorFieldUf_
 Face velocity field. More...
 
autoPtr< volScalarFielddivU_
 Dilatation rate. More...
 
autoPtr< phaseCompressible::momentumTransportModelmomentumTransport_
 Turbulence 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 momentumTransport 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 momentumTransport 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()

~MovingPhaseModel
virtual

Destructor.

Definition at line 216 of file MovingPhaseModel.C.

Member Function Documentation

◆ correct()

void correct
virtual

Correct the phase properties other than the thermo.

and momentumTransport

Definition at line 234 of file MovingPhaseModel.C.

References Foam::MULES::correct().

Here is the call graph for this function:

◆ correctContinuityError()

void correctContinuityError ( const volScalarField source)
virtual

Correct the continuity error.

Definition at line 223 of file MovingPhaseModel.C.

References Foam::fvc::ddt(), Foam::fvc::div(), and rho.

Here is the call graph for this function:

◆ correctKinematics()

void correctKinematics
virtual

Correct the kinematics.

Definition at line 241 of file MovingPhaseModel.C.

References Foam::magSqr(), and U.

Referenced by MovingPhaseModel< BasePhaseModel >::MovingPhaseModel().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ predictMomentumTransport()

void predictMomentumTransport
virtual

Predict the momentumTransport.

Definition at line 253 of file MovingPhaseModel.C.

◆ correctMomentumTransport()

void correctMomentumTransport
virtual

Correct the momentumTransport.

Definition at line 261 of file MovingPhaseModel.C.

◆ correctUf()

void correctUf
virtual

Correct the face velocity for moving meshes.

Definition at line 269 of file MovingPhaseModel.C.

References Foam::fvc::absolute(), Foam::fvc::interpolate(), fvMesh::magSf(), fvMesh::mesh(), n, and fvMesh::Sf().

Here is the call graph for this function:

◆ stationary()

bool stationary
virtual

Return whether the phase is stationary.

Definition at line 288 of file MovingPhaseModel.C.

◆ UEqn()

Return the momentum equation.

Definition at line 296 of file MovingPhaseModel.C.

References alpha(), Foam::fvc::DDt(), Foam::fvm::ddt(), Foam::fvm::div(), MRF(), rho, and Foam::fvm::SuSp().

Here is the call graph for this function:

◆ UfEqn()

Foam::tmp< Foam::fvVectorMatrix > UfEqn
virtual

Return the momentum equation for the face-based algorithm.

Definition at line 314 of file MovingPhaseModel.C.

References alpha(), Foam::fvc::ddt(), Foam::fvc::DDt(), Foam::fvm::div(), MRF(), rho, and Foam::fvm::SuSp().

Here is the call graph for this function:

◆ U()

Return the velocity.

Definition at line 333 of file MovingPhaseModel.C.

◆ URef() [1/2]

Foam::volVectorField & URef
virtual

Access the velocity.

Definition at line 341 of file MovingPhaseModel.C.

◆ URef() [2/2]

const Foam::volVectorField & URef
virtual

Access the velocity.

Definition at line 349 of file MovingPhaseModel.C.

◆ phi()

Return the volumetric flux.

Definition at line 357 of file MovingPhaseModel.C.

◆ phiRef() [1/2]

Foam::surfaceScalarField & phiRef
virtual

Access the volumetric flux.

Definition at line 365 of file MovingPhaseModel.C.

◆ phiRef() [2/2]

const Foam::surfaceScalarField & phiRef
virtual

Access the volumetric flux.

Definition at line 373 of file MovingPhaseModel.C.

◆ Uf()

const Foam::autoPtr< Foam::surfaceVectorField > & Uf
virtual

Return the face velocity.

Required for moving mesh cases

Definition at line 381 of file MovingPhaseModel.C.

◆ UfRef() [1/2]

Foam::surfaceVectorField & UfRef
virtual

Access the face velocity.

Required for moving mesh cases

Definition at line 389 of file MovingPhaseModel.C.

References Foam::exit(), Foam::FatalError, FatalErrorInFunction, and GeometricField< Type, PatchField, GeoMesh >::null().

Here is the call graph for this function:

◆ UfRef() [2/2]

const Foam::surfaceVectorField & UfRef
virtual

Access the face velocity.

Required for moving mesh cases

Definition at line 408 of file MovingPhaseModel.C.

References Foam::exit(), Foam::FatalError, FatalErrorInFunction, and GeometricField< Type, PatchField, GeoMesh >::null().

Here is the call graph for this function:

◆ alphaPhi()

Foam::tmp< Foam::surfaceScalarField > alphaPhi
virtual

Return the volumetric flux of the phase.

Definition at line 427 of file MovingPhaseModel.C.

◆ alphaPhiRef() [1/2]

Foam::surfaceScalarField & alphaPhiRef
virtual

Access the volumetric flux of the phase.

Definition at line 435 of file MovingPhaseModel.C.

◆ alphaPhiRef() [2/2]

const Foam::surfaceScalarField & alphaPhiRef
virtual

Access the volumetric flux of the phase.

Definition at line 443 of file MovingPhaseModel.C.

◆ alphaRhoPhi()

Foam::tmp< Foam::surfaceScalarField > alphaRhoPhi
virtual

Return the mass flux of the phase.

Definition at line 451 of file MovingPhaseModel.C.

◆ alphaRhoPhiRef() [1/2]

Foam::surfaceScalarField & alphaRhoPhiRef
virtual

Access the mass flux of the phase.

Definition at line 459 of file MovingPhaseModel.C.

◆ alphaRhoPhiRef() [2/2]

const Foam::surfaceScalarField & alphaRhoPhiRef
virtual

Access the mass flux of the phase.

Definition at line 467 of file MovingPhaseModel.C.

◆ UgradU()

Foam::tmp< Foam::fvVectorMatrix > UgradU
virtual

Return the velocity transport matrix.

Definition at line 475 of file MovingPhaseModel.C.

References Foam::fvc::absolute(), Foam::fvc::div(), Foam::fvm::div(), and Foam::fvm::Sp().

Here is the call graph for this function:

◆ DUDt()

Return the substantive acceleration matrix.

Definition at line 488 of file MovingPhaseModel.C.

References Foam::fvm::ddt().

Here is the call graph for this function:

◆ continuityError()

Foam::tmp< Foam::volScalarField > continuityError
virtual

Return the continuity error.

Definition at line 496 of file MovingPhaseModel.C.

◆ K()

Return the phase kinetic energy.

Definition at line 504 of file MovingPhaseModel.C.

References IOobject::groupName(), Foam::magSqr(), Foam::name(), and U.

Here is the call graph for this function:

◆ divU() [1/2]

const Foam::autoPtr< Foam::volScalarField > & divU
virtual

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

Definition at line 522 of file MovingPhaseModel.C.

◆ divU() [2/2]

void divU ( tmp< volScalarField divU)
virtual

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

Definition at line 529 of file MovingPhaseModel.C.

References IOobject::groupName(), Foam::name(), and tmp< T >::ptr().

Here is the call graph for this function:

◆ k()

Return the turbulent kinetic energy.

Definition at line 546 of file MovingPhaseModel.C.

◆ pPrimef()

Foam::tmp< Foam::surfaceScalarField > pPrimef
virtual

Return the face-phase-pressure'.

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

Definition at line 554 of file MovingPhaseModel.C.

Member Data Documentation

◆ U_

volVectorField U_
protected

Velocity field.

Definition at line 69 of file MovingPhaseModel.H.

Referenced by MovingPhaseModel< BasePhaseModel >::MovingPhaseModel().

◆ phi_

surfaceScalarField phi_
protected

◆ 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.

◆ Uf_

autoPtr<surfaceVectorField> Uf_
protected

Face velocity field.

Definition at line 81 of file MovingPhaseModel.H.

Referenced by MovingPhaseModel< BasePhaseModel >::MovingPhaseModel().

◆ divU_

autoPtr<volScalarField> divU_
protected

Dilatation rate.

Definition at line 84 of file MovingPhaseModel.H.

◆ momentumTransport_

autoPtr<phaseCompressible::momentumTransportModel> momentumTransport_
protected

Turbulence model.

Definition at line 87 of file MovingPhaseModel.H.

◆ continuityError_

volScalarField continuityError_
protected

Continuity error.

Definition at line 90 of file MovingPhaseModel.H.

◆ K_

tmp<volScalarField> K_
mutableprotected

Kinetic Energy.

Definition at line 93 of file MovingPhaseModel.H.


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