effectivenessHeatExchanger Class Reference

Heat exchanger model, based on an effectiveness. More...

Inheritance diagram for effectivenessHeatExchanger:
Collaboration diagram for effectivenessHeatExchanger:

Public Member Functions

 TypeName ("effectivenessHeatExchanger")
 Runtime type information. More...
 
 effectivenessHeatExchanger (const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
 Construct from components. More...
 
 effectivenessHeatExchanger (const effectivenessHeatExchanger &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~effectivenessHeatExchanger ()
 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, const volScalarField &he, fvMatrix< scalar > &eqn) 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 effectivenessHeatExchanger &)=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...
 
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 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 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

Heat exchanger model, based on an effectiveness.

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:
effectivenessHeatExchanger1
{
    type                    effectivenessHeatExchanger;

    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 158 of file effectivenessHeatExchanger.H.

Constructor & Destructor Documentation

◆ effectivenessHeatExchanger() [1/2]

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

Construct from components.

Definition at line 186 of file effectivenessHeatExchanger.C.

◆ effectivenessHeatExchanger() [2/2]

Disallow default bitwise copy construction.

◆ ~effectivenessHeatExchanger()

virtual ~effectivenessHeatExchanger ( )
inlinevirtual

Destructor.

Definition at line 247 of file effectivenessHeatExchanger.H.

Member Function Documentation

◆ TypeName()

TypeName ( "effectivenessHeatExchanger"  )

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 217 of file effectivenessHeatExchanger.C.

References thermo.

◆ addSup()

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

◆ movePoints()

bool movePoints ( )
virtual

Update for mesh motion.

Implements fvModel.

Definition at line 330 of file effectivenessHeatExchanger.C.

◆ topoChange()

void topoChange ( const polyTopoChangeMap map)
virtual

Update topology using the given map.

Implements fvModel.

Definition at line 337 of file effectivenessHeatExchanger.C.

◆ mapMesh()

void mapMesh ( const polyMeshMap map)
virtual

Update from another mesh using the given map.

Implements fvModel.

Definition at line 346 of file effectivenessHeatExchanger.C.

◆ distribute()

void distribute ( const polyDistributionMap map)
virtual

Redistribute or update using the given distribution map.

Implements fvModel.

Definition at line 352 of file effectivenessHeatExchanger.C.

◆ read()

bool read ( const dictionary dict)
virtual

Read dictionary.

Reimplemented from fvModel.

Definition at line 361 of file effectivenessHeatExchanger.C.

References dict, and fvModel::read().

Here is the call graph for this function:

◆ operator=()

void operator= ( const effectivenessHeatExchanger )
delete

Disallow default bitwise assignment.


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