39namespace diameterModels
57 tmp<volScalarField> tInvDsm
71 const sizeGroup& fi = sizeGroups_[i];
81Foam::diameterModels::velocityGroup::fSum()
const
83 tmp<volScalarField> tsumSizeGroups
97 sumSizeGroups += sizeGroups_[i];
100 return tsumSizeGroups;
104void Foam::diameterModels::velocityGroup::renormalize()
106 Info<<
"Renormalizing sizeGroups for velocityGroup "
113 sizeGroups_[i] *=
pos(sizeGroups_[i]);
118 sizeGroups_[i] /= fSum_;
124Foam::diameterModels::velocityGroup::mvconvection()
const
132 phase_.alphaRhoPhi(),
133 phase_.mesh().divScheme
135 "div(" + phase_.alphaRhoPhi()().name() +
",f)"
153 popBalName_(diameterProperties.
lookup(
"populationBalance")),
174 formFactor_(
"formFactor",
dimless, diameterProperties),
177 diameterProperties.
lookup(
"sizeGroups"),
227 "renormalizeAtRestart",
246 ||
mag(1 -
max(fSum_).value()) >= 1
e-5
247 ||
mag(1 -
min(fSum_).value()) >= 1
e-5
251 <<
" Initial values of the sizeGroups belonging to velocityGroup "
253 <<
" must add to" <<
nl <<
" unity. This condition might be"
254 <<
" violated due to wrong entries in the" <<
nl
255 <<
" velocityGroupCoeffs subdictionary or bad initial conditions in"
256 <<
" the startTime" <<
nl
257 <<
" directory. The sizeGroups can be renormalized at every"
258 <<
" timestep or at restart" <<
nl
259 <<
" only by setting the corresponding switch renormalize or"
260 <<
" renormalizeAtRestart" <<
nl
261 <<
" in the fvSolution subdictionary " << popBalName_ <<
"."
262 <<
" Note that boundary conditions are not" <<
nl <<
"renormalized."
268 fields_.add(sizeGroups_[i]);
286 mvConvection_ = mvconvection();
294 Info<< this->
phase().
name() <<
" Sauter mean diameter, min, max = "
295 << d_.weightedAverage(d_.mesh().V()).value()
296 <<
' ' <<
min(d_).value()
297 <<
' ' <<
max(d_).value()
302 Info<< phase_.name() <<
" sizeGroups-sum volume fraction, min, max = "
303 << fSum_.weightedAverage(phase_.mesh().V()).value()
304 <<
' ' <<
min(fSum_).value()
305 <<
' ' <<
max(fSum_).value()
310 phase_.mesh().solverDict(popBalName_).getOrDefault<
Switch>
tmp< fv::convectionScheme< scalar > > mvConvection(fv::convectionScheme< scalar >::New(mesh, fields, phi, mesh.divScheme("div(phi,Yi_h)")))
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &weightField) const
Calculate and return weighted average.
const Mesh & mesh() const
Return mesh.
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
virtual bool read()
Re-read model coefficients if they have changed.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Abstract base-class for dispersed-phase particle diameter models.
const phaseModel & phase() const
Return the phase.
const phaseModel & phase_
This class represents a single sizeGroup belonging to a velocityGroup. The main property of a sizeGro...
This diameterModel is intended for use with a populationBalanceModel in order to simulate polydispers...
void postSolve()
Corrections after populationBalanceModel::solve()
virtual tmp< volScalarField > d() const
Return the Sauter-mean diameter.
virtual ~velocityGroup()
Destructor.
void preSolve()
Corrections before populationBalanceModel::solve()
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
const word & name() const
Helper class to manage multi-specie phase properties.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
const word & name() const
Lookup type of boundary radiation properties.
const dictionary & solverDict(const word &name) const
Return the solver controls dictionary for the given field.
A class for managing temporary objects.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar pos(const dimensionedScalar &ds)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
const dimensionSet dimless
Dimensionless.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
constexpr char nl
The newline '\n' character (0x0a)
#define forAll(list, i)
Loop across all elements in list.