Go to the documentation of this file.
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];
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>
149 template<
class CloudType>
156 template<
class CloudType>
163 template<
class CloudType>
166 return this->SOI_ + duration_;
170 template<
class CloudType>
182 scalar gR = growthRate_.value(time1);
184 scalar dT = time1 - time0;
188 forAll(inflationCells_, iCI)
190 label cI = inflationCells_[iCI];
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)
240 <<
"Maximum particle split iterations ("
241 << maxIterations <<
") exceeded" <<
endl;
246 label cI = generationCells_
251 generationCells_.size() - 1
261 if (selfSeed_ && !cellCentresUsed.
found(cI))
263 scalar dNew = sizeDistribution_().sample();
274 volumeAccumulator_ -= CloudType::parcelType::volume(dNew);
276 cellCentresUsed.
insert(cI);
289 scalar pD = pPtr->d();
292 if ((pD/pPtr->dTarget()) < rnd.
sample01<scalar>())
297 const point& pP = pPtr->position();
298 const vector& pU = pPtr->U();
321 scalar
x = a/
sqrt(3.0);
322 scalar r = a/(2.0*
sqrt(6.0));
324 scalar d = a/(2.0*
sqrt(3.0));
326 scalar dNew = sizeDistribution_().sample();
327 scalar volNew = CloudType::parcelType::volume(dNew);
337 volumeAccumulator_ -= volNew;
339 dNew = sizeDistribution_().sample();
348 volumeAccumulator_ -= volNew;
350 dNew = sizeDistribution_().sample();
359 volumeAccumulator_ -= volNew;
361 dNew = sizeDistribution_().sample();
370 volumeAccumulator_ -= volNew;
374 volumeAccumulator_ += CloudType::parcelType::volume
405 gatheredNewParticles,
412 newParticles_ = combinedNewParticles;
418 return newParticles_.size();
422 template<
class CloudType>
429 if ((time0 >= 0.0) && (time0 < duration_))
431 return fraction_*flowRateProfile_.integrate(time0, time1);
440 template<
class CloudType>
452 position = newParticles_[parcelI].first().first();
454 this->findCellAtPosition
465 template<
class CloudType>
474 parcel.U() = newParticles_[parcelI].first().second();
476 parcel.d() = newParticles_[parcelI].second().first();
478 parcel.dTarget() = newParticles_[parcelI].second().second();
482 template<
class CloudType>
489 template<
class CloudType>
word inflationSetName_
Name of cellSet for inflating new particles.
Type sample01()
Return a sample whose components lie in the range 0-1.
A class for handling words, derived from string.
InflationInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
List< Key > toc() const
Return the table of contents.
#define forAll(list, i)
Loop across all elements in list.
const autoPtr< distributionModels::distributionModel > sizeDistribution_
Parcel size distribution model.
AccessType combine(const List< T > &, AccessOp aop=accessOp< T >())
Combines sublists into one big list.
static label nProcs(const label communicator=0)
Number of processes in parallel run.
static bool & parRun()
Is this a parallel run?
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.
static autoPtr< distributionModel > New(const dictionary &dict, cachedRandom &rndGen)
Selector.
Ostream & endl(Ostream &os)
Add newline and flush stream.
TimeDataEntry< scalar > flowRateProfile_
Flow rate profile relative to SOI [m3/s].
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
Inflation injection - creates new particles by splitting existing particles within in a set of genera...
Mesh consisting of general polyhedral cells.
Random & rndGen()
Return refernce to the random object.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
#define R(A, B, C, D, E, F, K, M)
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 List< DynamicList< molecule * > > & cellOccupancy
labelList generationCells_
Set of cells to generate particles in.
const fvMesh & mesh() const
Return refernce to the mesh.
Type position(const Type &start, const Type &end)
Return a sample between start and end.
void deleteParticle(ParticleType &)
Remove particle from cloud and delete.
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
Templated base class for dsmc cloud.
const List< DynamicList< ParcelType * > > & cellOccupancy() const
Return the cell occupancy addressing.
virtual ~InflationInjection()
Destructor.
A list of keyword definitions, which are a keyword followed by any number of values (e....
bool found(const Key &) const
Return true if hashedEntry is found in table.
const scalarField & cellVolumes() const
A collection of cell labels.
Vector< scalar > vector
A scalar version of the templated Vector.
scalar duration_
Injection duration [s].
Switch selfSeed_
Switch to control whether or not the injector is allowed.
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
static bool master(const label communicator=0)
Am I the master process.
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
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.
scalar dSeed_
Diameter with which to create new seed particles.
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
An ordered pair of two objects of type <T> with first() and second() elements.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
dimensionedScalar sqrt(const dimensionedScalar &ds)
scalar volumeAccumulator_
Accumulation variable to carry over volume from one injection.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
TimeDataEntry< scalar > growthRate_
Growth rate of particle diameters towards target [m/s].
word generationSetName_
Name of cellSet for generating new particles.
bool insert(const Key &key)
Insert a new entry.
scalar timeEnd() const
Return the end-of-injection time.
scalar fraction_
Fraction of injection controlled by this processor.
ParcelType parcelType
Type of parcel the cloud was instantiated for.
DynamicList< vectorPairScalarPair > newParticles_
Positions, velocities, diameters and target diameters of.
A 2-tuple for storing two objects of different types.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
#define WarningInFunction
Report a warning using Foam::Warning.
labelList inflationCells_
Set of cells to inflate particles in, includes all.
stressControl lookup("compactNormalStress") >> compactNormalStress
virtual void updateMesh()
Set injector locations when mesh is updated.