StationaryPhaseModel< BasePhaseModel > Class Template Reference

Class which represents a stationary (and therefore probably solid) phase. Generates, but does not store, zero velocity and flux field and turbulent qauantities. Throws an error when non-const access is requested to the motion fields or when the momentum equation is requested. Usage must must protect against such calls. More...

Inheritance diagram for StationaryPhaseModel< BasePhaseModel >:
[legend]
Collaboration diagram for StationaryPhaseModel< BasePhaseModel >:
[legend]

Public Member Functions

 StationaryPhaseModel (const phaseSystem &fluid, const word &phaseName, const label index)
 
virtual ~StationaryPhaseModel ()
 Destructor. 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< volScalarFieldcontinuityErrorFlow () const
 Return the continuity error due to the flow field. More...
 
virtual tmp< volScalarFieldcontinuityErrorSources () const
 Return the continuity error due to any sources. 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< volScalarFieldmut () const
 Return the turbulent dynamic viscosity. More...
 
virtual tmp< volScalarFieldmuEff () const
 Return the effective dynamic viscosity. More...
 
virtual tmp< volScalarFieldnut () const
 Return the turbulent kinematic viscosity. More...
 
virtual tmp< volScalarFieldnuEff () const
 Return the effective kinematic viscosity. More...
 
virtual tmp< volScalarFieldkappaEff () const
 Return the effective thermal conductivity. More...
 
virtual tmp< scalarFieldkappaEff (const label patchi) const
 Return the effective thermal conductivity on a patch. More...
 
virtual tmp< volScalarFieldalphaEff () const
 Return the effective thermal diffusivity. More...
 
