interRegionHeatTransfer Class Reference

Model for inter-region heat exchange. Requires specification of a model for the heat transfer coefficient (htc) and the area per unit volume (Av). These are then used to apply the following source to the energy equation: More...

Inheritance diagram for interRegionHeatTransfer:
Collaboration diagram for interRegionHeatTransfer:

Public Member Functions

 TypeName ("interRegionHeatTransfer")
 Runtime type information. More...
 
 interRegionHeatTransfer (const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
 Construct from dictionary. More...
 
virtual ~interRegionHeatTransfer ()
 Destructor. More...
 
const volScalarFieldhtc () const
 Return the heat transfer coefficient. More...
 
virtual wordList addSupFields () const
 Return the list of fields for which the fvModel adds source term. More...
 
virtual void addSup (const volScalarField &he, fvMatrix< scalar > &eqn) const
 Source term to energy equation. More...
 
virtual void addSup (const volScalarField &rho, const volScalarField &he, fvMatrix< scalar > &eqn) const
 Source term to compressible energy equation. More...
 
virtual void correct ()
 Correct the model. 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 dictionary. More...
 
- Public Member Functions inherited from interRegionModel
 TypeName ("interRegionModel")
 Runtime type information. More...
 
 interRegionModel (const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
 Construct from dictionary. More...
 
virtual ~interRegionModel ()
 Destructor. More...
 
bool master () const
 Return whether the master region. More...
 
const wordnbrRegionName () const
 Return const access to the neighbour region name. More...
 
const fvMeshnbrMesh () const
 Return const access to the neighbour mesh. More...
 
const cellsToCellsinterpolation () const
 Return const access to the interpolation engine. More...
 
template<class Type >
tmp< Field< Type > > interpolate (const Field< Type > &field) const
 Interpolate field. More...
 
template<class Type >
void interpolate (const Field< Type > &field, Field< Type > &result) const
 Interpolate field. More...
 
template<class Type >
Foam::tmp< Foam::Field< Type > > interpolate (const Field< Type > &field) const
 
- 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 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 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 interRegionModel
const interRegionModelnbrModel () const
 Get the neighbour interRegionModel. 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

Model for inter-region heat exchange. Requires specification of a model for the heat transfer coefficient (htc) and the area per unit volume (Av). These are then used to apply the following source to the energy equation:

\[ -htc*Av*(T_{nbr,mapped} - T) \]

If the semiImplicit option is set, then this becomes:

\[ -htc*Av*(T_{nbr,mapped} - T) + htc*Av/Cp*h - Sp(htc*Av/Cp, h); \]

Usage
Example usage:
interRegionHeatTransfer
{
    type            interRegionHeatTransfer;

    libs            ("libinterRegionFvModels.so");

    nbrRegion       other;

    interpolationMethod intersection;
    master          true;

    semiImplicit    no;

    Av              200;

    heatTransferCoefficientModel constant;

    htc             10;
}
See also
Foam::fv::heatTransferCoefficientModel
Source files

Definition at line 91 of file interRegionHeatTransfer.H.

Constructor & Destructor Documentation

◆ interRegionHeatTransfer()

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

Construct from dictionary.

Definition at line 79 of file interRegionHeatTransfer.C.

◆ ~interRegionHeatTransfer()

Destructor.

Definition at line 100 of file interRegionHeatTransfer.C.

Member Function Documentation

◆ TypeName()

TypeName ( "interRegionHeatTransfer"  )

Runtime type information.

◆ htc()

const volScalarField& htc ( ) const
inline

Return the heat transfer coefficient.

◆ addSupFields()

Foam::wordList addSupFields ( ) const
virtual

Return the list of fields for which the fvModel adds source term.

to the transport equation

Reimplemented from fvModel.

Definition at line 106 of file interRegionHeatTransfer.C.

References thermo.

◆ addSup() [1/2]

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

◆ addSup() [2/2]

void addSup ( const volScalarField rho,
const volScalarField he,
fvMatrix< scalar > &  eqn 
) const
virtual

Source term to compressible energy equation.

Definition at line 197 of file interRegionHeatTransfer.C.

References he().

Here is the call graph for this function:

◆ correct()

void correct ( )
virtual

Correct the model.

Reimplemented from fvModel.

Definition at line 208 of file interRegionHeatTransfer.C.

◆ movePoints()

bool movePoints ( )
virtual

Update for mesh motion.

Implements fvModel.

Definition at line 217 of file interRegionHeatTransfer.C.

References NotImplemented.

◆ topoChange()

void topoChange ( const polyTopoChangeMap )
virtual

Update topology using the given map.

Implements fvModel.

Definition at line 224 of file interRegionHeatTransfer.C.

References NotImplemented.

◆ mapMesh()

void mapMesh ( const polyMeshMap )
virtual

Update from another mesh using the given map.

Implements fvModel.

Definition at line 230 of file interRegionHeatTransfer.C.

References NotImplemented.

◆ distribute()

void distribute ( const polyDistributionMap )
virtual

Redistribute or update using the given distribution map.

Implements fvModel.

Definition at line 236 of file interRegionHeatTransfer.C.

References NotImplemented.

◆ read()

bool read ( const dictionary dict)
virtual

Read dictionary.

Reimplemented from interRegionModel.

Definition at line 245 of file interRegionHeatTransfer.C.

References dict, and interRegionModel::read().

Here is the call graph for this function:

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