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...
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 volScalarField & | htc () 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 (fvMatrix< scalar > &eqn, const word &fieldName) const |
Source term to energy equation. More... | |
virtual void | addSup (const volScalarField &rho, fvMatrix< scalar > &eqn, const word &fieldName) 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... | |
![]() | |
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 word & | nbrRegionName () const |
Return const access to the neighbour region name. More... | |
const fvMesh & | nbrMesh () const |
Return const access to the neighbour mesh. More... | |
const cellsToCells & | interpolation () 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 |
![]() | |
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< fvModel > | clone () const |
Return clone. More... | |
virtual | ~fvModel () |
Destructor. More... | |
const word & | name () const |
Return const access to the source name. More... | |
const fvMesh & | mesh () const |
Return const access to the mesh database. More... | |
const dictionary & | coeffs () 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... | |
FOR_ALL_FIELD_TYPES (DEFINE_FV_MODEL_ADD_SUP) | |
Add a source term to an equation. More... | |
FOR_ALL_FIELD_TYPES (DEFINE_FV_MODEL_ADD_RHO_SUP) | |
Add a source term to a compressible equation. More... | |
FOR_ALL_FIELD_TYPES (DEFINE_FV_MODEL_ADD_ALPHA_RHO_SUP) | |
Add a source term to a phase 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 > > | source (const VolField< Type > &field, const word &fieldName) const |
Return source for an equation with a specified name. 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 > > | source (const volScalarField &rho, const VolField< Type > &field, const word &fieldName) const |
Return source for a compressible equation with a specified name. 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 > > | source (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field, const word &fieldName) const |
Return source for a phase equation with a specified name. 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... | |
template<class Type > | |
tmp< fvMatrix< Type > > | d2dt2 (const VolField< Type > &field, const word &fieldName) const |
Return source for an equation with a second time derivative. More... | |
virtual void | preUpdateMesh () |
Prepare for mesh update. More... | |
template<class Type , class AlphaRhoFieldType , class ... AlphaRhoFieldTypes> | |
Foam::dimensionSet | sourceDims (const VolField< Type > &field, const dimensionSet &ds, const AlphaRhoFieldType &alphaRho, const AlphaRhoFieldTypes &... alphaRhos) |
template<class Type > | |
Foam::dimensionSet | sourceDims (const VolField< Type > &field, const dimensionSet &ds) |
template<class Type , class ... AlphaRhoFieldTypes> | |
Foam::tmp< Foam::fvMatrix< Type > > | source (const VolField< Type > &field, const word &fieldName, const dimensionSet &ds, const AlphaRhoFieldTypes &... alphaRhos) const |
template<class Type > | |
Foam::tmp< Foam::fvMatrix< Type > > | source (const VolField< Type > &field) const |
template<class Type > | |
Foam::tmp< Foam::fvMatrix< Type > > | source (const VolField< Type > &field, const word &fieldName) 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 > > | source (const volScalarField &rho, const VolField< Type > &field, const word &fieldName) 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 > > | source (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field, const word &fieldName) 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 |
template<class Type > | |
Foam::tmp< Foam::fvMatrix< Type > > | d2dt2 (const VolField< Type > &field, const word &fieldName) const |
Additional Inherited Members | |
![]() | |
template<class Type , class AlphaRhoFieldType , class ... AlphaRhoFieldTypes> | |
static dimensionSet | sourceDims (const VolField< Type > &field, const dimensionSet &ds, const AlphaRhoFieldType &alphaRho, const AlphaRhoFieldTypes &... alphaRhos) |
Return the dimensions of the matrix of a source term. More... | |
template<class Type > | |
static dimensionSet | sourceDims (const VolField< Type > &field, const dimensionSet &ds) |
Return the dimensions of the matrix of a source term (base. More... | |
static autoPtr< fvModel > | New (const word &name, const fvMesh &mesh, const dictionary &dict) |
Return a reference to the selected fvModel. More... | |
![]() | |
const interRegionModel & | nbrModel () const |
Get the neighbour interRegionModel. More... | |
![]() | |
template<class Type > | |
void | addSupType (fvMatrix< Type > &eqn, const word &fieldName) const |
Add a source term to an equation. More... | |
template<class Type > | |
void | addSupType (const volScalarField &rho, fvMatrix< Type > &eqn, const word &fieldName) const |
Add a source term to a compressible equation. More... | |
template<class Type > | |
void | addSupType (const volScalarField &alpha, const volScalarField &rho, fvMatrix< Type > &eqn, const word &fieldName) const |
Add a source term to a phase equation. More... | |
template<class Type , class ... AlphaRhoFieldTypes> | |
tmp< fvMatrix< Type > > | source (const VolField< Type > &field, const word &fieldName, const dimensionSet &ds, const AlphaRhoFieldTypes &... alphaRhos) const |
Return source for equation with specified name and dimensions. More... | |
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:
If the semiImplicit option is set, then this becomes:
interRegionHeatTransfer { type interRegionHeatTransfer; nbrRegion other; interpolationMethod intersection; master true; semiImplicit no; Av 200; heatTransferCoefficientModel constant; htc 10; }
Definition at line 89 of file interRegionHeatTransfer.H.
interRegionHeatTransfer | ( | const word & | name, |
const word & | modelType, | ||
const fvMesh & | mesh, | ||
const dictionary & | dict | ||
) |
Construct from dictionary.
Definition at line 79 of file interRegionHeatTransfer.C.
|
virtual |
Destructor.
Definition at line 100 of file interRegionHeatTransfer.C.
TypeName | ( | "interRegionHeatTransfer" | ) |
Runtime type information.
|
inline |
Return the heat transfer coefficient.
|
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.
Source term to energy equation.
Definition at line 115 of file interRegionHeatTransfer.C.
References Foam::dimEnergy, Foam::dimless, Foam::dimMass, Foam::dimTemperature, he(), Foam::fvc::interpolate(), Foam::name(), GeometricField< Type, PatchField, GeoMesh >::New(), fvMatrix< Type >::psi(), tmp< T >::ref(), Foam::fvm::Sp(), Foam::T(), and thermo.
|
virtual |
Source term to compressible energy equation.
Definition at line 199 of file interRegionHeatTransfer.C.
|
virtual |
Correct the model.
Reimplemented from fvModel.
Definition at line 210 of file interRegionHeatTransfer.C.
|
virtual |
Update for mesh motion.
Implements fvModel.
Definition at line 219 of file interRegionHeatTransfer.C.
References NotImplemented.
|
virtual |
Update topology using the given map.
Implements fvModel.
Definition at line 226 of file interRegionHeatTransfer.C.
References NotImplemented.
|
virtual |
Update from another mesh using the given map.
Implements fvModel.
Definition at line 232 of file interRegionHeatTransfer.C.
References NotImplemented.
|
virtual |
Redistribute or update using the given distribution map.
Implements fvModel.
Definition at line 238 of file interRegionHeatTransfer.C.
References NotImplemented.
|
virtual |
Read dictionary.
Reimplemented from interRegionModel.
Definition at line 247 of file interRegionHeatTransfer.C.
References dict, and interRegionModel::read().