virtual tmp< scalarFieldalphaEff (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...
 
template<class Type , template< class > class PatchField, class GeoMesh >
Foam::tmp< Foam::GeometricField< Type, PatchField, GeoMesh > > zeroField (const word &name, const dimensionSet &dims, const bool cache) const
 
template<class Type >
Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > zeroVolField (const word &name, const dimensionSet &dims, const bool cache) const
 
template<class Type >
Foam::tmp< Foam::GeometricField< Type, Foam::fvsPatchField, Foam::surfaceMesh > > zeroSurfaceField (const word &name, const dimensionSet &dims, const bool cache) const
 

Detailed Description

template<class BasePhaseModel>
class Foam::StationaryPhaseModel< BasePhaseModel >

Class which represents a stationary (and therefore probably solid) phase. Generates, but does not store, zero velocity and flux field and turbulent qauantities. Throws an error when non-const access is requested to the motion fields or when the momentum equation is requested. Usage must must protect against such calls.

See also
MovingPhaseModel
Source files

Definition at line 58 of file StationaryPhaseModel.H.

Constructor & Destructor Documentation

◆ StationaryPhaseModel()

StationaryPhaseModel ( const phaseSystem fluid,
const word phaseName,
const label  index 
)

Definition at line 90 of file StationaryPhaseModel.C.

◆ ~StationaryPhaseModel()

Destructor.

Definition at line 104 of file StationaryPhaseModel.C.

Member Function Documentation

◆ stationary()

bool stationary
virtual

Return whether the phase is stationary.

Definition at line 111 of file StationaryPhaseModel.C.

◆ UEqn()

Return the momentum equation.

Definition at line 119 of file StationaryPhaseModel.C.

References Foam::exit(), Foam::FatalError, and FatalErrorInFunction.

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 131 of file StationaryPhaseModel.C.

References Foam::exit(), Foam::FatalError, and FatalErrorInFunction.

Here is the call graph for this function:

◆ U()

Return the velocity.

Definition at line 143 of file StationaryPhaseModel.C.

References Foam::dimVelocity.

◆ URef()

Foam::volVectorField & URef
virtual

Access the velocity.

Definition at line 151 of file StationaryPhaseModel.C.

References Foam::exit(), Foam::FatalError, FatalErrorInFunction, and GeometricField< vector, fvPatchField, volMesh >::null().

Here is the call graph for this function:

◆ phi()

Return the volumetric flux.

Definition at line 163 of file StationaryPhaseModel.C.

References Foam::dimTime, and Foam::dimVolume.

◆ phiRef()

Foam::surfaceScalarField & phiRef
virtual

Access the volumetric flux.

Definition at line 171 of file StationaryPhaseModel.C.

References Foam::exit(), Foam::FatalError, FatalErrorInFunction, and GeometricField< scalar, fvsPatchField, surfaceMesh >::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 183 of file StationaryPhaseModel.C.

References Foam::dimTime, and Foam::dimVolume.

◆ alphaPhiRef()

Foam::surfaceScalarField & alphaPhiRef
virtual

Access the volumetric flux of the phase.

Definition at line 191 of file StationaryPhaseModel.C.

References Foam::exit(), Foam::FatalError, FatalErrorInFunction, and GeometricField< scalar, fvsPatchField, surfaceMesh >::null().

Here is the call graph for this function:

◆ alphaRhoPhi()

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

Return the mass flux of the phase.

Definition at line 203 of file StationaryPhaseModel.C.

References Foam::dimMass, and Foam::dimTime.

◆ alphaRhoPhiRef()

Foam::surfaceScalarField & alphaRhoPhiRef
virtual

Access the mass flux of the phase.

Definition at line 211 of file StationaryPhaseModel.C.

References Foam::exit(), Foam::FatalError, FatalErrorInFunction, and GeometricField< scalar, fvsPatchField, surfaceMesh >::null().

Here is the call graph for this function:

◆ DUDt()

Return the substantive acceleration.

Definition at line 223 of file StationaryPhaseModel.C.

References Foam::dimTime, and Foam::dimVelocity.

◆ DUDtf()

Return the substantive acceleration on the faces.

Definition at line 231 of file StationaryPhaseModel.C.

References Foam::dimArea, Foam::dimTime, and Foam::dimVelocity.

◆ continuityError()

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

Return the continuity error.

Definition at line 239 of file StationaryPhaseModel.C.

References Foam::dimDensity, and Foam::dimTime.

◆ continuityErrorFlow()

Foam::tmp< Foam::volScalarField > continuityErrorFlow
virtual

Return the continuity error due to the flow field.

Definition at line 247 of file StationaryPhaseModel.C.

References Foam::dimDensity, and Foam::dimTime.

◆ continuityErrorSources()

Foam::tmp< Foam::volScalarField > continuityErrorSources
virtual

Return the continuity error due to any sources.

Definition at line 255 of file StationaryPhaseModel.C.

References Foam::dimDensity, and Foam::dimTime.

◆ K()

Return the phase kinetic energy.

Definition at line 263 of file StationaryPhaseModel.C.

References Foam::dimVelocity, and Foam::sqr().

Here is the call graph for this function:

◆ divU() [1/2]

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

Definition at line 271 of file StationaryPhaseModel.C.

◆ divU() [2/2]

void divU ( tmp< volScalarField divU)
virtual

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

Definition at line 278 of file StationaryPhaseModel.C.

References Foam::exit(), Foam::FatalError, and FatalErrorInFunction.

Here is the call graph for this function:

◆ mut()

Return the turbulent dynamic viscosity.

Definition at line 291 of file StationaryPhaseModel.C.

References Foam::dimDynamicViscosity.

◆ muEff()

Foam::tmp< Foam::volScalarField > muEff
virtual

Return the effective dynamic viscosity.

Definition at line 299 of file StationaryPhaseModel.C.

◆ nut()

Return the turbulent kinematic viscosity.

Definition at line 307 of file StationaryPhaseModel.C.

References Foam::dimViscosity.

◆ nuEff()

Foam::tmp< Foam::volScalarField > nuEff
virtual

Return the effective kinematic viscosity.

Definition at line 315 of file StationaryPhaseModel.C.

◆ kappaEff() [1/2]

Foam::tmp< Foam::volScalarField > kappaEff
virtual

Return the effective thermal conductivity.

Definition at line 323 of file StationaryPhaseModel.C.

◆ kappaEff() [2/2]

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

Return the effective thermal conductivity on a patch.

Definition at line 331 of file StationaryPhaseModel.C.

◆ alphaEff() [1/2]

Foam::tmp< Foam::volScalarField > alphaEff
virtual

Return the effective thermal diffusivity.

Definition at line 339 of file StationaryPhaseModel.C.

◆ alphaEff() [2/2]

Foam::tmp< Foam::scalarField > alphaEff ( const label  patchi) const
virtual

Return the effective thermal conductivity on a patch.

Definition at line 347 of file StationaryPhaseModel.C.

◆ k()

Return the turbulent kinetic energy.

Definition at line 355 of file StationaryPhaseModel.C.

References Foam::dimVelocity, and Foam::sqr().

Here is the call graph for this function:

◆ pPrime()

Foam::tmp< Foam::volScalarField > pPrime
virtual

Return the phase-pressure'.

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

Definition at line 363 of file StationaryPhaseModel.C.

References Foam::dimPressure.

◆ zeroField()

Foam::tmp< Foam::GeometricField< Type, PatchField, GeoMesh > > zeroField ( const word name,
const dimensionSet dims,
const bool  cache 
) const

Definition at line 35 of file StationaryPhaseModel.C.

References mesh, Foam::name(), and timeName.

Here is the call graph for this function:

◆ zeroVolField()

Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > zeroVolField ( const word name,
const dimensionSet dims,
const bool  cache 
) const

Definition at line 62 of file StationaryPhaseModel.C.

References Foam::name().

Here is the call graph for this function:

◆ zeroSurfaceField()

Foam::tmp< Foam::GeometricField< Type, Foam::fvsPatchField, Foam::surfaceMesh > > zeroSurfaceField ( const word name,
const dimensionSet dims,
const bool  cache 
) const

Definition at line 76 of file StationaryPhaseModel.C.

References Foam::name().

Here is the call graph for this function:

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