Go to the documentation of this file.
35 template<
class Thermo,
class OtherThermo>
49 dict.getCheck<scalar>(
"C", scalarMinMax::ge(0))
57 dict.getCheck<scalar>(
"Tliquidus", scalarMinMax::ge(0))
65 dict.getCheck<scalar>(
"Tsolidus", scalarMinMax::ge(0))
73 dict.getCheck<scalar>(
"oxideCrit", scalarMinMax::ge(0))
81 this->mesh_.time().timeName(),
89 isoAlpha_(
dict.getOrDefault<scalar>(
"isoAlpha", 0.5))
95 template<
class Thermo,
class OtherThermo>
119 const label status =
cutCell.calcSubCell(celli, isoAlpha_);
122 Salpha[celli] = scalar(1);
128 max((oxideCrit_.value() - to)/oxideCrit_.value(), scalar(0));
135 - scalar(1)/
max((
T - Tsolidus_)/(Tliquidus_ - Tsolidus_),scalar(1
e-6))
139 mDotOxide_ = C_*tSalpha*tSoxide*tST;
145 if (isA<timeVaryingMassSorptionFvPatchScalarField>(alphab[patchi]))
148 refCast<const timeVaryingMassSorptionFvPatchScalarField>
160 mesh.time().timeName(),
168 rhoto = this->pair().
to().
rho();
172 const label cellI = fc[faceI];
173 const scalar rhoI = rhoto[cellI];
174 mDotOxide_[cellI] += rhoI*tsb()[faceI];
183 template<
class Thermo,
class OtherThermo>
195 template<
class Thermo,
class OtherThermo>
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Base class for interface composition models, templated on the two thermodynamic models either side of...
virtual tmp< volScalarField > KSu(label modelVariable, const volScalarField &field)
Explicit mass transfer coefficient.
A class for managing temporary objects.
static constexpr const zero Zero
Global zero (0)
const dimensionSet dimDensity
virtual const phaseModel & to() const
To phase.
interfaceOxideRate(const dictionary &dict, const phasePair &pair)
Construct from components.
dimensionedScalar exp(const dimensionedScalar &ds)
#define forAll(list, i)
Loop across all elements in list.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
virtual tmp< volScalarField > Kexp(const volScalarField &field)
Explicit total mass transfer coefficient.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
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.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Class for cutting a cell, celli, of an fvMesh, mesh_, at its intersection with an isosurface defined ...
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
virtual const phaseModel & from() const
From phase.
const dimensionedScalar e
Elementary charge.
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
virtual tmp< volScalarField > KSp(label modelVariable, const volScalarField &field)
Implicit mass transfer coefficient.
const dimensionedScalar & rho() const
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
const Boundary & boundaryField() const
Return const-reference to the boundary field.