BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / lagrangian / basic / InteractionLists / InteractionLists.C
blob6a7dfb1f96398fad1ab416ac7ec52705899b9489
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "InteractionLists.H"
27 #include "globalIndexAndTransform.H"
28 #include "indexedOctree.H"
29 #include "treeDataFace.H"
30 #include "treeDataCell.H"
31 #include "volFields.H"
32 #include "meshTools.H"
34 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
36 template<class ParticleType>
37 void Foam::InteractionLists<ParticleType>::buildInteractionLists()
39     Info<< "Building InteractionLists with interaction distance "
40         << maxDistance_ << endl;
42     Random rndGen(419715);
44     const vector interactionVec = maxDistance_*vector::one;
46     treeBoundBox procBb(treeBoundBox(mesh_.points()));
48     treeBoundBox extendedProcBb
49     (
50         procBb.min() - interactionVec,
51         procBb.max() + interactionVec
52     );
54     treeBoundBoxList allExtendedProcBbs(Pstream::nProcs());
56     allExtendedProcBbs[Pstream::myProcNo()] = extendedProcBb;
58     Pstream::gatherList(allExtendedProcBbs);
60     Pstream::scatterList(allExtendedProcBbs);
62     List<treeBoundBox> extendedProcBbsInRange;
63     List<label> extendedProcBbsTransformIndex;
64     List<label> extendedProcBbsOrigProc;
66     findExtendedProcBbsInRange
67     (
68         procBb,
69         allExtendedProcBbs,
70         mesh_.globalData().globalTransforms(),
71         extendedProcBbsInRange,
72         extendedProcBbsTransformIndex,
73         extendedProcBbsOrigProc
74     );
76     treeBoundBoxList cellBbs(mesh_.nCells());
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     }
90     const globalIndexAndTransform& globalTransforms =
91         mesh_.globalData().globalTransforms();
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);
98     // IAndT: index (=local cell index) and transform (from
99     // globalIndexAndTransform)
100     DynamicList<labelPair> cellIAndTToExchange;
102     DynamicList<treeBoundBox> cellBbsToExchange;
104     DynamicList<label> procToDistributeCellTo;
106     forAll(extendedProcBbsInRange, ePBIRI)
107     {
108         const treeBoundBox& otherExtendedProcBb =
109             extendedProcBbsInRange[ePBIRI];
111         label transformIndex = extendedProcBbsTransformIndex[ePBIRI];
113         label origProc = extendedProcBbsOrigProc[ePBIRI];
115         forAll(cellBbs, cellI)
116         {
117             const treeBoundBox& cellBb = cellBbs[cellI];
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
124                 cellInRangeOfCoupledPatch[cellI] = true;
126                 cellIAndTToExchange.append
127                 (
128                     globalTransforms.encode(cellI, transformIndex)
129                 );
131                 cellBbsToExchange.append(cellBb);
133                 procToDistributeCellTo.append(origProc);
134             }
135         }
136     }
138     buildMap(cellMapPtr_, procToDistributeCellTo);
140     // Needed for reverseDistribute
141     label preDistributionCellMapSize = procToDistributeCellTo.size();
143     cellMap().distribute(cellBbsToExchange);
145     cellMap().distribute(cellIAndTToExchange);
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;
151     forAll(cellInRangeOfCoupledPatch, cellI)
152     {
153         if (cellInRangeOfCoupledPatch[cellI])
154         {
155             coupledPatchRangeCells.append(cellI);
156         }
157     }
159     treeBoundBox procBbRndExt
160     (
161         treeBoundBox(mesh_.points()).extend(rndGen, 1e-4)
162     );
164     indexedOctree<treeDataCell> coupledPatchRangeTree
165     (
166         treeDataCell(true, mesh_, coupledPatchRangeCells),
167         procBbRndExt,
168         8,              // maxLevel,
169         10,             // leafSize,
170         100.0           // duplicity
171     );
173     ril_.setSize(cellBbsToExchange.size());
175     // This needs to be a boolList, not PackedBoolList if
176     // reverseDistribute is called.
177     boolList cellBbRequiredByAnyCell(cellBbsToExchange.size(), false);
179     Info<< "    Building referred interaction lists" << endl;
181     forAll(cellBbsToExchange, bbI)
182     {
183         const labelPair& ciat = cellIAndTToExchange[bbI];
185         const vectorTensorTransform& transform = globalTransforms.transform
186         (
187             globalTransforms.transformIndex(ciat)
188         );
190         treeBoundBox tempTransformedBb
191         (
192             transform.invTransformPosition(cellBbsToExchange[bbI].points())
193         );
195         treeBoundBox extendedBb
196         (
197             tempTransformedBb.min() - interactionVec,
198             tempTransformedBb.max() + interactionVec
199         );
201         // Find all elements intersecting box.
202         labelList interactingElems
203         (
204             coupledPatchRangeTree.findBox(extendedBb)
205         );
207         if (!interactingElems.empty())
208         {
209             cellBbRequiredByAnyCell[bbI] = true;
210         }
212         ril_[bbI].setSize(interactingElems.size(), -1);
214         forAll(interactingElems, i)
215         {
216             label elemI = interactingElems[i];
218             // Here, a more detailed geometric test could be applied,
219             // i.e. a more accurate bounding volume like a OBB or
220             // convex hull or an exact geometrical test.
222             label c = coupledPatchRangeTree.shapes().cellLabels()[elemI];
224             ril_[bbI][i] = c;
225         }
226     }
228     // Perform subset of ril_, to remove any referred cells that do
229     // not interact.  They will not be sent from originating
230     // processors.  This assumes that the disappearance of the cell
231     // from the sending list of the source processor, simply removes
232     // the referred cell from the ril_, all of the subsequent indices
233     // shuffle down one, but the order and structure is preserved,
234     // i.e. it, is as if the cell had never been referred in the first
235     // place.
237     inplaceSubset(cellBbRequiredByAnyCell, ril_);
239     // Send information about which cells are actually required back
240     // to originating processors.
242     // At this point, cellBbsToExchange does not need to be maintained
243     // or distributed as it is not longer needed.
245     cellBbsToExchange.setSize(0);
247     cellMap().reverseDistribute
248     (
249         preDistributionCellMapSize,
250         cellBbRequiredByAnyCell
251     );
253     cellMap().reverseDistribute
254     (
255         preDistributionCellMapSize,
256         cellIAndTToExchange
257     );
259     // Perform ordering of cells to send, this invalidates the
260     // previous value of preDistributionCellMapSize.
262     preDistributionCellMapSize = -1;
264     // Move all of the used cells to refer to the start of the list
265     // and truncate it
267     inplaceSubset(cellBbRequiredByAnyCell, cellIAndTToExchange);
269     inplaceSubset(cellBbRequiredByAnyCell, procToDistributeCellTo);
271     preDistributionCellMapSize = procToDistributeCellTo.size();
273     // Rebuild mapDistribute with only required referred cells
274     buildMap(cellMapPtr_, procToDistributeCellTo);
276     // Store cellIndexAndTransformToDistribute
277     cellIndexAndTransformToDistribute_.transfer(cellIAndTToExchange);
279     // Determine inverse addressing for referred cells
281     rilInverse_.setSize(mesh_.nCells());
283     // Temporary Dynamic lists for accumulation
284     List<DynamicList<label> > rilInverseTemp(rilInverse_.size());
286     // Loop over all referred cells
287     forAll(ril_, refCellI)
288     {
289         const labelList& realCells = ril_[refCellI];
291         // Loop over all real cells in that the referred cell is to
292         // supply interactions to and record the index of this
293         // referred cell in the real cells entry in rilInverse
295         forAll(realCells, realCellI)
296         {
297             rilInverseTemp[realCells[realCellI]].append(refCellI);
298         }
299     }
301     forAll(rilInverse_, cellI)
302     {
303         rilInverse_[cellI].transfer(rilInverseTemp[cellI]);
304     }
306     // Determine which wall faces to refer
308     // The referring of wall patch data relies on patches having the same
309     // index on each processor.
310     mesh_.boundaryMesh().checkParallelSync(true);
312     // Determine the index of all of the wall faces on this processor
313     DynamicList<label> localWallFaces;
315     forAll(mesh_.boundaryMesh(), patchI)
316     {
317         const polyPatch& patch = mesh_.boundaryMesh()[patchI];
319         if (isA<wallPolyPatch>(patch))
320         {
321             localWallFaces.append(identity(patch.size()) + patch.start());
322         }
323     }
325     treeBoundBoxList wallFaceBbs(localWallFaces.size());
327     forAll(wallFaceBbs, i)
328     {
329         wallFaceBbs[i] = treeBoundBox
330         (
331             mesh_.faces()[localWallFaces[i]].points(mesh_.points())
332         );
333     }
335     // IAndT: index and transform
336     DynamicList<labelPair> wallFaceIAndTToExchange;
338     DynamicList<treeBoundBox> wallFaceBbsToExchange;
340     DynamicList<label> procToDistributeWallFaceTo;
342     forAll(extendedProcBbsInRange, ePBIRI)
343     {
344         const treeBoundBox& otherExtendedProcBb =
345             extendedProcBbsInRange[ePBIRI];
347         label transformIndex = extendedProcBbsTransformIndex[ePBIRI];
349         label origProc = extendedProcBbsOrigProc[ePBIRI];
351         forAll(wallFaceBbs, i)
352         {
353             const treeBoundBox& wallFaceBb = wallFaceBbs[i];
355             if (wallFaceBb.overlaps(otherExtendedProcBb))
356             {
357                 // This wall face is in range of the Bb of the other
358                 // processor Bb, and so needs to be referred to it
360                 label wallFaceI = localWallFaces[i];
362                 wallFaceIAndTToExchange.append
363                 (
364                     globalTransforms.encode(wallFaceI, transformIndex)
365                 );
367                 wallFaceBbsToExchange.append(wallFaceBb);
369                 procToDistributeWallFaceTo.append(origProc);
370             }
371         }
372     }
374     buildMap(wallFaceMapPtr_, procToDistributeWallFaceTo);
376     // Needed for reverseDistribute
377     label preDistributionWallFaceMapSize = procToDistributeWallFaceTo.size();
379     wallFaceMap().distribute(wallFaceBbsToExchange);
381     wallFaceMap().distribute(wallFaceIAndTToExchange);
383     indexedOctree<treeDataCell> allCellsTree
384     (
385         treeDataCell(true, mesh_),
386         procBbRndExt,
387         8,              // maxLevel,
388         10,             // leafSize,
389         100.0           // duplicity
390     );
392     rwfil_.setSize(wallFaceBbsToExchange.size());
394     // This needs to be a boolList, not PackedBoolList if
395     // reverseDistribute is called.
396     boolList wallFaceBbRequiredByAnyCell(wallFaceBbsToExchange.size(), false);
398     forAll(wallFaceBbsToExchange, bbI)
399     {
400         const labelPair& wfiat = wallFaceIAndTToExchange[bbI];
402         const vectorTensorTransform& transform = globalTransforms.transform
403         (
404             globalTransforms.transformIndex(wfiat)
405         );
407         treeBoundBox tempTransformedBb
408         (
409             transform.invTransformPosition(wallFaceBbsToExchange[bbI].points())
410         );
412         treeBoundBox extendedBb
413         (
414             tempTransformedBb.min() - interactionVec,
415             tempTransformedBb.max() + interactionVec
416         );
418         // Find all elements intersecting box.
419         labelList interactingElems
420         (
421             coupledPatchRangeTree.findBox(extendedBb)
422         );
424         if (!interactingElems.empty())
425         {
426             wallFaceBbRequiredByAnyCell[bbI] = true;
427         }
429         rwfil_[bbI].setSize(interactingElems.size(), -1);
431         forAll(interactingElems, i)
432         {
433             label elemI = interactingElems[i];
435             // Here, a more detailed geometric test could be applied,
436             // i.e. a more accurate bounding volume like a OBB or
437             // convex hull or an exact geometrical test.
439             label c = coupledPatchRangeTree.shapes().cellLabels()[elemI];
441             rwfil_[bbI][i] = c;
442         }
443     }
445     // Perform subset of rwfil_, to remove any referred wallFaces that do
446     // not interact.  They will not be sent from originating
447     // processors.  This assumes that the disappearance of the wallFace
448     // from the sending list of the source processor, simply removes
449     // the referred wallFace from the rwfil_, all of the subsequent indices
450     // shuffle down one, but the order and structure is preserved,
451     // i.e. it, is as if the wallFace had never been referred in the first
452     // place.
454     inplaceSubset(wallFaceBbRequiredByAnyCell, rwfil_);
456     // Send information about which wallFaces are actually required
457     // back to originating processors.
459     // At this point, wallFaceBbsToExchange does not need to be
460     // maintained or distributed as it is not longer needed.
462     wallFaceBbsToExchange.setSize(0);
464     wallFaceMap().reverseDistribute
465     (
466         preDistributionWallFaceMapSize,
467         wallFaceBbRequiredByAnyCell
468     );
470     wallFaceMap().reverseDistribute
471     (
472         preDistributionWallFaceMapSize,
473         wallFaceIAndTToExchange
474     );
476     // Perform ordering of wallFaces to send, this invalidates the
477     // previous value of preDistributionWallFaceMapSize.
479     preDistributionWallFaceMapSize = -1;
481     // Move all of the used wallFaces to refer to the start of the
482     // list and truncate it
484     inplaceSubset(wallFaceBbRequiredByAnyCell, wallFaceIAndTToExchange);
486     inplaceSubset(wallFaceBbRequiredByAnyCell, procToDistributeWallFaceTo);
488     preDistributionWallFaceMapSize = procToDistributeWallFaceTo.size();
490     // Rebuild mapDistribute with only required referred wallFaces
491     buildMap(wallFaceMapPtr_, procToDistributeWallFaceTo);
493     // Store wallFaceIndexAndTransformToDistribute
494     wallFaceIndexAndTransformToDistribute_.transfer(wallFaceIAndTToExchange);
496     // Determine inverse addressing for referred wallFaces
498     rwfilInverse_.setSize(mesh_.nCells());
500     // Temporary Dynamic lists for accumulation
501     List<DynamicList<label> > rwfilInverseTemp(rwfilInverse_.size());
503     // Loop over all referred wall faces
504     forAll(rwfil_, refWallFaceI)
505     {
506         const labelList& realCells = rwfil_[refWallFaceI];
508         // Loop over all real cells in that the referred wall face is
509         // to supply interactions to and record the index of this
510         // referred wall face in the real cells entry in rwfilInverse
512         forAll(realCells, realCellI)
513         {
514             rwfilInverseTemp[realCells[realCellI]].append(refWallFaceI);
515         }
516     }
518     forAll(rwfilInverse_, cellI)
519     {
520         rwfilInverse_[cellI].transfer(rwfilInverseTemp[cellI]);
521     }
523     // Refer wall faces to the appropriate processor
524     referredWallFaces_.setSize(wallFaceIndexAndTransformToDistribute_.size());
526     forAll(referredWallFaces_, rWFI)
527     {
528         const labelPair& wfiat = wallFaceIndexAndTransformToDistribute_[rWFI];
530         label wallFaceIndex = globalTransforms.index(wfiat);
532         const vectorTensorTransform& transform = globalTransforms.transform
533         (
534             globalTransforms.transformIndex(wfiat)
535         );
537         const face& f = mesh_.faces()[wallFaceIndex];
539         label patchI = mesh_.boundaryMesh().patchID()
540         [
541             wallFaceIndex - mesh_.nInternalFaces()
542         ];
544         referredWallFaces_[rWFI] = referredWallFace
545         (
546             face(identity(f.size())),
547             transform.invTransformPosition(f.points(mesh_.points())),
548             patchI
549         );
550     }
552     wallFaceMap().distribute(referredWallFaces_);
554     if (writeCloud_)
555     {
556         writeReferredWallFaces();
557     }
559     // Direct interaction list and direct wall faces
561     Info<< "    Building direct interaction lists" << endl;
563     indexedOctree<treeDataFace> wallFacesTree
564     (
565         treeDataFace(true, mesh_, localWallFaces),
566         procBbRndExt,
567         8,              // maxLevel,
568         10,             // leafSize,
569         100.0
570     );
572     dil_.setSize(mesh_.nCells());
574     dwfil_.setSize(mesh_.nCells());
576     forAll(cellBbs, cellI)
577     {
578         const treeBoundBox& cellBb = cellBbs[cellI];
580         treeBoundBox extendedBb
581         (
582             cellBb.min() - interactionVec,
583             cellBb.max() + interactionVec
584         );
586         // Find all cells intersecting extendedBb
587         labelList interactingElems(allCellsTree.findBox(extendedBb));
589         // Reserve space to avoid multiple resizing
590         DynamicList<label> cellDIL(interactingElems.size());
592         forAll(interactingElems, i)
593         {
594             label elemI = interactingElems[i];
596             label c = allCellsTree.shapes().cellLabels()[elemI];
598             // Here, a more detailed geometric test could be applied,
599             // i.e. a more accurate bounding volume like a OBB or
600             // convex hull, or an exact geometrical test.
602             // The higher index cell is added to the lower index
603             // cell's DIL.  A cell is not added to its own DIL.
604             if (c > cellI)
605             {
606                 cellDIL.append(c);
607             }
608         }
610         dil_[cellI].transfer(cellDIL);
612         // Find all wall faces intersecting extendedBb
613         interactingElems = wallFacesTree.findBox(extendedBb);
615         dwfil_[cellI].setSize(interactingElems.size(), -1);
617         forAll(interactingElems, i)
618         {
619             label elemI = interactingElems[i];
621             label f = wallFacesTree.shapes().faceLabels()[elemI];
623             dwfil_[cellI][i] = f;
624         }
625     }
629 template<class ParticleType>
630 void Foam::InteractionLists<ParticleType>::findExtendedProcBbsInRange
632     const treeBoundBox& procBb,
633     const List<treeBoundBox>& allExtendedProcBbs,
634     const globalIndexAndTransform& globalTransforms,
635     List<treeBoundBox>& extendedProcBbsInRange,
636     List<label>& extendedProcBbsTransformIndex,
637     List<label>& extendedProcBbsOrigProc
640     extendedProcBbsInRange.setSize(0);
641     extendedProcBbsTransformIndex.setSize(0);
642     extendedProcBbsOrigProc.setSize(0);
644     DynamicList<treeBoundBox> tmpExtendedProcBbsInRange;
645     DynamicList<label> tmpExtendedProcBbsTransformIndex;
646     DynamicList<label> tmpExtendedProcBbsOrigProc;
648     label nTrans = globalTransforms.nIndependentTransforms();
650     forAll(allExtendedProcBbs, procI)
651     {
652         List<label> permutationIndices(nTrans, 0);
654         if (nTrans == 0 && procI != Pstream::myProcNo())
655         {
656             treeBoundBox extendedReferredProcBb = allExtendedProcBbs[procI];
658             if (procBb.overlaps(extendedReferredProcBb))
659             {
660                 tmpExtendedProcBbsInRange.append(extendedReferredProcBb);
662                 // Dummy index, there are no transforms, so there will
663                 // be no resultant transform when this is decoded.
664                 tmpExtendedProcBbsTransformIndex.append(0);
666                 tmpExtendedProcBbsOrigProc.append(procI);
667             }
668         }
669         else if (nTrans == 3)
670         {
671             label& i = permutationIndices[0];
672             label& j = permutationIndices[1];
673             label& k = permutationIndices[2];
675             for (i = -1; i <= 1; i++)
676             {
677                 for (j = -1; j <= 1; j++)
678                 {
679                     for (k = -1; k <= 1; k++)
680                     {
681                         if
682                         (
683                             i == 0
684                          && j == 0
685                          && k == 0
686                          && procI == Pstream::myProcNo()
687                         )
688                         {
689                             // Skip this processor's extended boundBox
690                             // when it has no transformation
691                             continue;
692                         }
694                         label transI = globalTransforms.encodeTransformIndex
695                         (
696                             permutationIndices
697                         );
699                         const vectorTensorTransform& transform =
700                             globalTransforms.transform(transI);
702                         treeBoundBox extendedReferredProcBb
703                         (
704                             transform.transformPosition
705                             (
706                                 allExtendedProcBbs[procI].points()
707                             )
708                         );
710                         if (procBb.overlaps(extendedReferredProcBb))
711                         {
712                             tmpExtendedProcBbsInRange.append
713                             (
714                                 extendedReferredProcBb
715                             );
717                             tmpExtendedProcBbsTransformIndex.append(transI);
719                             tmpExtendedProcBbsOrigProc.append(procI);
720                         }
721                     }
722                 }
723             }
724         }
725         else if (nTrans == 2)
726         {
727             label& i = permutationIndices[0];
728             label& j = permutationIndices[1];
730             for (i = -1; i <= 1; i++)
731             {
732                 for (j = -1; j <= 1; j++)
733                 {
734                     if (i == 0 && j == 0 && procI == Pstream::myProcNo())
735                     {
736                         // Skip this processor's extended boundBox
737                         // when it has no transformation
738                         continue;
739                     }
741                     label transI = globalTransforms.encodeTransformIndex
742                     (
743                         permutationIndices
744                     );
746                     const vectorTensorTransform& transform =
747                         globalTransforms.transform(transI);
749                     treeBoundBox extendedReferredProcBb
750                     (
751                         transform.transformPosition
752                         (
753                             allExtendedProcBbs[procI].points()
754                         )
755                     );
757                     if (procBb.overlaps(extendedReferredProcBb))
758                     {
759                         tmpExtendedProcBbsInRange.append
760                         (
761                             extendedReferredProcBb
762                         );
764                         tmpExtendedProcBbsTransformIndex.append(transI);
766                         tmpExtendedProcBbsOrigProc.append(procI);
767                     }
768                 }
769             }
770         }
771         else if (nTrans == 1)
772         {
773             label& i = permutationIndices[0];
775             for (i = -1; i <= 1; i++)
776             {
777                 if (i == 0 && procI == Pstream::myProcNo())
778                 {
779                     // Skip this processor's extended boundBox when it
780                     // has no transformation
781                     continue;
782                 }
784                 label transI = globalTransforms.encodeTransformIndex
785                 (
786                     permutationIndices
787                 );
789                 const vectorTensorTransform& transform =
790                     globalTransforms.transform(transI);
792                 treeBoundBox extendedReferredProcBb
793                 (
794                     transform.transformPosition
795                     (
796                         allExtendedProcBbs[procI].points()
797                     )
798                 );
800                 if (procBb.overlaps(extendedReferredProcBb))
801                 {
802                     tmpExtendedProcBbsInRange.append
803                     (
804                         extendedReferredProcBb
805                     );
807                     tmpExtendedProcBbsTransformIndex.append(transI);
809                     tmpExtendedProcBbsOrigProc.append(procI);
810                 }
811             }
812         }
813     }
815     extendedProcBbsInRange = tmpExtendedProcBbsInRange.xfer();
816     extendedProcBbsTransformIndex = tmpExtendedProcBbsTransformIndex.xfer();
817     extendedProcBbsOrigProc = tmpExtendedProcBbsOrigProc.xfer();
821 template<class ParticleType>
822 void Foam::InteractionLists<ParticleType>::buildMap
824     autoPtr<mapDistribute>& mapPtr,
825     const List<label>& toProc
828     // Determine send map
829     // ~~~~~~~~~~~~~~~~~~
831     // 1. Count
832     labelList nSend(Pstream::nProcs(), 0);
834     forAll(toProc, i)
835     {
836         label procI = toProc[i];
838         nSend[procI]++;
839     }
841     // Send over how many I need to receive
842     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
844     labelListList sendSizes(Pstream::nProcs());
846     sendSizes[Pstream::myProcNo()] = nSend;
848     combineReduce(sendSizes, UPstream::listEq());
850     // 2. Size sendMap
851     labelListList sendMap(Pstream::nProcs());
853     forAll(nSend, procI)
854     {
855         sendMap[procI].setSize(nSend[procI]);
857         nSend[procI] = 0;
858     }
860     // 3. Fill sendMap
861     forAll(toProc, i)
862     {
863         label procI = toProc[i];
865         sendMap[procI][nSend[procI]++] = i;
866     }
868     // Determine receive map
869     // ~~~~~~~~~~~~~~~~~~~~~
871     labelListList constructMap(Pstream::nProcs());
873     // Local transfers first
874     constructMap[Pstream::myProcNo()] = identity
875     (
876         sendMap[Pstream::myProcNo()].size()
877     );
879     label constructSize = constructMap[Pstream::myProcNo()].size();
881     forAll(constructMap, procI)
882     {
883         if (procI != Pstream::myProcNo())
884         {
885             label nRecv = sendSizes[procI][Pstream::myProcNo()];
887             constructMap[procI].setSize(nRecv);
889             for (label i = 0; i < nRecv; i++)
890             {
891                 constructMap[procI][i] = constructSize++;
892             }
893         }
894     }
896     mapPtr.reset
897     (
898         new mapDistribute
899         (
900             constructSize,
901             sendMap.xfer(),
902             constructMap.xfer()
903         )
904     );
908 template<class ParticleType>
909 void Foam::InteractionLists<ParticleType>::prepareParticlesToRefer
911     const List<DynamicList<ParticleType*> >& cellOccupancy
914     const globalIndexAndTransform& globalTransforms =
915         mesh_.globalData().globalTransforms();
917     referredParticles_.setSize(cellIndexAndTransformToDistribute_.size());
919     // Clear all existing referred particles
921     forAll(referredParticles_, i)
922     {
923         referredParticles_[i].clear();
924     }
926     // Clear all particles that may have been populated into the cloud
927     cloud_.clear();
929     forAll(cellIndexAndTransformToDistribute_, i)
930     {
931         const labelPair ciat = cellIndexAndTransformToDistribute_[i];
933         label cellIndex = globalTransforms.index(ciat);
935         List<ParticleType*> realParticles = cellOccupancy[cellIndex];
937         IDLList<ParticleType>& particlesToRefer = referredParticles_[i];
939         forAll(realParticles, rM)
940         {
941             const ParticleType& particle = *realParticles[rM];
943             particlesToRefer.append(particle.clone().ptr());
945             prepareParticleToBeReferred(particlesToRefer.last(), ciat);
946         }
947     }
951 template<class ParticleType>
952 void Foam::InteractionLists<ParticleType>::prepareParticleToBeReferred
954     ParticleType* particle,
955     labelPair ciat
958     const globalIndexAndTransform& globalTransforms =
959         mesh_.globalData().globalTransforms();
961     const vectorTensorTransform& transform = globalTransforms.transform
962     (
963         globalTransforms.transformIndex(ciat)
964     );
966     particle->position() = transform.invTransformPosition(particle->position());
968     particle->transformProperties(-transform.t());
970     if (transform.hasR())
971     {
972         particle->transformProperties(transform.R().T());
973     }
977 template<class ParticleType>
978 void Foam::InteractionLists<ParticleType>::fillReferredParticleCloud()
980     if (writeCloud_)
981     {
982         forAll(referredParticles_, refCellI)
983         {
984             const IDLList<ParticleType>& refCell =
985                 referredParticles_[refCellI];
987             forAllConstIter(typename IDLList<ParticleType>, refCell, iter)
988             {
989                 cloud_.addParticle
990                 (
991                     static_cast<ParticleType*>(iter().clone().ptr())
992                 );
993             }
994         }
995     }
999 template<class ParticleType>
1000 void Foam::InteractionLists<ParticleType>::prepareWallDataToRefer()
1002     const globalIndexAndTransform& globalTransforms =
1003         mesh_.globalData().globalTransforms();
1005     referredWallData_.setSize
1006     (
1007         wallFaceIndexAndTransformToDistribute_.size()
1008     );
1010     const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
1012     forAll(referredWallData_, rWVI)
1013     {
1014         const labelPair& wfiat = wallFaceIndexAndTransformToDistribute_[rWVI];
1016         label wallFaceIndex = globalTransforms.index(wfiat);
1018         const vectorTensorTransform& transform = globalTransforms.transform
1019         (
1020             globalTransforms.transformIndex(wfiat)
1021         );
1023         label patchI = mesh_.boundaryMesh().patchID()
1024         [
1025             wallFaceIndex - mesh_.nInternalFaces()
1026         ];
1028         label patchFaceI =
1029             wallFaceIndex
1030           - mesh_.boundaryMesh()[patchI].start();
1032         // Need to transform velocity when tensor transforms are
1033         // supported
1034         referredWallData_[rWVI] = U.boundaryField()[patchI][patchFaceI];
1036         if (transform.hasR())
1037         {
1038             referredWallData_[rWVI] =
1039                 transform.R().T() & referredWallData_[rWVI];
1040         }
1041     }
1045 template<class ParticleType>
1046 void Foam::InteractionLists<ParticleType>::writeReferredWallFaces() const
1048     if (referredWallFaces_.empty())
1049     {
1050         return;
1051     }
1053     fileName objDir = mesh_.time().timePath()/cloud::prefix;
1055     mkDir(objDir);
1057     fileName objFileName = "referredWallFaces.obj";
1059     OFstream str(objDir/objFileName);
1061     Info<< "    Writing "
1062         << mesh_.time().timeName()/cloud::prefix/objFileName
1063         << endl;
1065     label offset = 1;
1067     forAll(referredWallFaces_, rWFI)
1068     {
1069         const referredWallFace& rwf = referredWallFaces_[rWFI];
1071         forAll(rwf, fPtI)
1072         {
1073             meshTools::writeOBJ(str, rwf.points()[rwf[fPtI]]);
1074         }
1076         str<< 'f';
1078         forAll(rwf, fPtI)
1079         {
1080             str<< ' ' << fPtI + offset;
1081         }
1083         str<< nl;
1085         offset += rwf.size();
1086     }
1090 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
1092 template<class ParticleType>
1093 Foam::InteractionLists<ParticleType>::InteractionLists(const polyMesh& mesh)
1095     mesh_(mesh),
1096     cloud_(mesh_, "NULL_Cloud", IDLList<ParticleType>()),
1097     writeCloud_(false),
1098     cellMapPtr_(),
1099     wallFaceMapPtr_(),
1100     maxDistance_(0.0),
1101     dil_(),
1102     dwfil_(),
1103     ril_(),
1104     rilInverse_(),
1105     cellIndexAndTransformToDistribute_(),
1106     wallFaceIndexAndTransformToDistribute_(),
1107     referredWallFaces_(),
1108     UName_("unknown_UName"),
1109     referredWallData_(),
1110     referredParticles_()
1114 template<class ParticleType>
1115 Foam::InteractionLists<ParticleType>::InteractionLists
1117     const polyMesh& mesh,
1118     scalar maxDistance,
1119     Switch writeCloud,
1120     const word& UName
1123     mesh_(mesh),
1124     cloud_(mesh_, "referredParticleCloud", IDLList<ParticleType>()),
1125     writeCloud_(writeCloud),
1126     cellMapPtr_(),
1127     wallFaceMapPtr_(),
1128     maxDistance_(maxDistance),
1129     dil_(),
1130     dwfil_(),
1131     ril_(),
1132     rilInverse_(),
1133     cellIndexAndTransformToDistribute_(),
1134     wallFaceIndexAndTransformToDistribute_(),
1135     referredWallFaces_(),
1136     UName_(UName),
1137     referredWallData_(),
1138     referredParticles_()
1140     buildInteractionLists();
1144 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
1146 template<class ParticleType>
1147 Foam::InteractionLists<ParticleType>::~InteractionLists()
1151 // * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * * //
1153 template<class ParticleType>
1154 void Foam::InteractionLists<ParticleType>::sendReferredData
1156     const List<DynamicList<ParticleType*> >& cellOccupancy,
1157     PstreamBuffers& pBufs
1160     if (mesh_.changing())
1161     {
1162         WarningIn
1163         (
1164             "void Foam::InteractionLists<ParticleType>::sendReferredData"
1165             "("
1166                 "const List<DynamicList<ParticleType*> >& cellOccupancy,"
1167                 "PstreamBuffers& pBufs"
1168             ")"
1169         )
1170             << "Mesh changing, rebuilding InteractionLists form scratch."
1171             << endl;
1173         buildInteractionLists();
1174     }
1176     prepareWallDataToRefer();
1178     prepareParticlesToRefer(cellOccupancy);
1180     for (label domain = 0; domain < Pstream::nProcs(); domain++)
1181     {
1182         const labelList& subMap = cellMap().subMap()[domain];
1184         if (subMap.size())
1185         {
1186             UOPstream toDomain(domain, pBufs);
1188             UIndirectList<IDLList<ParticleType> > subMappedParticles
1189             (
1190                 referredParticles_,
1191                 subMap
1192             );
1194             forAll(subMappedParticles, i)
1195             {
1196                 toDomain << subMappedParticles[i];
1197             }
1198         }
1199     }
1201     // Using the mapDistribute to start sending and receiving the
1202     // buffer but not block, i.e. it is calling
1203     //     pBufs.finishedSends(false);
1204     wallFaceMap().send(pBufs, referredWallData_);
1208 template<class ParticleType>
1209 void Foam::InteractionLists<ParticleType>::receiveReferredData
1211     PstreamBuffers& pBufs
1214     Pstream::waitRequests();
1216     referredParticles_.setSize(cellMap().constructSize());
1218     for (label domain = 0; domain < Pstream::nProcs(); domain++)
1219     {
1220         const labelList& constructMap = cellMap().constructMap()[domain];
1222         if (constructMap.size())
1223         {
1224             UIPstream str(domain, pBufs);
1226             forAll(constructMap, i)
1227             {
1228                 referredParticles_[constructMap[i]] = IDLList<ParticleType>
1229                 (
1230                     str,
1231                     typename ParticleType::iNew(mesh_)
1232                 );
1233             }
1234         }
1235     }
1237     fillReferredParticleCloud();
1239     wallFaceMap().receive(pBufs, referredWallData_);
1243 // ************************************************************************* //