MovingPhaseModel< BasePhaseModel > Class Template Reference

Class which represents a moving fluid phase. Holds the velocity, fluxes and turbulence model. Provides access to the turbulent quantities. More...

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

Public Member Functions

 MovingPhaseModel (const multiphaseInterSystem &fluid, const word &phaseName)
 
virtual ~MovingPhaseModel ()=default
 Destructor. More...
 
virtual void correct ()
 Correct the phase properties other than the thermo and turbulence. More...
 
virtual tmp< surfaceScalarFieldphi () const
 Constant access the volumetric flux. More...
 
virtual const surfaceScalarFieldphi ()
 Access the volumetric flux. More...
 
virtual tmp< surfaceScalarFieldalphaPhi () const
 Constant access the volumetric flux of the phase. More...
 
virtual surfaceScalarFieldalphaPhi ()
 Access the volumetric flux of the phase. More...
 
virtual tmp< volVectorFieldU () const
 Access const reference to U. More...
 
virtual tmp< surfaceScalarFielddiffNo () const
 Diffusion number. More...
 
 MovingPhaseModel (const phaseSystem &fluid, const word &phaseName, const label index)
 Construct from phase system and phase name. More...
 
virtual ~MovingPhaseModel ()=default
 Destructor. More...
 
virtual void correct ()
 Correct the phase properties other than the thermo and turbulence. 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< 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...
 

Protected Attributes

volVectorField U_
 Velocity field. More...
 
surfaceScalarField phi_
 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< phaseCompressibleTurbulenceModelturbulence_
 Turbulence model. More...
 
volScalarField continuityErrorFlow_
 Continuity error due to the flow. More...
 
volScalarField continuityErrorSources_
 Continuity error due to any sources. 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. Provides access to the turbulent quantities.

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 fuctionality into another layer.

Source files

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

See also
StationaryPhaseModel
Source files

Definition at line 55 of file MovingPhaseModel.H.

Constructor & Destructor Documentation

◆ MovingPhaseModel() [1/2]

MovingPhaseModel ( const multiphaseInterSystem fluid,
const word phaseName 
)

Definition at line 47 of file MovingPhaseModel.C.

References Foam::name().

Here is the call graph for this function:

◆ ~MovingPhaseModel() [1/2]

virtual ~MovingPhaseModel ( )
virtualdefault

Destructor.

◆ MovingPhaseModel() [2/2]

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

Construct from phase system and phase name.

Definition at line 127 of file MovingPhaseModel.C.

References IOobject::AUTO_WRITE, MovingPhaseModel< BasePhaseModel >::correctKinematics(), Foam::dimDensity, Foam::dimTime, Foam::name(), Foam::New(), and IOobject::writeOpt().

Here is the call graph for this function:

◆ ~MovingPhaseModel() [2/2]

virtual ~MovingPhaseModel ( )
virtualdefault

Destructor.

Member Function Documentation

◆ correct() [1/2]

void correct
virtual

Correct the phase properties other than the thermo and turbulence.

Definition at line 77 of file MovingPhaseModel.C.

◆ phi() [1/3]

Constant access the volumetric flux.

Definition at line 85 of file MovingPhaseModel.C.

◆ phi() [2/3]

const Foam::surfaceScalarField & phi
virtual

Access the volumetric flux.

Definition at line 93 of file MovingPhaseModel.C.

◆ alphaPhi() [1/3]

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

Constant access the volumetric flux of the phase.

Definition at line 101 of file MovingPhaseModel.C.

◆ alphaPhi() [2/3]

Foam::surfaceScalarField & alphaPhi
virtual

Access the volumetric flux of the phase.

Definition at line 109 of file MovingPhaseModel.C.

◆ U() [1/2]

Access const reference to U.

Definition at line 117 of file MovingPhaseModel.C.

◆ diffNo()

Diffusion number.

Definition at line 124 of file MovingPhaseModel.C.

References Foam::dimless, IOobject::groupName(), phaseModel::name(), Time::New(), and Foam::Zero.

Here is the call graph for this function:

◆ correct() [2/2]

virtual void correct ( )
virtual

Correct the phase properties other than the thermo and turbulence.

◆ correctKinematics()

void correctKinematics
virtual

Correct the kinematics.

Definition at line 233 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:

◆ correctTurbulence()

void correctTurbulence
virtual

Correct the turbulence.

Definition at line 257 of file MovingPhaseModel.C.

◆ correctEnergyTransport()

void correctEnergyTransport
virtual

Correct the energy transport e.g. alphat.

Definition at line 266 of file MovingPhaseModel.C.

◆ stationary()

bool stationary
virtual

Return whether the phase is stationary.

Definition at line 275 of file MovingPhaseModel.C.

◆ UEqn()

Return the momentum equation.

Definition at line 283 of file MovingPhaseModel.C.

References alpha, Foam::fvm::ddt(), Foam::fvm::div(), fluid, 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 301 of file MovingPhaseModel.C.

