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 | Copyright (C) 2011-2016 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  Random rndGen(419715);
43 
44  const vector interactionVec = maxDistance_*vector::one;
45 
46  treeBoundBox procBb(treeBoundBox(mesh_.points()));
47 
48  treeBoundBox extendedProcBb
49  (
50  procBb.min() - interactionVec,
51  procBb.max() + interactionVec
52  );
53 
54  treeBoundBoxList allExtendedProcBbs(Pstream::nProcs());
55 
56  allExtendedProcBbs[Pstream::myProcNo()] = extendedProcBb;
57 
58  Pstream::gatherList(allExtendedProcBbs);
59 
60  Pstream::scatterList(allExtendedProcBbs);
61 
62  List<treeBoundBox> extendedProcBbsInRange;
63  List<label> extendedProcBbsTransformIndex;
64  List<label> extendedProcBbsOrigProc;
65 
66  findExtendedProcBbsInRange
67  (
68  procBb,
69  allExtendedProcBbs,
70  mesh_.globalData().globalTransforms(),
71  extendedProcBbsInRange,
72  extendedProcBbsTransformIndex,
73  extendedProcBbsOrigProc
74  );
75 
76  treeBoundBoxList cellBbs(mesh_.nCells());
77 
78  forAll(cellBbs, celli)
79  {
80  cellBbs[celli] = treeBoundBox
81  (
82  mesh_.cells()[celli].points
83  (
84  mesh_.faces(),
85  mesh_.points()
86  )
87  );
88  }
89 
90  const globalIndexAndTransform& globalTransforms =
91  mesh_.globalData().globalTransforms();
92 
93  // Recording which cells are in range of an extended boundBox, as
94  // only these cells will need to be tested to determine which
95  // referred cells that they interact with.
96  PackedBoolList cellInRangeOfCoupledPatch(mesh_.nCells(), false);
97 
98  // IAndT: index (=local cell index) and transform (from
99  // globalIndexAndTransform)
100  DynamicList<labelPair> cellIAndTToExchange;
101 
102  DynamicList<treeBoundBox> cellBbsToExchange;
103 
104  DynamicList<label> procToDistributeCellTo;
105 
106  forAll(extendedProcBbsInRange, ePBIRI)
107  {
108  const treeBoundBox& otherExtendedProcBb =
109  extendedProcBbsInRange[ePBIRI];
110 
111  label transformIndex = extendedProcBbsTransformIndex[ePBIRI];
112 
113  label origProc = extendedProcBbsOrigProc[ePBIRI];
114 
115  forAll(cellBbs, celli)
116  {
117  const treeBoundBox& cellBb = cellBbs[celli];
118 
119  if (cellBb.overlaps(otherExtendedProcBb))
120  {
121  // This cell is in range of the Bb of the other
122  // processor Bb, and so needs to be referred to it
123 
124  cellInRangeOfCoupledPatch[celli] = true;
125 
126  cellIAndTToExchange.append
127  (
128  globalTransforms.encode(celli, transformIndex)
129  );
130 
131  cellBbsToExchange.append(cellBb);
132 
133  procToDistributeCellTo.append(origProc);
134  }
135  }
136  }
137 
138  buildMap(cellMapPtr_, procToDistributeCellTo);
139 
140  // Needed for reverseDistribute
141  label preDistributionCellMapSize = procToDistributeCellTo.size();
142 
143  cellMap().distribute(cellBbsToExchange);
144 
145  cellMap().distribute(cellIAndTToExchange);
146 
147  // Determine labelList specifying only cells that are in range of
148  // a coupled boundary to build an octree limited to these cells.
149  DynamicList<label> coupledPatchRangeCells;
150 
151  forAll(cellInRangeOfCoupledPatch, celli)
152  {
153  if (cellInRangeOfCoupledPatch[celli])
154  {
155  coupledPatchRangeCells.append(celli);
156  }
157  }
158 
159  treeBoundBox procBbRndExt
160  (
161  treeBoundBox(mesh_.points()).extend(rndGen, 1e-4)
162  );
163 
164  indexedOctree<treeDataCell> coupledPatchRangeTree
165  (
166  treeDataCell
167  (
168  true, // Cache cell bb
169  mesh_,
170  coupledPatchRangeCells, // Subset of mesh
171  polyMesh::CELL_TETS // Consistent with tracking
172  ),
173  procBbRndExt,
174  8, // maxLevel,
175  10, // leafSize,
176  100.0 // duplicity
177  );
178 
179  ril_.setSize(cellBbsToExchange.size());
180 
181  // This needs to be a boolList, not PackedBoolList if
182  // reverseDistribute is called.
183  boolList cellBbRequiredByAnyCell(cellBbsToExchange.size(), false);
184 
185  Info<< " Building referred interaction lists" << endl;
186 
187  forAll(cellBbsToExchange, bbI)
188  {
189  const labelPair& ciat = cellIAndTToExchange[bbI];
190 
191  const vectorTensorTransform& transform = globalTransforms.transform
192  (
193  globalTransforms.transformIndex(ciat)
194  );
195 
196  treeBoundBox tempTransformedBb
197  (
198  transform.invTransformPosition(cellBbsToExchange[bbI].points())
199  );
200 
201  treeBoundBox extendedBb
202  (
203  tempTransformedBb.min() - interactionVec,
204  tempTransformedBb.max() + interactionVec
205  );
206 
207  // Find all elements intersecting box.
208  labelList interactingElems
209  (
210  coupledPatchRangeTree.findBox(extendedBb)
211  );
212 
213  if (!interactingElems.empty())
214  {
215  cellBbRequiredByAnyCell[bbI] = true;
216  }
217 
218  ril_[bbI].setSize(interactingElems.size(), -1);
219 
220  forAll(interactingElems, i)
221  {
222  label elemI = interactingElems[i];
223 
224  // Here, a more detailed geometric test could be applied,
225  // i.e. a more accurate bounding volume like a OBB or
226  // convex hull or an exact geometrical test.
227 
228  label c = coupledPatchRangeTree.shapes().cellLabels()[elemI];
229 
230  ril_[bbI][i] = c;
231  }
232  }
233 
234  // Perform subset of ril_, to remove any referred cells that do
235  // not interact. They will not be sent from originating
236  // processors. This assumes that the disappearance of the cell
237  // from the sending list of the source processor, simply removes
238  // the referred cell from the ril_, all of the subsequent indices
239  // shuffle down one, but the order and structure is preserved,
240  // i.e. it, is as if the cell had never been referred in the first
241  // place.
242 
243  inplaceSubset(cellBbRequiredByAnyCell, ril_);
244 
245  // Send information about which cells are actually required back
246  // to originating processors.
247 
248  // At this point, cellBbsToExchange does not need to be maintained
249  // or distributed as it is not longer needed.
250 
251  cellBbsToExchange.setSize(0);
252 
253  cellMap().reverseDistribute
254  (
255  preDistributionCellMapSize,
256  cellBbRequiredByAnyCell
257  );
258 
259  cellMap().reverseDistribute
260  (
261  preDistributionCellMapSize,
262  cellIAndTToExchange
263  );
264 
265  // Perform ordering of cells to send, this invalidates the
266  // previous value of preDistributionCellMapSize.
267 
268  preDistributionCellMapSize = -1;
269 
270  // Move all of the used cells to refer to the start of the list
271  // and truncate it
272 
273  inplaceSubset(cellBbRequiredByAnyCell, cellIAndTToExchange);
274 
275  inplaceSubset(cellBbRequiredByAnyCell, procToDistributeCellTo);
276 
277  preDistributionCellMapSize = procToDistributeCellTo.size();
278 
279  // Rebuild mapDistribute with only required referred cells
280  buildMap(cellMapPtr_, procToDistributeCellTo);
281 
282  // Store cellIndexAndTransformToDistribute
283  cellIndexAndTransformToDistribute_.transfer(cellIAndTToExchange);
284 
285  // Determine inverse addressing for referred cells
286 
287  rilInverse_.setSize(mesh_.nCells());
288 
289  // Temporary Dynamic lists for accumulation
290  List<DynamicList<label>> rilInverseTemp(rilInverse_.size());
291 
292  // Loop over all referred cells
293  forAll(ril_, refCelli)
294  {
295  const labelList& realCells = ril_[refCelli];
296 
297  // Loop over all real cells in that the referred cell is to
298  // supply interactions to and record the index of this
299  // referred cell in the real cells entry in rilInverse
300 
301  forAll(realCells, realCelli)
302  {
303  rilInverseTemp[realCells[realCelli]].append(refCelli);
304  }
305  }
306 
307  forAll(rilInverse_, celli)
308  {
309  rilInverse_[celli].transfer(rilInverseTemp[celli]);
310  }
311 
312  // Determine which wall faces to refer
313 
314  // The referring of wall patch data relies on patches having the same
315  // index on each processor.
316  mesh_.boundaryMesh().checkParallelSync(true);
317 
318  // Determine the index of all of the wall faces on this processor
319  DynamicList<label> localWallFaces;
320 
321  forAll(mesh_.boundaryMesh(), patchi)
322  {
323  const polyPatch& patch = mesh_.boundaryMesh()[patchi];
324 
325  if (isA<wallPolyPatch>(patch))
326  {
327  localWallFaces.append(identity(patch.size()) + patch.start());
328  }
329  }
330 
331  treeBoundBoxList wallFaceBbs(localWallFaces.size());
332 
333  forAll(wallFaceBbs, i)
334  {
335  wallFaceBbs[i] = treeBoundBox
336  (
337  mesh_.faces()[localWallFaces[i]].points(mesh_.points())
338  );
339  }
340 
341  // IAndT: index and transform
342  DynamicList<labelPair> wallFaceIAndTToExchange;
343 
344  DynamicList<treeBoundBox> wallFaceBbsToExchange;
345 
346  DynamicList<label> procToDistributeWallFaceTo;
347 
348  forAll(extendedProcBbsInRange, ePBIRI)
349  {
350  const treeBoundBox& otherExtendedProcBb =
351  extendedProcBbsInRange[ePBIRI];
352 
353  label transformIndex = extendedProcBbsTransformIndex[ePBIRI];
354 
355  label origProc = extendedProcBbsOrigProc[ePBIRI];
356 
357  forAll(wallFaceBbs, i)
358  {
359  const treeBoundBox& wallFaceBb = wallFaceBbs[i];
360 
361  if (wallFaceBb.overlaps(otherExtendedProcBb))
362  {
363  // This wall face is in range of the Bb of the other
364  // processor Bb, and so needs to be referred to it
365 
366  label wallFacei = localWallFaces[i];
367 
368  wallFaceIAndTToExchange.append
369  (
370  globalTransforms.encode(wallFacei, transformIndex)
371  );
372 
373  wallFaceBbsToExchange.append(wallFaceBb);
374 
375  procToDistributeWallFaceTo.append(origProc);
376  }
377  }
378  }
379 
380  buildMap(wallFaceMapPtr_, procToDistributeWallFaceTo);
381 
382  // Needed for reverseDistribute
383  label preDistributionWallFaceMapSize = procToDistributeWallFaceTo.size();
384 
385  wallFaceMap().distribute(wallFaceBbsToExchange);
386 
387  wallFaceMap().distribute(wallFaceIAndTToExchange);
388 
389  indexedOctree<treeDataCell> allCellsTree
390  (
391  treeDataCell(true, mesh_, polyMesh::CELL_TETS),
392  procBbRndExt,
393  8, // maxLevel,
394  10, // leafSize,
395  100.0 // duplicity
396  );
397 
398  rwfil_.setSize(wallFaceBbsToExchange.size());
399 
400  // This needs to be a boolList, not PackedBoolList if
401  // reverseDistribute is called.
402  boolList wallFaceBbRequiredByAnyCell(wallFaceBbsToExchange.size(), false);
403 
404  forAll(wallFaceBbsToExchange, bbI)
405  {
406  const labelPair& wfiat = wallFaceIAndTToExchange[bbI];
407 
408  const vectorTensorTransform& transform = globalTransforms.transform
409  (
410  globalTransforms.transformIndex(wfiat)
411  );
412 
413  treeBoundBox tempTransformedBb
414  (
415  transform.invTransformPosition(wallFaceBbsToExchange[bbI].points())
416  );
417 
418  treeBoundBox extendedBb
419  (
420  tempTransformedBb.min() - interactionVec,
421  tempTransformedBb.max() + interactionVec
422  );
423 
424  // Find all elements intersecting box.
425  labelList interactingElems
426  (
427  coupledPatchRangeTree.findBox(extendedBb)
428  );
429 
430  if (!interactingElems.empty())
431  {
432  wallFaceBbRequiredByAnyCell[bbI] = true;
433  }
434 
435  rwfil_[bbI].setSize(interactingElems.size(), -1);
436 
437  forAll(interactingElems, i)
438  {
439  label elemI = interactingElems[i];
440 
441  // Here, a more detailed geometric test could be applied,
442  // i.e. a more accurate bounding volume like a OBB or
443  // convex hull or an exact geometrical test.
444 
445  label c = coupledPatchRangeTree.shapes().cellLabels()[elemI];
446 
447  rwfil_[bbI][i] = c;
448  }
449  }
450 
451  // Perform subset of rwfil_, to remove any referred wallFaces that do
452  // not interact. They will not be sent from originating
453  // processors. This assumes that the disappearance of the wallFace
454  // from the sending list of the source processor, simply removes
455  // the referred wallFace from the rwfil_, all of the subsequent indices
456  // shuffle down one, but the order and structure is preserved,
457  // i.e. it, is as if the wallFace had never been referred in the first
458  // place.
459 
460  inplaceSubset(wallFaceBbRequiredByAnyCell, rwfil_);
461 
462  // Send information about which wallFaces are actually required
463  // back to originating processors.
464 
465  // At this point, wallFaceBbsToExchange does not need to be
466  // maintained or distributed as it is not longer needed.
467 
468  wallFaceBbsToExchange.setSize(0);
469 
470  wallFaceMap().reverseDistribute
471  (
472  preDistributionWallFaceMapSize,
473  wallFaceBbRequiredByAnyCell
474  );
475 
476  wallFaceMap().reverseDistribute
477  (
478  preDistributionWallFaceMapSize,
479  wallFaceIAndTToExchange
480  );
481 
482  // Perform ordering of wallFaces to send, this invalidates the
483  // previous value of preDistributionWallFaceMapSize.
484 
485  preDistributionWallFaceMapSize = -1;
486 
487  // Move all of the used wallFaces to refer to the start of the
488  // list and truncate it
489 
490  inplaceSubset(wallFaceBbRequiredByAnyCell, wallFaceIAndTToExchange);
491 
492  inplaceSubset(wallFaceBbRequiredByAnyCell, procToDistributeWallFaceTo);
493 
494  preDistributionWallFaceMapSize = procToDistributeWallFaceTo.size();
495 
496  // Rebuild mapDistribute with only required referred wallFaces
497  buildMap(wallFaceMapPtr_, procToDistributeWallFaceTo);
498 
499  // Store wallFaceIndexAndTransformToDistribute
500  wallFaceIndexAndTransformToDistribute_.transfer(wallFaceIAndTToExchange);
501 
502  // Determine inverse addressing for referred wallFaces
503 
504  rwfilInverse_.setSize(mesh_.nCells());
505 
506  // Temporary Dynamic lists for accumulation
507  List<DynamicList<label>> rwfilInverseTemp(rwfilInverse_.size());
508 
509  // Loop over all referred wall faces
510  forAll(rwfil_, refWallFacei)
511  {
512  const labelList& realCells = rwfil_[refWallFacei];
513 
514  // Loop over all real cells in that the referred wall face is
515  // to supply interactions to and record the index of this
516  // referred wall face in the real cells entry in rwfilInverse
517 
518  forAll(realCells, realCelli)
519  {
520  rwfilInverseTemp[realCells[realCelli]].append(refWallFacei);
521  }
522  }
523 
524  forAll(rwfilInverse_, celli)
525  {
526  rwfilInverse_[celli].transfer(rwfilInverseTemp[celli]);
527  }
528 
529  // Refer wall faces to the appropriate processor
530  referredWallFaces_.setSize(wallFaceIndexAndTransformToDistribute_.size());
531 
532  forAll(referredWallFaces_, rWFI)
533  {
534  const labelPair& wfiat = wallFaceIndexAndTransformToDistribute_[rWFI];
535 
536  label wallFaceIndex = globalTransforms.index(wfiat);
537 
538  const vectorTensorTransform& transform = globalTransforms.transform
539  (
540  globalTransforms.transformIndex(wfiat)
541  );
542 
543  const face& f = mesh_.faces()[wallFaceIndex];
544 
545  label patchi = mesh_.boundaryMesh().patchID()
546  [
547  wallFaceIndex - mesh_.nInternalFaces()
548  ];
549 
550  referredWallFaces_[rWFI] = referredWallFace
551  (
552  face(identity(f.size())),
553  transform.invTransformPosition(f.points(mesh_.points())),
554  patchi
555  );
556  }
557 
558  wallFaceMap().distribute(referredWallFaces_);
559 
560  if (writeCloud_)
561  {
562  writeReferredWallFaces();
563  }
564 
565  // Direct interaction list and direct wall faces
566 
567  Info<< " Building direct interaction lists" << endl;
568 
569  indexedOctree<treeDataFace> wallFacesTree
570  (
571  treeDataFace(true, mesh_, localWallFaces),
572  procBbRndExt,
573  8, // maxLevel,
574  10, // leafSize,
575  100.0
576  );
577 
578  dil_.setSize(mesh_.nCells());
579 
580  dwfil_.setSize(mesh_.nCells());
581 
582  forAll(cellBbs, celli)
583  {
584  const treeBoundBox& cellBb = cellBbs[celli];
585 
586  treeBoundBox extendedBb
587  (
588  cellBb.min() - interactionVec,
589  cellBb.max() + interactionVec
590  );
591 
592  // Find all cells intersecting extendedBb
593  labelList interactingElems(allCellsTree.findBox(extendedBb));
594 
595  // Reserve space to avoid multiple resizing
596  DynamicList<label> cellDIL(interactingElems.size());
597 
598  forAll(interactingElems, i)
599  {
600  label elemI = interactingElems[i];
601 
602  label c = allCellsTree.shapes().cellLabels()[elemI];
603 
604  // Here, a more detailed geometric test could be applied,
605  // i.e. a more accurate bounding volume like a OBB or
606  // convex hull, or an exact geometrical test.
607 
608  // The higher index cell is added to the lower index
609  // cell's DIL. A cell is not added to its own DIL.
610  if (c > celli)
611  {
612  cellDIL.append(c);
613  }
614  }
615 
616  dil_[celli].transfer(cellDIL);
617 
618  // Find all wall faces intersecting extendedBb
619  interactingElems = wallFacesTree.findBox(extendedBb);
620 
621  dwfil_[celli].setSize(interactingElems.size(), -1);
622 
623  forAll(interactingElems, i)
624  {
625  label elemI = interactingElems[i];
626 
627  label f = wallFacesTree.shapes().faceLabels()[elemI];
628 
629  dwfil_[celli][i] = f;
630  }
631  }
632 }
633 
634 
635 template<class ParticleType>
637 (
638  const treeBoundBox& procBb,
639  const List<treeBoundBox>& allExtendedProcBbs,
640  const globalIndexAndTransform& globalTransforms,
641  List<treeBoundBox>& extendedProcBbsInRange,
642  List<label>& extendedProcBbsTransformIndex,
643  List<label>& extendedProcBbsOrigProc
644 )
645 {
646  extendedProcBbsInRange.setSize(0);
647  extendedProcBbsTransformIndex.setSize(0);
648  extendedProcBbsOrigProc.setSize(0);
649 
650  DynamicList<treeBoundBox> tmpExtendedProcBbsInRange;
651  DynamicList<label> tmpExtendedProcBbsTransformIndex;
652  DynamicList<label> tmpExtendedProcBbsOrigProc;
653 
654  label nTrans = globalTransforms.nIndependentTransforms();
655 
656  forAll(allExtendedProcBbs, proci)
657  {
658  List<label> permutationIndices(nTrans, 0);
659 
660  if (nTrans == 0 && proci != Pstream::myProcNo())
661  {
662  treeBoundBox extendedReferredProcBb = allExtendedProcBbs[proci];
663 
664  if (procBb.overlaps(extendedReferredProcBb))
665  {
666  tmpExtendedProcBbsInRange.append(extendedReferredProcBb);
667 
668  // Dummy index, there are no transforms, so there will
669  // be no resultant transform when this is decoded.
670  tmpExtendedProcBbsTransformIndex.append(0);
671 
672  tmpExtendedProcBbsOrigProc.append(proci);
673  }
674  }
675  else if (nTrans == 3)
676  {
677  label& i = permutationIndices[0];
678  label& j = permutationIndices[1];
679  label& k = permutationIndices[2];
680 
681  for (i = -1; i <= 1; i++)
682  {
683  for (j = -1; j <= 1; j++)
684  {
685  for (k = -1; k <= 1; k++)
686  {
687  if
688  (
689  i == 0
690  && j == 0
691  && k == 0
692  && proci == Pstream::myProcNo()
693  )
694  {
695  // Skip this processor's extended boundBox
696  // when it has no transformation
697  continue;
698  }
699 
700  label transI = globalTransforms.encodeTransformIndex
701  (
702  permutationIndices
703  );
704 
705  const vectorTensorTransform& transform =
706  globalTransforms.transform(transI);
707 
708  treeBoundBox extendedReferredProcBb
709  (
710  transform.transformPosition
711  (
712  allExtendedProcBbs[proci].points()
713  )
714  );
715 
716  if (procBb.overlaps(extendedReferredProcBb))
717  {
718  tmpExtendedProcBbsInRange.append
719  (
720  extendedReferredProcBb
721  );
722 
723  tmpExtendedProcBbsTransformIndex.append(transI);
724 
725  tmpExtendedProcBbsOrigProc.append(proci);
726  }
727  }
728  }
729  }
730  }
731  else if (nTrans == 2)
732  {
733  label& i = permutationIndices[0];
734  label& j = permutationIndices[1];
735 
736  for (i = -1; i <= 1; i++)
737  {
738  for (j = -1; j <= 1; j++)
739  {
740  if (i == 0 && j == 0 && proci == Pstream::myProcNo())
741  {
742  // Skip this processor's extended boundBox
743  // when it has no transformation
744  continue;
745  }
746 
747  label transI = globalTransforms.encodeTransformIndex
748  (
749  permutationIndices
750  );
751 
752  const vectorTensorTransform& transform =
753  globalTransforms.transform(transI);
754 
755  treeBoundBox extendedReferredProcBb
756  (
757  transform.transformPosition
758  (
759  allExtendedProcBbs[proci].points()
760  )
761  );
762 
763  if (procBb.overlaps(extendedReferredProcBb))
764  {
765  tmpExtendedProcBbsInRange.append
766  (
767  extendedReferredProcBb
768  );
769 
770  tmpExtendedProcBbsTransformIndex.append(transI);
771 
772  tmpExtendedProcBbsOrigProc.append(proci);
773  }
774  }
775  }
776  }
777  else if (nTrans == 1)
778  {
779  label& i = permutationIndices[0];
780 
781  for (i = -1; i <= 1; i++)
782  {
783  if (i == 0 && proci == Pstream::myProcNo())
784  {
785  // Skip this processor's extended boundBox when it
786  // has no transformation
787  continue;
788  }
789 
790  label transI = globalTransforms.encodeTransformIndex
791  (
792  permutationIndices
793  );
794 
795  const vectorTensorTransform& transform =
796  globalTransforms.transform(transI);
797 
798  treeBoundBox extendedReferredProcBb
799  (
800  transform.transformPosition
801  (
802  allExtendedProcBbs[proci].points()
803  )
804  );
805 
806  if (procBb.overlaps(extendedReferredProcBb))
807  {
808  tmpExtendedProcBbsInRange.append
809  (
810  extendedReferredProcBb
811  );
812 
813  tmpExtendedProcBbsTransformIndex.append(transI);
814 
815  tmpExtendedProcBbsOrigProc.append(proci);
816  }
817  }
818  }
819  }
820 
821  extendedProcBbsInRange = tmpExtendedProcBbsInRange.xfer();
822  extendedProcBbsTransformIndex = tmpExtendedProcBbsTransformIndex.xfer();
823  extendedProcBbsOrigProc = tmpExtendedProcBbsOrigProc.xfer();
824 }
825 
826 
827 template<class ParticleType>
829 (
830  autoPtr<mapDistribute>& mapPtr,
831  const List<label>& toProc
832 )
833 {
834  // Determine send map
835  // ~~~~~~~~~~~~~~~~~~
836 
837  // 1. Count
838  labelList nSend(Pstream::nProcs(), 0);
839 
840  forAll(toProc, i)
841  {
842  label proci = toProc[i];
843 
844  nSend[proci]++;
845  }
846 
847  // 2. Size sendMap
848  labelListList sendMap(Pstream::nProcs());
849 
850  forAll(nSend, proci)
851  {
852  sendMap[proci].setSize(nSend[proci]);
853 
854  nSend[proci] = 0;
855  }
856 
857  // 3. Fill sendMap
858  forAll(toProc, i)
859  {
860  label proci = toProc[i];
861 
862  sendMap[proci][nSend[proci]++] = i;
863  }
864 
865  // 4. Send over how many I need to receive
866  labelList recvSizes;
867  Pstream::exchangeSizes(sendMap, recvSizes);
868 
869 
870  // Determine receive map
871  // ~~~~~~~~~~~~~~~~~~~~~
872 
873  labelListList constructMap(Pstream::nProcs());
874 
875  // Local transfers first
876  constructMap[Pstream::myProcNo()] = identity
877  (
878  sendMap[Pstream::myProcNo()].size()
879  );
880 
881  label constructSize = constructMap[Pstream::myProcNo()].size();
882 
883  forAll(constructMap, proci)
884  {
885  if (proci != Pstream::myProcNo())
886  {
887  label nRecv = recvSizes[proci];
888 
889  constructMap[proci].setSize(nRecv);
890 
891  for (label i = 0; i < nRecv; i++)
892  {
893  constructMap[proci][i] = constructSize++;
894  }
895  }
896  }
897 
898  mapPtr.reset
899  (
900  new mapDistribute
901  (
902  constructSize,
903  sendMap.xfer(),
904  constructMap.xfer()
905  )
906  );
907 }
908 
909 
910 template<class ParticleType>
912 (
913  const List<DynamicList<ParticleType*>>& cellOccupancy
914 )
915 {
916  const globalIndexAndTransform& globalTransforms =
917  mesh_.globalData().globalTransforms();
918 
919  referredParticles_.setSize(cellIndexAndTransformToDistribute_.size());
920 
921  // Clear all existing referred particles
922 
923  forAll(referredParticles_, i)
924  {
925  referredParticles_[i].clear();
926  }
927 
928  // Clear all particles that may have been populated into the cloud
929  cloud_.clear();
930 
931  forAll(cellIndexAndTransformToDistribute_, i)
932  {
933  const labelPair ciat = cellIndexAndTransformToDistribute_[i];
934 
935  label cellIndex = globalTransforms.index(ciat);
936 
937  List<ParticleType*> realParticles = cellOccupancy[cellIndex];
938 
939  IDLList<ParticleType>& particlesToRefer = referredParticles_[i];
940 
941  forAll(realParticles, rM)
942  {
943  const ParticleType& particle = *realParticles[rM];
944 
945  particlesToRefer.append(particle.clone().ptr());
946 
947  prepareParticleToBeReferred(particlesToRefer.last(), ciat);
948  }
949  }
950 }
951 
952 
953 template<class ParticleType>
955 (
956  ParticleType* particle,
957  labelPair ciat
958 )
959 {
960  const globalIndexAndTransform& globalTransforms =
961  mesh_.globalData().globalTransforms();
962 
963  const vectorTensorTransform& transform = globalTransforms.transform
964  (
965  globalTransforms.transformIndex(ciat)
966  );
967 
968  particle->position() = transform.invTransformPosition(particle->position());
969 
970  particle->transformProperties(-transform.t());
971 
972  if (transform.hasR())
973  {
974  particle->transformProperties(transform.R().T());
975  }
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  forAllConstIter(typename IDLList<ParticleType>, refCell, iter)
990  {
991  cloud_.addParticle
992  (
993  static_cast<ParticleType*>(iter().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_, "NULL_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  Switch 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 form scratch."
1166  << endl;
1167 
1168  buildInteractionLists();
1169  }
1170 
1171  prepareWallDataToRefer();
1172 
1173  prepareParticlesToRefer(cellOccupancy);
1174 
1175  for (label domain = 0; domain < Pstream::nProcs(); domain++)
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 (label domain = 0; domain < Pstream::nProcs(); domain++)
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  fillReferredParticleCloud();
1234 
1235  wallFaceMap().receive(pBufs, referredWallData_);
1236 }
1237 
1238 
1239 // ************************************************************************* //
void inplaceSubset(const UList< T > &select, const T &value, ListType &)
Inplace extract elements of List when select is a certain value.
cachedRandom rndGen(label(0),-1)
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 ...
const labelListList & constructMap() const
From subsetted data to new reconstructed data.
void sendReferredData(const List< DynamicList< ParticleType * >> &cellOccupancy, PstreamBuffers &pBufs)
Prepare and send referred particles and wall data,.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
const labelListList & subMap() const
From subsetted data back to original data.
const mapDistribute & cellMap() const
Return access to the cellMap.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
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: doubleFloat.H:78
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
bool changing() const
Is mesh changing (topology changing and/or moving)
Definition: polyMesh.H:521
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:262
labelList f(nPoints)
Buffers for inter-processor communications streams (UOPstream, UIPstream).
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
Definition: UPstream.C:117
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:295
Intrusive doubly-linked list.
Definition: IDLList.H:47
void setSize(const label)
Reset size of List.
Definition: List.C:295
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:399
label patchi
const List< DynamicList< molecule * > > & cellOccupancy
#define WarningInFunction
Report a warning using Foam::Warning.
void send(PstreamBuffers &, const List< T > &) const
Do all sends using PstreamBuffers.
A List with indirect addressing.
Definition: fvMatrix.H:106
const dimensionedScalar c
Speed of light in a vacuum.
messageStream Info
const mapDistribute & wallFaceMap() const
Return access to the wallFaceMap.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
void receive(PstreamBuffers &, List< T > &) const
Do all receives using PstreamBuffers.