33 namespace functionObjects
53 >::names[] = {
"integerMoment",
"mean",
"variance",
"stdDev"};
72 >::names[] = {
"volume",
"area",
"diameter"};
93 "numberConcentration",
94 "volumeConcentration",
111 const char* NamedEnum
115 >::names[] = {
"arithmetic",
"geometric",
"notApplicable"};
130 Foam::functionObjects::populationBalanceMoments::coordinateTypeSymbolicName()
134 switch (coordinateType_)
138 coordinateTypeSymbolicName =
"v";
144 coordinateTypeSymbolicName =
"a";
150 coordinateTypeSymbolicName =
"d";
156 return coordinateTypeSymbolicName;
161 Foam::functionObjects::populationBalanceMoments::weightTypeSymbolicName()
167 case weightType::numberConcentration:
169 weightTypeSymbolicName =
"N";
173 case weightType::volumeConcentration:
175 weightTypeSymbolicName =
"V";
179 case weightType::areaConcentration:
181 weightTypeSymbolicName =
"A";
187 return weightTypeSymbolicName;
191 Foam::word Foam::functionObjects::populationBalanceMoments::defaultFldName()
195 meanType_ == meanType::geometric
196 ? word(meanTypeNames_[meanType_]).capitalise()
207 + word(momentTypeNames_[momentType_]).capitalise()
209 + weightTypeSymbolicName()
211 + coordinateTypeSymbolicName()
220 Foam::functionObjects::populationBalanceMoments::integerMomentFldName()
227 word(momentTypeNames_[momentType_])
230 + weightTypeSymbolicName()
232 + coordinateTypeSymbolicName()
240 void Foam::functionObjects::populationBalanceMoments::setDimensions
248 case momentType::integerMoment:
250 switch (coordinateType_)
252 case coordinateType::volume:
254 fld.dimensions().reset
261 case coordinateType::area:
263 fld.dimensions().reset
270 case coordinateType::diameter:
272 fld.dimensions().reset
283 case weightType::volumeConcentration:
289 case weightType::areaConcentration:
303 case momentType::mean:
305 switch (coordinateType_)
307 case coordinateType::volume:
313 case coordinateType::area:
319 case coordinateType::diameter:
329 case momentType::variance:
331 switch (coordinateType_)
333 case coordinateType::volume:
339 case coordinateType::area:
345 case coordinateType::diameter:
353 if (meanType_ == meanType::geometric)
360 case momentType::stdDev:
362 switch (coordinateType_)
364 case coordinateType::volume:
370 case coordinateType::area:
376 case coordinateType::diameter:
384 if (meanType_ == meanType::geometric)
396 Foam::functionObjects::populationBalanceMoments::totalConcentration()
398 tmp<volScalarField> tTotalConcentration
402 "totalConcentration",
412 case weightType::volumeConcentration:
414 totalConcentration.dimensions().reset
416 totalConcentration.dimensions()*
dimVolume
421 case weightType::areaConcentration:
423 totalConcentration.dimensions().reset
425 totalConcentration.dimensions()*
dimArea
436 forAll(popBal_.sizeGroups(), i)
442 case weightType::numberConcentration:
444 totalConcentration += fi*fi.
phase()/fi.
x();
448 case weightType::volumeConcentration:
450 totalConcentration += fi*fi.
phase();
454 case weightType::areaConcentration:
456 totalConcentration += fi.
a()*fi*fi.
phase()/fi.
x();
463 return tTotalConcentration;
468 Foam::functionObjects::populationBalanceMoments::mean()
470 tmp<volScalarField> tMean
482 setDimensions(mean, momentType::mean);
486 forAll(popBal_.sizeGroups(), i)
494 case weightType::volumeConcentration:
496 concentration *= fi.
x();
500 case weightType::areaConcentration:
502 concentration *= fi.
a();
514 case meanType::geometric:
516 mean.dimensions().reset(
dimless);
518 switch (coordinateType_)
520 case coordinateType::volume:
526 *concentration/totalConcentration;
530 case coordinateType::area:
536 *concentration/totalConcentration;
540 case coordinateType::diameter:
546 *concentration/totalConcentration;
556 switch (coordinateType_)
558 case coordinateType::volume:
560 mean += fi.
x()*concentration/totalConcentration;
564 case coordinateType::area:
566 mean += fi.
a()*concentration/totalConcentration;
570 case coordinateType::diameter:
572 mean += fi.
d()*concentration/totalConcentration;
583 if (meanType_ == meanType::geometric)
587 setDimensions(mean, momentType::mean);
595 Foam::functionObjects::populationBalanceMoments::variance()
597 tmp<volScalarField> tVariance
609 setDimensions(variance, momentType::variance);
614 forAll(popBal_.sizeGroups(), i)
622 case weightType::volumeConcentration:
624 concentration *= fi.
x();
628 case weightType::areaConcentration:
630 concentration *= fi.
a();
642 case meanType::geometric:
644 switch (coordinateType_)
646 case coordinateType::volume:
650 *concentration/totalConcentration;
654 case coordinateType::area:
658 *concentration/totalConcentration;
662 case coordinateType::diameter:
666 *concentration/totalConcentration;
676 switch (coordinateType_)
678 case coordinateType::volume:
681 sqr(fi.
x() - mean)*concentration/totalConcentration;
685 case coordinateType::area:
688 sqr(fi.
a() - mean)*concentration/totalConcentration;
692 case coordinateType::diameter:
695 sqr(fi.
d() - mean)*concentration/totalConcentration;
711 Foam::functionObjects::populationBalanceMoments::stdDev()
715 case meanType::geometric:
717 return exp(
sqrt(this->variance()));
721 return sqrt(this->variance());
739 obr_.lookupObject<
Foam::diameterModels::populationBalanceModel>
741 dict.lookup(
"populationBalance")
744 momentType_(momentTypeNames_.
read(
dict.lookup(
"momentType"))),
745 coordinateType_(coordinateTypeNames_.
read(
dict.lookup(
"coordinateType"))),
749 ? weightTypeNames_.
read(
dict.lookup(
"weightType"))
775 case momentType::integerMoment:
777 order_ =
dict.lookup<
int>(
"order");
784 dict.found(
"meanType")
785 ? meanTypeNames_.read(
dict.lookup(
"meanType"))
786 : meanType::arithmetic;
794 case momentType::integerMoment:
802 this->integerMomentFldName(),
815 setDimensions(integerMoment, momentType::integerMoment);
819 case momentType::mean:
827 this->defaultFldName(),
839 case momentType::variance:
847 this->defaultFldName(),
859 case momentType::stdDev:
867 this->defaultFldName(),
889 case momentType::integerMoment:
893 integerMoment =
Zero;
895 forAll(popBal_.sizeGroups(), i)
898 popBal_.sizeGroups()[i];
904 case weightType::volumeConcentration:
906 concentration *= fi.
x();
910 case weightType::areaConcentration:
912 concentration *= fi.
a();
922 switch (coordinateType_)
924 case coordinateType::volume:
927 pow(fi.
x(), order_)*concentration;
931 case coordinateType::area:
934 pow(fi.
a(), order_)*concentration;
938 case coordinateType::diameter:
941 pow(fi.
d(), order_)*concentration;
950 case momentType::mean:
952 fldPtr_() = this->mean();
956 case momentType::variance:
958 fldPtr_() = this->variance();
962 case momentType::stdDev:
964 fldPtr_() =
sqrt(this->variance());
976 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, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, 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.
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 keyword definitions, which are a keyword followed by any number of values (e....
Abstract base-class for Time/database functionObjects.
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
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.
virtual bool read(const dictionary &)
Read optional controls.
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(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(), 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 &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
word name(const bool)
Return a word representation of a bool.
const dimensionSet dimless
const dimensionSet dimLength
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimVolume
VolField< scalar > volScalarField
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
const dimensionSet dimArea
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.