43 VoFSolidificationMeltingSource,
52 void Foam::fv::VoFSolidificationMeltingSource::readCoeffs()
56 relax_ = coeffs().lookupOrDefault<scalar>(
"relax", 0.9);
57 Cu_ = coeffs().lookupOrDefault<scalar>(
"Cu", 100000);
58 q_ = coeffs().lookupOrDefault<scalar>(
"q", 0.001);
62 Foam::word Foam::fv::VoFSolidificationMeltingSource::alphaSolidName()
const 64 const compressibleTwoPhaseMixture&
thermo 66 mesh().lookupObject<compressibleTwoPhaseMixture>
83 const word& modelType,
84 const dictionary& dict,
88 fvModel(name, modelType, dict, mesh),
102 IOobject::READ_IF_PRESENT,
107 zeroGradientFvPatchScalarField::typeName
125 fvMatrix<scalar>& eqn,
126 const word& fieldName
131 Info<<
type() <<
": applying source to " << eqn.psi().name() <<
endl;
134 const compressibleTwoPhaseMixture&
thermo 136 mesh().lookupObject<compressibleTwoPhaseMixture>
146 eqn += L_/CpVoF*(
fvc::ddt(rho, alphaSolid_));
150 eqn += L_*(
fvc::ddt(rho, alphaSolid_));
158 fvMatrix<vector>& eqn,
159 const word& fieldName
164 Info<<
type() <<
": applying source to " << eqn.psi().name() <<
endl;
174 const label celli = cells[i];
175 const scalar Vc = V[celli];
176 const scalar alphaFluid = 1 - alphaSolid_[celli];
178 const scalar
S = Cu_*
sqr(1 - alphaFluid)/(
pow3(alphaFluid) + q_);
180 Sp[celli] -= Vc*rho[celli]*
S;
190 <<
" - updating solid phase fraction" <<
endl;
193 alphaSolid_.oldTime();
195 const compressibleTwoPhaseMixture&
thermo 197 mesh().lookupObject<compressibleTwoPhaseMixture>
211 const label celli = cells[i];
213 alphaSolid_[celli] =
min 215 relax_*alphaVoF[celli]*alphaSolidT_->value(TVoF[celli])
216 + (1 - relax_)*alphaSolid_[celli],
221 alphaSolid_.correctBoundaryConditions();
227 const polyTopoChangeMap& map
230 set_.topoChange(map);
242 const polyDistributionMap& map
245 set_.distribute(map);
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
#define forAll(list, i)
Loop across all elements in list.
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
fluidReactionThermo & thermo
const word & name() const
Return const access to the source name.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual bool read(const dictionary &dict)
Read source dictionary.
Ostream & endl(Ostream &os)
Add newline and flush stream.
virtual bool movePoints()
Update for mesh motion.
virtual void correct()
Correct the model.
const dimensionSet dimless
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const fvMesh & mesh() const
Return const access to the mesh database.
Macros for easy insertion into run-time selection tables.
VoFSolidificationMeltingSource(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Calculate the first temporal derivative.
A class for handling words, derived from string.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const GeometricField< Type, fvPatchField, volMesh > &)
static word groupName(Name name, const word &group)
const dictionary & coeffs() const
Return dictionary.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
List< label > labelList
A List of labels.
virtual bool read(const dictionary &dict)
Read source dictionary.
const dimensionSet dimEnergy
const dimensionSet dimMass
word name(const complex &)
Return a string representation of a complex.
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
List< word > wordList
A List of words.
dimensionedScalar pow3(const dimensionedScalar &ds)
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.
virtual void addSup(const volScalarField &rho, fvMatrix< scalar > &eqn, const word &fieldName) const
Add explicit contribution to compressible enthalpy equation.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
const dimensionSet dimTemperature
static autoPtr< Function1< scalar > > New(const word &name, const dictionary &dict)
Selector.