37template<
class CloudType>
52 tag[particlei] =
p.tag();
53 position[particlei] =
p.position();
55 soi[particlei] =
p.soi();
112 maxTag =
max(maxTag, tag[particlei]);
115 label nInjectors = maxTag + 1;
125 const label tagi = tag[i];
126 const scalar t = soi[i];
127 injStartTime[tagi] =
min(t, injStartTime[tagi]);
128 injEndTime[tagi] =
max(t, injEndTime[tagi]);
129 injPosition[tagi].
append(position[i]);
131 injDiameter[tagi].
append(d[i]);
136 scalar sumVolume = 0;
141 U_.setSize(nInjectors);
144 scalar minTime = GREAT;
152 const label nParticle = diameters.
size();
153 const scalar dTime = injEndTime[i] - injStartTime[i];
155 if ((nParticle > 1) && (dTime > ROOTVSMALL))
157 minTime =
min(minTime, injStartTime[i]);
160 endTime_[injectori] = injEndTime[i];
172 Ui[samplei] = injU[i][posi];
177 forAll(diameters, particlei)
179 sumPow3 +=
pow3(diameters[particlei]);
206 U_.setSize(injectori);
229template<
class CloudType>
235 const word& modelName
239 cloudName_(this->coeffDict().
lookup(
"cloud")),
240 startTime_(this->template getModelProperty<
scalarList>(
"startTime")),
241 endTime_(this->template getModelProperty<
scalarList>(
"endTime")),
242 position_(this->template getModelProperty<
List<
vectorList>>(
"position")),
243 positionOffset_(this->coeffDict().
lookup(
"positionOffset")),
246 this->template getModelProperty<
scalarList>(
"volumeFlowRate")
249 binWidth_(this->coeffDict().getScalar(
"binWidth")),
253 ceil(this->coeffDict().getScalar(
"parcelsPerInjector"))
257 this->coeffDict().getOrDefault(
"resampleSize", label(100))
259 applyDistributionMassTotal_
261 this->coeffDict().getBool(
"applyDistributionMassTotal")
265 this->coeffDict().getOrDefault(
"ignoreOutOfBounds", false)
267 nParcelsInjected_(this->parcelsAddedTotal()),
268 nParcelsInjected0_(0),
269 currentInjectori_(0),
275 sizeDistribution_.setSize(startTime_.size());
276 forAll(sizeDistribution_, i)
278 const word dictName(
"distribution" + Foam::name(i));
280 this->getModelDict(dictName, dict);
282 sizeDistribution_.set
285 new distributionModels::general(dict, this->owner().rndGen())
296 if (applyDistributionMassTotal_)
298 this->massTotal_ = this->volumeTotal_*this->owner().constProps().rho0();
299 Info<<
" Set mass to inject from distribution: "
305template<
class CloudType>
313 cloudName_(im.cloudName_),
314 startTime_(im.startTime_),
315 endTime_(im.endTime_),
316 position_(im.position_),
317 positionOffset_(im.positionOffset_),
318 volumeFlowRate_(im.volumeFlowRate_),
320 binWidth_(im.binWidth_),
321 sizeDistribution_(im.sizeDistribution_.size()),
322 parcelsPerInjector_(im.parcelsPerInjector_),
323 resampleSize_(im.resampleSize_),
324 applyDistributionMassTotal_(im.applyDistributionMassTotal_),
325 ignoreOutOfBounds_(im.ignoreOutOfBounds_),
326 nParcelsInjected_(im.nParcelsInjected_),
327 nParcelsInjected0_(im.nParcelsInjected0_),
328 currentInjectori_(0),
350template<
class CloudType>
358template<
class CloudType>
363template<
class CloudType>
371template<
class CloudType>
381 nParcelsInjected0_ = 0;
383 if (startTime_.empty() || this->volumeTotal_ < ROOTVSMALL)
388 scalar targetVolume = 0;
389 forAll(startTime_, injectori)
391 if (time1 > startTime_[injectori])
393 scalar totalDuration =
394 min(time1, endTime_[injectori]) - startTime_[injectori];
396 targetVolume += volumeFlowRate_[injectori]*totalDuration;
400 const label targetParcels =
403 scalar(startTime_.size()*parcelsPerInjector_)
404 *targetVolume/this->volumeTotal_
407 const label nParcels = targetParcels - nParcelsInjected_;
413template<
class CloudType>
422 forAll(startTime_, injectori)
424 if ((time1 > startTime_[injectori]) && (time1 <= endTime_[injectori]))
426 scalar duration =
min(time1, endTime_[injectori]) - time0;
427 volume += volumeFlowRate_[injectori]*duration;
435template<
class CloudType>
439 const label nParcels,
447 Random& rnd = this->owner().rndGen();
448 currentInjectori_ = rnd.
globalPosition<label>(0, position_.size() - 1);
449 currentSamplei_ = rnd.
globalPosition<label>(0, resampleSize_ - 1);
451 position = position_[currentInjectori_][currentSamplei_];
454 this->findCellAtPosition
464template<
class CloudType>
474 parcel.
U() = U_[currentInjectori_][currentSamplei_];
477 parcel.d() = sizeDistribution_[currentInjectori_].sample();
485template<
class CloudType>
493template<
class CloudType>
503template<
class CloudType>
508 if (this->writeTime())
510 this->setModelProperty(
"startTime", startTime_);
511 this->setModelProperty(
"endTime", endTime_);
512 this->setModelProperty(
"position", position_);
513 this->setModelProperty(
"volumeFlowRate", volumeFlowRate_);
514 this->setModelProperty(
"U", U_);
515 forAll(sizeDistribution_, i)
const CloudType & owner() const
Return const access to the owner cloud.
virtual bool writeTime() const
Flag to indicate when to write a property.
const vector & U() const
Return const access to velocity.
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...
InjectedParticleDistributionInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
label nParcelsInjected0_
Number of parcels injected in last step (local proc only).
void initialise()
Initialise injectors.
bool applyDistributionMassTotal_
Flag to apply mass calculated from distribution instead of.
scalarList endTime_
List of end time per injector.
const word cloudName_
Name of cloud used to seed the new particles.
label currentSamplei_
Current sample.
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.
label nParcelsInjected_
Running total of number of parcels injected.
scalarList startTime_
List of start time per injector.
vector positionOffset_
Position offset to apply to input positions.
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
scalar binWidth_
Bin width when generating particle distributions.
scalarList volumeFlowRate_
List of volume flow rate per injector [m3/s].
List< vectorList > position_
List of position per injector.
List< vectorList > U_
List of parcel velocity per injector.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
Switch ignoreOutOfBounds_
Flag to suppress errors if particle injection site is out-of-bounds.
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.
void info()
Write injection info.
scalar parcelsPerInjector_
Target number of parcels to inject per injector.
label currentInjectori_
Current injector.
virtual void updateMesh()
Set injector locations when mesh is updated.
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
label resampleSize_
Resample size.
scalar timeEnd() const
Return the end-of-injection time.
Templated injection model class.
virtual void info()
Write injection info.
label parcelsAddedTotal() const
Return the total number parcels added.
scalar volumeTotal_
Total volume of particles introduced by this injector [m^3] Note: scaled to ensure massTotal is achie...
InjectionModel(CloudType &owner)
Construct null from owner.
scalar massTotal_
Total mass to inject [kg].
virtual bool findCellAtPosition(label &celli, label &tetFacei, label &tetPti, vector &position, bool errorOnNotFound=true)
Find the cell that contains the supplied position.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
void append(const T &val)
Append an element at the end of the list.
static void allGatherList(UList< T > &values, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Gather data, but keep individual values separate. Uses MPI_Allgather or manual communication.
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.
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run.
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.
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 getModelProperty(const word &entryName, Type &value) const
Retrieve generic property from the sub-model.
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
const word & modelName() const
Return const access to the name of the sub-model.
void setModelProperty(const word &entryName, const Type &value)
Add generic property to the sub-model.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
const word dictName("faMeshDefinition")
AccessType combine(const UList< T > &lists, AccessOp aop=accessOp< T >())
Combines sub-lists into a single list.
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.
DSMCCloud< dsmcParcel > CloudType
dimensionedScalar pow3(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< vector > vectorList
List of vector.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
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 expressions::valueTypeCode::INVALID.
List< scalar > scalarList
List of scalar.
#define forAll(list, i)
Loop across all elements in list.
Object access operator or list access operator (default is pass-through).