32template<
class CloudType>
45 this->
writeTabbed(
os, outPatchName +
"_massRemoved_" + suffix);
54 this->
writeTabbed(
os, inPatchName +
"_massInjected_" + suffix);
62template<
class CloudType>
69 PatchInteractionModel<CloudType>(
dict, cloud, typeName),
71 recyclePatches_(this->coeffDict().lookup(
"recyclePatches")),
72 recyclePatchesIds_(recyclePatches_.size()),
73 recycledParcels_(recyclePatches_.size()),
74 nRemoved_(recyclePatches_.size()),
75 massRemoved_(nRemoved_.size()),
76 nInjected_(nRemoved_.size()),
77 massInjected_(nRemoved_.size()),
78 injectionPatchPtr_(nRemoved_.size()),
81 this->coeffDict().template getCheck<scalar>
84 scalarMinMax::zero_one()
89 this->coeffDict().getOrDefault(
"outputByInjectorId", false)
96 for (
const auto& inj : cloud.injectors())
132template<
class CloudType>
135 const RecycleInteraction<CloudType>& pim
140 recyclePatches_(pim.recyclePatches_),
141 recyclePatchesIds_(pim.recyclePatchesIds_),
142 nRemoved_(pim.nRemoved_),
143 massRemoved_(pim.massRemoved_),
144 nInjected_(pim.nInjected_),
145 massInjected_(pim.massInjected_),
146 injIdToIndex_(pim.injIdToIndex_),
147 injectionPatchPtr_(),
148 recycleFraction_(pim.recycleFraction_),
155template<
class CloudType>
167 ? injIdToIndex_.lookup(
p.typeId(), 0)
173 forAll(recyclePatchesIds_, i)
175 if (recyclePatchesIds_[i].first() ==
pp.index())
190 keepParticle =
false;
191 recycledParcels_[addri].append
193 static_cast<parcelType*
>(
p.
clone().ptr())
196 ++nRemoved_[addri][idx];
203template<
class CloudType>
216 auto& rnd = this->owner().rndGen();
218 forAll(recycledParcels_, addri)
220 auto& patchParcels = recycledParcels_[addri];
221 auto& injectionPatch = injectionPatchPtr_[addri];
223 for (parcelType&
p : patchParcels)
226 const scalar fraction01 = rnd.template sample01<scalar>();
229 const label toProci = injectionPatch.whichProc(fraction01);
232 auto* osptr = UOPstreamPtrs.
get(toProci);
236 UOPstreamPtrs.
set(toProci, osptr);
240 (*osptr) << addri << fraction01 <<
p;
243 delete(patchParcels.remove(&
p));
258 for (
const int proci : pBufs.
allProcs())
269 auto* newp =
new parcelType(this->owner().
mesh(), is);
275 injectionPatchPtr_[addri].setPositionAndCell
285 newp->relocate(newPosition, cellOwner);
286 newp->nParticle() *= recycleFraction_;
289 newp->U() = this->owner().U()[cellOwner];
295 ? injIdToIndex_.lookup(newp->typeId(), 0)
299 ++nInjected_[addri][idx];
300 massInjected_[addri][idx] += newp->nParticle()*newp->mass();
302 this->owner().addParticle(newp);
309 forAll(recycledParcels_, addri)
311 for (parcelType&
p : recycledParcels_[addri])
313 parcelType* newp = recycledParcels_[addri].remove(&
p);
319 injectionPatchPtr_[addri].setPositionAndCell
330 newp->relocate(newPosition, cellOwner);
331 newp->nParticle() *= recycleFraction_;
334 newp->U() = this->owner().U()[cellOwner];
340 ? injIdToIndex_.lookup(newp->typeId(), 0)
343 ++nInjected_[addri][idx];
344 massInjected_[addri][idx] += newp->nParticle()*newp->mass();
346 this->owner().addParticle(newp);
353template<
class CloudType>
365 const label lsd = nRemoved_[patchi].size();
366 npr0[patchi].setSize(lsd,
Zero);
367 mpr0[patchi].setSize(lsd,
Zero);
368 npi0[patchi].setSize(lsd,
Zero);
369 mpi0[patchi].setSize(lsd,
Zero);
372 this->getModelProperty(
"nRemoved", npr0);
373 this->getModelProperty(
"massRemoved", mpr0);
374 this->getModelProperty(
"nInjected", npi0);
375 this->getModelProperty(
"massInjected", mpi0);
383 npr[i] = npr[i] + npr0[i];
390 mpr[i] = mpr[i] + mpr0[i];
397 npi[i] = npi[i] + npi0[i];
404 mpi[i] = mpi[i] + mpi0[i];
407 if (injIdToIndex_.size())
411 labelList indexToInjector(injIdToIndex_.size());
414 indexToInjector[iter.val()] = iter.key();
419 const word& outPatchName = recyclePatches_[i].first();
421 Log_<<
" Parcel fate: patch " << outPatchName
422 <<
" (number, mass)" <<
nl;
426 Log_<<
" - removed (injector " << indexToInjector[indexi]
427 <<
") = " << npr[i][indexi]
428 <<
", " << mpr[i][indexi] <<
nl;
431 <<
tab << npr[i][indexi] <<
tab << mpr[i][indexi];
434 const word& inPatchName = recyclePatches_[i].second();
436 Log_<<
" Parcel fate: patch " << inPatchName
437 <<
" (number, mass)" <<
nl;
441 Log_<<
" - injected (injector " << indexToInjector[indexi]
442 <<
") = " << npi[i][indexi]
443 <<
", " << mpi[i][indexi] <<
nl;
445 <<
tab << npi[i][indexi] <<
tab << mpi[i][indexi];
449 this->file() <<
endl;
455 const word& outPatchName = recyclePatches_[i].first();
457 Log_<<
" Parcel fate: patch " << outPatchName
458 <<
" (number, mass)" <<
nl
459 <<
" - removed = " << npr[i][0] <<
", " << mpr[i][0]
463 <<
tab << npr[i][0] <<
tab << mpr[i][0];
468 const word& inPatchName = recyclePatches_[i].second();
470 Log_<<
" Parcel fate: patch " << inPatchName
471 <<
" (number, mass)" <<
nl
472 <<
" - injected = " << npi[i][0] <<
", " << mpi[i][0]
476 <<
tab << npi[i][0] <<
tab << mpi[i][0];
479 this->file() <<
endl;
482 if (this->writeTime())
484 this->setModelProperty(
"nRemoved", npr);
485 this->setModelProperty(
"massRemoved", mpr);
486 this->setModelProperty(
"nInjected", npi);
487 this->setModelProperty(
"massInjected", mpi);
492 massInjected_ =
Zero;
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const CloudType & owner() const
Return const access to the owner cloud.
virtual bool writeTime() const
Flag to indicate when to write a property.
label typeId() const
Return type id.
const vector & U() const
Return const access to velocity.
tmp< GeometricField< Type, PatchField, GeoMesh > > clone() const
Clone.
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
bool eof() const noexcept
True if end of input seen.
void setSize(label n)
Alias for resize().
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
PatchInteractionModel(CloudType &owner)
Construct null from owner.
virtual void info()
Write patch interaction info.
virtual void writeFileHeader(Ostream &os)
Output file header information.
Buffers for inter-processor communications streams (UOPstream, UIPstream).
UPstream::rangeType allProcs() const noexcept
Range of ranks indices associated with PstreamBuffers.
label recvDataCount(const label proci) const
Number of unconsumed receive bytes for the specified processor. Must call finishedSends() or other fi...
void finishedSends(const bool wait=true)
Mark the send phase as being finished.
static void listGather(UList< T > &values, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Gather (reduce) list elements, applying bop to each list element.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
virtual void postEvolve()
Post-evolve hook.
const fvMesh & mesh_
Reference to mesh.
RecycleInteraction(const dictionary &dict, CloudType &cloud)
Construct from dictionary.
List< Pair< label > > recyclePatchesIds_
Patch IDs of recyclePatches.
Map< label > injIdToIndex_
Injector ID to local index map.
List< List< scalar > > massInjected_
Mass of parcels injected.
virtual bool correct(typename CloudType::parcelType &p, const polyPatch &pp, bool &keepParticle)
Apply velocity correction.
List< List< label > > nRemoved_
Number of parcels removed.
List< Pair< word > > recyclePatches_
Outlet-inlet patch pair to apply parcel recycling.
List< List< label > > nInjected_
Number of parcels injected.
virtual void info()
Write patch interaction info.
virtual void writeFileHeader(Ostream &os)
Output file header information.
List< List< scalar > > massRemoved_
Mass of parcels removed.
const scalar recycleFraction_
Parcel fraction to be recycled from outlet to inlet.
List< IDLList< parcelType > > recycledParcels_
Parcel IDs of recycled parcels.
CloudType::parcelType parcelType
bool outputByInjectorId_
Flag to output escaped/mass particles sorted by injectorID.
PtrList< patchInjectionBase > injectionPatchPtr_
Injection patch pointer.
Input inter-processor communications stream using MPI send/recv etc. - operating on external buffer.
Output inter-processor communications stream using MPI send/recv etc. - operating on external buffer.
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.
const T * get(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
A cloud is a registry collection of lagrangian particles.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
virtual OFstream & file()
Return access to the file (if only 1).
A traits class, which is primarily used for primitives and vector-space.
void relocate(const point &position, const label celli=-1)
Set the addressing based on the provided position.
A patch is a list of labels that address the faces in the global face list.
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.
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.
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
OBJstream os(runTime.globalPath()/outputName)
#define Log_
Report write to Foam::Info if the class log switch is true.
List< scalarList > scalarListList
List of scalarList.
DSMCCloud< dsmcParcel > CloudType
List< labelList > labelListList
List of labelList.
List< label > labelList
A List of labels.
MinMax< scalar > scalarMinMax
A scalar min/max range.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
static constexpr const zero Zero
Global zero (0).
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
constexpr char nl
The newline '\n' character (0x0a).
constexpr char tab
The tab '\t' character(0x09).
#define forAll(list, i)
Loop across all elements in list.
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.