36 template<
class CloudType>
45 generationSetName_(this->coeffDict().lookup(
"generationCellSet")),
46 inflationSetName_(this->coeffDict().lookup(
"inflationCellSet")),
69 volumeAccumulator_(0.0),
71 selfSeed_(this->coeffDict().lookupOrDefault(
"selfSeed",
false)),
77 this->coeffDict().subDict(
"sizeDistribution"),
82 duration_ = owner.db().time().userTimeToTime(duration_);
89 cellSet generationCells(this->owner().
mesh(), generationSetName_);
91 generationCells_ = generationCells.
toc();
93 cellSet inflationCells(this->owner().
mesh(), inflationSetName_);
96 inflationCells |= generationCells;
98 inflationCells_ = inflationCells.
toc();
102 scalar generationVolume = 0.0;
104 forAll(generationCells_, gCI)
106 label cI = generationCells_[gCI];
108 generationVolume += this->owner().mesh().cellVolumes()[cI];
111 scalar totalGenerationVolume = generationVolume;
115 fraction_ = generationVolume/totalGenerationVolume;
119 this->volumeTotal_ = fraction_*flowRateProfile_.integrate(0.0, duration_);
120 this->massTotal_ *= fraction_;
124 template<
class CloudType>
131 generationSetName_(im.generationSetName_),
132 inflationSetName_(im.inflationSetName_),
133 generationCells_(im.generationCells_),
134 inflationCells_(im.inflationCells_),
135 duration_(im.duration_),
136 flowRateProfile_(im.flowRateProfile_),
137 growthRate_(im.growthRate_),
138 newParticles_(im.newParticles_),
139 volumeAccumulator_(im.volumeAccumulator_),
140 fraction_(im.fraction_),
141 selfSeed_(im.selfSeed_),
143 sizeDistribution_(im.sizeDistribution_().clone().ptr())
149 template<
class CloudType>
156 template<
class CloudType>
163 template<
class CloudType>
166 return this->SOI_ + duration_;
170 template<
class CloudType>
180 this->owner().cellOccupancy();
182 scalar gR = growthRate_.value(time1);
184 scalar dT = time1 - time0;
188 forAll(inflationCells_, iCI)
190 label cI = inflationCells_[iCI];
194 forAll(cellOccupancy[cI], cPI)
196 pPtr = cellOccupancy[cI][cPI];
198 scalar dTarget = pPtr->dTarget();
200 pPtr->d() =
min(dTarget, pPtr->d() + gR*dT);
206 newParticles_.clear();
214 if ((time0 >= 0.0) && (time0 < duration_))
216 volumeAccumulator_ +=
217 fraction_*flowRateProfile_.integrate(time0, time1);
226 (10*volumeAccumulator_)
227 /CloudType::parcelType::volume(sizeDistribution_().
minValue())
230 label iterationNo = 0;
235 while (!generationCells_.empty() && volumeAccumulator_ > 0)
237 if (iterationNo > maxIterations)
242 "Foam::InflationInjection<CloudType>::parcelsToInject" 248 <<
"Maximum particle split iterations (" 249 << maxIterations <<
") exceeded" <<
endl;
254 label cI = generationCells_
259 generationCells_.size() - 1
267 if (cellOccupancy[cI].empty())
269 if (selfSeed_ && !cellCentresUsed.
found(cI))
271 scalar dNew = sizeDistribution_().sample();
282 volumeAccumulator_ -= CloudType::parcelType::volume(dNew);
284 cellCentresUsed.
insert(cI);
297 scalar pD = pPtr->d();
300 if ((pD/pPtr->dTarget()) < rnd.
sample01<scalar>())
305 const point& pP = pPtr->position();
306 const vector& pU = pPtr->U();
329 scalar
x = a/
sqrt(3.0);
330 scalar r = a/(2.0*
sqrt(6.0));
332 scalar d = a/(2.0*
sqrt(3.0));
334 scalar dNew = sizeDistribution_().sample();
335 scalar volNew = CloudType::parcelType::volume(dNew);
345 volumeAccumulator_ -= volNew;
347 dNew = sizeDistribution_().sample();
356 volumeAccumulator_ -= volNew;
358 dNew = sizeDistribution_().sample();
367 volumeAccumulator_ -= volNew;
369 dNew = sizeDistribution_().sample();
378 volumeAccumulator_ -= volNew;
382 volumeAccumulator_ += CloudType::parcelType::volume
387 this->owner().deleteParticle(*pPtr);
413 gatheredNewParticles,
420 newParticles_ = combinedNewParticles;
426 return newParticles_.size();
430 template<
class CloudType>
437 if ((time0 >= 0.0) && (time0 < duration_))
439 return fraction_*flowRateProfile_.integrate(time0, time1);
448 template<
class CloudType>
460 position = newParticles_[parcelI].first().first();
462 this->findCellAtPosition
473 template<
class CloudType>
482 parcel.U() = newParticles_[parcelI].first().second();
484 parcel.d() = newParticles_[parcelI].second().first();
486 parcel.dTarget() = newParticles_[parcelI].second().second();
490 template<
class CloudType>
497 template<
class CloudType>
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensionedScalar sqrt(const dimensionedScalar &ds)
static bool & parRun()
Is this a parallel run?
A 2-tuple for storing two objects of different types.
Inflation injection - creates new particles by splitting existing particles within in a set of genera...
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
An ordered pair of two objects of type <T> with first() and second() elements.
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 fullyDescribed() const
Flag to identify whether model fully describes the parcel.
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
#define R(A, B, C, D, E, F, K, M)
A collection of cell labels.
A class for handling words, derived from string.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const vectorField & cellCentres() const
A list of keyword definitions, which are a keyword followed by any number of values (e...
ParcelType parcelType
Type of parcel the cloud was instantiated for.
AccessType combine(const List< T > &, AccessOp aop=accessOp< T >())
Combines sublists into one big list.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Ostream & endl(Ostream &os)
Add newline and flush stream.
stressControl lookup("compactNormalStress") >> compactNormalStress
#define WarningIn(functionName)
Report a warning using Foam::Warning.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
virtual void updateMesh()
Set injector locations when mesh is updated.
Templated injection model class.
static void scatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Distribute without modification. Reverse of gather.
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
static autoPtr< distributionModel > New(const dictionary &dict, cachedRandom &rndGen)
Selector.
List< Key > toc() const
Return the table of contents.
scalar timeEnd() const
Return the end-of-injection time.
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
Mesh consisting of general polyhedral cells.
const List< DynamicList< molecule * > > & cellOccupancy
Type sample01()
Return a sample whose components lie in the range 0-1.
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
Vector< scalar > vector
A scalar version of the templated Vector.
static bool master(const label communicator=0)
Am I the master process.
Type position(const Type &start, const Type &end)
Return a sample between start and end.
Templated base class for dsmc cloud.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
InflationInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
virtual ~InflationInjection()
Destructor.
bool found(const Key &) const
Return true if hashedEntry is found in table.
bool insert(const Key &key)
Insert a new entry.