Go to the documentation of this file.
45 VoFSolidificationMeltingSource,
54 void Foam::fv::VoFSolidificationMeltingSource::update()
64 <<
" - updating solid phase fraction" <<
endl;
69 const twoPhaseMixtureThermo&
thermo
83 const label celli =
cells_[i];
85 alphaSolid_[celli] =
min
87 relax_*alphaVoF[celli]*alphaSolidT_->value(TVoF[celli])
88 + (1 - relax_)*alphaSolid_[celli],
99 Foam::word Foam::fv::VoFSolidificationMeltingSource::alphaSolidName()
const
101 const twoPhaseMixtureThermo&
thermo
103 mesh_.lookupObject<twoPhaseMixtureThermo>
117 Foam::fv::VoFSolidificationMeltingSource::VoFSolidificationMeltingSource
119 const word& sourceName,
120 const word& modelType,
121 const dictionary&
dict,
125 fv::cellSetOption(sourceName, modelType,
dict,
mesh),
126 alphaSolidT_(Function1<scalar>::
New(
"alphaSolidT", coeffs_, &
mesh)),
128 relax_(coeffs_.getOrDefault(
"relax", 0.9)),
129 Cu_(coeffs_.getOrDefault<scalar>(
"Cu", 100000)),
130 q_(coeffs_.getOrDefault<scalar>(
"q", 0.001)),
138 IOobject::READ_IF_PRESENT,
143 zeroGradientFvPatchScalarField::typeName
147 fieldNames_.resize(2);
148 fieldNames_[0] =
"U";
149 fieldNames_[1] =
"T";
159 fvMatrix<scalar>& eqn,
163 apply(geometricOneField(), eqn);
170 fvMatrix<scalar>& eqn,
180 fvMatrix<vector>& eqn,
186 Info<<
type() <<
": applying source to " << eqn.psi().name() <<
endl;
196 const label celli = cells_[i];
197 const scalar Vc = V[celli];
198 const scalar alphaFluid = 1 - alphaSolid_[celli];
200 const scalar S = Cu_*
sqr(1 - alphaFluid)/(
pow3(alphaFluid) + q_);
210 fvMatrix<vector>& eqn,
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldi)
A class for handling words, derived from Foam::string.
static constexpr const zero Zero
const dimensionSet dimEnergy
static const word dictName
psiReactionThermo & thermo
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Ostream & endl(Ostream &os)
label min(const labelHashSet &set, label minValue=labelMax)
label timeIndex() const noexcept
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const Type & lookupObject(const word &name, const bool recursive=false) const
Macros for easy insertion into run-time selection tables.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
void correctBoundaryConditions()
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static void apply(bitSet &selection, const Detail::parcelSelection::actionType action, const Predicate &accept, const UList< Type > &list, const AccessOp &aop)
fileName::Type type(const fileName &name, const bool followLink=true)
const Time & time() const
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
static word groupName(StringType base, const word &group)
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
const dimensionSet dimless