fvModel Class Referenceabstract

Finite volume model abstract base class. More...

Inheritance diagram for fvModel:

Classes

class  iNew
 Return pointer to new fvModel object created. More...
 

Public Member Functions

 TypeName ("fvModel")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, fvModel, dictionary,(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict),(name, modelType, mesh, dict))
 
 fvModel (const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
 Construct from components. More...
 
autoPtr< fvModelclone () const
 Return clone. More...
 
virtual ~fvModel ()
 Destructor. More...
 
const wordname () const
 Return const access to the source name. More...
 
const fvMeshmesh () const
 Return const access to the mesh database. More...
 
const dictionarycoeffs () const
 Return dictionary. More...
 
virtual wordList addSupFields () const
 Return the list of fields for which the fvModel adds source term. More...
 
virtual bool addsSupToField (const word &fieldName) const
 Return true if the fvModel adds a source term to the given. More...
 
virtual scalar maxDeltaT () const
 Return the maximum time-step for stable operation. More...
 
virtual void addSup (fvMatrix< scalar > &eqn) const
 Add a source term to a field-less proxy equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > sourceProxy (const VolField< Type > &eqnField) const
 Add a source term to an equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const VolField< Type > &field) const
 Return source for an equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > sourceProxy (const VolField< Type > &field, const VolField< Type > &eqnField) const
 Return source for an equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const volScalarField &rho, const VolField< Type > &field) const
 Return source for a compressible equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > sourceProxy (const volScalarField &rho, const VolField< Type > &field, const VolField< Type > &eqnField) const
 Return source for a compressible equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field) const
 Return source for a phase equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > sourceProxy (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field, const VolField< Type > &eqnField) const
 Return source for a phase equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const volScalarField &alpha, const geometricOneField &rho, const VolField< Type > &field) const
 Return source for a phase equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const geometricOneField &alpha, const volScalarField &rho, const VolField< Type > &field) const
 Return source for a phase equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const geometricOneField &alpha, const geometricOneField &rho, const VolField< Type > &field) const
 Return source for a phase equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > d2dt2 (const VolField< Type > &field) const
 Return source for an equation with a second time derivative. More...
 
virtual void preUpdateMesh ()
 Prepare for mesh update. More...
 
virtual bool movePoints ()=0
 Update for mesh motion. More...
 
virtual void topoChange (const polyTopoChangeMap &)=0
 Update topology using the given map. More...
 
virtual void mapMesh (const polyMeshMap &)=0
 Update from another mesh using the given map. More...
 
virtual void distribute (const polyDistributionMap &)=0
 Redistribute or update using the given distribution map. More...
 
virtual void correct ()
 Correct the fvModel. More...
 
virtual bool read (const dictionary &dict)
 Read source dictionary. More...
 
virtual bool write (const bool write=true) const
 Write fvModel data. More...
 
template<class AlphaRhoFieldType , class ... AlphaRhoFieldTypes>
Foam::dimensionSet sourceDims (const dimensionSet &ds, const AlphaRhoFieldType &alphaRhoField, const AlphaRhoFieldTypes &... alphaRhoFields)
 
template<class AlphaRhoFieldType , class ... AlphaRhoFieldTypes>
const Foam::wordfieldName (const AlphaRhoFieldType &alphaRhoField, const AlphaRhoFieldTypes &... alphaRhoFields)
 
template<class AlphaRhoFieldType >
const Foam::wordfieldName (const AlphaRhoFieldType &alphaRhoField)
 
template<class Type , class ... AlphaRhoFieldTypes>
Foam::tmp< Foam::fvMatrix< Type > > sourceTerm (const VolField< Type > &eqnField, const dimensionSet &ds, const AlphaRhoFieldTypes &... alphaRhoFields) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > sourceProxy (const VolField< Type > &eqnField) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > sourceProxy (const VolField< Type > &field, const VolField< Type > &eqnField) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const volScalarField &rho, const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > sourceProxy (const volScalarField &rho, const VolField< Type > &field, const VolField< Type > &eqnField) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > sourceProxy (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field, const VolField< Type > &eqnField) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const geometricOneField &alpha, const geometricOneField &rho, const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const volScalarField &alpha, const geometricOneField &rho, const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const geometricOneField &alpha, const volScalarField &rho, const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > d2dt2 (const VolField< Type > &field) const
 

