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