References alpha, Foam::fvc::div(), Foam::fvm::div(), fluid, MRF, rho, Foam::fvm::Sp(), and Foam::fvm::SuSp().

Here is the call graph for this function:

◆ U() [2/2]

virtual tmp< volVectorField > U ( ) const
virtual

Return the velocity.

◆ URef()

Foam::volVectorField & URef
virtual

Access the velocity.

Definition at line 329 of file MovingPhaseModel.C.

◆ phi() [3/3]

virtual tmp< surfaceScalarField > phi ( ) const
virtual

Return the volumetric flux.

◆ phiRef()

Foam::surfaceScalarField & phiRef
virtual

Access the volumetric flux.

Definition at line 345 of file MovingPhaseModel.C.

◆ alphaPhi() [3/3]

virtual tmp< surfaceScalarField > alphaPhi ( ) const
virtual

Return the volumetric flux of the phase.

◆ alphaPhiRef()

Foam::surfaceScalarField & alphaPhiRef
virtual

Access the volumetric flux of the phase.

Definition at line 361 of file MovingPhaseModel.C.

◆ alphaRhoPhi()

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

Return the mass flux of the phase.

Definition at line 369 of file MovingPhaseModel.C.

◆ alphaRhoPhiRef()

Foam::surfaceScalarField & alphaRhoPhiRef
virtual

Access the mass flux of the phase.

Definition at line 377 of file MovingPhaseModel.C.

◆ DUDt()

Return the substantive acceleration.

Definition at line 385 of file MovingPhaseModel.C.

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

Here is the call graph for this function:

◆ DUDtf()

Return the substantive acceleration on the faces.

Definition at line 398 of file MovingPhaseModel.C.

References Foam::byDt().

Here is the call graph for this function:

◆ continuityError()

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

Return the continuity error.

Definition at line 411 of file MovingPhaseModel.C.

◆ continuityErrorFlow()

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

Return the continuity error due to the flow field.

Definition at line 419 of file MovingPhaseModel.C.

◆ continuityErrorSources()

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

Return the continuity error due to any sources.

Definition at line 427 of file MovingPhaseModel.C.

◆ K()

Return the phase kinetic energy.

Definition at line 435 of file MovingPhaseModel.C.

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

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 452 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 459 of file MovingPhaseModel.C.

References divU.

◆ mut()

Return the turbulent dynamic viscosity.

Definition at line 467 of file MovingPhaseModel.C.

◆ muEff()

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

Return the effective dynamic viscosity.

Definition at line 475 of file MovingPhaseModel.C.

◆ nut()

Return the turbulent kinematic viscosity.

Definition at line 483 of file MovingPhaseModel.C.

◆ nuEff()

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

Return the effective kinematic viscosity.

Definition at line 491 of file MovingPhaseModel.C.

◆ kappaEff() [1/2]

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

Return the effective thermal conductivity.

Definition at line 499 of file MovingPhaseModel.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 507 of file MovingPhaseModel.C.

◆ alphaEff() [1/2]

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

Return the effective thermal diffusivity.

Definition at line 515 of file MovingPhaseModel.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 523 of file MovingPhaseModel.C.

◆ k()

Return the turbulent kinetic energy.

Definition at line 531 of file MovingPhaseModel.C.

◆ pPrime()

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

Return the phase-pressure'.

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

Definition at line 539 of file MovingPhaseModel.C.

Member Data Documentation

◆ U_

volVectorField U_
protected

Velocity field.

Definition at line 71 of file MovingPhaseModel.H.

◆ phi_

surfaceScalarField phi_
protected

Flux.

Definition at line 74 of file MovingPhaseModel.H.

◆ alphaRhoPhi_

surfaceScalarField alphaRhoPhi_
protected

Mass flux.

Definition at line 80 of file MovingPhaseModel.H.

◆ DUDt_

tmp<volVectorField> DUDt_
mutableprotected

Lagrangian acceleration field (needed for virtual-mass)

Definition at line 83 of file MovingPhaseModel.H.

◆ DUDtf_

tmp<surfaceScalarField> DUDtf_
mutableprotected

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

Definition at line 86 of file MovingPhaseModel.H.

◆ divU_

tmp<volScalarField> divU_
protected

Dilatation rate.

Definition at line 89 of file MovingPhaseModel.H.

◆ turbulence_

autoPtr<phaseCompressibleTurbulenceModel> turbulence_
protected

Turbulence model.

Definition at line 92 of file MovingPhaseModel.H.

◆ continuityErrorFlow_

volScalarField continuityErrorFlow_
protected

Continuity error due to the flow.

Definition at line 95 of file MovingPhaseModel.H.

◆ continuityErrorSources_

volScalarField continuityErrorSources_
protected

Continuity error due to any sources.

Definition at line 98 of file MovingPhaseModel.H.

◆ K_

tmp<volScalarField> K_
mutableprotected

Kinetic Energy.

Definition at line 101 of file MovingPhaseModel.H.


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