92#ifndef meltingEvaporationModels_Lee_H
93#define meltingEvaporationModels_Lee_H
95#include "InterfaceCompositionModel.H"
101namespace meltingEvaporationModels
108template<
class Thermo,
class OtherThermo>
111 public InterfaceCompositionModel<Thermo, OtherThermo>
142 virtual ~Lee() =
default;
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,...
Mass transfer Lee model. Simple model driven by field value difference as:
TypeName("Lee")
Runtime type information.
virtual ~Lee()=default
Destructor.
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 const dimensionedScalar & Tactivate() const noexcept
Return T transition between phases.
virtual bool includeDivU() const noexcept
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.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.