33 namespace functionObjects
51 {
"integerMoment",
"mean",
"variance",
"stdDev"};
59 {
"volume",
"area",
"diameter"};
68 "numberConcentration",
69 "volumeConcentration",
79 {
"arithmetic",
"geometric",
"notApplicable"};
85 Foam::functionObjects::populationBalanceMoments::coordinateTypeSymbolicName()
89 switch (coordinateType_)
93 coordinateTypeSymbolicName =
"v";
99 coordinateTypeSymbolicName =
"a";
105 coordinateTypeSymbolicName =
"d";
111 return coordinateTypeSymbolicName;
116 Foam::functionObjects::populationBalanceMoments::weightTypeSymbolicName()
122 case weightType::numberConcentration:
124 weightTypeSymbolicName =
"N";
128 case weightType::volumeConcentration:
130 weightTypeSymbolicName =
"V";
134 case weightType::areaConcentration:
136 weightTypeSymbolicName =
"A";
142 return weightTypeSymbolicName;
146 Foam::word Foam::functionObjects::populationBalanceMoments::defaultFldName()
150 meanType_ == meanType::geometric
151 ? word(meanTypeNames_[meanType_]).capitalise()
162 + word(momentTypeNames_[momentType_]).capitalise()
164 + weightTypeSymbolicName()
166 + coordinateTypeSymbolicName()
175 Foam::functionObjects::populationBalanceMoments::integerMomentFldName()
182 word(momentTypeNames_[momentType_])
185 + weightTypeSymbolicName()
187 + coordinateTypeSymbolicName()
195 void Foam::functionObjects::populationBalanceMoments::setDimensions
203 case momentType::integerMoment:
205 switch (coordinateType_)
207 case coordinateType::volume:
209 fld.dimensions().reset
216 case coordinateType::area:
218 fld.dimensions().reset
225 case coordinateType::diameter:
227 fld.dimensions().reset
238 case weightType::volumeConcentration:
244 case weightType::areaConcentration:
258 case momentType::mean:
260 switch (coordinateType_)
262 case coordinateType::volume:
268 case coordinateType::area:
274 case coordinateType::diameter:
284 case momentType::variance:
286 switch (coordinateType_)
288 case coordinateType::volume:
294 case coordinateType::area:
300 case coordinateType::diameter:
308 if (meanType_ == meanType::geometric)
315 case momentType::stdDev:
317 switch (coordinateType_)
319 case coordinateType::volume:
325 case coordinateType::area:
331 case coordinateType::diameter:
339 if (meanType_ == meanType::geometric)
351 Foam::functionObjects::populationBalanceMoments::totalConcentration
353 const diameterModels::populationBalanceModel& popBal
356 tmp<volScalarField> tTotalConcentration
360 "totalConcentration",
370 case weightType::volumeConcentration:
372 totalConcentration.dimensions().reset
374 totalConcentration.dimensions()*
dimVolume
379 case weightType::areaConcentration:
381 totalConcentration.dimensions().reset
383 totalConcentration.dimensions()*
dimArea
394 forAll(popBal.sizeGroups(), i)
400 case weightType::numberConcentration:
402 totalConcentration += fi*fi.
phase()/fi.
x();
406 case weightType::volumeConcentration:
408 totalConcentration += fi*fi.
phase();
412 case weightType::areaConcentration:
414 totalConcentration += fi.
a()*fi*fi.
phase()/fi.
x();
421 return tTotalConcentration;
426 Foam::functionObjects::populationBalanceMoments::mean
428 const diameterModels::populationBalanceModel& popBal
431 tmp<volScalarField> tMean
443 setDimensions(mean, momentType::mean);
445 volScalarField totalConcentration(this->totalConcentration(popBal));
447 forAll(popBal.sizeGroups(), i)
455 case weightType::volumeConcentration:
457 concentration *= fi.
x();
461 case weightType::areaConcentration:
463 concentration *= fi.
a();
475 case meanType::geometric:
477 mean.dimensions().reset(
dimless);
479 switch (coordinateType_)
481 case coordinateType::volume:
487 *concentration/totalConcentration;
491 case coordinateType::area:
497 *concentration/totalConcentration;
501 case coordinateType::diameter:
507 *concentration/totalConcentration;
517 switch (coordinateType_)
519 case coordinateType::volume:
521 mean += fi.
x()*concentration/totalConcentration;
525 case coordinateType::area:
527 mean += fi.
a()*concentration/totalConcentration;
531 case coordinateType::diameter:
533 mean += fi.
d()*concentration/totalConcentration;
544 if (meanType_ == meanType::geometric)
548 setDimensions(mean, momentType::mean);
556 Foam::functionObjects::populationBalanceMoments::variance
558 const diameterModels::populationBalanceModel& popBal
561 tmp<volScalarField> tVariance
573 setDimensions(variance, momentType::variance);
575 volScalarField totalConcentration(this->totalConcentration(popBal));
578 forAll(popBal.sizeGroups(), i)
586 case weightType::volumeConcentration:
588 concentration *= fi.
x();
592 case weightType::areaConcentration:
594 concentration *= fi.
a();
606 case meanType::geometric:
608 switch (coordinateType_)
610 case coordinateType::volume:
614 *concentration/totalConcentration;
618 case coordinateType::area:
622 *concentration/totalConcentration;
626 case coordinateType::diameter:
630 *concentration/totalConcentration;
640 switch (coordinateType_)
642 case coordinateType::volume:
645 sqr(fi.
x() - mean)*concentration/totalConcentration;
649 case coordinateType::area:
652 sqr(fi.
a() - mean)*concentration/totalConcentration;
656 case coordinateType::diameter:
659 sqr(fi.
d() - mean)*concentration/totalConcentration;
675 Foam::functionObjects::populationBalanceMoments::stdDev
677 const diameterModels::populationBalanceModel& popBal
682 case meanType::geometric:
684 return exp(
sqrt(this->variance(popBal)));
688 return sqrt(this->variance(popBal));
704 popBalName_(
dict.lookup(
"populationBalance")),
705 momentType_(momentTypeNames_.
read(
dict.lookup(
"momentType"))),
706 coordinateType_(coordinateTypeNames_.
read(
dict.lookup(
"coordinateType"))),
710 ? weightTypeNames_.
read(
dict.lookup(
"weightType"))
736 case momentType::integerMoment:
738 order_ =
dict.lookup<
int>(
"order");
745 dict.found(
"meanType")
746 ? meanTypeNames_.read(
dict.lookup(
"meanType"))
747 : meanType::arithmetic;
755 case momentType::integerMoment:
763 this->integerMomentFldName(),
776 setDimensions(integerMoment, momentType::integerMoment);
780 case momentType::mean:
788 this->defaultFldName(),
799 setDimensions(fldPtr_(), momentType::mean);
803 case momentType::variance:
811 this->defaultFldName(),
822 setDimensions(fldPtr_(), momentType::variance);
826 case momentType::stdDev:
834 this->defaultFldName(),
845 setDimensions(fldPtr_(), momentType::stdDev);
865 case momentType::integerMoment:
869 integerMoment =
Zero;
880 case weightType::volumeConcentration:
882 concentration *= fi.
x();
886 case weightType::areaConcentration:
888 concentration *= fi.
a();
898 switch (coordinateType_)
900 case coordinateType::volume:
903 pow(fi.
x(), order_)*concentration;
907 case coordinateType::area:
910 pow(fi.
a(), order_)*concentration;
914 case coordinateType::diameter:
917 pow(fi.
d(), order_)*concentration;
926 case momentType::mean:
928 fldPtr_() = this->mean(popBal);
932 case momentType::variance:
934 fldPtr_() = this->variance(popBal);
938 case momentType::stdDev:
940 fldPtr_() =
sqrt(this->variance(popBal));
952 writeObject(fldPtr_->name());
#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 word groupName(Name name, const word &group)
Initialise the NamedEnum HashTable from the static list of names.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Model for tracking the evolution of a dispersed phase size distribution due to coalescence (synonymou...
const UPtrList< sizeGroup > & sizeGroups() const
Return the size groups belonging to this populationBalance.
Single size class fraction field representing a fixed particle volume as defined by the user through ...
const dimensionedScalar & x() const
Return representative volume of the sizeGroup.
const tmp< volScalarField > a() const
Return representative surface area of the sizeGroup.
const phaseModel & phase() const
Return const-reference to the phase.
const tmp< volScalarField > d() const
Return representative diameter of the sizeGroup.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Abstract base-class for Time/database functionObjects.
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
virtual bool read(const dictionary &)
Read optional controls.
Calculates and writes out integral (integer moments) or mean properties (mean, variance,...
static const NamedEnum< coordinateType, 3 > coordinateTypeNames_
Names of the coordinate types.
weightType
Enumeration for the weight types.
static const NamedEnum< meanType, 3 > meanTypeNames_
Names of the mean types.
static const NamedEnum< momentType, 4 > momentTypeNames_
Names of the moment types.
static const NamedEnum< weightType, 3 > weightTypeNames_
Names of the weight types.
coordinateType
Enumeration for the coordinate types.
populationBalanceMoments(const word &name, const Time &runTime, const dictionary &)
Construct from Time and dictionary.
virtual ~populationBalanceMoments()
Destructor.
meanType
Enumeration for the mean types.
virtual bool execute()
Calculate the moment fields.
virtual bool write()
Write the moment fields.
virtual bool read(const dictionary &)
Read the data.
momentType
Enumeration for the moment types.
A class for managing temporary objects.
A class for handling words, derived from string.
static const word null
An empty word.
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
dimensionedScalar exp(const dimensionedScalar &ds)
bool read(const char *, int32_t &)
const dimensionSet dimless
const dimensionSet dimLength
dimensionedScalar log(const dimensionedScalar &ds)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimVolume
VolField< scalar > volScalarField
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 dimArea
void inv(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.