74#ifndef meltingEvaporationModels_interfaceHeatResistance_H
75#define meltingEvaporationModels_interfaceHeatResistance_H
77#include "InterfaceCompositionModel.H"
87namespace meltingEvaporationModels
94template<
class Thermo,
class OtherThermo>
97 public InterfaceCompositionModel<Thermo, OtherThermo>
134 TypeName(
"interfaceHeatResistance");
196# include "interfaceHeatResistance.C"
Base class for interface composition models, templated on the two thermodynamic models either side of...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Interface Heat Resistance type of condensation/saturation model using spread source distribution foll...
virtual tmp< volScalarField > Kexp(const volScalarField &field)
Explicit total mass transfer coefficient.
virtual tmp< volScalarField > KSu(label modelVariable, const volScalarField &field)
Explicit mass transfer coefficient.
virtual ~interfaceHeatResistance()=default
Destructor.
virtual const dimensionedScalar & Tactivate() const noexcept
Return Tactivate.
virtual bool includeDivU() const noexcept
TypeName("interfaceHeatResistance")
Runtime type information.
virtual tmp< volScalarField > KSp(label modelVariable, const volScalarField &field)
Implicit mass transfer coefficient.
const phasePair & pair() const
The phase pair.
modelVariable
Enumeration for variable based mass transfer models.
Interface Heat Resistance type of condensation/saturation model using spread source distribution foll...
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
A class for managing temporary objects.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.