Loading...
Searching...
No Matches
InteractionLists.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2019-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "InteractionLists.H"
31#include "indexedOctree.H"
32#include "treeDataCell.H"
33#include "treeDataFace.H"
34#include "volFields.H"
35#include "meshTools.H"
36
37// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
38
39template<class ParticleType>
40void Foam::InteractionLists<ParticleType>::buildInteractionLists()
41{
42 Info<< "Building InteractionLists with interaction distance "
43 << maxDistance_ << endl;
44
45 Random rndGen(419715);
46
47 treeBoundBox procBb(mesh_.points());
48
49 treeBoundBoxList allExtendedProcBbs(Pstream::nProcs());
50
51 allExtendedProcBbs[Pstream::myProcNo()] = procBb;
52 allExtendedProcBbs[Pstream::myProcNo()].grow(maxDistance_);
53
54 Pstream::allGatherList(allExtendedProcBbs);
55
56 List<treeBoundBox> extendedProcBbsInRange;
57 List<label> extendedProcBbsTransformIndex;
58 List<label> extendedProcBbsOrigProc;
59
60 findExtendedProcBbsInRange
61 (
62 procBb,
63 allExtendedProcBbs,
64 mesh_.globalData().globalTransforms(),
65 extendedProcBbsInRange,
66 extendedProcBbsTransformIndex,
67 extendedProcBbsOrigProc
68 );
69
70 treeBoundBoxList cellBbs(treeDataCell::boxes(mesh_));
71
72 const globalIndexAndTransform& globalTransforms =
73 mesh_.globalData().globalTransforms();
74
75 // Recording which cells are in range of an extended boundBox, as
76 // only these cells will need to be tested to determine which
77 // referred cells that they interact with.
78 bitSet cellInRangeOfCoupledPatch(mesh_.nCells(), false);
79
80 // IAndT: index (=local cell index) and transform (from
81 // globalIndexAndTransform)
82 DynamicList<labelPair> cellIAndTToExchange;
83
84 DynamicList<treeBoundBox> cellBbsToExchange;
85
86 DynamicList<label> procToDistributeCellTo;
87
88 forAll(extendedProcBbsInRange, ePBIRI)
89 {
90 const treeBoundBox& otherExtendedProcBb =
91 extendedProcBbsInRange[ePBIRI];
92
93 label transformIndex = extendedProcBbsTransformIndex[ePBIRI];
94
95 label origProc = extendedProcBbsOrigProc[ePBIRI];
96
97 forAll(cellBbs, celli)
98 {
99 const treeBoundBox& cellBb = cellBbs[celli];
100
101 if (cellBb.overlaps(otherExtendedProcBb))
102 {
103 // This cell is in range of the Bb of the other
104 // processor Bb, and so needs to be referred to it
105
106 cellInRangeOfCoupledPatch.set(celli);
107
108 cellIAndTToExchange.append
109 (
110 globalTransforms.encode(celli, transformIndex)
111 );
112
113 cellBbsToExchange.append(cellBb);
114
115 procToDistributeCellTo.append(origProc);
116 }
117 }
118 }
119
120 buildMap(cellMapPtr_, procToDistributeCellTo);
121
122 // Needed for reverseDistribute
123 label preDistributionCellMapSize = procToDistributeCellTo.size();
124
125 cellMap().distribute(cellBbsToExchange);
126
127 cellMap().distribute(cellIAndTToExchange);
128
129 // Determine labelList specifying only cells that are in range of
130 // a coupled boundary to build an octree limited to these cells.
131 DynamicList<label> coupledPatchRangeCells;
132
133 forAll(cellInRangeOfCoupledPatch, celli)
134 {
135 if (cellInRangeOfCoupledPatch[celli])
136 {
137 coupledPatchRangeCells.append(celli);
138 }
139 }
140
141 treeBoundBox procBbRndExt
142 (
143 treeBoundBox(mesh_.points()).extend(rndGen, 1e-4)
144 );
145
146 indexedOctree<treeDataCell> coupledPatchRangeTree
147 (
148 treeDataCell
149 (
150 true, // Cache cell bb
151 mesh_,
152 coupledPatchRangeCells, // Subset of mesh
153 polyMesh::CELL_TETS // Consistent with tracking
154 ),
155 procBbRndExt,
156 8, // maxLevel,
157 10, // leafSize,
158 100.0 // duplicity
159 );
160
161 ril_.setSize(cellBbsToExchange.size());
162
163 // This needs to be a boolList, not bitSet if
164 // reverseDistribute is called.
165 boolList cellBbRequiredByAnyCell(cellBbsToExchange.size(), false);
166
167 Info<< " Building referred interaction lists" << endl;
168
169 forAll(cellBbsToExchange, bbI)
170 {
171 const labelPair& ciat = cellIAndTToExchange[bbI];
172
173 const vectorTensorTransform& transform = globalTransforms.transform
174 (
175 globalTransforms.transformIndex(ciat)
176 );
177
178 treeBoundBox extendedBb
179 (
180 transform.invTransformPosition(cellBbsToExchange[bbI].points())
181 );
182 extendedBb.grow(maxDistance_);
183
184 // Find all elements intersecting box.
185 labelList interactingElems
186 (
187 coupledPatchRangeTree.findBox(extendedBb)
188 );
189
190 if (!interactingElems.empty())
191 {
192 cellBbRequiredByAnyCell[bbI] = true;
193 }
194
195 ril_[bbI].setSize(interactingElems.size(), -1);
196
197 forAll(interactingElems, i)
198 {
199 label elemI = interactingElems[i];
200
201 // Here, a more detailed geometric test could be applied,
202 // i.e. a more accurate bounding volume like a OBB or
203 // convex hull or an exact geometrical test.
204
205 label c = coupledPatchRangeTree.shapes().objectIndex(elemI);
206
207 ril_[bbI][i] = c;
208 }
209 }
210
211 // Perform subset of ril_, to remove any referred cells that do
212 // not interact. They will not be sent from originating
213 // processors. This assumes that the disappearance of the cell
214 // from the sending list of the source processor, simply removes
215 // the referred cell from the ril_, all of the subsequent indices
216 // shuffle down one, but the order and structure is preserved,
217 // i.e. it, is as if the cell had never been referred in the first
218 // place.
219
220 inplaceSubset(cellBbRequiredByAnyCell, ril_);
221
222 // Send information about which cells are actually required back
223 // to originating processors.
224
225 // At this point, cellBbsToExchange does not need to be maintained
226 // or distributed as it is not longer needed.
227
228 cellBbsToExchange.clearStorage();
229
230 cellMap().reverseDistribute
231 (
232 preDistributionCellMapSize,
233 cellBbRequiredByAnyCell
234 );
235
236 cellMap().reverseDistribute
237 (
238 preDistributionCellMapSize,
239 cellIAndTToExchange
240 );
241
242 // Perform ordering of cells to send, this invalidates the
243 // previous value of preDistributionCellMapSize.
244
245 preDistributionCellMapSize = -1;
246
247 // Move all of the used cells to refer to the start of the list
248 // and truncate it
249
250 inplaceSubset(cellBbRequiredByAnyCell, cellIAndTToExchange);
251
252 inplaceSubset(cellBbRequiredByAnyCell, procToDistributeCellTo);
253
254 preDistributionCellMapSize = procToDistributeCellTo.size();
255
256 // Rebuild mapDistribute with only required referred cells
257 buildMap(cellMapPtr_, procToDistributeCellTo);
258
259 // Store cellIndexAndTransformToDistribute
260 cellIndexAndTransformToDistribute_.transfer(cellIAndTToExchange);
261
262 // Determine inverse addressing for referred cells
263
264 rilInverse_.setSize(mesh_.nCells());
265
266 // Temporary Dynamic lists for accumulation
267 List<DynamicList<label>> rilInverseTemp(rilInverse_.size());
268
269 // Loop over all referred cells
270 forAll(ril_, refCelli)
271 {
272 const labelList& realCells = ril_[refCelli];
273
274 // Loop over all real cells in that the referred cell is to
275 // supply interactions to and record the index of this
276 // referred cell in the real cells entry in rilInverse
277
278 forAll(realCells, realCelli)
279 {
280 rilInverseTemp[realCells[realCelli]].append(refCelli);
281 }
282 }
283
284 forAll(rilInverse_, celli)
285 {
286 rilInverse_[celli].transfer(rilInverseTemp[celli]);
287 }
288
289 // Determine which wall faces to refer
290
291 // The referring of wall patch data relies on patches having the same
292 // index on each processor.
293 mesh_.boundaryMesh().checkParallelSync(true);
294
295 // Determine the index of all of the wall faces on this processor
296 DynamicList<label> localWallFaces;
298 for (const polyPatch& patch : mesh_.boundaryMesh())
299 {
300 if (isA<wallPolyPatch>(patch))
301 {
302 const auto areaFraction(patch.areaFraction());
303
304 forAll(patch, facei)
306 if (!areaFraction || areaFraction()[facei] > 0.5)
307 {
308 localWallFaces.append(facei + patch.start());
309 }
310 }
311 }
312 }
313
314 treeBoundBoxList wallFaceBbs(localWallFaces.size());
315
316 forAll(wallFaceBbs, i)
317 {
318 wallFaceBbs[i] =
320 (
321 mesh_.points(),
322 mesh_.faces()[localWallFaces[i]]
323 );
324 }
325
326 // IAndT: index and transform
327 DynamicList<labelPair> wallFaceIAndTToExchange;
328
329 DynamicList<treeBoundBox> wallFaceBbsToExchange;
330
331 DynamicList<label> procToDistributeWallFaceTo;
332
333 forAll(extendedProcBbsInRange, ePBIRI)
334 {
335 const treeBoundBox& otherExtendedProcBb =
336 extendedProcBbsInRange[ePBIRI];
337
338 label transformIndex = extendedProcBbsTransformIndex[ePBIRI];
339
340 label origProc = extendedProcBbsOrigProc[ePBIRI];
341
342 forAll(wallFaceBbs, i)
343 {
344 const treeBoundBox& wallFaceBb = wallFaceBbs[i];
345
346 if (wallFaceBb.overlaps(otherExtendedProcBb))
347 {
348 // This wall face is in range of the Bb of the other
349 // processor Bb, and so needs to be referred to it
350
351 label wallFacei = localWallFaces[i];
352
353 wallFaceIAndTToExchange.append
354 (
355 globalTransforms.encode(wallFacei, transformIndex)
356 );
357
358 wallFaceBbsToExchange.append(wallFaceBb);
359
360 procToDistributeWallFaceTo.append(origProc);
361 }
362 }
363 }
364
365 buildMap(wallFaceMapPtr_, procToDistributeWallFaceTo);
366
367 // Needed for reverseDistribute
368 label preDistributionWallFaceMapSize = procToDistributeWallFaceTo.size();
369
370 wallFaceMap().distribute(wallFaceBbsToExchange);
371
372 wallFaceMap().distribute(wallFaceIAndTToExchange);
373
374 indexedOctree<treeDataCell> allCellsTree
375 (
376 treeDataCell(true, mesh_, polyMesh::CELL_TETS),
377 procBbRndExt,
378 8, // maxLevel,
379 10, // leafSize,
380 100.0 // duplicity
381 );
382
383 rwfil_.setSize(wallFaceBbsToExchange.size());
384
385 // This needs to be a boolList, not bitSet if
386 // reverseDistribute is called.
387 boolList wallFaceBbRequiredByAnyCell(wallFaceBbsToExchange.size(), false);
388
389 forAll(wallFaceBbsToExchange, bbI)
390 {
391 const labelPair& wfiat = wallFaceIAndTToExchange[bbI];
392
393 const vectorTensorTransform& transform = globalTransforms.transform
394 (
395 globalTransforms.transformIndex(wfiat)
396 );
397
398 treeBoundBox extendedBb
399 (
400 transform.invTransformPosition(wallFaceBbsToExchange[bbI].points())
401 );
402 extendedBb.grow(maxDistance_);
403
404 // Find all elements intersecting box.
405 labelList interactingElems
406 (
407 coupledPatchRangeTree.findBox(extendedBb)
408 );
409
410 if (!interactingElems.empty())
411 {
412 wallFaceBbRequiredByAnyCell[bbI] = true;
413 }
414
415 rwfil_[bbI].setSize(interactingElems.size(), -1);
416
417 forAll(interactingElems, i)
418 {
419 label elemI = interactingElems[i];
420
421 // Here, a more detailed geometric test could be applied,
422 // i.e. a more accurate bounding volume like a OBB or
423 // convex hull or an exact geometrical test.
424
425 label c = coupledPatchRangeTree.shapes().objectIndex(elemI);
426
427 rwfil_[bbI][i] = c;
428 }
429 }
430
431 // Perform subset of rwfil_, to remove any referred wallFaces that do
432 // not interact. They will not be sent from originating
433 // processors. This assumes that the disappearance of the wallFace
434 // from the sending list of the source processor, simply removes
435 // the referred wallFace from the rwfil_, all of the subsequent indices
436 // shuffle down one, but the order and structure is preserved,
437 // i.e. it, is as if the wallFace had never been referred in the first
438 // place.
439
440 inplaceSubset(wallFaceBbRequiredByAnyCell, rwfil_);
441
442 // Send information about which wallFaces are actually required
443 // back to originating processors.
444
445 // At this point, wallFaceBbsToExchange does not need to be
446 // maintained or distributed as it is not longer needed.
447
448 wallFaceBbsToExchange.clear();
449
450 wallFaceMap().reverseDistribute
451 (
452 preDistributionWallFaceMapSize,
453 wallFaceBbRequiredByAnyCell
454 );
455
456 wallFaceMap().reverseDistribute
457 (
458 preDistributionWallFaceMapSize,
459 wallFaceIAndTToExchange
460 );
461
462 // Perform ordering of wallFaces to send, this invalidates the
463 // previous value of preDistributionWallFaceMapSize.
464
465 preDistributionWallFaceMapSize = -1;
466
467 // Move all of the used wallFaces to refer to the start of the
468 // list and truncate it
469
470 inplaceSubset(wallFaceBbRequiredByAnyCell, wallFaceIAndTToExchange);
471
472 inplaceSubset(wallFaceBbRequiredByAnyCell, procToDistributeWallFaceTo);
473
474 preDistributionWallFaceMapSize = procToDistributeWallFaceTo.size();
475
476 // Rebuild mapDistribute with only required referred wallFaces
477 buildMap(wallFaceMapPtr_, procToDistributeWallFaceTo);
478
479 // Store wallFaceIndexAndTransformToDistribute
480 wallFaceIndexAndTransformToDistribute_.transfer(wallFaceIAndTToExchange);
481
482 // Determine inverse addressing for referred wallFaces
483
484 rwfilInverse_.setSize(mesh_.nCells());
485
486 // Temporary Dynamic lists for accumulation
487 List<DynamicList<label>> rwfilInverseTemp(rwfilInverse_.size());
488
489 // Loop over all referred wall faces
490 forAll(rwfil_, refWallFacei)
491 {
492 const labelList& realCells = rwfil_[refWallFacei];
493
494 // Loop over all real cells in that the referred wall face is
495 // to supply interactions to and record the index of this
496 // referred wall face in the real cells entry in rwfilInverse
497
498 forAll(realCells, realCelli)
499 {
500 rwfilInverseTemp[realCells[realCelli]].append(refWallFacei);
501 }
502 }
503
504 forAll(rwfilInverse_, celli)
505 {
506 rwfilInverse_[celli].transfer(rwfilInverseTemp[celli]);
507 }
508
509 // Refer wall faces to the appropriate processor
510 referredWallFaces_.setSize(wallFaceIndexAndTransformToDistribute_.size());
511
512 forAll(referredWallFaces_, rWFI)
513 {
514 const labelPair& wfiat = wallFaceIndexAndTransformToDistribute_[rWFI];
515
516 label wallFaceIndex = globalTransforms.index(wfiat);
517
518 const vectorTensorTransform& transform = globalTransforms.transform
519 (
520 globalTransforms.transformIndex(wfiat)
521 );
522
523 const face& f = mesh_.faces()[wallFaceIndex];
524
525 const label patchi = mesh_.boundaryMesh().patchID(wallFaceIndex);
526
527 referredWallFaces_[rWFI] = referredWallFace
528 (
529 face(identity(f.size())),
530 transform.invTransformPosition(f.points(mesh_.points())),
531 patchi
532 );
533 }
534
535 wallFaceMap().distribute(referredWallFaces_);
536
537 if (writeCloud_)
538 {
539 writeReferredWallFaces();
540 }
541
542 // Direct interaction list and direct wall faces
543
544 Info<< " Building direct interaction lists" << endl;
545
546 indexedOctree<treeDataFace> wallFacesTree
547 (
548 treeDataFace(true, mesh_, localWallFaces),
549 procBbRndExt,
550 8, // maxLevel,
551 10, // leafSize,
552 100.0
553 );
554
555 dil_.setSize(mesh_.nCells());
556
557 dwfil_.setSize(mesh_.nCells());
558
559 forAll(cellBbs, celli)
560 {
561 const treeBoundBox& cellBb = cellBbs[celli];
562
563 treeBoundBox extendedBb(cellBb);
564 extendedBb.grow(maxDistance_);
565
566 // Find all cells intersecting extendedBb
567 labelList interactingElems(allCellsTree.findBox(extendedBb));
568
569 // Reserve space to avoid multiple resizing
570 DynamicList<label> cellDIL(interactingElems.size());
571
572 for (const label elemi : interactingElems)
573 {
574 const label c = allCellsTree.shapes().objectIndex(elemi);
575
576 // Here, a more detailed geometric test could be applied,
577 // i.e. a more accurate bounding volume like a OBB or
578 // convex hull, or an exact geometrical test.
579
580 // The higher index cell is added to the lower index
581 // cell's DIL. A cell is not added to its own DIL.
582 if (c > celli)
583 {
584 cellDIL.append(c);
585 }
586 }
587
588 dil_[celli].transfer(cellDIL);
589
590 // Find all wall faces intersecting extendedBb
591 interactingElems = wallFacesTree.findBox(extendedBb);
592
593 dwfil_[celli].setSize(interactingElems.size(), -1);
594
595 forAll(interactingElems, i)
596 {
597 const label elemi = interactingElems[i];
598
599 const label f = wallFacesTree.shapes().objectIndex(elemi);
600
601 dwfil_[celli][i] = f;
602 }
603 }
604}
605
606
607template<class ParticleType>
608void Foam::InteractionLists<ParticleType>::findExtendedProcBbsInRange
609(
610 const treeBoundBox& procBb,
611 const List<treeBoundBox>& allExtendedProcBbs,
612 const globalIndexAndTransform& globalTransforms,
613 List<treeBoundBox>& extendedProcBbsInRange,
614 List<label>& extendedProcBbsTransformIndex,
615 List<label>& extendedProcBbsOrigProc
616)
617{
618 extendedProcBbsInRange.clear();
619 extendedProcBbsTransformIndex.clear();
620 extendedProcBbsOrigProc.clear();
621
622 DynamicList<treeBoundBox> tmpExtendedProcBbsInRange;
623 DynamicList<label> tmpExtendedProcBbsTransformIndex;
624 DynamicList<label> tmpExtendedProcBbsOrigProc;
625
626 label nTrans = globalTransforms.nIndependentTransforms();
627
628 forAll(allExtendedProcBbs, proci)
629 {
630 labelList permutationIndices(nTrans, Zero);
631
632 if (nTrans == 0 && proci != Pstream::myProcNo())
633 {
634 treeBoundBox extendedReferredProcBb = allExtendedProcBbs[proci];
635
636 if (procBb.overlaps(extendedReferredProcBb))
637 {
638 tmpExtendedProcBbsInRange.append(extendedReferredProcBb);
639
640 // Dummy index, there are no transforms, so there will
641 // be no resultant transform when this is decoded.
642 tmpExtendedProcBbsTransformIndex.append(0);
643
644 tmpExtendedProcBbsOrigProc.append(proci);
645 }
646 }
647 else if (nTrans == 3)
648 {
649 label& i = permutationIndices[0];
650 label& j = permutationIndices[1];
651 label& k = permutationIndices[2];
652
653 for (i = -1; i <= 1; i++)
654 {
655 for (j = -1; j <= 1; j++)
656 {
657 for (k = -1; k <= 1; k++)
658 {
659 if
660 (
661 i == 0
662 && j == 0
663 && k == 0
664 && proci == Pstream::myProcNo()
665 )
666 {
667 // Skip this processor's extended boundBox
668 // when it has no transformation
669 continue;
670 }
671
672 label transI = globalTransforms.encodeTransformIndex
673 (
674 permutationIndices
675 );
676
677 const vectorTensorTransform& transform =
678 globalTransforms.transform(transI);
679
680 treeBoundBox extendedReferredProcBb
681 (
682 transform.transformPosition
683 (
684 allExtendedProcBbs[proci].points()
685 )
686 );
687
688 if (procBb.overlaps(extendedReferredProcBb))
689 {
690 tmpExtendedProcBbsInRange.append
691 (
692 extendedReferredProcBb
693 );
694
695 tmpExtendedProcBbsTransformIndex.append(transI);
696
697 tmpExtendedProcBbsOrigProc.append(proci);
698 }
699 }
700 }
701 }
702 }
703 else if (nTrans == 2)
704 {
705 label& i = permutationIndices[0];
706 label& j = permutationIndices[1];
707
708 for (i = -1; i <= 1; i++)
709 {
710 for (j = -1; j <= 1; j++)
711 {
712 if (i == 0 && j == 0 && proci == Pstream::myProcNo())
713 {
714 // Skip this processor's extended boundBox
715 // when it has no transformation
716 continue;
717 }
718
719 label transI = globalTransforms.encodeTransformIndex
720 (
721 permutationIndices
722 );
723
724 const vectorTensorTransform& transform =
725 globalTransforms.transform(transI);
726
727 treeBoundBox extendedReferredProcBb
728 (
729 transform.transformPosition
730 (
731 allExtendedProcBbs[proci].points()
732 )
733 );
734
735 if (procBb.overlaps(extendedReferredProcBb))
736 {
737 tmpExtendedProcBbsInRange.append
738 (
739 extendedReferredProcBb
740 );
741
742 tmpExtendedProcBbsTransformIndex.append(transI);
743
744 tmpExtendedProcBbsOrigProc.append(proci);
745 }
746 }
747 }
748 }
749 else if (nTrans == 1)
750 {
751 label& i = permutationIndices[0];
752
753 for (i = -1; i <= 1; i++)
754 {
755 if (i == 0 && proci == Pstream::myProcNo())
756 {
757 // Skip this processor's extended boundBox when it
758 // has no transformation
759 continue;
760 }
761
762 label transI = globalTransforms.encodeTransformIndex
763 (
764 permutationIndices
765 );
766
767 const vectorTensorTransform& transform =
768 globalTransforms.transform(transI);
769
770 treeBoundBox extendedReferredProcBb
771 (
772 transform.transformPosition
773 (
774 allExtendedProcBbs[proci].points()
775 )
776 );
777
778 if (procBb.overlaps(extendedReferredProcBb))
779 {
780 tmpExtendedProcBbsInRange.append
781 (
782 extendedReferredProcBb
783 );
784
785 tmpExtendedProcBbsTransformIndex.append(transI);
786
787 tmpExtendedProcBbsOrigProc.append(proci);
788 }
789 }
790 }
791 }
792
793 extendedProcBbsInRange.transfer(tmpExtendedProcBbsInRange);
794 extendedProcBbsTransformIndex.transfer(tmpExtendedProcBbsTransformIndex);
795 extendedProcBbsOrigProc.transfer(tmpExtendedProcBbsOrigProc);
796}
797
798
799template<class ParticleType>
800void Foam::InteractionLists<ParticleType>::buildMap
801(
802 autoPtr<mapDistribute>& mapPtr,
803 const List<label>& toProc
804)
805{
806 // Determine send map
807 // ~~~~~~~~~~~~~~~~~~
808
809 // 1. Count
810 labelList nSend(Pstream::nProcs(), Zero);
811
812 for (const label proci : toProc)
813 {
814 nSend[proci]++;
815 }
816
817 // 2. Size sendMap
818 labelListList sendMap(Pstream::nProcs());
819
820 forAll(nSend, proci)
821 {
822 sendMap[proci].resize_nocopy(nSend[proci]);
823 nSend[proci] = 0;
824 }
825
826 // 3. Fill sendMap
827 forAll(toProc, i)
828 {
829 label proci = toProc[i];
830 sendMap[proci][nSend[proci]++] = i;
831 }
832
833 mapPtr.reset(new mapDistribute(std::move(sendMap)));
834}
835
836
837template<class ParticleType>
838void Foam::InteractionLists<ParticleType>::prepareParticlesToRefer
839(
840 const List<DynamicList<ParticleType*>>& cellOccupancy
841)
842{
843 const globalIndexAndTransform& globalTransforms =
844 mesh_.globalData().globalTransforms();
845
846 referredParticles_.setSize(cellIndexAndTransformToDistribute_.size());
847
848 // Clear all existing referred particles
849
850 forAll(referredParticles_, i)
851 {
852 referredParticles_[i].clear();
853 }
854
855 // Clear all particles that may have been populated into the cloud
856 cloud_.clear();
857
858 forAll(cellIndexAndTransformToDistribute_, i)
859 {
860 const labelPair ciat = cellIndexAndTransformToDistribute_[i];
861
862 label cellIndex = globalTransforms.index(ciat);
863
864 List<ParticleType*> realParticles = cellOccupancy[cellIndex];
865
866 IDLList<ParticleType>& particlesToRefer = referredParticles_[i];
867
868 forAll(realParticles, rM)
869 {
870 const ParticleType& particle = *realParticles[rM];
871
872 particlesToRefer.append(particle.clone().ptr());
873
874 prepareParticleToBeReferred(particlesToRefer.last(), ciat);
875 }
876 }
877}
878
879
880template<class ParticleType>
881void Foam::InteractionLists<ParticleType>::prepareParticleToBeReferred
882(
883 ParticleType* particle,
884 labelPair ciat
885)
886{
887 const globalIndexAndTransform& globalTransforms =
888 mesh_.globalData().globalTransforms();
889
890 const vectorTensorTransform& transform = globalTransforms.transform
891 (
892 globalTransforms.transformIndex(ciat)
893 );
894
895 particle->prepareForInteractionListReferral(transform);
896}
897
898
899template<class ParticleType>
900void Foam::InteractionLists<ParticleType>::fillReferredParticleCloud()
901{
902 if (writeCloud_)
903 {
904 forAll(referredParticles_, refCelli)
905 {
906 const IDLList<ParticleType>& refCell =
907 referredParticles_[refCelli];
908
909 for (const ParticleType& p : refCell)
910 {
911 cloud_.addParticle
912 (
913 static_cast<ParticleType*>(p.clone().ptr())
914 );
915 }
916 }
917 }
918}
919
920
921template<class ParticleType>
922void Foam::InteractionLists<ParticleType>::prepareWallDataToRefer()
923{
924 const globalIndexAndTransform& globalTransforms =
925 mesh_.globalData().globalTransforms();
926
927 referredWallData_.setSize
928 (
929 wallFaceIndexAndTransformToDistribute_.size()
930 );
931
932 const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
933
934 forAll(referredWallData_, rWVI)
935 {
936 const labelPair& wfiat = wallFaceIndexAndTransformToDistribute_[rWVI];
937
938 label wallFaceIndex = globalTransforms.index(wfiat);
939
940 const vectorTensorTransform& transform = globalTransforms.transform
941 (
942 globalTransforms.transformIndex(wfiat)
943 );
944
945 const label patchi = mesh_.boundaryMesh().patchID(wallFaceIndex);
946
947 const label patchFacei =
948 mesh_.boundaryMesh()[patchi].whichFace(wallFaceIndex);
949
950 // Need to transform velocity when tensor transforms are
951 // supported
952 referredWallData_[rWVI] = U.boundaryField()[patchi][patchFacei];
953
954 if (transform.hasR())
955 {
956 referredWallData_[rWVI] =
957 transform.R().T() & referredWallData_[rWVI];
958 }
959 }
960}
961
962
963template<class ParticleType>
964void Foam::InteractionLists<ParticleType>::writeReferredWallFaces() const
965{
966 if (referredWallFaces_.empty())
967 {
968 return;
969 }
970
971 fileName objDir = mesh_.time().timePath()/cloud::prefix;
972
973 mkDir(objDir);
974
975 fileName objFileName = "referredWallFaces.obj";
976
977 OFstream str(objDir/objFileName);
978
979 Info<< " Writing "
980 << mesh_.time().timeName()/cloud::prefix/objFileName
981 << endl;
982
983 label offset = 1;
984
985 forAll(referredWallFaces_, rWFI)
986 {
987 const referredWallFace& rwf = referredWallFaces_[rWFI];
988
989 forAll(rwf, fPtI)
990 {
991 meshTools::writeOBJ(str, rwf.points()[rwf[fPtI]]);
992 }
993
994 str<< 'f';
995
996 forAll(rwf, fPtI)
997 {
998 str<< ' ' << fPtI + offset;
999 }
1000
1001 str<< nl;
1002
1003 offset += rwf.size();
1005}
1006
1007
1008// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1009
1010template<class ParticleType>
1011Foam::InteractionLists<ParticleType>::InteractionLists(const polyMesh& mesh)
1012:
1013 mesh_(mesh),
1014 cloud_(mesh_, "nullptr_Cloud", IDLList<ParticleType>()),
1015 writeCloud_(false),
1016 cellMapPtr_(),
1017 wallFaceMapPtr_(),
1018 maxDistance_(0.0),
1019 dil_(),
1020 dwfil_(),
1021 ril_(),
1022 rilInverse_(),
1023 cellIndexAndTransformToDistribute_(),
1024 wallFaceIndexAndTransformToDistribute_(),
1025 referredWallFaces_(),
1026 UName_("unknown_U"),
1027 referredWallData_(),
1028 referredParticles_()
1029{}
1030
1031
1032template<class ParticleType>
1033Foam::InteractionLists<ParticleType>::InteractionLists
1034(
1035 const polyMesh& mesh,
1036 scalar maxDistance,
1037 bool writeCloud,
1038 const word& UName
1039)
1040:
1041 mesh_(mesh),
1042 cloud_(mesh_, "referredParticleCloud", IDLList<ParticleType>()),
1043 writeCloud_(writeCloud),
1044 cellMapPtr_(),
1045 wallFaceMapPtr_(),
1046 maxDistance_(maxDistance),
1047 dil_(),
1048 dwfil_(),
1049 ril_(),
1050 rilInverse_(),
1051 cellIndexAndTransformToDistribute_(),
1052 wallFaceIndexAndTransformToDistribute_(),
1053 referredWallFaces_(),
1054 UName_(UName),
1055 referredWallData_(),
1056 referredParticles_()
1057{
1058 buildInteractionLists();
1059}
1060
1061
1062// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1063
1064template<class ParticleType>
1066{}
1067
1068
1069// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
1070
1071template<class ParticleType>
1073(
1075 PstreamBuffers& pBufs
1076)
1077{
1078 if (mesh_.changing())
1079 {
1081 << "Mesh changing, rebuilding InteractionLists from scratch."
1082 << endl;
1083
1084 buildInteractionLists();
1085 }
1086
1087 prepareWallDataToRefer();
1088
1089 prepareParticlesToRefer(cellOccupancy);
1090
1091 for (const int domain : Pstream::allProcs())
1092 {
1093 const labelList& subMap = cellMap().subMap()[domain];
1094
1095 if (subMap.size())
1096 {
1097 UOPstream toDomain(domain, pBufs);
1098
1099 UIndirectList<IDLList<ParticleType>> subMappedParticles
1100 (
1101 referredParticles_,
1102 subMap
1103 );
1104
1105 forAll(subMappedParticles, i)
1106 {
1107 toDomain << subMappedParticles[i];
1108 }
1109 }
1110 }
1111
1112 // Using the mapDistribute to start sending and receiving the
1113 // buffer but not block, i.e. it is calling
1114 // pBufs.finishedSends(false);
1115 wallFaceMap().send(pBufs, referredWallData_);
1116}
1117
1118
1119template<class ParticleType>
1121(
1122 PstreamBuffers& pBufs,
1123 const label startOfRequests
1124)
1125{
1126 Pstream::waitRequests(startOfRequests);
1127
1128 referredParticles_.setSize(cellMap().constructSize());
1129
1130 for (const int domain : Pstream::allProcs())
1131 {
1132 const labelList& constructMap = cellMap().constructMap()[domain];
1133
1134 if (constructMap.size())
1135 {
1136 UIPstream str(domain, pBufs);
1137
1138 forAll(constructMap, i)
1139 {
1140 referredParticles_[constructMap[i]] = IDLList<ParticleType>
1141 (
1142 str,
1143 typename ParticleType::iNew(mesh_)
1144 );
1145 }
1146 }
1147 }
1148
1149 forAll(referredParticles_, refCelli)
1150 {
1151 IDLList<ParticleType>& refCell = referredParticles_[refCelli];
1152 for (ParticleType& p : refCell)
1153 {
1154 p.correctAfterInteractionListReferral(ril_[refCelli][0]);
1155 }
1156 }
1157
1158 fillReferredParticleCloud();
1159
1160 wallFaceMap().receive(pBufs, referredWallData_);
1161}
1162
1163
1164// ************************************************************************* //
label k
const word UName(propsDict.getOrDefault< word >("U", "U"))
const List< DynamicList< molecule * > > & cellOccupancy
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
void append(const T &val)
Copy append an element to the end of this list.
const mapDistribute & cellMap() const
Return access to the cellMap.
const polyMesh & mesh() const
Return access to the mesh.
const mapDistribute & wallFaceMap() const
Return access to the wallFaceMap.
void sendReferredData(const List< DynamicList< ParticleType * > > &cellOccupancy, PstreamBuffers &pBufs)
Prepare and send referred particles and wall data (nonBlocking).
const word & UName() const
Return the name of the velocity field.
void receiveReferredData(PstreamBuffers &pBufs, const label startReq=0)
Receive referred data.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void append(const T &val)
Append an element at the end of the list.
Definition List.H:497
Buffers for inter-processor communications streams (UOPstream, UIPstream).
Input inter-processor communications stream using MPI send/recv etc. - operating on external buffer.
Definition UIPstream.H:313
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Output inter-processor communications stream using MPI send/recv etc. - operating on external buffer.
Definition UOPstream.H:408
static rangeType allProcs(const label communicator=worldComm)
Range of process indices for all processes.
Definition UPstream.H:1857
static void waitRequests()
Wait for all requests to finish.
Definition UPstream.H:2497
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundaryMesh reference.
Definition polyPatch.C:295
void transform(const polyMesh &, const tensor &, TrackingData &td)
Apply rotation matrix to any coordinates.
Standard boundBox with extra functionality for use in octree.
bool overlaps(const boundBox &bb) const
Overlaps with other bounding box, sphere etc?
Definition boundBoxI.H:439
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
volScalarField & p
dynamicFvMesh & mesh
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar c
Speed of light in a vacuum.
const std::string patch
OpenFOAM patch number as a std::string.
Pair< label > labelPair
A pair of labels.
Definition Pair.H:54
ILList< DLListBase, T > IDLList
Definition IDLList.H:39
void inplaceSubset(const BoolListType &select, ListType &input, const bool invert=false)
Inplace extract elements of the input list when select is true.
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
List< treeBoundBox > treeBoundBoxList
A List of treeBoundBox.
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
List< bool > boolList
A List of bools.
Definition List.H:60
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
labelList f(nPoints)
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
mkDir(pdfPath)
Random rndGen