effectivenessHeatExchangerSource Class Reference

Heat exchanger source model, in which the heat exchanger is defined as a selection of cells. More...

Inheritance diagram for effectivenessHeatExchangerSource:
Collaboration diagram for effectivenessHeatExchangerSource:

Public Member Functions

 TypeName ("effectivenessHeatExchangerSource")
 Runtime type information. More...
 
 effectivenessHeatExchangerSource (const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
 Construct from components. More...
 
 effectivenessHeatExchangerSource (const effectivenessHeatExchangerSource &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~effectivenessHeatExchangerSource ()
 Destructor. More...
 
virtual wordList addSupFields () const
 Return the list of fields for which the fvModel adds source term. More...
 
virtual void addSup (const volScalarField &rho, fvMatrix< scalar > &eqn, const word &fieldName) const
 Explicit and implicit source for compressible equation. 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...
 
void operator= (const effectivenessHeatExchangerSource &)=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 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...
 
virtual void correct ()
 Correct the fvModel. 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

- Static Public Member Functions inherited from fvModel
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< fvModelNew (const word &name, const fvMesh &mesh, const dictionary &dict)
 Return a reference to the selected fvModel. More...
 
- Protected Member Functions inherited from fvModel
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...
 

Detailed Description

Heat exchanger source model, in which the heat exchanger is defined as a selection of cells.

The total heat exchange source is given by:

\[ Q_t = e(\phi, \dot{m}_2) (T_2 - T_1) \phi c_p \]

where:

$ Q_t $ = total heat source
$ e(\phi,\dot{m}_2) $ = effectiveness table
$ \phi $ = net mass flux entering heat exchanger [kg/s]
$ \dot{m}_2 $ = secondary mass flow rate [kg/s]
$ T_1 $ = primary inlet temperature [K]
$ T_2 $ = secondary inlet temperature [K]
$ c_p $ = specific heat capacity [J/kg/K]

The distribution inside the hear exchanger is given by:

\[ Q_c = \frac{V_c |U_c| (T_c - T_{ref})}{\sum(V_c |U_c| (T_c - T_{ref}))} \]

where:

$ Q_c $ = source for cell
$ V_c $ = volume of the cell [m^3]
$ U_c $ = local cell velocity [m/s]
$ T_c $ = local call temperature [K]
$ T_{ref} $ = min or max(T) in cell zone depending on the sign of Q_t [K]
Usage
Example usage:
effectivenessHeatExchangerSource1
{
    type                    effectivenessHeatExchangerSource;

    select                  cellZone;
    cellZone                porosity;

    secondaryMassFlowRate   1.0;
    secondaryInletT         336;
    primaryInletT           293;

    faceZone                facesZoneInletOriented;

    effectiveness           <function2>;
}

Note:

  • The effectiveness Function2 is described in terms of the primary and secondary mass flow rates and has the same units as the secondary mass flow rate and kg/s for phi.
  • faceZone is the faces at the inlet of the cellzone, it needs to be created with flip map flags. It is used to integrate the net mass flow rate into the heat exchanger
Source files

Definition at line 160 of file effectivenessHeatExchangerSource.H.

Constructor & Destructor Documentation

◆ effectivenessHeatExchangerSource() [1/2]

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

Construct from components.

Definition at line 170 of file effectivenessHeatExchangerSource.C.

◆ effectivenessHeatExchangerSource() [2/2]

Disallow default bitwise copy construction.

◆ ~effectivenessHeatExchangerSource()

virtual ~effectivenessHeatExchangerSource ( )
inlinevirtual

Destructor.

Definition at line 249 of file effectivenessHeatExchangerSource.H.

Member Function Documentation

◆ TypeName()

TypeName ( "effectivenessHeatExchangerSource"  )

Runtime type information.

◆ 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 201 of file effectivenessHeatExchangerSource.C.

References thermo.

◆ addSup()

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

◆ movePoints()

bool movePoints ( )
virtual

Update for mesh motion.

Implements fvModel.

Definition at line 314 of file effectivenessHeatExchangerSource.C.

◆ topoChange()

void topoChange ( const polyTopoChangeMap map)
virtual

Update topology using the given map.

Implements fvModel.

Definition at line 321 of file effectivenessHeatExchangerSource.C.

◆ mapMesh()

void mapMesh ( const polyMeshMap map)
virtual

Update from another mesh using the given map.

Implements fvModel.

Definition at line 330 of file effectivenessHeatExchangerSource.C.

◆ distribute()

void distribute ( const polyDistributionMap map)
virtual

Redistribute or update using the given distribution map.

Implements fvModel.

Definition at line 336 of file effectivenessHeatExchangerSource.C.

◆ read()

bool read ( const dictionary dict)
virtual

Read dictionary.

Reimplemented from fvModel.

Definition at line 345 of file effectivenessHeatExchangerSource.C.

References dict, and fvModel::read().

Here is the call graph for this function:

◆ operator=()

void operator= ( const effectivenessHeatExchangerSource )
delete

Disallow default bitwise assignment.


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