TwoResistanceHeatTransferPhaseSystem< BasePhaseSystem > Class Template Referenceabstract

Class which models interfacial heat transfer between a number of phases. Two heat transfer models are stored at each interface, one for each phase. This permits definition of an interface temperature with which heat transfer occurs. It also allows derived systems to define other thermodynamic properties at the interface and therefore represent phase changes. More...

Inheritance diagram for TwoResistanceHeatTransferPhaseSystem< BasePhaseSystem >:
Collaboration diagram for TwoResistanceHeatTransferPhaseSystem< BasePhaseSystem >:

Public Member Functions

 TwoResistanceHeatTransferPhaseSystem (const fvMesh &)
 Construct from fvMesh. More...
 
virtual ~TwoResistanceHeatTransferPhaseSystem ()
 Destructor. More...
 
virtual autoPtr< phaseSystem::heatTransferTableheatTransfer () const
 Return the heat transfer matrices. More...
 
virtual void predictThermophysicalTransport ()
 Predict the energy transport e.g. alphat. More...
 
virtual void correctThermophysicalTransport ()
 Correct the energy transport e.g. alphat. More...
 
virtual void correctInterfaceThermo ()=0
 Correct the interface thermodynamics. More...
 
virtual bool read ()
 Read base phaseProperties dictionary. More...
 
- Public Member Functions inherited from HeatTransferPhaseSystem< BasePhaseSystem >
 HeatTransferPhaseSystem (const fvMesh &)
 Construct from fvMesh. More...
 
virtual ~HeatTransferPhaseSystem ()
 Destructor. More...
 
virtual tmp< volScalarFieldL (const phaseInterface &interface, const volScalarField &dmdtf, const volScalarField &Tf, const latentHeatScheme scheme) const
 Return the latent heat for a given interface, mass transfer rate. More...
 
virtual tmp< scalarFieldL (const phaseInterface &interface, const scalarField &dmdtf, const scalarField &Tf, const labelUList &cells, const latentHeatScheme scheme) const
 As above, but for a cell-set. More...
 
virtual tmp< volScalarFieldLi (const phaseInterface &interface, const word &member, const volScalarField &dmidtf, const volScalarField &Tf, const latentHeatScheme scheme) const
 Return the latent heat for a given interface, specie, mass transfer. More...
 
virtual tmp< scalarFieldLi (const phaseInterface &interface, const word &member, const scalarField &dmidtf, const scalarField &Tf, const labelUList &cells, const latentHeatScheme scheme) const
 As above, but for a cell-set. More...
 
- Public Member Functions inherited from heatTransferPhaseSystem
 TypeName ("heatTransferPhaseSystem")
 Runtime type information. More...
 
 heatTransferPhaseSystem ()
 Default constructor. More...
 
virtual ~heatTransferPhaseSystem ()
 Destructor. More...
 

Protected Types

typedef HashTable< autoPtr< sidedBlendedHeatTransferModel >, phaseInterfaceKey, phaseInterfaceKey::hashheatTransferModelTable
 
typedef HeatTransferPhaseSystem< BasePhaseSystem >::latentHeatScheme latentHeatScheme
 

Protected Member Functions

void addDmdtHefs (const phaseSystem::dmdtfTable &dmdtfs, const phaseSystem::dmdtfTable &Tfs, const latentHeatScheme scheme, const latentHeatTransfer transfer, phaseSystem::heatTransferTable &eqns) const
 Add energy transfer terms which result from bulk phase changes. More...
 
void addDmidtHefs (const phaseSystem::dmidtfTable &dmidtfs, const phaseSystem::dmdtfTable &Tfs, const latentHeatScheme scheme, const latentHeatTransfer transfer, phaseSystem::heatTransferTable &eqns) const
 Add energy transfer terms which result from specie phase changes. More...
 
- Protected Member Functions inherited from HeatTransferPhaseSystem< BasePhaseSystem >
void addDmdtHefs (const phaseSystem::dmdtfTable &dmdtfs, phaseSystem::heatTransferTable &eqns) const
 Add energy transfer terms which result from bulk mass transfers. More...
 
void addDmidtHefs (const phaseSystem::dmidtfTable &dmidtfs, phaseSystem::heatTransferTable &eqns) const
 Add energy transfer terms which result from specie mass transfers. More...
 
void addDmdtHefsWithoutL (const phaseSystem::dmdtfTable &dmdtfs, const phaseSystem::dmdtfTable &Tfs, const latentHeatScheme scheme, phaseSystem::heatTransferTable &eqns) const
 Add energy transfer terms which result from bulk phase changes,. More...
 
void addDmdtL (const phaseSystem::dmdtfTable &dmdtfs, const phaseSystem::dmdtfTable &Tfs, const scalar weight, const latentHeatScheme scheme, phaseSystem::heatTransferTable &eqns) const
 Add latent heat terms which result from bulk phase changes. More...
 
void addDmdtHefs (const phaseSystem::dmdtfTable &dmdtfs, const phaseSystem::dmdtfTable &Tfs, const scalar weight, const latentHeatScheme scheme, phaseSystem::heatTransferTable &eqns) const
 Add energy transfer terms which result from bulk phase changes. More...
 
void addDmidtHefsWithoutL (const phaseSystem::dmidtfTable &dmidtfs, const phaseSystem::dmdtfTable &Tfs, const latentHeatScheme scheme, phaseSystem::heatTransferTable &eqns) const
 Add energy transfer terms which result from specie phase changes,. More...
 
