massSource Class Reference

This fvModel applies a mass source to the continuity equation and to all field equations. It can be applied to compressible solvers, such as fluid, isothermalFluid, compressibleVoF and multiphaseEuler. For incompressible solvers, use the volumeSource model instead. More...

Inheritance diagram for massSource:
Collaboration diagram for massSource:

Public Member Functions

 TypeName ("massSource")
 Runtime type information. More...
 
 massSource (const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
 Construct from explicit source name and mesh. More...
 
virtual labelUList cells () const
 Return the cells that the source applies to. More...
 
virtual label nCells () const
 Return the number of cells that the source applies to. More...
 
virtual scalar V () const
 Return the volume of cells that the source applies to. More...
 
virtual dimensionedScalar S () const
 Return the source value. More...
 
virtual bool movePoints ()
 Update for mesh motion. More...
 
virtual void topoChange (const polyTopoChangeMap &)
 Update topology using the given map. More...
 
virtual void mapMesh (const polyMeshMap &)
 Update from another mesh using the given map. More...
 
virtual void distribute (const polyDistributionMap &)
 Redistribute or update using the given distribution map. More...
 
virtual bool read (const dictionary &dict)
 Read source dictionary. More...
 
- Public Member Functions inherited from massSourceBase
 TypeName ("massSourceBase")
 Runtime type information. More...
 
 massSourceBase (const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
 Construct from explicit source name and mesh. More...
 
virtual void addSup (fvMatrix< scalar > &eqn) const
 Add a source term to a field-less proxy equation. More...
 
 FOR_ALL_FIELD_TYPES (DEFINE_FV_MODEL_ADD_FIELD_SUP)
 Add a source term to an equation. More...
 
 FOR_ALL_FIELD_TYPES (DEFINE_FV_MODEL_ADD_RHO_FIELD_SUP)
 Add a source term to a compressible equation. More...
 
 FOR_ALL_FIELD_TYPES (DEFINE_FV_MODEL_ADD_ALPHA_RHO_FIELD_SUP)
 Add a source term to a phase equation. More...
 
- Public Member Functions inherited from fvTotalSource
 TypeName ("fvTotalSource")
 Runtime type information. More...
 
 fvTotalSource (const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
 Construct from explicit source name and mesh. More...
 
virtual ~fvTotalSource ()
 Destructor. More...
 
virtual bool addsSupToField (const word &fieldName) const
 Return true if the fvModel adds a source term to the given. More...
 
const wordphaseName () const
 Return the phase name. More...
 
virtual tmp< scalarFieldsource (const word &fieldName) const
 Return the source value. More...
 
- Public Member Functions inherited from fvSource
 TypeName ("fvSource")
 Runtime type information. More...
 
 fvSource (const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
 Construct from explicit source name and mesh. More...
 
 fvSource (const fvSource &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~fvSource ()
 Destructor. More...
 
virtual wordList addSupFields () const
 Return the list of fields for which the fvModel adds source term. More...
 
void operator= (const fvSource &)=delete
 Disallow default bitwise assignment. More...
 
- Public Member Functions inherited from fvModel
 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 scalar maxDeltaT () const
 Return the maximum time-step for stable operation. 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 void correct ()
 Correct the fvModel. 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
 

Additional Inherited Members

- Static Public Member Functions inherited from fvModel
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 inherited from fvTotalSource
void addSource (fvMatrix< scalar > &eqn) const
 Add a source term to a field-less proxy equation. More...
 
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...
 
- Protected Member Functions inherited from fvModel
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

This fvModel applies a mass source to the continuity equation and to all field equations. It can be applied to compressible solvers, such as fluid, isothermalFluid, compressibleVoF and multiphaseEuler. For incompressible solvers, use the volumeSource model instead.

This model requires a corresponding field source to be specified for all solved-for fields.

Usage
Example usage for a constant mass flow rate applied to a cell set:
massSource
{
    type            massSource;

    select          cellSet;
    cellSet         massSource;

    massFlowRate    1e-4;
}

Example usage for a pulsing flow rate applied at a point:

massSource
{
    type            massSource;

    select          points;
    points          ((2.75 0.5 0));

    massFlowRate
    {
        type            scale;
        scale           squarePulse;
        start           0.2;
        duration        2;
        value           1e-4;
    }
}
Source files
See also
Foam::fv::volumeSource

Definition at line 94 of file massSource.H.

Constructor & Destructor Documentation

◆ massSource()

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

Construct from explicit source name and mesh.

Definition at line 63 of file massSource.C.

Member Function Documentation

◆ TypeName()

TypeName ( "massSource"  )

Runtime type information.

◆ cells()

Foam::labelUList cells ( ) const
virtual

Return the cells that the source applies to.

Implements fvSource.

Definition at line 81 of file massSource.C.

◆ nCells()

Foam::label nCells ( ) const
virtual

Return the number of cells that the source applies to.

Reimplemented from fvSource.

Definition at line 87 of file massSource.C.

◆ V()

Foam::scalar V ( ) const
virtual

Return the volume of cells that the source applies to.

Implements fvTotalSource.

Definition at line 93 of file massSource.C.

◆ S()

Foam::dimensionedScalar S ( ) const
virtual

Return the source value.

Implements fvTotalSource.

Definition at line 99 of file massSource.C.

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

◆ movePoints()

bool movePoints ( )
virtual

Update for mesh motion.

Implements fvModel.

Definition at line 110 of file massSource.C.

◆ topoChange()

void topoChange ( const polyTopoChangeMap map)
virtual

Update topology using the given map.

Implements fvModel.

Definition at line 117 of file massSource.C.

◆ mapMesh()

void mapMesh ( const polyMeshMap map)
virtual

Update from another mesh using the given map.

Implements fvModel.

Definition at line 123 of file massSource.C.

◆ distribute()

void distribute ( const polyDistributionMap map)
virtual

Redistribute or update using the given distribution map.

Implements fvModel.

Definition at line 129 of file massSource.C.

◆ read()

bool read ( const dictionary dict)
virtual

Read source dictionary.

Reimplemented from massSourceBase.

Definition at line 135 of file massSource.C.

References dict, and massSourceBase::read().

Here is the call graph for this function:

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