Go to the documentation of this file.
55 { thermoMode::mdThermo,
"thermo" },
56 { thermoMode::mdLookup,
"lookup" },
63 Foam::fv::solidificationMeltingSource::Cp()
const
77 if (CpName_ ==
"CpRef")
79 const scalar CpRef =
coeffs_.
get<scalar>(
"CpRef");
98 extrapolatedCalculatedFvPatchScalarField::typeName
120 void Foam::fv::solidificationMeltingSource::update(
const volScalarField&
Cp)
122 if (curTimeIndex_ == mesh_.time().timeIndex())
129 Info<<
type() <<
": " << name_ <<
" - updating phase indicator" <<
endl;
139 label celli = cells_[i];
141 scalar Tc =
T[celli];
142 scalar Cpc =
Cp[celli];
143 scalar alpha1New = alpha1_[celli] + relax_*Cpc*(Tc - Tmelt_)/L_;
145 alpha1_[celli] =
max(0,
min(alpha1New, 1));
146 deltaT_[i] = Tc - Tmelt_;
149 alpha1_.correctBoundaryConditions();
151 curTimeIndex_ = mesh_.time().timeIndex();
159 const word& sourceName,
160 const word& modelType,
166 Tmelt_(coeffs_.get<scalar>(
"Tmelt")),
167 L_(coeffs_.get<scalar>(
"L")),
168 relax_(coeffs_.getOrDefault<scalar>(
"relax", 0.9)),
169 mode_(thermoModeTypeNames_.get(
"thermoMode", coeffs_)),
170 rhoRef_(coeffs_.get<scalar>(
"rhoRef")),
171 TName_(coeffs_.getOrDefault<
word>(
"T",
"T")),
172 CpName_(coeffs_.getOrDefault<
word>(
"Cp",
"Cp")),
173 UName_(coeffs_.getOrDefault<
word>(
"U",
"U")),
174 phiName_(coeffs_.getOrDefault<
word>(
"phi",
"phi")),
175 Cu_(coeffs_.getOrDefault<scalar>(
"Cu", 100000)),
176 q_(coeffs_.getOrDefault<scalar>(
"q", 0.001)),
177 beta_(coeffs_.get<scalar>(
"beta")),
183 mesh.time().timeName(),
190 zeroGradientFvPatchScalarField::typeName
193 deltaT_(cells_.size(), 0)
195 fieldNames_.setSize(2);
196 fieldNames_[0] = UName_;
205 fieldNames_[1] =
thermo.he().name();
210 fieldNames_[1] = TName_;
216 <<
"Unhandled thermo mode: " << thermoModeTypeNames_[mode_]
221 applied_.setSize(fieldNames_.size(),
false);
271 label celli = cells_[i];
273 scalar Vc = V[celli];
274 scalar alpha1c = alpha1_[celli];
276 scalar S = -Cu_*
sqr(1.0 - alpha1c)/(
pow3(alpha1c) + q_);
277 vector Sb(rhoRef_*
g*beta_*deltaT_[i]);
301 coeffs_.readEntry(
"Tmelt", Tmelt_);
302 coeffs_.readEntry(
"L", L_);
304 coeffs_.readIfPresent(
"relax", relax_);
306 thermoModeTypeNames_.readEntry(
"thermoMode", coeffs_, mode_);
308 coeffs_.readEntry(
"rhoRef", rhoRef_);
309 coeffs_.readIfPresent(
"T", TName_);
310 coeffs_.readIfPresent(
"U", UName_);
312 coeffs_.readIfPresent(
"Cu", Cu_);
313 coeffs_.readIfPresent(
"q", q_);
315 coeffs_.readEntry(
"beta", beta_);
int debug
Static debugging option.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
A class for handling words, derived from Foam::string.
Intermediate abstract class for handling cell-set options for the derived fvOptions.
A class for managing temporary objects.
static constexpr const zero Zero
Global zero (0)
const dimensionSet dimEnergy
solidificationMeltingSource(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
static word timeName(const scalar t, const int precision=precision_)
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Abstract base-class for fluid and solid thermodynamic properties.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const word name_
Source name.
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
const GeometricField< Type, fvPatchField, volMesh > & psi() const
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
thermoMode
Options for the thermo mode specification.
const fvMesh & mesh_
Reference to the mesh database.
#define forAll(list, i)
Loop across all elements in list.
A special matrix type and solver, designed for finite volume solutions of scalar equations.
dimensionedScalar pow3(const dimensionedScalar &ds)
messageStream Info
Information stream (uses stdout - output is on the master only)
dictionary coeffs_
Dictionary containing source coefficients.
virtual bool read(const dictionary &dict)
Read source dictionary.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldi)
Add explicit contribution to enthalpy equation.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const Type & lookupObject(const word &name, const bool recursive=false) const
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
word dictName() const
The local dictionary name (final part of scoped name)
static const Enum< thermoMode > thermoModeTypeNames_
Names for thermoMode.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Macros for easy insertion into run-time selection tables.
Mesh data needed to do the Finite Volume discretisation.
const uniformDimensionedVectorField & g
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
errorManip< error > abort(error &err)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual bool read(const dictionary &dict)
Read source dictionary.
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)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
const volScalarField & Cp
A special matrix type and solver, designed for finite volume solutions of scalar equations....
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
const Time & time() const
Return the top-level database.
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
static const gravity & New(const Time &runTime)
Construct on Time.