38template<
class Thermo,
class OtherThermo>
39void Foam::meltingEvaporationModels::
40diffusionGasEvaporation<Thermo, OtherThermo>::updateInterface
56 forAll(interfaceArea_, celli)
58 const label status =
cutCell.calcSubCell(celli, isoAlpha_);
59 interfaceArea_[celli] = 0;
70template<
class Thermo,
class OtherThermo>
83 dict.subDict(
"saturationPressure"),
87 isoAlpha_(
dict.getOrDefault<scalar>(
"isoAlpha", 0.5)),
121template<
class Thermo,
class OtherThermo>
133 const typename OtherThermo::thermoType& vapourThermo =
154 auto&
rhog = tRhog.ref();
165 auto& Dvg = tDvg.ref();
167 Dvg = this->
Dto(speciesName);
178 Wv/(XvSat*Wv + (1-XvSat)*this->
toThermo_.W())
192 *C_*
rhog*Dvg*gradYgm*interfaceArea_
217template<
class Thermo,
class OtherThermo>
230template<
class Thermo,
class OtherThermo>
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,...
static word member(const word &name)
Return member (name without the extension).
Base class for interface composition models, templated on the two thermodynamic models either side of...
virtual tmp< volScalarField > Dto(const word &speciesName) const
Specie mass diffusivity for specie in a multicomponent.
InterfaceCompositionModel(const dictionary &dict, const phasePair &pair)
Construct from components.
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,...
Mesh data needed to do the Finite Volume discretisation.
Gas diffusion based evaporation/condensation mass transfer model.
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 multiphaseInterSystem & fluid() const
Return the multiphaseInterSystem this interface belongs to.
static autoPtr< interfaceCompositionModel > New(const dictionary &dict, const phasePair &pair)
const word & variable() const
Returns the variable on which the model is based.
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
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.
Calculate the gradient of the given field.
Different types of constants.
Namespace for handling debugging switches.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
dimensionedScalar pos(const dimensionedScalar &ds)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
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)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0).
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
dimensionedScalar rhog("rhog", dimDensity, transportProperties)
#define forAll(list, i)
Loop across all elements in list.