40template<
class Thermo,
class OtherThermo>
55 forAll(interfaceArea_, celli)
57 label status =
cutCell.calcSubCell(celli, isoAlpha_);
58 interfaceArea_[celli] = 0;
61 interfaceArea_[celli] =
70template<
class Thermo,
class OtherThermo>
121 isoAlpha_(
dict.getOrDefault<scalar>(
"isoAlpha", 0.5))
126 const typename OtherThermo::thermoType& toThermo =
134 Mv_.
value() = toThermo.W()*1
e-3;
136 if (Mv_.
value() == -1)
139 <<
" Please provide the molar weight (Mv) of vapour [g/mol] "
147template<
class Thermo,
class OtherThermo>
177 auto& rhov = tRhov.ref();
186 auto& deltaT = tdeltaT.ref();
190 if (
sign(C_.value()) > 0)
192 rhov = this->
pair().to().rho();
193 deltaT =
max(
T - Tactivate_,
T0);
197 rhov = this->
pair().from().rho();
198 deltaT =
max(Tactivate_ -
T,
T0);
201 htc_ = 2*
mag(C_)/(2-
mag(C_))*(
L()*rhov/HerztKnudsConst);
203 mDotc_ = htc_*deltaT*interfaceArea_;
209template<
class Thermo,
class OtherThermo>
217 if (this->modelVariable_ == variable)
221 if (
sign(C_.value()) > 0)
223 return(coeff*
pos(refValue - Tactivate_));
227 return(coeff*
pos(Tactivate_ - refValue));
235template<
class Thermo,
class OtherThermo>
243 if (this->modelVariable_ == variable)
247 if (
sign(C_.value()) > 0)
249 return(-coeff*
pos(refValue - Tactivate_));
253 return(coeff*
pos(Tactivate_ - refValue));
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
Defines the attributes of an object for which implicit objectRegistry management is supported,...
word member() const
Return member (name without the extension).
Base class for interface composition models, templated on the two thermodynamic models either side of...
InterfaceCompositionModel(const dictionary &dict, const phasePair &pair)
Construct from components.
virtual tmp< volScalarField > L(const word &speciesName, const volScalarField &Tf) const
Latent heat (to - from)(thermo - otherThermo).
const OtherThermo & toThermo_
Other Thermo (to).
const pureMixture< ThermoType >::thermoType & getLocalThermo(const word &speciesName, const pureMixture< ThermoType > &globalThermo) const
Get a reference to the local thermo for a pure mixture.
static FOAM_NO_DANGLING_REFERENCE const volPointInterpolation & New(const fvMesh &mesh, Args &&... args)
Class for cutting a cell, celli, of an fvMesh, mesh_, at its intersection with an isosurface defined ...
Service routines for cutting a cell, celli, of an fvMesh, mesh_, at its intersection with a surface.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
const Type & value() const noexcept
Return const reference to value.
Mesh data needed to do the Finite Volume discretisation.
Considering the Hertz Knudsen formula, which gives the evaporation-condensation flux based on the kin...
virtual tmp< volScalarField > KSu(label modelVariable, const volScalarField &field)
Explicit mass transfer coefficient.
virtual tmp< volScalarField > KSp(label modelVariable, const volScalarField &field)
Implicit mass transfer coefficient.
const fvMesh & mesh_
Reference to mesh.
const word transferSpecie() const
Return the transferring species name.
const phasePair & pair() const
The phase pair.
const word & variable() const
Returns the variable on which the model is based.
modelVariable modelVariable_
Enumeration for the model variable.
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
A class for managing temporary objects.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
A class for handling words, derived from Foam::string.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
constexpr scalar pi(M_PI)
const dimensionedScalar R
Universal gas constant: default SI units: [J/mol/K].
Different types of constants.
dimensionedScalar pos(const dimensionedScalar &ds)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
const dimensionSet dimless
Dimensionless.
dimensionedScalar sign(const dimensionedScalar &ds)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
const dimensionSet dimArea(sqr(dimLength))
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
#define forAll(list, i)
Loop across all elements in list.
const vector L(dict.get< vector >("L"))