RecycleInteraction.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2020-2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "RecycleInteraction.H"
29
30// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
31
32template<class CloudType>
34{
36
37 forAll(nRemoved_, i)
38 {
39 const word& outPatchName = recyclePatches_[i].first();
40
41 forAll(nRemoved_[i], injectori)
42 {
43 const word suffix = Foam::name(injectori);
44 this->writeTabbed(os, outPatchName + "_nRemoved_" + suffix);
45 this->writeTabbed(os, outPatchName + "_massRemoved_" + suffix);
46 }
47
48 const word& inPatchName = recyclePatches_[i].second();
49
50 forAll(nInjected_[i], injectori)
51 {
52 const word suffix = Foam::name(injectori);
53 this->writeTabbed(os, inPatchName + "_nInjected_" + suffix);
54 this->writeTabbed(os, inPatchName + "_massInjected_" + suffix);
55 }
56 }
57}
58
59
60// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
61
62template<class CloudType>
64(
65 const dictionary& dict,
67)
68:
70 mesh_(cloud.mesh()),
71 recyclePatches_(this->coeffDict().lookup("recyclePatches")),
72 recyclePatchesIds_(recyclePatches_.size()),
73 recycledParcels_(recyclePatches_.size()),
74 nRemoved_(recyclePatches_.size()), // per patch the no. of parcels
75 massRemoved_(nRemoved_.size()),
76 nInjected_(nRemoved_.size()),
77 massInjected_(nRemoved_.size()),
78 injectionPatchPtr_(nRemoved_.size()),
79 recycleFraction_
80 (
81 this->coeffDict().template getCheck<scalar>
82 (
83 "recycleFraction",
84 scalarMinMax::zero_one()
85 )
86 ),
87 outputByInjectorId_
88 (
89 this->coeffDict().getOrDefault("outputByInjectorId", false)
90 )
91{
92 // Determine the number of injectors and the injector mapping
93 label nInjectors = 0;
95 {
96 for (const auto& inj : cloud.injectors())
97 {
98 injIdToIndex_.insert(inj.injectorID(), ++nInjectors);
99 }
100 }
101
102 // The normal case, and safety if injector mapping was somehow null
103 if (injIdToIndex_.empty())
104 {
105 nInjectors = 1;
106 }
107
108 forAll(nRemoved_, i)
109 {
110 // Create injection helper for each inflow patch
112 (
113 i,
115 );
116
117 // Mappings from patch names to patch IDs
118 recyclePatchesIds_[i].first() =
120 recyclePatchesIds_[i].second() =
122
123 // Storage for reporting
124 nRemoved_[i].setSize(nInjectors, Zero);
125 massRemoved_[i].setSize(nInjectors, Zero);
126 nInjected_[i].setSize(nInjectors, Zero);
127 massInjected_[i].setSize(nInjectors, Zero);
128 }
129}
130
131
132template<class CloudType>
134(
136)
137:
139 mesh_(pim.mesh_),
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_),
149 outputByInjectorId_(pim.outputByInjectorId_)
150{}
151
152
153// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
154
155template<class CloudType>
157(
158 typename CloudType::parcelType& p,
159 const polyPatch& pp,
160 bool& keepParticle
161)
162{
163 // Injector ID
164 const label idx =
165 (
166 injIdToIndex_.size()
167 ? injIdToIndex_.lookup(p.typeId(), 0)
168 : 0
169 );
170
171 // See if this patch is designated an outflow patch
172 label addri = -1;
173 forAll(recyclePatchesIds_, i)
174 {
175 if (recyclePatchesIds_[i].first() == pp.index())
176 {
177 addri = i;
178 break;
179 }
180 }
181
182 if (addri == -1)
183 {
184 // Not processing this outflow patch
185 keepParticle = true;
186 return false;
187 }
188
189 // Flag to remove current parcel and copy to local storage
190 keepParticle = false;
191 recycledParcels_[addri].append
192 (
193 static_cast<parcelType*>(p.clone().ptr())
194 );
195
196 ++nRemoved_[addri][idx];
197 massRemoved_[addri][idx] += p.nParticle()*p.mass();
198
199 return true;
200}
201
202
203template<class CloudType>
205{
206 if (Pstream::parRun())
207 {
208 // See comments in Cloud::move() about transfer particles handling
209
210 // Allocate transfer buffers
212
213 // Cache of opened UOPstream wrappers
214 PtrList<UOPstream> UOPstreamPtrs(Pstream::nProcs());
215
216 auto& rnd = this->owner().rndGen();
217
218 forAll(recycledParcels_, addri)
219 {
220 auto& patchParcels = recycledParcels_[addri];
221 auto& injectionPatch = injectionPatchPtr_[addri];
222
223 for (parcelType& p : patchParcels)
224 {
225 // Choose a random location to insert the parcel
226 const scalar fraction01 = rnd.template sample01<scalar>();
227
228 // Identify the processor that owns the location
229 const label toProci = injectionPatch.whichProc(fraction01);
230
231 // Get/create output stream
232 auto* osptr = UOPstreamPtrs.get(toProci);
233 if (!osptr)
234 {
235 osptr = new UOPstream(toProci, pBufs);
236 UOPstreamPtrs.set(toProci, osptr);
237 }
238
239 // Tuple: (address fraction particle)
240 (*osptr) << addri << fraction01 << p;
241
242 // Can now remove from list and delete
243 delete(patchParcels.remove(&p));
244 }
245 }
246
247 pBufs.finishedSends();
248
249 if (!returnReduce(pBufs.hasRecvData(), orOp<bool>()))
250 {
251 // No parcels to recycle
252 return;
253 }
254
255 // Retrieve from receive buffers
256 for (const int proci : pBufs.allProcs())
257 {
258 if (pBufs.recvDataCount(proci))
259 {
260 UIPstream is(proci, pBufs);
261
262 // Read out each (address fraction particle) tuple
263 while (!is.eof())
264 {
265 const label addri = pTraits<label>(is);
266 const scalar fraction01 = pTraits<scalar>(is);
267 auto* newp = new parcelType(this->owner().mesh(), is);
268
269 // Parcel to be recycled
270 vector newPosition;
271 label cellOwner;
272 label dummy;
273 injectionPatchPtr_[addri].setPositionAndCell
274 (
275 mesh_,
276 fraction01,
277 this->owner().rndGen(),
278 newPosition,
279 cellOwner,
280 dummy,
281 dummy
282 );
283 newp->relocate(newPosition, cellOwner);
284 newp->nParticle() *= recycleFraction_;
285
286 // Assume parcel velocity is same as the carrier velocity
287 newp->U() = this->owner().U()[cellOwner];
288
289 // Injector ID
290 const label idx =
291 (
292 injIdToIndex_.size()
293 ? injIdToIndex_.lookup(newp->typeId(), 0)
294 : 0
295 );
296
297 ++nInjected_[addri][idx];
298 massInjected_[addri][idx] += newp->nParticle()*newp->mass();
299
300 this->owner().addParticle(newp);
301 }
302 }
303 }
304 }
305 else
306 {
307 forAll(recycledParcels_, addri)
308 {
309 for (parcelType& p : recycledParcels_[addri])
310 {
311 parcelType* newp = recycledParcels_[addri].remove(&p);
312
313 // Parcel to be recycled
314 vector newPosition;
315 label cellOwner;
316 label dummy;
317 injectionPatchPtr_[addri].setPositionAndCell
318 (
319 mesh_,
320 this->owner().rndGen(),
321 newPosition,
322 cellOwner,
323 dummy,
324 dummy
325 );
326
327 // Update parcel properties
328 newp->relocate(newPosition, cellOwner);
329 newp->nParticle() *= recycleFraction_;
330
331 // Assume parcel velocity is same as the carrier velocity
332 newp->U() = this->owner().U()[cellOwner];
333
334 // Injector ID
335 const label idx =
336 (
337 injIdToIndex_.size()
338 ? injIdToIndex_.lookup(newp->typeId(), 0)
339 : 0
340 );
341 ++nInjected_[addri][idx];
342 massInjected_[addri][idx] += newp->nParticle()*newp->mass();
343
344 this->owner().addParticle(newp);
345 }
346 }
347 }
348}
349
350
351template<class CloudType>
353{
355
356 labelListList npr0(nRemoved_.size());
357 scalarListList mpr0(massRemoved_.size());
358 labelListList npi0(nInjected_.size());
359 scalarListList mpi0(massInjected_.size());
360
361 forAll(nRemoved_, patchi)
362 {
363 const label lsd = nRemoved_[patchi].size();
364 npr0[patchi].setSize(lsd, Zero);
365 mpr0[patchi].setSize(lsd, Zero);
366 npi0[patchi].setSize(lsd, Zero);
367 mpi0[patchi].setSize(lsd, Zero);
368 }
369
370 this->getModelProperty("nRemoved", npr0);
371 this->getModelProperty("massRemoved", mpr0);
372 this->getModelProperty("nInjected", npi0);
373 this->getModelProperty("massInjected", mpi0);
374
375 // Accumulate current data
376 labelListList npr(nRemoved_);
377
378 forAll(npr, i)
379 {
381 npr[i] = npr[i] + npr0[i];
382 }
383
384 scalarListList mpr(massRemoved_);
385 forAll(mpr, i)
386 {
388 mpr[i] = mpr[i] + mpr0[i];
389 }
390
391 labelListList npi(nInjected_);
392 forAll(npi, i)
393 {
395 npi[i] = npi[i] + npi0[i];
396 }
397
398 scalarListList mpi(massInjected_);
399 forAll(mpi, i)
400 {
402 mpi[i] = mpi[i] + mpi0[i];
403 }
404
405 if (injIdToIndex_.size())
406 {
407 // Since injIdToIndex_ is a one-to-one mapping (starting as zero),
408 // can simply invert it.
409 labelList indexToInjector(injIdToIndex_.size());
410 forAllConstIters(injIdToIndex_, iter)
411 {
412 indexToInjector[iter.val()] = iter.key();
413 }
414
415 forAll(npr, i)
416 {
417 const word& outPatchName = recyclePatches_[i].first();
418
419 os << " Parcel fate: patch " << outPatchName
420 << " (number, mass)" << nl;
421
422 forAll(mpr[i], indexi)
423 {
424 os << " - removed (injector " << indexToInjector[indexi]
425 << ") = " << npr[i][indexi]
426 << ", " << mpr[i][indexi] << nl;
427
428 this->file()
429 << tab << npr[i][indexi] << tab << mpr[i][indexi];
430 }
431
432 const word& inPatchName = recyclePatches_[i].second();
433
434 os << " Parcel fate: patch " << inPatchName
435 << " (number, mass)" << nl;
436
437 forAll(mpi[i], indexi)
438 {
439 os << " - injected (injector " << indexToInjector[indexi]
440 << ") = " << npi[i][indexi]
441 << ", " << mpi[i][indexi] << nl;
442 this->file()
443 << tab << npi[i][indexi] << tab << mpi[i][indexi];
444 }
445 }
446
447 this->file() << endl;
448 }
449 else
450 {
451 forAll(npr, i)
452 {
453 const word& outPatchName = recyclePatches_[i].first();
454
455 os << " Parcel fate: patch " << outPatchName
456 << " (number, mass)" << nl
457 << " - removed = " << npr[i][0] << ", " << mpr[i][0]
458 << nl;
459
460 this->file()
461 << tab << npr[i][0] << tab << mpr[i][0];
462 }
463
464 forAll(npi, i)
465 {
466 const word& inPatchName = recyclePatches_[i].second();
467
468 os << " Parcel fate: patch " << inPatchName
469 << " (number, mass)" << nl
470 << " - injected = " << npi[i][0] << ", " << mpi[i][0]
471 << nl;
472
473 this->file()
474 << tab << npi[i][0] << tab << mpi[i][0];
475 }
476
477 this->file() << endl;
478 }
479
480 if (this->writeTime())
481 {
482 this->setModelProperty("nRemoved", npr);
483 this->setModelProperty("massRemoved", mpr);
484 this->setModelProperty("nInjected", npi);
485 this->setModelProperty("massInjected", mpi);
486
487 nRemoved_ = Zero;
488 massRemoved_ = Zero;
489 nInjected_ = Zero;
490 massInjected_ = Zero;
491 }
492}
493
494
495// ************************************************************************* //
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
tmp< GeometricField< Type, PatchField, GeoMesh > > clone() const
Clone.
bool empty() const noexcept
True if the hash table is empty.
Definition: HashTableI.H:59
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:180
bool eof() const noexcept
True if end of input seen.
Definition: IOstream.H:239
void setSize(const label n)
Alias for resize()
Definition: List.H:218
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Templated patch interaction model class.
Buffers for inter-processor communications streams (UOPstream, UIPstream).
UPstream::rangeType allProcs() const noexcept
Range of ranks indices associated with PstreamBuffers.
label nProcs() const noexcept
Number of ranks associated with PstreamBuffers.
bool hasRecvData() const
label recvDataCount(const label proci) const
void finishedSends(const bool wait=true)
Mark sends as done.
static void listCombineGather(const List< commsStruct > &comms, List< T > &values, const CombineOp &cop, const int tag, const label comm)
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
const T * set(const label i) const
Definition: PtrList.H:138
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Patch interaction model to perform parcel transfer and recycle from one patch to another.
virtual void postEvolve()
Post-evolve hook.
const fvMesh & mesh_
Reference to mesh.
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.
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 writeFileHeader(Ostream &os)
Output file header information.
List< List< scalar > > massRemoved_
Mass of parcels removed.
CloudType::parcelType parcelType
bool outputByInjectorId_
Flag to output escaped/mass particles sorted by injectorID.
PtrList< patchInjectionBase > injectionPatchPtr_
Injection patch pointer.
@ nonBlocking
"nonBlocking"
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
T * get(const label i)
Definition: UPtrListI.H:127
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:60
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
InfoProxy< ensightCells > info() const
Return info proxy.
Definition: ensightCells.H:254
A traits class, which is primarily used for primitives.
Definition: pTraits.H:59
label index() const noexcept
The index of this patch in the boundaryMesh.
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
Lookup type of boundary radiation properties.
Definition: lookup.H:66
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
constexpr char tab
The tab '\t' character(0x09)
Definition: Ostream.H:52
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278
Random rndGen
Definition: createFields.H:23