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