37template<
class CloudType>
52 tag[particlei] =
p.tag();
53 position[particlei] =
p.position();
55 soi[particlei] =
p.soi();
61 if (Pstream::parRun())
64 procTag[Pstream::myProcNo()].
transfer(tag);
65 Pstream::gatherList(procTag);
66 Pstream::scatterList(procTag);
68 ListListOps::combine<List<label>>
74 procPosition[Pstream::myProcNo()].
transfer(position);
75 Pstream::gatherList(procPosition);
76 Pstream::scatterList(procPosition);
78 ListListOps::combine<List<point>>
85 Pstream::gatherList(procU);
86 Pstream::scatterList(procU);
88 ListListOps::combine<List<vector>>
94 procSOI[Pstream::myProcNo()].
transfer(soi);
95 Pstream::gatherList(procSOI);
96 Pstream::scatterList(procSOI);
98 ListListOps::combine<List<scalar>>
104 procD[Pstream::myProcNo()].
transfer(d);
105 Pstream::gatherList(procD);
106 Pstream::scatterList(procD);
108 ListListOps::combine<List<scalar>>
117 maxTag =
max(maxTag, tag[particlei]);
120 label nInjectors = maxTag + 1;
130 const label tagi = tag[i];
131 const scalar t = soi[i];
132 injStartTime[tagi] =
min(t, injStartTime[tagi]);
133 injEndTime[tagi] =
max(t, injEndTime[tagi]);
134 injPosition[tagi].
append(position[i]);
136 injDiameter[tagi].
append(d[i]);
141 scalar sumVolume = 0;
142 startTime_.setSize(nInjectors, 0);
143 endTime_.setSize(nInjectors, 0);
144 sizeDistribution_.setSize(nInjectors);
145 position_.setSize(nInjectors);
146 U_.setSize(nInjectors);
147 volumeFlowRate_.setSize(nInjectors, 0);
149 scalar minTime = GREAT;
152 Random& rnd = this->owner().rndGen();
157 const label nParticle = diameters.
size();
158 const scalar dTime = injEndTime[i] - injStartTime[i];
160 if ((nParticle > 1) && (dTime > ROOTVSMALL))
162 minTime =
min(minTime, injStartTime[i]);
164 startTime_[injectori] = injStartTime[i];
165 endTime_[injectori] = injEndTime[i];
168 position_[injectori].
setSize(resampleSize_);
169 U_[injectori].setSize(resampleSize_);
173 for (label samplei = 0; samplei < resampleSize_; ++samplei)
176 positioni[samplei] = injPosition[i][posi] + positionOffset_;
177 Ui[samplei] = injU[i][posi];
182 forAll(diameters, particlei)
184 sumPow3 +=
pow3(diameters[particlei]);
189 volumeFlowRate_[injectori] = volume/dTime;
192 sizeDistribution_.
set
208 startTime_.setSize(injectori);
209 endTime_.setSize(injectori);
210 position_.setSize(injectori);
211 U_.setSize(injectori);
212 volumeFlowRate_.setSize(injectori);
213 sizeDistribution_.setSize(injectori);
216 forAll(startTime_, injectori)
218 startTime_[injectori] -= minTime;
219 endTime_[injectori] -= minTime;
224 this->volumeTotal_ = sumVolume;
227 Info<<
" Read " << position_.size() <<
" injectors with "
228 << tag.
size() <<
" total particles" <<
endl;
234template<
class CloudType>
240 const word& modelName
244 cloudName_(this->coeffDict().
lookup(
"cloud")),
245 startTime_(this->template getModelProperty<
scalarList>(
"startTime")),
246 endTime_(this->template getModelProperty<
scalarList>(
"endTime")),
247 position_(this->template getModelProperty<
List<
vectorList>>(
"position")),
248 positionOffset_(this->coeffDict().
lookup(
"positionOffset")),
251 this->template getModelProperty<
scalarList>(
"volumeFlowRate")
254 binWidth_(this->coeffDict().getScalar(
"binWidth")),
258 ceil(this->coeffDict().getScalar(
"parcelsPerInjector"))
262 this->coeffDict().getOrDefault(
"resampleSize", label(100))
264 applyDistributionMassTotal_
266 this->coeffDict().getBool(
"applyDistributionMassTotal")
270 this->coeffDict().getOrDefault(
"ignoreOutOfBounds", false)
272 nParcelsInjected_(this->parcelsAddedTotal()),
273 nParcelsInjected0_(0),
274 currentInjectori_(0),
304 Info<<
" Set mass to inject from distribution: "
310template<
class CloudType>
318 cloudName_(im.cloudName_),
319 startTime_(im.startTime_),
320 endTime_(im.endTime_),
321 position_(im.position_),
322 positionOffset_(im.positionOffset_),
323 volumeFlowRate_(im.volumeFlowRate_),
325 binWidth_(im.binWidth_),
326 sizeDistribution_(im.sizeDistribution_.size()),
327 parcelsPerInjector_(im.parcelsPerInjector_),
328 resampleSize_(im.resampleSize_),
329 applyDistributionMassTotal_(im.applyDistributionMassTotal_),
330 ignoreOutOfBounds_(im.ignoreOutOfBounds_),
331 nParcelsInjected_(im.nParcelsInjected_),
332 nParcelsInjected0_(im.nParcelsInjected0_),
333 currentInjectori_(0),
355template<
class CloudType>
363template<
class CloudType>
368template<
class CloudType>
372 return max(endTime_);
376template<
class CloudType>
386 nParcelsInjected0_ = 0;
388 if (startTime_.empty() || this->volumeTotal_ < ROOTVSMALL)
393 scalar targetVolume = 0;
394 forAll(startTime_, injectori)
396 if (time1 > startTime_[injectori])
398 scalar totalDuration =
399 min(time1, endTime_[injectori]) - startTime_[injectori];
401 targetVolume += volumeFlowRate_[injectori]*totalDuration;
405 const label targetParcels =
408 scalar(startTime_.size()*parcelsPerInjector_)
409 *targetVolume/this->volumeTotal_
412 const label nParcels = targetParcels - nParcelsInjected_;
418template<
class CloudType>
427 forAll(startTime_, injectori)
429 if ((time1 > startTime_[injectori]) && (time1 <= endTime_[injectori]))
431 scalar duration =
min(time1, endTime_[injectori]) - time0;
432 volume += volumeFlowRate_[injectori]*duration;
440template<
class CloudType>
444 const label nParcels,
452 Random& rnd = this->owner().rndGen();
453 currentInjectori_ = rnd.
globalPosition<label>(0, position_.size() - 1);
454 currentSamplei_ = rnd.
globalPosition<label>(0, resampleSize_ - 1);
456 position = position_[currentInjectori_][currentSamplei_];
459 this->findCellAtPosition
469template<
class CloudType>
479 parcel.U() = U_[currentInjectori_][currentSamplei_];
482 parcel.d() = sizeDistribution_[currentInjectori_].sample();
486 nParcelsInjected0_++;
490template<
class CloudType>
498template<
class CloudType>
508template<
class CloudType>
513 if (this->writeTime())
515 this->setModelProperty(
"startTime", startTime_);
516 this->setModelProperty(
"endTime", endTime_);
517 this->setModelProperty(
"position", position_);
518 this->setModelProperty(
"volumeFlowRate", volumeFlowRate_);
519 this->setModelProperty(
"U", U_);
520 forAll(sizeDistribution_, i)
const CloudType & owner() const
Return const access to the owner cloud.
Templated base class for dsmc cloud.
const List< typename ParcelType::constantProperties > & constProps() const
Return all of the constant properties.
ParcelType parcelType
Type of parcel the cloud was instantiated for.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
label size() const noexcept
The number of elements in table.
Interrogates an injectedParticleCloud to convert the raw particle data into a set of 'binned' injecto...
void initialise()
Initialise injectors.
bool applyDistributionMassTotal_
Flag to apply mass calculated from distribution instead of.
PtrList< distributionModels::general > sizeDistribution_
List of size distribution model per injector.
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.
scalarList startTime_
List of start time per injector.
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 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 ~InjectedParticleDistributionInjection()
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...
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
void transfer(List< T > &list)
void setSize(const label n)
Alias for resize()
void append(const T &val)
Append an element at the end of the list.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Type globalPosition(const Type &start, const Type &end)
Return a sample on the interval [start,end].
void size(const label n)
Older name for setAddressableSize.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Particle-size distribution model wherein random samples are drawn from a given arbitrary probability ...
virtual dictionary writeDict(const word &dictName) const
Write data in dictionary format.
InfoProxy< ensightCells > info() const
Return info proxy.
Primarily stores particle properties so that it can be injected at a later time. Note that this store...
Lookup type of boundary radiation properties.
bool getModelDict(const word &entryName, dictionary &dict) const
Retrieve dictionary, return true if set.
const dictionary & dict() const
Return const access to the cloud dictionary.
A class for handling words, derived from Foam::string.
OBJstream os(runTime.globalPath()/outputName)
const word dictName("faMeshDefinition")
constexpr scalar pi(M_PI)
Different types of constants.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
dimensionedScalar pow3(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
#define forAll(list, i)
Loop across all elements in list.