50 solidificationMeltingSource,
51 "solidificationMeltingSource"
67 void Foam::fv::solidificationMelting::readCoeffs(
const dictionary&
dict)
69 Tsol_ =
dict.lookup<scalar>(
"Tsol");
70 Tliq_ =
dict.lookupOrDefault<scalar>(
"Tliq", Tsol_);
71 alpha1e_ =
dict.lookupOrDefault<scalar>(
"alpha1e", 0.0);
72 L_ =
dict.lookup<scalar>(
"L");
74 relax_ =
dict.lookupOrDefault(
"relax", 0.9);
78 rhoRef_ =
dict.lookup<scalar>(
"rhoRef");
79 TName_ =
dict.lookupOrDefault<word>(
"T",
"T");
80 CpName_ =
dict.lookupOrDefault<word>(
"Cp",
"Cp");
81 UName_ =
dict.lookupOrDefault<word>(
"U",
"U");
82 phiName_ =
dict.lookupOrDefault<word>(
"phi",
"phi");
84 Cu_ =
dict.lookupOrDefault<scalar>(
"Cu", 100000);
85 q_ =
dict.lookupOrDefault(
"q", 0.001);
87 beta_ =
dict.lookup<scalar>(
"beta");
91 CpRef_ =
dict.lookup<scalar>(
"CpRef");
96 g_ =
dict.lookup(
"g");
102 Foam::fv::solidificationMelting::Cp()
const
108 const basicThermo&
thermo =
114 case thermoMode::lookup:
116 if (CpName_ ==
"CpRef")
127 extrapolatedCalculatedFvPatchScalarField::typeName
140 <<
"Unhandled thermo mode: " << thermoModeTypeNames_[mode_]
145 return tmp<volScalarField>(
nullptr);
151 if (
mesh().foundObject<uniformDimensionedVectorField>(
"g"))
155 return value.
value();
164 void Foam::fv::solidificationMelting::update
177 <<
" - updating phase indicator" <<
endl;
191 const scalar Tc =
T[celli];
192 const scalar Cpc =
Cp[celli];
193 const scalar alpha1New =
202 + (Tliq_ - Tsol_)*(alpha1_[celli] - alpha1e_)/(1 - alpha1e_)
206 alpha1_[celli] =
max(0,
min(alpha1New, 1));
213 + (Tliq_ - Tsol_)*(alpha1_[celli] - alpha1e_)/(1 - alpha1e_)
217 alpha1_.correctBoundaryConditions();
223 template<
class RhoFieldType>
224 void Foam::fv::solidificationMelting::apply
226 const RhoFieldType&
rho,
227 fvMatrix<scalar>& eqn
232 Info<<
type() <<
": applying source to " << eqn.psi().name() <<
endl;
258 const word& modelType,
275 phiName_(
word::null),
283 this->
name() +
":alpha1",
291 zeroGradientFvPatchScalarField::typeName
294 deltaT_(zone_.nCells(), 0)
313 case thermoMode::lookup:
359 const vector g = this->g();
371 const scalar Vc = V[celli];
372 const scalar alpha1c = alpha1_[celli];
374 const scalar
S = -Cu_*
sqr(1.0 - alpha1c)/(
pow3(alpha1c) + q_);
375 const vector Sb = rhoRef_*g*beta_*deltaT_[i];
406 zone_.topoChange(map);
421 zone_.distribute(map);
429 zone_.read(coeffs(
dict));
430 readCoeffs(coeffs(
dict));
scalar Cp(const scalar p, const scalar T) const
#define forAll(list, i)
Loop across all elements in list.
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
static const List< word > & null()
Return a null List.
Initialise the NamedEnum HashTable from the static list of names.
Enum read(Istream &) const
Read a word from Istream and return the corresponding.
label timeIndex() const
Return current time index.
Base-class for fluid and solid thermodynamic properties.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Mesh data needed to do the Finite Volume discretisation.
const Time & time() const
Return the top-level database.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Finite volume model abstract base class.
static const dictionary & coeffs(const word &modelType, const dictionary &)
Return the coefficients sub-dictionary for a given model type.
virtual bool read(const dictionary &dict)
Read source dictionary.
const fvMesh & mesh() const
Return const access to the mesh database.
This source is designed to model the effect of solidification and melting processes,...
virtual bool movePoints()
Update for mesh motion.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
virtual void addSup(const volScalarField &he, fvMatrix< scalar > &eqn) const
Add explicit contribution to enthalpy equation.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual bool read(const dictionary &dict)
Read source dictionary.
static const NamedEnum< thermoMode, 2 > thermoModeTypeNames_
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
solidificationMelting(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
bool foundObject(const word &name) const
Is the named Type in registry.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
A class for handling words, derived from string.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Calculate the first temporal derivative.
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
addBackwardCompatibleToRunTimeSelectionTable(fvConstraint, fixedTemperature, dictionary, fixedTemperatureConstraint, "fixedTemperatureConstraint")
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
List< word > wordList
A List of words.
List< label > labelList
A List of labels.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const dimensionSet dimEnergy
Ostream & endl(Ostream &os)
Add newline and flush stream.
errorManip< error > abort(error &err)
const dimensionSet dimless
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimTemperature
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
UniformDimensionedField< vector > uniformDimensionedVectorField
VolField< scalar > volScalarField
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void pow3(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimMass
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
fluidMulticomponentThermo & thermo