115#ifndef meltingEvaporationModels_kineticGasEvaporation_H
116#define meltingEvaporationModels_kineticGasEvaporation_H
118#include "InterfaceCompositionModel.H"
128namespace meltingEvaporationModels
135template<
class Thermo,
class OtherThermo>
136class kineticGasEvaporation
138 public InterfaceCompositionModel<Thermo, OtherThermo>
181 const dictionary&
dict,
182 const phasePair&
pair
193 virtual tmp<volScalarField>
Kexp
199 virtual tmp<volScalarField>
KSp
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,...
Considering the Hertz Knudsen formula, which gives the evaporation-condensation flux based on the kin...
virtual tmp< volScalarField > Kexp(const volScalarField &field)
Explicit total mass transfer coefficient.
virtual ~kineticGasEvaporation()=default
Destructor.
virtual tmp< volScalarField > KSu(label modelVariable, const volScalarField &field)
Explicit mass transfer coefficient.
virtual const dimensionedScalar & Tactivate() const noexcept
Return Tactivate.
virtual bool includeDivU() const noexcept
kineticGasEvaporation(const dictionary &dict, const phasePair &pair)
Construct from components.
TypeName("kineticGasEvaporation")
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.
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.