Go to the documentation of this file.
30 #include "phaseSystem.H"
37 namespace diameterModels
39 namespace nucleationModels
65 refCast<const velocityGroup>
83 turbulenceModel::propertiesName,
84 popBal_.continuousPhase().name()
96 const volScalarField::Boundary& alphatBf = talphat().boundaryField();
99 alphatWallBoilingWallFunction;
105 isA<alphatWallBoilingWallFunction>(alphatBf[patchi])
108 const alphatWallBoilingWallFunction& alphatw =
109 refCast<const alphatWallBoilingWallFunction>(alphatBf[patchi]);
116 <<
"Minimum departure diameter " <<
min(dDep)
117 <<
" m outside of range ["
118 << velGroup_.
sizeGroups().first().d().value() <<
", "
119 << velGroup_.
sizeGroups().last().d().value() <<
"] m"
120 <<
" at patch " << alphatw.patch().name()
122 <<
" The nucleation rate in populationBalance "
124 <<
" Adjust discretization over property space to"
125 <<
" suppress this warning."
128 else if (
max(dDep) > velGroup_.
sizeGroups().last().d().value())
131 <<
"Maximum departure diameter " <<
max(dDep)
132 <<
" m outside of range ["
133 << velGroup_.
sizeGroups().first().d().value() <<
", "
134 << velGroup_.
sizeGroups().last().d().value() <<
"] m"
135 <<
" at patch " << alphatw.patch().name()
137 <<
" The nucleation rate in populationBalance "
139 <<
" Adjust discretization over property space to"
140 <<
" suppress this warning."
155 const sizeGroup& fi = popBal_.sizeGroups()[i];
159 const volScalarField::Boundary& alphatBf = talphat().boundaryField();
162 alphatWallBoilingWallFunction;
168 isA<alphatWallBoilingWallFunction>(alphatBf[patchi])
171 const alphatWallBoilingWallFunction& alphatw =
172 refCast<const alphatWallBoilingWallFunction>(alphatBf[patchi]);
183 if (dmdt[facei] > SMALL)
185 const label faceCelli =
faceCells[facei];
187 nucleationRate[faceCelli] +=
191 velGroup_.formFactor()*
pow3(dDep[facei]*unitLength)
193 *dmdt[facei]/
rho[faceCelli]/fi.
x().
value();
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
const phaseModel & phase() const
Return const-reference to the phase.
const word & name() const
Return name.
A class for handling words, derived from Foam::string.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
virtual void correct()
Correct diameter independent expressions.
Ostream & endl(Ostream &os)
Add newline and flush stream.
virtual void addToNucleationRate(volScalarField &nucleationRate, const label i)
Add to nucleationRate.
const Type & value() const
Return const reference to value.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
A thermal wall function for simulation of boiling wall.
Base class for nucleation models.
#define forAll(list, i)
Loop across all elements in list.
const fvMesh & mesh() const
Return reference to the mesh.
wallBoiling(const populationBalanceModel &popBal, const dictionary &dict)
dimensionedScalar pow3(const dimensionedScalar &ds)
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.
This class represents a single sizeGroup belonging to a velocityGroup. The main property of a sizeGro...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Class that solves the univariate population balance equation by means of a class method (also called ...
Macros for easy insertion into run-time selection tables.
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
virtual tmp< volScalarField > alphat() const
Return the turbulent thermal diffusivity for enthalpy [kg/m/s].
defineTypeNameAndDebug(constantNucleation, 0)
addToRunTimeSelectionTable(nucleationModel, constantNucleation, dictionary)
const PtrList< sizeGroup > & sizeGroups() const
Return sizeGroups belonging to this velocityGroup.
const dimensionedScalar & x() const
Return representative volume of the sizeGroup.
Smooth ATC in cells next to a set of patches supplied by type.
const dimensionedScalar & rho() const
Return const-access to phase1 density.