39 #ifndef interfaceCompositionModel_H 40 #define interfaceCompositionModel_H 56 namespace multiphaseInter
107 TypeName(
"interfaceCompositionModel");
162 const word& speciesName,
169 const word& speciesName,
176 const word& speciesName
182 const word& speciesName
188 const word& speciesName,
static const Enum< modelVariable > modelVariableNames_
Selection names for the modelVariable.
const fvMesh & mesh_
Reference to mesh.
word speciesName_
Names of the transferring specie.
modelVariable
Enumeration for variable based mass transfer models.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
virtual tmp< volScalarField > KSp(label modelVariable, const volScalarField &field)=0
Implicit mass transfer.
virtual bool includeDivU() const noexcept
Add/subtract alpha*div(U) as a source term for alpha, substituting div(U) = mDot(1/rho1 - 1/rho2) ...
const word & variable() const
Returns the variable on which the model is based.
bool includeVolChange()
Add volume change in pEq.
const multiphaseInterSystem & fluid() const
Return the multiphaseInterSystem this interface belongs to.
virtual tmp< volScalarField > KSu(label modelVariable, const volScalarField &field)=0
Explicit mass transfer.
virtual tmp< volScalarField > Dto(const word &speciesName) const =0
Specie mass diffusivity for specie in a multicomponent.
virtual const dimensionedScalar & Tactivate() const noexcept=0
Reference value.
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const =0
Interface mass fraction.
const word transferSpecie() const
Return the transferring species name.
A class for handling words, derived from Foam::string.
bool includeVolChange_
Add volume change in pEq.
virtual tmp< volScalarField > Kexp(const volScalarField &field)=0
Explicit full mass transfer.
interfaceCompositionModel(const dictionary &dict, const phasePair &pair)
Construct from a dictionary and a phase pair.
const phasePair & pair_
Phase pair.
modelVariable modelVariable_
Enumeration for the model variable.
const phasePair & pair() const
The phase pair.
virtual tmp< volScalarField > Dfrom(const word &speciesName) const =0
Specie mass diffusivity for pure mixture.
virtual tmp< volScalarField > dY(const word &speciesName, const volScalarField &Tf) const =0
Mass fraction difference between the interface and the field.
TypeName("interfaceCompositionModel")
Runtime type information.
static autoPtr< interfaceCompositionModel > New(const dictionary &dict, const phasePair &pair)
virtual ~interfaceCompositionModel()=default
Destructor.
Mesh data needed to do the Finite Volume discretisation.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Macros to ease declaration of run-time selection tables.
A class for managing temporary objects.
declareRunTimeSelectionTable(autoPtr, interfaceCompositionModel, dictionary,(const dictionary &dict, const phasePair &pair),(dict, pair))
virtual tmp< volScalarField > L(const word &speciesName, const volScalarField &Tf) const =0
Latent heat (delta Hc)
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...