void addDmidtL (const phaseSystem::dmidtfTable &dmidtfs, const phaseSystem::dmdtfTable &Tfs, const scalar weight, const latentHeatScheme scheme, phaseSystem::heatTransferTable &eqns) const
 Add latent heat terms which result from specie phase changes. More...
 
void addDmidtHefs (const phaseSystem::dmidtfTable &dmidtfs, const phaseSystem::dmdtfTable &Tfs, const scalar weight, const latentHeatScheme scheme, phaseSystem::heatTransferTable &eqns) const
 Add energy transfer terms which result from specie phase changes. More...
 

Protected Attributes

heatTransferModelTable heatTransferModels_
 Heat transfer models. More...
 
- Protected Attributes inherited from HeatTransferPhaseSystem< BasePhaseSystem >
scalar residualY_
 Residual mass fraction used for linearisation of heat transfer. More...
 

Additional Inherited Members

- Public Types inherited from twoResistanceHeatTransferPhaseSystem
enum class  latentHeatTransfer { heat , mass }
 Enumeration for the form of the latent heat transfer. More...
 
- Public Types inherited from heatTransferPhaseSystem
enum class  latentHeatScheme { symmetric , upwind }
 Enumeration for the latent heat scheme. More...
 
- Static Public Attributes inherited from twoResistanceHeatTransferPhaseSystem
static const NamedEnum< latentHeatTransfer, 2 > latentHeatTransferNames_
 Names of the forms of the latent heat transfer. More...
 
- Static Public Attributes inherited from heatTransferPhaseSystem
static const NamedEnum< latentHeatScheme, 2 > latentHeatSchemeNames_
 Names of the latent heat schemes. More...
 

Detailed Description

template<class BasePhaseSystem>
class Foam::TwoResistanceHeatTransferPhaseSystem< BasePhaseSystem >

Class which models interfacial heat transfer between a number of phases. Two heat transfer models are stored at each interface, one for each phase. This permits definition of an interface temperature with which heat transfer occurs. It also allows derived systems to define other thermodynamic properties at the interface and therefore represent phase changes.

See also
OneResistanceHeatTransferPhaseSystem
Source files

Definition at line 59 of file TwoResistanceHeatTransferPhaseSystem.H.

Member Typedef Documentation

◆ heatTransferModelTable

◆ latentHeatScheme

typedef HeatTransferPhaseSystem<BasePhaseSystem>::latentHeatScheme latentHeatScheme
protected

Definition at line 77 of file TwoResistanceHeatTransferPhaseSystem.H.

Constructor & Destructor Documentation

◆ TwoResistanceHeatTransferPhaseSystem()

◆ ~TwoResistanceHeatTransferPhaseSystem()

Destructor.

Definition at line 234 of file TwoResistanceHeatTransferPhaseSystem.C.

Member Function Documentation

◆ addDmdtHefs()

void addDmdtHefs ( const phaseSystem::dmdtfTable dmdtfs,
const phaseSystem::dmdtfTable Tfs,
const latentHeatScheme  scheme,
const latentHeatTransfer  transfer,
phaseSystem::heatTransferTable eqns 
) const
protected

Add energy transfer terms which result from bulk phase changes.

that are coupled to the two-resistance heat transfer model

Definition at line 34 of file TwoResistanceHeatTransferPhaseSystem.C.

References forAllConstIter, heatTransferModel::K(), phaseModel::name(), phaseInterface::phase1(), phaseInterface::phase2(), Foam::fvc::scheme(), basicThermo::T(), and phaseModel::thermo().

Here is the call graph for this function:

◆ addDmidtHefs()

void addDmidtHefs ( const phaseSystem::dmidtfTable dmidtfs,
const phaseSystem::dmdtfTable Tfs,
const latentHeatScheme  scheme,
const latentHeatTransfer  transfer,
phaseSystem::heatTransferTable eqns 
) const
protected

Add energy transfer terms which result from specie phase changes.

that are coupled to the two-resistance heat transfer model

Definition at line 101 of file TwoResistanceHeatTransferPhaseSystem.C.

References forAllConstIter, heatTransferModel::K(), phaseModel::name(), phaseInterface::phase1(), phaseInterface::phase2(), Foam::fvc::scheme(), basicThermo::T(), and phaseModel::thermo().

Here is the call graph for this function:

◆ heatTransfer()

◆ predictThermophysicalTransport()

void predictThermophysicalTransport
virtual

Predict the energy transport e.g. alphat.

and interface properties e.g. Tf

Definition at line 297 of file TwoResistanceHeatTransferPhaseSystem.C.

◆ correctThermophysicalTransport()

void correctThermophysicalTransport
virtual

Correct the energy transport e.g. alphat.

Definition at line 306 of file TwoResistanceHeatTransferPhaseSystem.C.

◆ correctInterfaceThermo()

virtual void correctInterfaceThermo ( )
pure virtual

Correct the interface thermodynamics.

◆ read()

bool read
virtual

Read base phaseProperties dictionary.

Reimplemented from HeatTransferPhaseSystem< BasePhaseSystem >.

Definition at line 314 of file TwoResistanceHeatTransferPhaseSystem.C.

References Foam::blockMeshTools::read().

Here is the call graph for this function:

Member Data Documentation

◆ heatTransferModels_


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