Static Public Member Functions

template<class AlphaRhoFieldType , class ... AlphaRhoFieldTypes>
static dimensionSet sourceDims (const dimensionSet &ds, const AlphaRhoFieldType &alphaRhoField, const AlphaRhoFieldTypes &... alphaRhoFields)
 Return the dimensions of the matrix of a source term. More...
 
static const dimensionSetsourceDims (const dimensionSet &ds)
 Return the dimensions of the matrix of a source term (base. More...
 
template<class AlphaRhoFieldType , class ... AlphaRhoFieldTypes>
static const wordfieldName (const AlphaRhoFieldType &alphaRhoField, const AlphaRhoFieldTypes &... alphaRhoFields)
 Return the name of the field associated with a source term. More...
 
template<class AlphaRhoFieldType >
static const wordfieldName (const AlphaRhoFieldType &alphaRhoField)
 Return the name of the field associated with a source term (base. More...
 
static const wordfieldName ()
 Return the name of the field associated with a source term. Special. More...
 
static autoPtr< fvModelNew (const word &name, const fvMesh &mesh, const dictionary &dict)
 Return a reference to the selected fvModel. More...
 

Protected Member Functions

template<class Type >
void addSupType (const VolField< Type > &field, fvMatrix< Type > &eqn) const
 Add a source term to an equation. More...
 
template<class Type >
void addSupType (const volScalarField &rho, const VolField< Type > &field, fvMatrix< Type > &eqn) const
 Add a source term to a compressible equation. More...
 
template<class Type >
void addSupType (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field, fvMatrix< Type > &eqn) const
 Add a source term to a phase equation. More...
 
template<class Type , class ... AlphaRhoFieldTypes>
tmp< fvMatrix< Type > > sourceTerm (const VolField< Type > &eqnField, const dimensionSet &ds, const AlphaRhoFieldTypes &... alphaRhoFields) const
 Return a source for an equation. More...
 

Detailed Description

Finite volume model abstract base class.

Source files

Definition at line 58 of file fvModel.H.

Constructor & Destructor Documentation

◆ fvModel()

fvModel ( const word name,
const word modelType,
const fvMesh mesh,
const dictionary dict 
)

Construct from components.

Definition at line 72 of file fvModel.C.

References Foam::decrIndent(), Foam::endl(), Foam::incrIndent(), Foam::indent(), and Foam::Info.

Here is the call graph for this function:

◆ ~fvModel()

~fvModel ( )
virtual

Destructor.

Definition at line 154 of file fvModel.C.

Member Function Documentation

◆ addSupType() [1/3]

void addSupType ( const VolField< Type > &  field,
fvMatrix< Type > &  eqn 
) const
protected

Add a source term to an equation.

Definition at line 41 of file fvModel.C.

◆ addSupType() [2/3]

void addSupType ( const volScalarField rho,
const VolField< Type > &  field,
fvMatrix< Type > &  eqn 
) const
protected

Add a source term to a compressible equation.

Definition at line 50 of file fvModel.C.

◆ addSupType() [3/3]

void addSupType ( const volScalarField alpha,
const volScalarField rho,
const VolField< Type > &  field,
fvMatrix< Type > &  eqn 
) const
protected

Add a source term to a phase equation.

Definition at line 60 of file fvModel.C.

◆ sourceTerm() [1/2]

tmp<fvMatrix<Type> > sourceTerm ( const VolField< Type > &  eqnField,
const dimensionSet ds,
const AlphaRhoFieldTypes &...  alphaRhoFields 
) const
protected

Return a source for an equation.

◆ TypeName()

TypeName ( "fvModel"  )

Runtime type information.

◆ declareRunTimeSelectionTable()

declareRunTimeSelectionTable ( autoPtr  ,
fvModel  ,
dictionary  ,
(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict ,
(name, modelType, mesh, dict  
)

◆ sourceDims() [1/3]

static dimensionSet sourceDims ( const dimensionSet ds,
const AlphaRhoFieldType &  alphaRhoField,
const AlphaRhoFieldTypes &...  alphaRhoFields 
)
static

Return the dimensions of the matrix of a source term.

Referenced by fvModel::sourceDims(), and fvModels::sourceTerm().

Here is the caller graph for this function:

◆ sourceDims() [2/3]

const Foam::dimensionSet & sourceDims ( const dimensionSet ds)
inlinestatic

Return the dimensions of the matrix of a source term (base.

condition for the above)

Definition at line 30 of file fvModelI.H.

◆ fieldName() [1/5]

static const word& fieldName ( const AlphaRhoFieldType &  alphaRhoField,
const AlphaRhoFieldTypes &...  alphaRhoFields 
)
static

Return the name of the field associated with a source term.

◆ fieldName() [2/5]

static const word& fieldName ( const AlphaRhoFieldType &  alphaRhoField)
static

Return the name of the field associated with a source term (base.

condition for the above)

◆ fieldName() [3/5]

const Foam::word & fieldName ( )
inlinestatic

Return the name of the field associated with a source term. Special.

overload for volume sources with no associated field.

Definition at line 39 of file fvModelI.H.

References word::null.

Referenced by fvModels::sourceTerm().

Here is the caller graph for this function:

◆ clone()

autoPtr<fvModel> clone ( ) const
inline

Return clone.

Definition at line 187 of file fvModel.H.

References NotImplemented.

◆ New()

Foam::autoPtr< Foam::fvModel > New ( const word name,
const fvMesh mesh,
const dictionary dict 
)
static

Return a reference to the selected fvModel.

Definition at line 93 of file fvModel.C.

References dict, Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, Foam::FatalIOError, FatalIOErrorInFunction, Foam::indent(), Foam::Info, Foam::libs, dictionary::lookup(), Foam::name(), Foam::nl, dlLibraryTable::open(), and string::remove().

Referenced by fvModels::fvModels(), and fvModel::iNew::operator()().

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

◆ name()

◆ mesh()

◆ coeffs()

const Foam::dictionary & coeffs ( ) const
inline

Return dictionary.

Definition at line 59 of file fvModelI.H.

Referenced by interRegionPorosityForce::interRegionPorosityForce(), propellerDisk::propellerDisk(), and forcing::readLambda().

Here is the caller graph for this function:

◆ addSupFields()

◆ addsSupToField()

bool addsSupToField ( const word fieldName) const
virtual

Return true if the fvModel adds a source term to the given.

field's transport equation

Reimplemented in fvSource, volumeBlockage, massTransfer, zeroDimensionalFixedPressureModel, and fvTotalSource.

Definition at line 166 of file fvModel.C.

References Foam::findIndex().

Referenced by fvModels::sourceTerm().

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

◆ maxDeltaT()

Foam::scalar maxDeltaT ( ) const
virtual

Return the maximum time-step for stable operation.

Definition at line 172 of file fvModel.C.

◆ addSup()

void addSup ( fvMatrix< scalar > &  eqn) const
virtual

Add a source term to a field-less proxy equation.

Reimplemented in volumeSource, massTransfer, and massSourceBase.

Definition at line 178 of file fvModel.C.

Referenced by fvModels::sourceTerm().

Here is the caller graph for this function:

◆ sourceProxy() [1/8]

tmp<fvMatrix<Type> > sourceProxy ( const VolField< Type > &  eqnField) const

Add a source term to an equation.

Add a source term to a compressible equation Add a source term to a phase equation Return source for an equation

◆ source() [1/12]

tmp<fvMatrix<Type> > source ( const VolField< Type > &  field) const

Return source for an equation.

◆ sourceProxy() [2/8]

tmp<fvMatrix<Type> > sourceProxy ( const VolField< Type > &  field,
const VolField< Type > &  eqnField 
) const

Return source for an equation.

◆ source() [2/12]

tmp<fvMatrix<Type> > source ( const volScalarField rho,
const VolField< Type > &  field 
) const

Return source for a compressible equation.

◆ sourceProxy() [3/8]

tmp<fvMatrix<Type> > sourceProxy ( const volScalarField rho,
const VolField< Type > &  field,
const VolField< Type > &  eqnField 
) const

Return source for a compressible equation.

◆ source() [3/12]

tmp<fvMatrix<Type> > source ( const volScalarField alpha,
const volScalarField rho,
const VolField< Type > &  field 
) const

Return source for a phase equation.

◆ sourceProxy() [4/8]

tmp<fvMatrix<Type> > sourceProxy ( const volScalarField alpha,
const volScalarField rho,
const VolField< Type > &  field,
const VolField< Type > &  eqnField 
) const

Return source for a phase equation.

◆ source() [4/12]

tmp<fvMatrix<Type> > source ( const volScalarField alpha,
const geometricOneField rho,
const VolField< Type > &  field 
) const

Return source for a phase equation.

◆ source() [5/12]

tmp<fvMatrix<Type> > source ( const geometricOneField alpha,
const volScalarField rho,
const VolField< Type > &  field 
) const

Return source for a phase equation.

◆ source() [6/12]

tmp<fvMatrix<Type> > source ( const geometricOneField alpha,
const geometricOneField rho,
const VolField< Type > &  field 
) const

Return source for a phase equation.

◆ d2dt2() [1/2]

tmp<fvMatrix<Type> > d2dt2 ( const VolField< Type > &  field) const

Return source for an equation with a second time derivative.

◆ preUpdateMesh()

void preUpdateMesh ( )
virtual

Prepare for mesh update.

Reimplemented in clouds, and VoFClouds.

Definition at line 191 of file fvModel.C.

◆ movePoints()

◆ topoChange()

◆ mapMesh()

◆ distribute()

◆ correct()

void correct ( )
virtual

◆ read()

bool read ( const dictionary dict)
virtual

Read source dictionary.

Reimplemented in waveForcing, verticalDamping, isotropicDamping, forcing, OUForce, rotorDisk, propellerDisk, interRegionPorosityForce, interRegionModel, interRegionHeatTransfer, heatTransfer, zeroDimensionalMassSource, volumeSource, volumeBlockage, viscousHeating, solidThermalEquilibrium, solidificationMelting, sixDoFAcceleration, semiImplicitSource, radialActuationDisk, porosityForce, phaseLimitStabilisation, singleComponentPhaseChange, multicomponentPhaseChange, coefficientPhaseChange, massTransfer, coefficientMassTransfer, massSourceBase, massSource, heatSource, effectivenessHeatExchanger, codedFvModel, buoyancyForce, buoyancyEnergy, actuationDisk, acceleration, fvTotalSource, homogeneousLiquidPhaseSeparation, homogeneousCondensation, and VoFSolidificationMelting.

Definition at line 199 of file fvModel.C.

References dict, and dictionary::optionalSubDict().

Referenced by fvModels::read(), VoFSolidificationMelting::read(), fvTotalSource::read(), acceleration::read(), actuationDisk::read(), buoyancyEnergy::read(), buoyancyForce::read(), codedFvModel::read(), effectivenessHeatExchanger::read(), heatSource::read(), massTransfer::read(), phaseLimitStabilisation::read(), porosityForce::read(), semiImplicitSource::read(), sixDoFAcceleration::read(), solidificationMelting::read(), solidThermalEquilibrium::read(), viscousHeating::read(), volumeBlockage::read(), heatTransfer::read(), interRegionModel::read(), propellerDisk::read(), rotorDisk::read(), OUForce::read(), and forcing::read().

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

◆ write()

bool write ( const bool  write = true) const
virtual

Write fvModel data.

Reimplemented in OUForce.

Definition at line 207 of file fvModel.C.

◆ sourceDims() [3/3]

Foam::dimensionSet sourceDims ( const dimensionSet ds,
const AlphaRhoFieldType &  alphaRhoField,
const AlphaRhoFieldTypes &...  alphaRhoFields 
)

Definition at line 31 of file fvModelTemplates.C.

References fvModel::sourceDims().

Here is the call graph for this function:

◆ fieldName() [4/5]

const Foam::word& fieldName ( const AlphaRhoFieldType &  alphaRhoField,
const AlphaRhoFieldTypes &...  alphaRhoFields 
)

Definition at line 43 of file fvModelTemplates.C.

◆ fieldName() [5/5]

const Foam::word& fieldName ( const AlphaRhoFieldType &  alphaRhoField)

Definition at line 54 of file fvModelTemplates.C.

◆ sourceTerm() [2/2]

Foam::tmp<Foam::fvMatrix<Type> > sourceTerm ( const VolField< Type > &  eqnField,
const dimensionSet ds,
const AlphaRhoFieldTypes &...  alphaRhoFields 
) const

Definition at line 66 of file fvModelTemplates.C.

References tmp< T >::ref().

Here is the call graph for this function:

◆ sourceProxy() [5/8]

Foam::tmp<Foam::fvMatrix<Type> > sourceProxy ( const VolField< Type > &  eqnField) const

Definition at line 95 of file fvModelTemplates.C.

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

◆ source() [7/12]

Foam::tmp<Foam::fvMatrix<Type> > source ( const VolField< Type > &  field) const

Definition at line 105 of file fvModelTemplates.C.

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

◆ sourceProxy() [6/8]

Foam::tmp<Foam::fvMatrix<Type> > sourceProxy ( const VolField< Type > &  field,
const VolField< Type > &  eqnField 
) const

Definition at line 115 of file fvModelTemplates.C.

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

◆ source() [8/12]

Foam::tmp<Foam::fvMatrix<Type> > source ( const volScalarField rho,
const VolField< Type > &  field 
) const

Definition at line 126 of file fvModelTemplates.C.

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

◆ sourceProxy() [7/8]

Foam::tmp<Foam::fvMatrix<Type> > sourceProxy ( const volScalarField rho,
const VolField< Type > &  field,
const VolField< Type > &  eqnField 
) const

Definition at line 137 of file fvModelTemplates.C.

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

◆ source() [9/12]

Foam::tmp<Foam::fvMatrix<Type> > source ( const volScalarField alpha,
const volScalarField rho,
const VolField< Type > &  field 
) const

Definition at line 149 of file fvModelTemplates.C.

References alpha(), Foam::dimTime, Foam::dimVolume, and rho.

Here is the call graph for this function:

◆ sourceProxy() [8/8]

Foam::tmp<Foam::fvMatrix<Type> > sourceProxy ( const volScalarField alpha,
const volScalarField rho,
const VolField< Type > &  field,
const VolField< Type > &  eqnField 
) const

Definition at line 161 of file fvModelTemplates.C.

References alpha(), Foam::dimTime, Foam::dimVolume, and rho.

Here is the call graph for this function:

◆ source() [10/12]

Foam::tmp<Foam::fvMatrix<Type> > source ( const geometricOneField alpha,
const geometricOneField rho,
const VolField< Type > &  field 
) const

Definition at line 174 of file fvModelTemplates.C.

◆ source() [11/12]

Foam::tmp<Foam::fvMatrix<Type> > source ( const volScalarField alpha,
const geometricOneField rho,
const VolField< Type > &  field 
) const

Definition at line 186 of file fvModelTemplates.C.

References alpha().

Here is the call graph for this function:

◆ source() [12/12]

Foam::tmp<Foam::fvMatrix<Type> > source ( const geometricOneField alpha,
const volScalarField rho,
const VolField< Type > &  field 
) const

Definition at line 198 of file fvModelTemplates.C.

References rho.

◆ d2dt2() [2/2]

Foam::tmp<Foam::fvMatrix<Type> > d2dt2 ( const VolField< Type > &  field) const

Definition at line 210 of file fvModelTemplates.C.

References Foam::dimTime, Foam::dimVolume, and Foam::sqr().

Here is the call graph for this function:

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