51 VoFSolidificationMeltingSource,
52 "VoFSolidificationMeltingSource"
60 void Foam::fv::VoFSolidificationMelting::readCoeffs(
const dictionary&
dict)
73 relax_ =
dict.lookupOrDefault<scalar>(
"relax",
dimless, 0.9);
75 q_ =
dict.lookupOrDefault<scalar>(
"q",
dimless, 0.001);
79 Foam::word Foam::fv::VoFSolidificationMelting::alphaSolidName()
const
81 const compressibleTwoPhaseVoFMixture&
thermo
83 mesh().lookupObject<compressibleTwoPhaseVoFMixture>
100 const word& modelType,
133 zeroGradientFvPatchScalarField::typeName
185 const scalar Vc = V[celli];
186 const scalar alphaFluid = 1 - alphaSolid_[celli];
188 const scalar
S = Cu_*
sqr(1 - alphaFluid)/(
pow3(alphaFluid) + q_);
190 Sp[celli] -= Vc*
rho[celli]*
S;
200 <<
" - updating solid phase fraction" <<
endl;
203 alphaSolid_.oldTime();
207 mesh().lookupObject<compressibleTwoPhaseVoFMixture>
222 alphaSolid_[celli] =
min
224 relax_*alphaVoF[celli]*alphaSolidT_->value(TVoF[celli])
225 + (1 - relax_)*alphaSolid_[celli],
230 alphaSolid_.correctBoundaryConditions();
239 zone_.topoChange(map);
254 zone_.distribute(map);
269 zone_.read(coeffs(
dict));
270 readCoeffs(coeffs(
dict));
#define forAll(list, i)
Loop across all elements in list.
Macros for easy insertion into run-time selection tables.
static autoPtr< Function1< scalar > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
static word groupName(Name name, const word &group)
Class to represent a mixture of two rhoFluidThermo-based phases.
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 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.
Solidification and melting model for VoF simulations.
virtual bool movePoints()
Update for mesh motion.
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
virtual void addSup(const volScalarField &alpha, const volScalarField &rho, const volScalarField &he, fvMatrix< scalar > &eqn) const
Add explicit contribution to phase energy equation.
virtual void correct()
Correct the model.
VoFSolidificationMelting(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
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.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
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 handling words, derived from string.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
Calculate the first temporal derivative.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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 > > 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.
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.
const dimensionSet dimless
const dimensionSet dimTemperature
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
VolField< scalar > volScalarField
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.
const unitConversion unitFraction
fluidMulticomponentThermo & thermo