33 namespace functionObjects
53 "numberConcentration",
55 "volumeConcentration",
78 "projectedAreaDiameter"
98 "numberConcentration",
99 "volumeConcentration",
119 Foam::functionObjects::populationBalanceSizeDistribution::
120 functionTypeSymbolicName()
124 switch (functionType_)
126 case functionType::numberConcentration:
128 functionTypeSymbolicName =
"N";
132 case functionType::numberDensity:
134 functionTypeSymbolicName =
"n";
138 case functionType::volumeConcentration:
140 functionTypeSymbolicName =
"V";
144 case functionType::volumeDensity:
146 functionTypeSymbolicName =
"v";
150 case functionType::areaConcentration:
152 functionTypeSymbolicName =
"A";
156 case functionType::areaDensity:
158 functionTypeSymbolicName =
"a";
164 return functionTypeSymbolicName;
169 Foam::functionObjects::populationBalanceSizeDistribution::
170 coordinateTypeSymbolicName
172 const coordinateType& cType
179 case coordinateType::volume:
181 coordinateTypeSymbolicName =
"v";
185 case coordinateType::area:
187 coordinateTypeSymbolicName =
"a";
191 case coordinateType::diameter:
193 coordinateTypeSymbolicName =
"d";
197 case coordinateType::projectedAreaDiameter:
199 coordinateTypeSymbolicName =
"dPa";
205 return coordinateTypeSymbolicName;
210 Foam::functionObjects::populationBalanceSizeDistribution::filterField
227 Foam::functionObjects::populationBalanceSizeDistribution::averageCoordinateValue
230 const coordinateType& cType
233 scalar averageCoordinateValue(
Zero);
237 case coordinateType::volume:
239 averageCoordinateValue = fi.
x().
value();
243 case coordinateType::area:
245 averageCoordinateValue =
246 weightedAverage(fi.
a(), fi);
250 case coordinateType::diameter:
252 averageCoordinateValue =
253 weightedAverage(fi.
d(), fi);
257 case coordinateType::projectedAreaDiameter:
259 averageCoordinateValue =
260 weightedAverage(
sqrt(fi.
a()/
pi), fi);
266 return averageCoordinateValue;
271 Foam::functionObjects::populationBalanceSizeDistribution::weightedAverage
277 scalar weightedAverage(
Zero);
281 case weightType::numberConcentration:
288 gSum(filterField(mesh_.V()*
fld))/this->V();
298 case weightType::volumeConcentration:
305 gSum(filterField(mesh_.V()*
fld))/this->V();
315 case weightType::areaConcentration:
322 gSum(filterField(mesh_.V()*
fld))/this->V();
332 case weightType::cellVolume:
335 gSum(filterField(mesh_.V()*
fld))/this->V();
341 return weightedAverage;
361 obr_.lookupObject<
Foam::diameterModels::populationBalanceModel>
363 dict.lookup(
"populationBalance")
366 functionType_(functionTypeNames_.
read(
dict.lookup(
"functionType"))),
367 coordinateType_(coordinateTypeNames_.
read(
dict.lookup(
"coordinateType"))),
370 dict.lookupOrDefault<
Switch>(
"allCoordinates", false)
372 normalise_(
dict.lookupOrDefault<
Switch>(
"normalise", false)),
375 dict.lookupOrDefaultBackwardsCompatible<
Switch>
377 {
"logTransform",
"geometric"},
383 dict.found(
"weightType")
384 ? weightTypeNames_.read(
dict.lookup(
"weightType"))
385 : weightType::numberConcentration
387 formatterPtr_(
nullptr)
428 popBal_.sizeGroups();
438 coordinateValues[i] = averageCoordinateValue(fi, coordinateType_);
443 functionType_ == functionType::numberDensity
444 || functionType_ == functionType::volumeDensity
445 || functionType_ == functionType::areaDensity
448 boundaryValues.
first() = coordinateValues.
first();
449 boundaryValues.
last() = coordinateValues.
last();
451 for (
label i = 1; i < boundaryValues.
size() - 1; i++)
454 0.5*(coordinateValues[i] + coordinateValues[i-1]);
459 boundaryValues =
Foam::log(boundaryValues);
463 switch (functionType_)
465 case functionType::numberConcentration:
472 gSum(filterField(mesh_.V()*fi*fi.
phase()/fi.
x()))/this->V();
475 if (normalise_ &&
sum(resultValues) != 0)
477 resultValues /=
sum(resultValues);
482 case functionType::numberDensity:
489 gSum(filterField(mesh_.V()*fi*fi.
phase()/fi.
x()))/this->V();
492 if (normalise_ &&
sum(resultValues) != 0)
494 resultValues /=
sum(resultValues);
499 resultValues[i] /= (boundaryValues[i+1] - boundaryValues[i]);
504 case functionType::volumeConcentration:
511 gSum(filterField(mesh_.V()*fi*fi.
phase()))/this->V();
514 if (normalise_ &&
sum(resultValues) != 0)
516 resultValues /=
sum(resultValues);
521 case functionType::volumeDensity:
528 gSum(filterField(mesh_.V()*fi*fi.
phase()))/this->V();
531 if (normalise_ &&
sum(resultValues) != 0)
533 resultValues /=
sum(resultValues);
538 resultValues[i] /= (boundaryValues[i+1] - boundaryValues[i]);
543 case functionType::areaConcentration:
552 filterField(mesh_.V()*fi.
a().ref()*fi*fi.
phase()/fi.
x())
557 if (normalise_ &&
sum(resultValues) != 0)
559 resultValues /=
sum(resultValues);
564 case functionType::areaDensity:
573 filterField(mesh_.V()*fi.
a().ref()*fi*fi.
phase()/fi.
x())
578 if (normalise_ &&
sum(resultValues) != 0)
580 resultValues /=
sum(resultValues);
585 resultValues[i] /= (boundaryValues[i+1] - boundaryValues[i]);
595 wordList otherCoordinateSymbolicNames(coordinateTypeNames_.size());
603 otherCoordinateSymbolicNames[cType] =
604 coordinateTypeSymbolicName(cType);
606 otherCoordinateValues.
set
616 otherCoordinateValues[cType][i] =
617 averageCoordinateValue(fi, cType);
630 coordinateTypeSymbolicName(coordinateType_),
633 functionTypeSymbolicName(),
635 otherCoordinateSymbolicNames,
636 otherCoordinateValues
651 coordinateTypeSymbolicName(coordinateType_),
654 functionTypeSymbolicName(),
#define forAll(list, i)
Loop across all elements in list.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Macros for easy insertion into run-time selection tables.
void size(const label)
Override size to be inconsistent with allocated storage.
Initialise the NamedEnum HashTable from the static list of names.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
bool set(const label) const
Is element set.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
T & first()
Return the first element of the list.
T & last()
Return the last element of the list.
static bool master(const label communicator=0)
Am I the master process.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
label size() const
Return the number of elements in the UPtrList.
Holds list of sampling positions.
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....
const Type & value() const
Return const reference to value.
Abstract base-class for Time/database functionObjects.
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Writes out the size distribution determined by a population balance model, either for the entire doma...
static const NamedEnum< coordinateType, 4 > coordinateTypeNames_
Coordinate type names.
weightType
Enumeration for the weight types.
populationBalanceSizeDistribution(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
functionType
Function type enumeration.
coordinateType
Coordinate type enumeration.
virtual ~populationBalanceSizeDistribution()
Destructor.
static const NamedEnum< weightType, 4 > weightTypeNames_
Names of the weight types.
static const NamedEnum< functionType, 6 > functionTypeNames_
Function type names.
virtual bool execute()
Execute, currently does nothing.
virtual bool write()
Calculate and write the size distribution.
virtual bool read(const dictionary &)
Read the populationBalanceSizeDistribution data.
virtual bool read(const dictionary &)
Read optional controls.
General run-time selected cell set selection class for fvMesh.
static autoPtr< setWriter > New(const word &writeType, const IOstream::streamFormat writeFormat=IOstream::ASCII, const IOstream::compressionType writeCompression=IOstream::UNCOMPRESSED)
Select given write options.
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))
#define Log
Report write to Foam::Info if the local log switch is true.
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
Type gSum(const FieldField< Field, Type > &f)
bool read(const char *, int32_t &)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
word name(const bool)
Return a word representation of a bool.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar sqrt(const dimensionedScalar &ds)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.