39template<
class CloudType>
48 generationSetName_(this->coeffDict().
lookup(
"generationCellSet")),
49 inflationSetName_(this->coeffDict().
lookup(
"inflationCellSet")),
52 duration_(this->coeffDict().getScalar(
"duration")),
72 volumeAccumulator_(0.0),
74 selfSeed_(this->coeffDict().getOrDefault(
"selfSeed", false)),
80 this->coeffDict().subDict(
"sizeDistribution"),
88 flowRateProfile_->userTimeToTime(time);
89 growthRate_->userTimeToTime(time);
98 generationCells_ = generationCells.
toc();
103 inflationCells |= generationCells;
105 inflationCells_ = inflationCells.
toc();
109 scalar generationVolume = 0.0;
111 forAll(generationCells_, gCI)
113 label cI = generationCells_[gCI];
118 scalar totalGenerationVolume = generationVolume;
122 fraction_ = generationVolume/totalGenerationVolume;
126 this->
volumeTotal_ = fraction_*flowRateProfile_->integrate(0.0, duration_);
131template<
class CloudType>
138 generationSetName_(im.generationSetName_),
139 inflationSetName_(im.inflationSetName_),
140 generationCells_(im.generationCells_),
141 inflationCells_(im.inflationCells_),
142 duration_(im.duration_),
143 flowRateProfile_(im.flowRateProfile_.clone()),
144 growthRate_(im.growthRate_.clone()),
145 newParticles_(im.newParticles_),
146 volumeAccumulator_(im.volumeAccumulator_),
147 fraction_(im.fraction_),
148 selfSeed_(im.selfSeed_),
150 sizeDistribution_(im.sizeDistribution_.clone())
156template<
class CloudType>
163template<
class CloudType>
168template<
class CloudType>
171 return this->SOI_ + duration_;
175template<
class CloudType>
185 this->owner().cellOccupancy();
187 scalar gR = growthRate_->value(time1);
189 scalar dT = time1 - time0;
193 forAll(inflationCells_, iCI)
195 label cI = inflationCells_[iCI];
203 scalar dTarget = pPtr->dTarget();
205 pPtr->d() =
min(dTarget, pPtr->d() + gR*dT);
211 newParticles_.clear();
213 Random& rnd = this->owner().rndGen();
219 if ((time0 >= 0.0) && (time0 < duration_))
221 volumeAccumulator_ +=
222 fraction_*flowRateProfile_->integrate(time0, time1);
228 label maxIterations =
max
231 (10*volumeAccumulator_)
235 label iterationNo = 0;
240 while (!generationCells_.empty() && volumeAccumulator_ > 0)
242 if (iterationNo > maxIterations)
245 <<
"Maximum particle split iterations ("
246 << maxIterations <<
") exceeded" <<
endl;
251 label cI = generationCells_
256 generationCells_.size() - 1
266 if (selfSeed_ && !cellCentresUsed.
found(cI))
268 scalar dNew = sizeDistribution_().sample();
281 cellCentresUsed.
insert(cI);
294 scalar pD = pPtr->d();
297 if ((pD/pPtr->dTarget()) < rnd.
sample01<scalar>())
302 const point pP(pPtr->position());
303 const vector& pU = pPtr->U();
326 scalar
x = a/
sqrt(3.0);
327 scalar r = a/(2.0*
sqrt(6.0));
329 scalar d = a/(2.0*
sqrt(3.0));
331 scalar dNew = sizeDistribution_().sample();
342 volumeAccumulator_ -= volNew;
344 dNew = sizeDistribution_().sample();
353 volumeAccumulator_ -= volNew;
355 dNew = sizeDistribution_().sample();
364 volumeAccumulator_ -= volNew;
366 dNew = sizeDistribution_().sample();
375 volumeAccumulator_ -= volNew;
384 this->owner().deleteParticle(*pPtr);
410 gatheredNewParticles,
417 newParticles_ = combinedNewParticles;
423 return newParticles_.size();
427template<
class CloudType>
434 if ((time0 >= 0.0) && (time0 < duration_))
436 return fraction_*flowRateProfile_->integrate(time0, time1);
443template<
class CloudType>
455 position = newParticles_[parcelI].first().first();
457 this->findCellAtPosition
468template<
class CloudType>
477 parcel.U() = newParticles_[parcelI].first().second();
479 parcel.d() = newParticles_[parcelI].second().first();
481 parcel.dTarget() = newParticles_[parcelI].second().second();
485template<
class CloudType>
492template<
class CloudType>
#define R(A, B, C, D, E, F, K, M)
const List< DynamicList< molecule * > > & cellOccupancy
const CloudType & owner() const
Return const access to the owner cloud.
Templated base class for dsmc cloud.
ParcelType parcelType
Type of parcel the cloud was instantiated for.
const fvMesh & mesh() const
Return reference to the mesh.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
List< Key > toc() const
The table of contents (the keys) in unsorted order.
bool found(const Key &key) const
Return true if hashed entry is found in table.
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Inflation injection - creates new particles by splitting existing particles within in a set of genera...
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
virtual void setPositionAndCell(const label parcelI, const label nParcels, const scalar time, vector &position, label &cellOwner, label &tetFacei, label &tetPti)
Set the injection position and owner cell, tetFace and tetPt.
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
virtual ~InflationInjection()
Destructor.
virtual void updateMesh()
Set injector locations when mesh is updated.
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
scalar timeEnd() const
Return the end-of-injection time.
Templated injection model class.
scalar massTotal_
Total mass to inject [kg].
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
An ordered pair of two objects of type <T> with first() and second() elements.
label nProcs() const noexcept
Number of ranks associated with PstreamBuffers.
static void gatherList(const List< commsStruct > &comms, List< T > &values, const int tag, const label comm)
static void scatter(const List< commsStruct > &comms, T &value, const int tag, const label comm)
Broadcast data: Distribute without modification.
Type position(const Type &start, const Type &end)
Return a sample on the interval [start,end].
Type sample01()
Return a sample whose components lie in the range [0,1].
virtual scalar userTimeToTime(const scalar theta) const
Convert the user-time (e.g. CA deg) to real-time (s).
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
static bool & parRun() noexcept
Test if this a parallel run.
A collection of cell labels.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
A library of runtime-selectable doubly-truncated probability distribution models. Returns random samp...
const Time & time() const noexcept
Return time registry.
Mesh consisting of general polyhedral cells.
const scalarField & cellVolumes() const
const vectorField & cellCentres() const
int myProcNo() const noexcept
Return processor number.
Lookup type of boundary radiation properties.
splitCell * master() const
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
A class for handling words, derived from Foam::string.
#define WarningInFunction
Report a warning using Foam::Warning.
AccessType combine(const UList< T > &lists, AccessOp aop=accessOp< T >())
Combines sub-lists into a single list.
const wordList volume
Standard volume field types (scalar, vector, tensor, etc)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar sqrt(const dimensionedScalar &ds)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
static constexpr const zero Zero
Global zero (0)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
#define forAll(list, i)
Loop across all elements in list.