55 { thermoMode::mdThermo,
"thermo" },
56 { thermoMode::mdLookup,
"lookup" },
63Foam::fv::solidificationMeltingSource::Cp()
const
77 if (CpName_ ==
"CpRef")
79 const scalar CpRef =
coeffs_.get<scalar>(
"CpRef");
115 if (curTimeIndex_ == mesh_.time().timeIndex())
122 Info<<
type() <<
": " << name_ <<
" - updating phase indicator" <<
endl;
125 if (mesh_.topoChanging())
127 deltaT_.resize(cells_.size());
137 label celli = cells_[i];
139 scalar Tc =
T[celli];
140 scalar Cpc =
Cp[celli];
141 scalar alpha1New = alpha1_[celli] + relax_*Cpc*(Tc - Tmelt_)/L_;
143 alpha1_[celli] =
clamp(alpha1New, zero_one{});
144 deltaT_[i] = Tc - Tmelt_;
147 alpha1_.correctBoundaryConditions();
149 curTimeIndex_ =
mesh_.time().timeIndex();
157 const word& sourceName,
158 const word& modelType,
164 Tmelt_(coeffs_.get<scalar>(
"Tmelt")),
165 L_(coeffs_.get<scalar>(
"L")),
166 relax_(coeffs_.getOrDefault<scalar>(
"relax", 0.9)),
167 mode_(thermoModeTypeNames_.get(
"thermoMode", coeffs_)),
168 rhoRef_(coeffs_.get<scalar>(
"rhoRef")),
169 TName_(coeffs_.getOrDefault<
word>(
"T",
"T")),
170 CpName_(coeffs_.getOrDefault<
word>(
"Cp",
"Cp")),
171 UName_(coeffs_.getOrDefault<
word>(
"U",
"U")),
172 phiName_(coeffs_.getOrDefault<
word>(
"phi",
"phi")),
173 Cu_(coeffs_.getOrDefault<scalar>(
"Cu", 100000)),
174 q_(coeffs_.getOrDefault<scalar>(
"q", 0.001)),
175 beta_(coeffs_.get<scalar>(
"beta")),
180 IOobject::scopedName(name_,
"alpha1"),
192 deltaT_(cells_.size(), 0)
270 label celli = cells_[i];
272 scalar Vc = V[celli];
273 scalar alpha1c = alpha1_[celli];
275 scalar S = -Cu_*
sqr(1.0 - alpha1c)/(
pow3(alpha1c) + q_);
276 vector Sb(rhoRef_*
g*beta_*deltaT_[i]);
300 coeffs_.readEntry(
"Tmelt", Tmelt_);
301 coeffs_.readEntry(
"L", L_);
303 coeffs_.readIfPresent(
"relax", relax_);
305 thermoModeTypeNames_.readEntry(
"thermoMode", coeffs_, mode_);
307 coeffs_.readEntry(
"rhoRef", rhoRef_);
308 coeffs_.readIfPresent(
"T", TName_);
309 coeffs_.readIfPresent(
"U", UName_);
311 coeffs_.readIfPresent(
"Cu", Cu_);
312 coeffs_.readIfPresent(
"q", q_);
314 coeffs_.readEntry(
"beta", beta_);
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const uniformDimensionedVectorField & g
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
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())
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const word & name() const noexcept
Return the object name.
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
void resize(const label len)
Adjust allocated size of list.
Abstract base-class for fluid and solid thermodynamic properties.
static const word dictName
The dictionary name ("thermophysicalProperties").
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.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Field< Type > & source() noexcept
Mesh data needed to do the Finite Volume discretisation.
Template invariant parts for fvPatchField.
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Intermediate abstract class for handling cell-set options for the derived fvOptions.
scalar V() const noexcept
Return const access to the total cell volume.
labelList cells_
Set of cells to apply source to.
virtual bool read(const dictionary &dict)
Read source dictionary.
cellSetOption(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
Base abstract class for handling finite volume options (i.e. fvOption).
const fvMesh & mesh_
Reference to the mesh database.
wordList fieldNames_
Field names to apply source to - populated by derived models.
dictionary coeffs_
Dictionary containing source coefficients.
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
const word name_
Source name.
This source is designed to model the effect of solidification and melting processes,...
solidificationMeltingSource(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
thermoMode
Options for the thermo mode specification.
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldi)
Add explicit contribution to enthalpy equation.
virtual bool read(const dictionary &dict)
Read source dictionary.
static const Enum< thermoMode > thermoModeTypeNames_
Names for thermoMode.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
const scalarField & diag() const
static const gravity & New(const word &name, const Time &runTime)
Return named gravity field cached or construct on Time.
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
const volScalarField & Cp
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Namespace for handling debugging switches.
Namespace for finite-volume.
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimEnergy
dimensionedScalar pow3(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
static void apply(bitSet &selection, const Detail::parcelSelection::actionType action, const Predicate &accept, const UList< Type > &list, const AccessOp &aop)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
errorManip< error > abort(error &err)
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
void Sp(fvMatrix< typename Expr::value_type > &m, const Expr2 &mult, const Expr &expression)
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.
void Su(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
psiReactionThermo & thermo
#define forAll(list, i)
Loop across all elements in list.