inverseDistanceCellCellStencil.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2017-2023 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify i
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
30 #include "OBJstream.H"
31 #include "Time.H"
32 #include "fvMeshSubset.H"
33 
34 #include "globalIndex.H"
35 #include "oversetFvPatch.H"
37 #include "syncTools.H"
38 #include "treeBoundBoxList.H"
39 #include "waveMethod.H"
40 
41 #include "regionSplit.H"
42 #include "oversetFvPatchFields.H"
43 #include "topoDistanceData.H"
44 #include "FaceCellWave.H"
45 
46 #include "OBJstream.H"
48 
49 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 namespace cellCellStencils
54 {
55  defineTypeNameAndDebug(inverseDistance, 0);
57 }
58 }
59 
60 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
61 
63 (
64  const labelVector& nDivs,
65  const labelVector& ijk
66 )
67 {
68  return (ijk[0]*nDivs[1] + ijk[1])*nDivs[2] + ijk[2];
69 }
70 
71 
73 (
74  const labelVector& nDivs,
75  const label boxI
76 )
77 {
78  label ij = boxI/nDivs[2];
79  label k = boxI-ij*nDivs[2];
80  label i = ij/nDivs[1];
81  label j = ij-i*nDivs[1];
82 
83  return labelVector(i, j, k);
84 }
85 
86 
88 (
89  const boundBox& bb,
90  const labelVector& nDivs,
91  const point& pt
92 )
93 {
94  const vector d(bb.span());
95  const point relPt(pt-bb.min());
96 
97  return labelVector
98  (
99  floor(relPt[0]/d[0]*nDivs[0]),
100  floor(relPt[1]/d[1]*nDivs[1]),
101  floor(relPt[2]/d[2]*nDivs[2])
102  );
103 }
104 
105 
107 (
108  const boundBox& bb,
109  const labelVector& nDivs,
110  const label boxI
111 )
112 {
113  // Return midpoint of box indicated by boxI
114  labelVector ids(index3(nDivs, boxI));
115 
116  const vector d(bb.span());
117  const vector sz(d[0]/nDivs[0], d[1]/nDivs[1], d[2]/nDivs[2]);
118 
119  return bb.min()+0.5*sz+vector(sz[0]*ids[0], sz[1]*ids[1], sz[2]*ids[2]);
120 }
121 
122 
124 (
125  PackedList<2>& elems,
126  const boundBox& bb,
127  const labelVector& nDivs,
128  const boundBox& subBb,
129  const unsigned int val
130 )
131 {
132  labelVector minIds(index3(bb, nDivs, subBb.min()));
133  labelVector maxIds(index3(bb, nDivs, subBb.max()));
134 
135  for (direction cmpt = 0; cmpt < 3; cmpt++)
136  {
137  if (maxIds[cmpt] < 0 || minIds[cmpt] > nDivs[cmpt])
138  {
139  return;
140  }
141  }
142 
143  labelVector maxIndex(labelVector(nDivs[0]-1, nDivs[1]-1, nDivs[2]-1));
144  minIds = max(labelVector::zero, minIds);
145  maxIds = min(maxIndex, maxIds);
146 
147  for (label i = minIds[0]; i <= maxIds[0]; i++)
148  {
149  for (label j = minIds[1]; j <= maxIds[1]; j++)
150  {
151  for (label k = minIds[2]; k <= maxIds[2]; k++)
152  {
153  label i1 = index(nDivs, labelVector(i, j, k));
154  elems[i1] = val;
155  }
156  }
157  }
158 }
159 
160 
162 (
163  const fvMesh& mesh,
164  const vector& smallVec,
165 
166  const boundBox& bb,
167  const labelVector& nDivs,
168  PackedList<2>& patchTypes,
169 
170  const labelList& cellMap,
171  labelList& patchCellTypes
172 )
173 {
174  // Mark all voxels that overlap the bounding box of any patch
175 
176  const fvBoundaryMesh& pbm = mesh.boundary();
177 
178  patchTypes = patchCellType::OTHER;
179 
180  // Mark wall boundaries
181  forAll(pbm, patchI)
182  {
183  const fvPatch& fvp = pbm[patchI];
184  const labelList& fc = fvp.faceCells();
185 
186  if (!fvPatch::constraintType(fvp.type()))
187  {
188  //Info<< "Marking cells on proper patch " << fvp.name()
189  // << " with type " << fvp.type() << endl;
190  const polyPatch& pp = fvp.patch();
191  forAll(pp, i)
192  {
193  // Mark in overall patch types
194  patchCellTypes[cellMap[fc[i]]] = patchCellType::PATCH;
195 
196  // Mark in voxel mesh
197  boundBox faceBb(pp.points(), pp[i], false);
198  faceBb.grow(smallVec);
199 
200  if (bb.overlaps(faceBb))
201  {
202  fill(patchTypes, bb, nDivs, faceBb, patchCellType::PATCH);
203  }
204  }
205  }
206  }
207 
208  // Override with overset boundaries
209  forAll(pbm, patchI)
210  {
211  const fvPatch& fvp = pbm[patchI];
212  const labelList& fc = fvp.faceCells();
213 
214  if (isA<oversetFvPatch>(fvp))
215  {
216  //Info<< "Marking cells on overset patch " << fvp.name() << endl;
217  const polyPatch& pp = fvp.patch();
218  forAll(pp, i)
219  {
220  // Mark in overall patch types
221  patchCellTypes[cellMap[fc[i]]] = patchCellType::OVERSET;
222 
223  // Mark in voxel mesh
224  boundBox faceBb(pp.points(), pp[i], false);
225  faceBb.grow(smallVec);
226 
227  if (bb.overlaps(faceBb))
228  {
229  fill(patchTypes, bb, nDivs, faceBb, patchCellType::OVERSET);
230  }
231  }
232  }
233  }
234 }
235 
236 
238 (
239  const boundBox& bb,
240  const labelVector& nDivs,
241  const PackedList<2>& vals,
242  const treeBoundBox& subBb,
243  const unsigned int val
244 )
245 {
246  // Checks if subBb overlaps any voxel set to val
247 
248  labelVector minIds(index3(bb, nDivs, subBb.min()));
249  labelVector maxIds(index3(bb, nDivs, subBb.max()));
250 
251  for (direction cmpt = 0; cmpt < 3; cmpt++)
252  {
253  if (maxIds[cmpt] < 0 || minIds[cmpt] > nDivs[cmpt])
254  {
255  return false;
256  }
257  }
258 
259  labelVector maxIndex(labelVector(nDivs[0]-1, nDivs[1]-1, nDivs[2]-1));
260  minIds = max(labelVector::zero, minIds);
261  maxIds = min(maxIndex, maxIds);
262 
263  for (label i = minIds[0]; i <= maxIds[0]; i++)
264  {
265  for (label j = minIds[1]; j <= maxIds[1]; j++)
266  {
267  for (label k = minIds[2]; k <= maxIds[2]; k++)
268  {
269  label i1 = index(nDivs, labelVector(i, j, k));
270  if (vals[i1] == patchCellType::PATCH)
271  {
272  return true;
273  }
274  }
275  }
276  }
277  return false;
278 }
279 
280 
282 (
283  PstreamBuffers& pBufs,
284 
285  const PtrList<fvMeshSubset>& meshParts,
286 
287  const List<treeBoundBoxList>& patchBb,
288  const List<labelVector>& patchDivisions,
289  const PtrList<PackedList<2>>& patchParts,
290 
291  const label srcI,
292  const label tgtI,
293  labelList& allCellTypes
294 ) const
295 {
296  const treeBoundBoxList& srcPatchBbs = patchBb[srcI];
297  const treeBoundBoxList& tgtPatchBbs = patchBb[tgtI];
298  const labelList& tgtCellMap = meshParts[tgtI].cellMap();
299 
300  // 1. do processor-local src-tgt patch overlap
301  {
302  const treeBoundBox& srcPatchBb = srcPatchBbs[Pstream::myProcNo()];
303  const treeBoundBox& tgtPatchBb = tgtPatchBbs[Pstream::myProcNo()];
304 
305  if (srcPatchBb.overlaps(tgtPatchBb))
306  {
307  const PackedList<2>& srcPatchTypes = patchParts[srcI];
308  const labelVector& zoneDivs = patchDivisions[srcI];
309 
310  forAll(tgtCellMap, tgtCelli)
311  {
312  label celli = tgtCellMap[tgtCelli];
313  treeBoundBox cBb(mesh_.cellBb(celli));
314  cBb.grow(smallVec_);
315 
316  if
317  (
318  overlaps
319  (
320  srcPatchBb,
321  zoneDivs,
322  srcPatchTypes,
323  cBb,
325  )
326  )
327  {
328  allCellTypes[celli] = HOLE;
329  }
330  }
331  }
332  }
333 
334 
335  // 2. Send over srcMesh bits that overlap tgt and do calculation
336  pBufs.clear();
337  for (const int procI : Pstream::allProcs())
338  {
339  if (procI != Pstream::myProcNo())
340  {
341  const treeBoundBox& srcPatchBb = srcPatchBbs[Pstream::myProcNo()];
342  const treeBoundBox& tgtPatchBb = tgtPatchBbs[procI];
343 
344  if (srcPatchBb.overlaps(tgtPatchBb))
345  {
346  // Send over complete patch voxel map. Tbd: could
347  // subset
348  UOPstream os(procI, pBufs);
349  os << srcPatchBb << patchDivisions[srcI] << patchParts[srcI];
350  }
351  }
352  }
353  pBufs.finishedSends();
354  for (const int procI : Pstream::allProcs())
355  {
356  if (procI != Pstream::myProcNo())
357  {
358  const treeBoundBox& srcPatchBb = srcPatchBbs[procI];
359  const treeBoundBox& tgtPatchBb = tgtPatchBbs[Pstream::myProcNo()];
360 
361  if (srcPatchBb.overlaps(tgtPatchBb))
362  {
363  UIPstream is(procI, pBufs);
364  const treeBoundBox receivedBb(is);
365  const labelVector zoneDivs(is);
366  const PackedList<2> srcPatchTypes(is);
367 
368  // Verify validity
369  if (srcPatchBb != receivedBb)
370  {
372  << "proc:" << procI
373  << " srcPatchBb:" << srcPatchBb
374  << " receivedBb:" << receivedBb
375  << exit(FatalError);
376  }
377 
378  forAll(tgtCellMap, tgtCelli)
379  {
380  label celli = tgtCellMap[tgtCelli];
381  treeBoundBox cBb(mesh_.cellBb(celli));
382  cBb.grow(smallVec_);
383 
384  if
385  (
386  overlaps
387  (
388  srcPatchBb,
389  zoneDivs,
390  srcPatchTypes,
391  cBb,
393  )
394  )
395  {
396  allCellTypes[celli] = HOLE;
397  }
398  }
399  }
400  }
401  }
402 }
403 
404 
406 (
407  const label destMesh,
408  const label currentDonorMesh,
409  const label newDonorMesh
410 ) const
411 {
412  // This determines for multiple overlapping meshes which one provides
413  // the best donors. Is very basic and only looks at indices of meshes:
414  // - 'nearest' mesh index wins, i.e. on mesh 0 it preferentially uses donors
415  // from mesh 1 over mesh 2 (if applicable)
416  // - if same 'distance' the highest mesh wins. So on mesh 1 it
417  // preferentially uses donors from mesh 2 over mesh 0. This particular
418  // rule helps to avoid some interpolation loops where mesh 1 uses donors
419  // from mesh 0 (usually the background) but mesh 0 then uses
420  // donors from 1.
421 
422  if (currentDonorMesh == -1)
423  {
424  return true;
425  }
426  else
427  {
428  const label currentDist = mag(currentDonorMesh-destMesh);
429  const label newDist = mag(newDonorMesh-destMesh);
430 
431  if (newDist < currentDist)
432  {
433  return true;
434  }
435  else if (newDist == currentDist && newDonorMesh > currentDonorMesh)
436  {
437  return true;
438  }
439  else
440  {
441  return false;
442  }
443  }
444 }
445 
446 
448 (
449  const globalIndex& globalCells,
450  PstreamBuffers& pBufs,
451  const PtrList<fvMeshSubset>& meshParts,
452  const List<treeBoundBoxList>& meshBb,
453 
454  const labelList& allCellTypes,
455 
456  const label srcI,
457  const label tgtI,
458  labelListList& allStencil,
459  labelList& allDonor
460 ) const
461 {
462  const treeBoundBoxList& srcBbs = meshBb[srcI];
463  const treeBoundBoxList& tgtBbs = meshBb[tgtI];
464 
465  const fvMesh& srcMesh = meshParts[srcI].subMesh();
466  const labelList& srcCellMap = meshParts[srcI].cellMap();
467  const fvMesh& tgtMesh = meshParts[tgtI].subMesh();
468  const pointField& tgtCc = tgtMesh.cellCentres();
469  const labelList& tgtCellMap = meshParts[tgtI].cellMap();
470 
471  // 1. do processor-local src/tgt overlap
472  {
473  labelList tgtToSrcAddr;
474  waveMethod::calculate(tgtMesh, srcMesh, tgtToSrcAddr);
475  forAll(tgtCellMap, tgtCelli)
476  {
477  const label srcCelli = tgtToSrcAddr[tgtCelli];
478  if
479  (
480  srcCelli != -1
481  && allCellTypes[tgtCellMap[tgtCelli]] != HOLE
482  )
483  {
484  label celli = tgtCellMap[tgtCelli];
485 
486  // TBD: check for multiple donors. Maybe better one? For
487  // now check 'nearer' mesh
488  if (betterDonor(tgtI, allDonor[celli], srcI))
489  {
490  label globalDonor =
491  globalCells.toGlobal(srcCellMap[srcCelli]);
492  allStencil[celli].setSize(1);
493  allStencil[celli][0] = globalDonor;
494  allDonor[celli] = srcI;
495  }
496  }
497  }
498  }
499 
500 
501  // 2. Send over tgtMesh bits that overlap src and do calculation on
502  // srcMesh.
503 
504 
505  // (remote) processors where the tgt overlaps my src
506  DynamicList<label> tgtOverlapProcs(Pstream::nProcs());
507  // (remote) processors where the src overlaps my tgt
508  DynamicList<label> srcOverlapProcs(Pstream::nProcs());
509  for (const int procI : Pstream::allProcs())
510  {
511  if (procI != Pstream::myProcNo())
512  {
513  if (tgtBbs[procI].overlaps(srcBbs[Pstream::myProcNo()]))
514  {
515  tgtOverlapProcs.append(procI);
516  }
517  if (srcBbs[procI].overlaps(tgtBbs[Pstream::myProcNo()]))
518  {
519  srcOverlapProcs.append(procI);
520  }
521  }
522  }
523 
524 
525  // Indices of tgtcells to send over to each processor
526  List<DynamicList<label>> tgtSendCells(Pstream::nProcs());
527  forAll(srcOverlapProcs, i)
528  {
529  label procI = srcOverlapProcs[i];
530  tgtSendCells[procI].reserve(tgtMesh.nCells()/srcOverlapProcs.size());
531  }
532 
533 
534  forAll(tgtCellMap, tgtCelli)
535  {
536  label celli = tgtCellMap[tgtCelli];
537  if (srcOverlapProcs.size())
538  {
539  treeBoundBox subBb(mesh_.cellBb(celli));
540  subBb.grow(smallVec_);
541 
542  forAll(srcOverlapProcs, i)
543  {
544  const label procI = srcOverlapProcs[i];
545  if (subBb.overlaps(srcBbs[procI]))
546  {
547  tgtSendCells[procI].append(tgtCelli);
548  }
549  }
550  }
551  }
552 
553  // Send target cell centres to overlapping processors
554  pBufs.clear();
555 
556  forAll(srcOverlapProcs, i)
557  {
558  label procI = srcOverlapProcs[i];
559  const labelList& cellIDs = tgtSendCells[procI];
560 
561  UOPstream os(procI, pBufs);
562  os << UIndirectList<point>(tgtCc, cellIDs);
563  }
564  pBufs.finishedSends();
565 
566  // Receive bits of target processors; find; send back
567  (void)srcMesh.tetBasePtIs();
568  forAll(tgtOverlapProcs, i)
569  {
570  label procI = tgtOverlapProcs[i];
571 
572  UIPstream is(procI, pBufs);
573  pointList samples(is);
574 
575  labelList donors(samples.size(), -1);
576  forAll(samples, sampleI)
577  {
578  const point& sample = samples[sampleI];
579  label srcCelli = srcMesh.findCell(sample, polyMesh::CELL_TETS);
580  if (srcCelli != -1)
581  {
582  donors[sampleI] = globalCells.toGlobal(srcCellMap[srcCelli]);
583  }
584  }
585 
586  // Use same pStreamBuffers to send back.
587  UOPstream os(procI, pBufs);
588  os << donors;
589  }
590  pBufs.finishedSends();
591 
592  forAll(srcOverlapProcs, i)
593  {
594  label procI = srcOverlapProcs[i];
595  const labelList& cellIDs = tgtSendCells[procI];
596 
597  UIPstream is(procI, pBufs);
598  labelList donors(is);
599 
600  if (donors.size() != cellIDs.size())
601  {
602  FatalErrorInFunction<< "problem : cellIDs:" << cellIDs.size()
603  << " donors:" << donors.size() << abort(FatalError);
604  }
605 
606  forAll(donors, donorI)
607  {
608  label globalDonor = donors[donorI];
609 
610  if (globalDonor != -1)
611  {
612  label celli = tgtCellMap[cellIDs[donorI]];
613  {
614  // TBD: check for multiple donors. Maybe better one? For
615  // now check 'nearer' mesh
616  if (betterDonor(tgtI, allDonor[celli], srcI))
617  {
618  allStencil[celli].setSize(1);
619  allStencil[celli][0] = globalDonor;
620  allDonor[celli] = srcI;
621  }
622  }
623  }
624  }
625  }
626 
627 }
628 
629 
630 //void Foam::cellCellStencils::inverseDistance::uncompactedRegionSplit
631 //(
632 // const fvMesh& mesh,
633 // const globalIndex& globalFaces,
634 // const label nZones,
635 // const labelList& zoneID,
636 // const labelList& cellTypes,
637 // const boolList& isBlockedFace,
638 // labelList& cellRegion
639 //) const
640 //{
641 // // Pass 1: locally seed 2 cells per zone (one unblocked, one blocked).
642 // // This avoids excessive numbers of front
643 //
644 // // Field on cells and faces.
645 // List<minData> cellData(mesh.nCells());
646 // List<minData> faceData(mesh.nFaces());
647 //
648 // // Take over blockedFaces by seeding a negative number
649 // // (so is always less than the decomposition)
650 //
651 // forAll(isBlockedFace, facei)
652 // {
653 // if (isBlockedFace[facei])
654 // {
655 // faceData[facei] = minData(-2);
656 // }
657 // }
658 //
659 //
660 // labelList seedFace(nZones, -1);
661 //
662 // const labelList& owner = mesh.faceOwner();
663 // const labelList& neighbour = mesh.faceNeighbour();
664 //
665 // forAll(owner, facei)
666 // {
667 // label own = owner[facei];
668 // if (seedFace[zoneID[own]] == -1)
669 // {
670 // if (cellTypes[own] != HOLE)
671 // {
672 // const cell& cFaces = mesh.cells()[own];
673 // forAll(cFaces, i)
674 // {
675 // if (!isBlockedFace[cFaces[i]])
676 // {
677 // seedFace[zoneID[own]] = cFaces[i];
678 // }
679 // }
680 // }
681 // }
682 // }
683 // forAll(neighbour, facei)
684 // {
685 // label nei = neighbour[facei];
686 // if (seedFace[zoneID[nei]] == -1)
687 // {
688 // if (cellTypes[nei] != HOLE)
689 // {
690 // const cell& cFaces = mesh.cells()[nei];
691 // forAll(cFaces, i)
692 // {
693 // if (!isBlockedFace[cFaces[i]])
694 // {
695 // seedFace[zoneID[nei]] = cFaces[i];
696 // }
697 // }
698 // }
699 // }
700 // }
701 //
702 // DynamicList<label> seedFaces(nZones);
703 // DynamicList<minData> seedData(seedFaces.size());
704 // forAll(seedFace, zonei)
705 // {
706 // if (seedFace[zonei] != -1)
707 // {
708 // seedFaces.append(seedFace[zonei]);
709 // seedData.append(minData(globalFaces.toGlobal(seedFace[zonei])));
710 // }
711 // }
712 //
713 // // Propagate information inwards
714 // FaceCellWave<minData> deltaCalc
715 // (
716 // mesh,
717 // List<labelPair>(),
718 // false, // disable walking through cyclicAMI for backwards
719 // // compatibility
720 // seedFaces,
721 // seedData,
722 // faceData,
723 // cellData,
724 // mesh.globalData().nTotalCells()+1
725 // );
726 //
727 // // Extract
728 // cellRegion.setSize(mesh.nCells());
729 // forAll(cellRegion, celli)
730 // {
731 // if (cellData[celli].valid(deltaCalc.data()))
732 // {
733 // cellRegion[celli] = cellData[celli].data();
734 // }
735 // else
736 // {
737 // // Unvisited cell -> only possible if surrounded by blocked faces.
738 // // If so make up region from any of the faces
739 // const cell& cFaces = mesh.cells()[celli];
740 // label facei = cFaces[0];
741 // cellRegion[celli] = globalFaces.toGlobal(facei);
742 // }
743 // }
744 //}
745 //Foam::autoPtr<Foam::globalIndex>
746 //Foam::cellCellStencils::inverseDistance::compactedRegionSplit
747 //(
748 // const fvMesh& mesh,
749 // const globalIndex& globalRegions,
750 // labelList& cellRegion
751 //) const
752 //{
753 // // Now our cellRegion will have
754 // // - non-local regions (i.e. originating from other processors)
755 // // - non-compact locally originating regions
756 // // so we'll need to compact
757 //
758 // // 4a: count per originating processor the number of regions
759 // labelList nOriginating(Pstream::nProcs(), Zero);
760 // {
761 // labelHashSet haveRegion(mesh.nCells()/8);
762 //
763 // forAll(cellRegion, celli)
764 // {
765 // label region = cellRegion[celli];
766 // label proci = globalRegions.whichProcID(region);
767 // if (haveRegion.insert(region))
768 // {
769 // nOriginating[proci]++;
770 // }
771 // }
772 // }
773 //
774 // if (debug)
775 // {
776 // Pout<< "Counted " << nOriginating[Pstream::myProcNo()]
777 // << " local regions." << endl;
778 // }
779 //
780 //
781 // // Global numbering for compacted local regions
782 // autoPtr<globalIndex> globalCompactPtr
783 // (
784 // new globalIndex(nOriginating[Pstream::myProcNo()])
785 // );
786 // const globalIndex& globalCompact = globalCompactPtr();
787 //
788 //
789 // // 4b: renumber
790 // // Renumber into compact indices. Note that since we've already made
791 // // all regions global we now need a Map to store the compacting
792 // // information
793 // // instead of a labelList - otherwise we could have used a straight
794 // // labelList.
795 //
796 // // Local compaction map
797 // Map<label> globalToCompact(2*nOriginating[Pstream::myProcNo()]);
798 // // Remote regions we want the compact number for
799 // List<labelHashSet> nonLocal(Pstream::nProcs());
800 // forAll(nonLocal, proci)
801 // {
802 // if (proci != Pstream::myProcNo())
803 // {
804 // nonLocal[proci].reserve(nOriginating[proci]);
805 // }
806 // }
807 //
808 // forAll(cellRegion, celli)
809 // {
810 // label region = cellRegion[celli];
811 // if (globalRegions.isLocal(region))
812 // {
813 // // Insert new compact region (if not yet present)
814 // globalToCompact.insert
815 // (
816 // region,
817 // globalCompact.toGlobal(globalToCompact.size())
818 // );
819 // }
820 // else
821 // {
822 // nonLocal[globalRegions.whichProcID(region)].insert(region);
823 // }
824 // }
825 //
826 //
827 // // Now we have all the local regions compacted. Now we need to get the
828 // // non-local ones from the processors to whom they are local.
829 // // Convert the nonLocal (labelHashSets) to labelLists.
830 //
831 // labelListList sendNonLocal(Pstream::nProcs());
832 // forAll(sendNonLocal, proci)
833 // {
834 // sendNonLocal[proci] = nonLocal[proci].toc();
835 // }
836 //
837 // if (debug)
838 // {
839 // forAll(sendNonLocal, proci)
840 // {
841 // Pout<< " from processor " << proci
842 // << " want " << sendNonLocal[proci].size()
843 // << " region numbers." << endl;
844 // }
845 // Pout<< endl;
846 // }
847 //
848 //
849 // // Get the wanted region labels into recvNonLocal
850 // labelListList recvNonLocal;
851 // Pstream::exchange<labelList, label>(sendNonLocal, recvNonLocal);
852 //
853 // // Now we have the wanted compact region labels that proci wants in
854 // // recvNonLocal[proci]. Construct corresponding list of compact
855 // // region labels to send back.
856 //
857 // labelListList sendWantedLocal(Pstream::nProcs());
858 // forAll(recvNonLocal, proci)
859 // {
860 // const labelList& nonLocal = recvNonLocal[proci];
861 // sendWantedLocal[proci].setSize(nonLocal.size());
862 //
863 // forAll(nonLocal, i)
864 // {
865 // sendWantedLocal[proci][i] = globalToCompact[nonLocal[i]];
866 // }
867 // }
868 //
869 //
870 // // Send back (into recvNonLocal)
871 // recvNonLocal.clear();
872 // Pstream::exchange<labelList, label>(sendWantedLocal, recvNonLocal);
873 // sendWantedLocal.clear();
874 //
875 // // Now recvNonLocal contains for every element in setNonLocal the
876 // // corresponding compact number. Insert these into the local compaction
877 // // map.
878 //
879 // forAll(recvNonLocal, proci)
880 // {
881 // const labelList& wantedRegions = sendNonLocal[proci];
882 // const labelList& compactRegions = recvNonLocal[proci];
883 //
884 // forAll(wantedRegions, i)
885 // {
886 // globalToCompact.insert(wantedRegions[i], compactRegions[i]);
887 // }
888 // }
889 //
890 // // Finally renumber the regions
891 // forAll(cellRegion, celli)
892 // {
893 // cellRegion[celli] = globalToCompact[cellRegion[celli]];
894 // }
895 //
896 // return globalCompactPtr;
897 //}
898 
899 
901 (
902  const globalIndex& globalCells,
903  const fvMesh& mesh,
904  const labelList& zoneID,
905  const labelListList& stencil,
907 ) const
908 {
909  const fvBoundaryMesh& pbm = mesh.boundary();
910  const labelList& own = mesh.faceOwner();
911  const labelList& nei = mesh.faceNeighbour();
912 
913 
914  // The input cellTypes will be
915  // - HOLE : cell part covered by other-mesh patch
916  // - INTERPOLATED : cell fully covered by other-mesh patch
917  // or next to 'overset' patch
918  // - CALCULATED : otherwise
919  //
920  // so we start a walk from our patches and any cell we cannot reach
921  // (because we walk is stopped by other-mesh patch) is a hole.
922 
923 
924  DebugInfo<< FUNCTION_NAME << " : Starting hole flood filling" << endl;
925 
926  DebugInfo<< FUNCTION_NAME << " : Starting hole cells : "
927  << findIndices(cellTypes, HOLE).size() << endl;
928 
929  boolList isBlockedFace(mesh.nFaces(), false);
930  label nBlocked = 0;
931 
932  for (label faceI = 0; faceI < mesh.nInternalFaces(); faceI++)
933  {
934  label ownType = cellTypes[own[faceI]];
935  label neiType = cellTypes[nei[faceI]];
936  if
937  (
938  (ownType == HOLE && neiType != HOLE)
939  || (ownType != HOLE && neiType == HOLE)
940  )
941  {
942  isBlockedFace[faceI] = true;
943  nBlocked++;
944  }
945  }
946  DebugInfo<< FUNCTION_NAME << " : Marked internal hole boundaries : "
947  << nBlocked << endl;
948 
949 
950  labelList nbrCellTypes;
952 
953  for (label faceI = mesh.nInternalFaces(); faceI < mesh.nFaces(); faceI++)
954  {
955  label ownType = cellTypes[own[faceI]];
956  label neiType = nbrCellTypes[faceI-mesh.nInternalFaces()];
957 
958  if
959  (
960  (ownType == HOLE && neiType != HOLE)
961  || (ownType != HOLE && neiType == HOLE)
962  )
963  {
964  isBlockedFace[faceI] = true;
965  nBlocked++;
966  }
967  }
968 
969  DebugInfo<< FUNCTION_NAME << " : Marked all hole boundaries : "
970  << nBlocked << endl;
971 
972  // Determine regions
973  regionSplit cellRegion(mesh, isBlockedFace);
974  const label nRegions = cellRegion.nRegions();
975 
976  //labelList cellRegion;
977  //label nRegions = -1;
978  //{
979  // const globalIndex globalFaces(mesh.nFaces());
980  // uncompactedRegionSplit
981  // (
982  // mesh,
983  // globalFaces,
984  // gMax(zoneID)+1,
985  // zoneID,
986  // cellTypes,
987  // isBlockedFace,
988  // cellRegion
989  // );
990  // autoPtr<globalIndex> globalRegions
991  // (
992  // compactedRegionSplit
993  // (
994  // mesh,
995  // globalFaces,
996  // cellRegion
997  // )
998  // );
999  // nRegions = globalRegions().size();
1000  //}
1001  DebugInfo<< FUNCTION_NAME << " : Determined regions : "
1002  << nRegions << endl;
1003 
1004  //Info<< typeName << " : detected " << nRegions
1005  // << " mesh regions after overset" << nl << endl;
1006 
1007 
1008 
1009  // Now we'll have a mesh split according to where there are cells
1010  // covered by the other-side patches. See what we can reach from our
1011  // real patches
1012 
1013  // 0 : region not yet determined
1014  // 1 : borders blockage so is not ok (but can be overridden by real
1015  // patch)
1016  // 2 : has real patch in it so is reachable
1017  labelList regionType(nRegions, Zero);
1018 
1019 
1020  // See if any regions borders blockage. Note: isBlockedFace is already
1021  // parallel synchronised.
1022  {
1023  for (label faceI = 0; faceI < mesh.nInternalFaces(); faceI++)
1024  {
1025  if (isBlockedFace[faceI])
1026  {
1027  label ownRegion = cellRegion[own[faceI]];
1028 
1029  if (cellTypes[own[faceI]] != HOLE)
1030  {
1031  if (regionType[ownRegion] == 0)
1032  {
1033  regionType[ownRegion] = 1;
1034  }
1035  }
1036 
1037  label neiRegion = cellRegion[nei[faceI]];
1038 
1039  if (cellTypes[nei[faceI]] != HOLE)
1040  {
1041  if (regionType[neiRegion] == 0)
1042  {
1043  regionType[neiRegion] = 1;
1044  }
1045  }
1046  }
1047  }
1048  for
1049  (
1050  label faceI = mesh.nInternalFaces();
1051  faceI < mesh.nFaces();
1052  faceI++
1053  )
1054  {
1055  if (isBlockedFace[faceI])
1056  {
1057  label ownRegion = cellRegion[own[faceI]];
1058 
1059  if (regionType[ownRegion] == 0)
1060  {
1061  regionType[ownRegion] = 1;
1062  }
1063  }
1064  }
1065  }
1066 
1067 
1068  // Override with real patches
1069  forAll(pbm, patchI)
1070  {
1071  const fvPatch& fvp = pbm[patchI];
1072 
1073  if (isA<oversetFvPatch>(fvp))
1074  {}
1075  else if (!fvPatch::constraintType(fvp.type()))
1076  {
1077  const labelList& fc = fvp.faceCells();
1078  forAll(fc, i)
1079  {
1080  label regionI = cellRegion[fc[i]];
1081 
1082  if (cellTypes[fc[i]] != HOLE && regionType[regionI] != 2)
1083  {
1084  regionType[regionI] = 2;
1085  }
1086  }
1087  }
1088  }
1089 
1090  DebugInfo<< FUNCTION_NAME << " : Done local analysis" << endl;
1091 
1092  // Now we've handled
1093  // - cells next to blocked cells
1094  // - coupled boundaries
1095  // Only thing to handle is the interpolation between regions
1096 
1097 
1098  labelListList compactStencil(stencil);
1099  List<Map<label>> compactMap;
1100  mapDistribute map(globalCells, compactStencil, compactMap);
1101 
1102  DebugInfo<< FUNCTION_NAME << " : Converted stencil into compact form"
1103  << endl;
1104 
1105 
1106  while (true)
1107  {
1108  // Synchronise region status on processors
1109  // (could instead swap status through processor patches)
1111 
1112  DebugInfo<< FUNCTION_NAME << " : Gathered region type" << endl;
1113 
1114  // Communicate region status through interpolative cells
1115  labelList cellRegionType(labelUIndList(regionType, cellRegion));
1116  map.distribute(cellRegionType);
1117 
1118  DebugInfo<< FUNCTION_NAME << " : Interpolated region type" << endl;
1119 
1120 
1121 
1122  label nChanged = 0;
1123  forAll(pbm, patchI)
1124  {
1125  const fvPatch& fvp = pbm[patchI];
1126 
1127  if (isA<oversetFvPatch>(fvp))
1128  {
1129  const labelUList& fc = fvp.faceCells();
1130  forAll(fc, i)
1131  {
1132  label cellI = fc[i];
1133  label regionI = cellRegion[cellI];
1134 
1135  if (regionType[regionI] != 2)
1136  {
1137  const labelList& slots = compactStencil[cellI];
1138  forAll(slots, i)
1139  {
1140  label otherType = cellRegionType[slots[i]];
1141 
1142  if (otherType == 2)
1143  {
1144  //Pout<< "Reachable through interpolation : "
1145  // << regionI << " at cell "
1146  // << mesh.cellCentres()[cellI] << endl;
1147  regionType[regionI] = 2;
1148  nChanged++;
1149  break;
1150  }
1151  }
1152  }
1153  }
1154  }
1155  }
1156 
1157  reduce(nChanged, sumOp<label>());
1158  DebugInfo<< FUNCTION_NAME << " : Determined regions changed : "
1159  << nChanged << endl;
1160 
1161  if (nChanged == 0)
1162  {
1163  break;
1164  }
1165  }
1166 
1167 
1168  // See which regions have not been visited (regionType == 1)
1169  // label count = 0;
1170  forAll(cellRegion, cellI)
1171  {
1172  label type = regionType[cellRegion[cellI]];
1173  if (type == 1 && cellTypes[cellI] != HOLE)
1174  {
1175  cellTypes[cellI] = HOLE;
1176  // ++count;
1177  }
1178  }
1179 }
1180 
1181 
1183 (
1184  const point& sample,
1185  const pointList& donorCcs,
1186  scalarList& weights
1187 ) const
1188 {
1189  // Inverse-distance weighting
1190 
1191  weights.setSize(donorCcs.size());
1192  scalar sum = 0;
1193  forAll(donorCcs, i)
1194  {
1195  const scalar d = mag(sample-donorCcs[i]);
1196 
1197  if (d > ROOTVSMALL)
1198  {
1199  weights[i] = scalar(1)/d;
1200  sum += weights[i];
1201  }
1202  else
1203  {
1204  // Short circuit
1205  weights = scalar(0);
1206  weights[i] = scalar(1);
1207  return;
1208  }
1209  }
1210  forAll(weights, i)
1211  {
1212  weights[i] /= sum;
1213  }
1214 }
1215 
1216 
1218 (
1219  const globalIndex& globalCells
1220 )
1221 {
1222  // Detects holes that are used for interpolation. If so
1223  // - make type interpolated
1224  // - add interpolation to nearest non-hole cell(s)
1225 
1226  const labelList& owner = mesh_.faceOwner();
1227  const labelList& neighbour = mesh_.faceNeighbour();
1228 
1229 
1230  // Get hole-cells that are used in an interpolation stencil
1231  boolList isDonor(cellInterpolationMap().constructSize(), false);
1232  label nHoleDonors = 0;
1233  {
1234  for (const label celli : interpolationCells_)
1235  {
1236  const labelList& slots = cellStencil_[celli];
1237  UIndirectList<bool>(isDonor, slots) = true;
1238  }
1239 
1240  mapDistributeBase::distribute<bool, orEqOp<bool>, flipOp>
1241  (
1243  List<labelPair>(),
1244  mesh_.nCells(),
1245  cellInterpolationMap().constructMap(),
1246  false,
1247  cellInterpolationMap().subMap(),
1248  false,
1249  isDonor,
1250  false,
1251  orEqOp<bool>(),
1252  flipOp() // negateOp
1253  );
1254 
1255 
1256  nHoleDonors = 0;
1257  forAll(cellTypes_, celli)
1258  {
1259  if (cellTypes_[celli] == cellCellStencil::HOLE && isDonor[celli])
1260  {
1261  nHoleDonors++;
1262  }
1263  }
1264  reduce(nHoleDonors, sumOp<label>());
1265  if (debug)
1266  {
1267  Pout<< "Detected " << nHoleDonors
1268  << " hole cells that are used as donors" << endl;
1269  }
1270  }
1271 
1272 
1273  if (nHoleDonors)
1274  {
1275  // Convert map and stencil back to global indices by 'interpolating'
1276  // global cells. Note: since we're only adding interpolations
1277  // (HOLE->SPECIAL) this could probably be done more efficiently.
1278 
1279  labelList globalCellIDs(identity(mesh_.nCells()));
1280  globalCells.inplaceToGlobal(Pstream::myProcNo(), globalCellIDs);
1281  cellInterpolationMap().distribute(globalCellIDs);
1282 
1283  forAll(interpolationCells_, i)
1284  {
1285  label cellI = interpolationCells_[i];
1286  const labelList& slots = cellStencil_[cellI];
1287  cellStencil_[cellI] = UIndirectList<label>(globalCellIDs, slots)();
1288  }
1289 
1290 
1291  labelList nbrTypes;
1292  syncTools::swapBoundaryCellList(mesh_, cellTypes_, nbrTypes);
1293 
1294  label nSpecialNear = 0;
1295  label nSpecialFar = 0;
1296 
1297 
1298  // 1. Find hole cells next to live cells
1299  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1300 
1301  {
1302  List<pointList> donorCcs(mesh_.nCells());
1303  labelListList donorCells(mesh_.nCells());
1304 
1305  for (label facei = 0; facei < mesh_.nInternalFaces(); ++facei)
1306  {
1307  const label own = owner[facei];
1308  const bool ownHole = (cellTypes_[own] == cellCellStencil::HOLE);
1309  const label nbr = neighbour[facei];
1310  const bool nbrHole = (cellTypes_[nbr] == cellCellStencil::HOLE);
1311 
1312  if (isDonor[own] && ownHole && !nbrHole)
1313  {
1314  donorCells[own].append(globalCells.toGlobal(nbr));
1315  donorCcs[own].append(mesh_.cellCentres()[nbr]);
1316  }
1317  else if (isDonor[nbr] && nbrHole && !ownHole)
1318  {
1319  donorCells[nbr].append(globalCells.toGlobal(own));
1320  donorCcs[nbr].append(mesh_.cellCentres()[own]);
1321  }
1322  }
1323  labelList globalCellIDs(identity(mesh_.nCells()));
1324  globalCells.inplaceToGlobal(Pstream::myProcNo(), globalCellIDs);
1325  labelList nbrCells;
1326  syncTools::swapBoundaryCellList(mesh_, globalCellIDs, nbrCells);
1327  pointField nbrCc;
1328  syncTools::swapBoundaryCellList(mesh_, mesh_.cellCentres(), nbrCc);
1329  forAll(nbrTypes, bFacei)
1330  {
1331  const label facei = bFacei+mesh_.nInternalFaces();
1332  const label own = owner[facei];
1333  const bool ownHole = (cellTypes_[own] == cellCellStencil::HOLE);
1334  const bool nbrHole = (nbrTypes[bFacei] == cellCellStencil::HOLE);
1335 
1336  if (isDonor[own] && ownHole && !nbrHole)
1337  {
1338  donorCells[own].append(nbrCells[bFacei]);
1339  donorCcs[own].append(nbrCc[bFacei]);
1340  }
1341  }
1342 
1343  forAll(donorCcs, celli)
1344  {
1345  if (donorCcs[celli].size())
1346  {
1347  cellStencil_[celli] = std::move(donorCells[celli]);
1348  stencilWeights
1349  (
1350  mesh_.cellCentres()[celli],
1351  donorCcs[celli],
1352  cellInterpolationWeights_[celli]
1353  );
1354  cellTypes_[celli] = SPECIAL;
1355  interpolationCells_.append(celli);
1356  cellInterpolationWeight_[celli] = scalar(1);
1357  nSpecialNear++;
1358  }
1359  }
1360  }
1361 
1362 
1363 
1364 
1365  // 2. Walk to find (topologically) nearest 'live' cell
1366  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1367  // This will do the remainder but will only find a single donor
1368 
1369  // Field on cells and faces.
1370  List<topoDistanceData<int>> cellData(mesh_.nCells());
1371  List<topoDistanceData<int>> faceData(mesh_.nFaces());
1372 
1373  int dummyTrackData(0);
1374 
1375  {
1376  // Mark all non-hole cells to disable any wave across them
1377  forAll(cellTypes_, celli)
1378  {
1379  label cType = cellTypes_[celli];
1380  if
1381  (
1383  || cType == cellCellStencil::INTERPOLATED
1384  )
1385  {
1386  cellData[celli] = topoDistanceData<int>(labelMin, 0);
1387  }
1388  }
1389 
1390  DynamicList<label> seedFaces(nHoleDonors);
1391  DynamicList<topoDistanceData<int>> seedData(nHoleDonors);
1392 
1393  for (label facei = 0; facei < mesh_.nInternalFaces(); ++facei)
1394  {
1395  const label own = owner[facei];
1396  const bool ownLive =
1397  (
1398  cellTypes_[own] == cellCellStencil::CALCULATED
1399  || cellTypes_[own] == cellCellStencil::INTERPOLATED
1400  );
1401 
1402  const label nbr = neighbour[facei];
1403  const bool nbrLive =
1404  (
1405  cellTypes_[nbr] == cellCellStencil::CALCULATED
1406  || cellTypes_[nbr] == cellCellStencil::INTERPOLATED
1407  );
1408 
1409  if (ownLive && !nbrLive)
1410  {
1411  seedFaces.append(facei);
1412  const label globalOwn = globalCells.toGlobal(own);
1413  seedData.append(topoDistanceData<int>(globalOwn, 0));
1414  }
1415  else if (!ownLive && nbrLive)
1416  {
1417  seedFaces.append(facei);
1418  const label globalNbr = globalCells.toGlobal(nbr);
1419  seedData.append(topoDistanceData<int>(globalNbr, 0));
1420  }
1421  }
1422 
1423  forAll(nbrTypes, bFacei)
1424  {
1425  const label facei = bFacei+mesh_.nInternalFaces();
1426  const label own = owner[facei];
1427  const bool ownLive =
1428  (
1429  cellTypes_[own] == cellCellStencil::CALCULATED
1430  || cellTypes_[own] == cellCellStencil::INTERPOLATED
1431  );
1432  const bool nbrLive =
1433  (
1434  nbrTypes[bFacei] == cellCellStencil::CALCULATED
1435  || nbrTypes[bFacei] == cellCellStencil::INTERPOLATED
1436  );
1437 
1438  if (ownLive && !nbrLive)
1439  {
1440  seedFaces.append(facei);
1441  const label globalOwn = globalCells.toGlobal(own);
1442  seedData.append(topoDistanceData<int>(globalOwn, 0));
1443  }
1444  }
1445 
1446  // Propagate information inwards
1447  FaceCellWave<topoDistanceData<int>> distanceCalc
1448  (
1449  mesh_,
1450  seedFaces,
1451  seedData,
1452  faceData,
1453  cellData,
1454  mesh_.globalData().nTotalCells()+1,
1455  dummyTrackData
1456  );
1457  }
1458 
1459 
1460  // Add the new donor ones
1461  forAll(cellData, celli)
1462  {
1463  if
1464  (
1465  cellTypes_[celli] == cellCellStencil::HOLE
1466  && isDonor[celli]
1467  && cellData[celli].valid(dummyTrackData)
1468  )
1469  {
1470  cellStencil_[celli].setSize(1);
1471  cellStencil_[celli] = cellData[celli].data(); // global cellID
1472  cellInterpolationWeights_[celli].setSize(1);
1473  cellInterpolationWeights_[celli] = 1.0;
1474  cellTypes_[celli] = SPECIAL; // handled later on
1475  interpolationCells_.append(celli);
1476  cellInterpolationWeight_[celli] = 1.0;
1477  nSpecialFar++;
1478  }
1479  }
1480 
1481  if (debug & 2)
1482  {
1483  reduce(nSpecialNear, sumOp<label>());
1484  reduce(nSpecialFar, sumOp<label>());
1485 
1486  Pout<< "Detected " << nHoleDonors
1487  << " hole cells that are used as donors of which" << nl
1488  << " next to live cells : " << nSpecialNear << nl
1489  << " other : " << nSpecialFar << nl
1490  << endl;
1491  }
1492 
1493 
1494  // Re-do the mapDistribute
1495  List<Map<label>> compactMap;
1496  cellInterpolationMap_.reset
1497  (
1498  new mapDistribute
1499  (
1500  globalCells,
1501  cellStencil_,
1502  compactMap
1503  )
1504  );
1505  }
1506 }
1507 
1508 
1510 (
1511  const globalIndex& globalCells,
1512  const bool allowHoleDonors
1513 )
1514 {
1515  // Send cell centre back to donor
1516  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1517  // The complication is that multiple acceptors need the same donor
1518  // (but with different weights obviously)
1519  // So we do multi-pass:
1520  // - send over cc of acceptor for which we want stencil.
1521  // Consistently choose the acceptor with smallest magSqr in case of
1522  // multiple acceptors for the containing cell/donor.
1523  // - find the cell-cells and weights for the donor
1524  // - send back together with the acceptor cc
1525  // - use the acceptor cc to see if it was 'me' that sent it. If so
1526  // mark me as complete so it doesn't get involved in the next loop.
1527  // - loop until all resolved.
1528 
1529  // Special value for unused points
1530  const vector greatPoint(GREAT, GREAT, GREAT);
1531 
1532  boolList isValidDonor(mesh_.nCells(), true);
1533  if (!allowHoleDonors)
1534  {
1535  forAll(cellTypes_, celli)
1536  {
1537  if (cellTypes_[celli] == HOLE)
1538  {
1539  isValidDonor[celli] = false;
1540  }
1541  }
1542  }
1543 
1544  // Has acceptor been handled already?
1545  bitSet doneAcceptor(interpolationCells_.size());
1546 
1547  while (true)
1548  {
1549  pointField samples(cellInterpolationMap().constructSize(), greatPoint);
1550 
1551  // Fill remote slots (override old content). We'll find out later
1552  // on which one has won and mark this one in doneAcceptor.
1553  label nSamples = 0;
1554  forAll(interpolationCells_, i)
1555  {
1556  if (!doneAcceptor[i])
1557  {
1558  label cellI = interpolationCells_[i];
1559  const point& cc = mesh_.cellCentres()[cellI];
1560  const labelList& slots = cellStencil_[cellI];
1561 
1562  //- Note: empty slots can happen for interpolated cells with
1563  // bad/insufficient donors. These are handled later on.
1564  if (slots.size() > 1)
1565  {
1566  FatalErrorInFunction<< "Problem:" << slots
1567  << abort(FatalError);
1568  }
1569  else if (slots.size())
1570  {
1571  forAll(slots, sloti)
1572  {
1573  const label elemi = slots[sloti];
1574  //Pout<< " acceptor:" << cellI
1575  // << " at:" << mesh_.cellCentres()[cellI]
1576  // << " global:" << globalCells.toGlobal(cellI)
1577  // << " found in donor:" << elemi << endl;
1578  minMagSqrEqOp<point>()(samples[elemi], cc);
1579  }
1580  nSamples++;
1581  }
1582  }
1583  }
1584 
1585 
1586  if (!returnReduceOr(nSamples))
1587  {
1588  break;
1589  }
1590 
1591  // Send back to donor. Make sure valid point takes priority
1592  mapDistributeBase::distribute<point, minMagSqrEqOp<point>, flipOp>
1593  (
1596  mesh_.nCells(),
1597  cellInterpolationMap().constructMap(),
1598  false,
1599  cellInterpolationMap().subMap(),
1600  false,
1601  samples,
1602  greatPoint, // nullValue
1603  minMagSqrEqOp<point>(),
1604  flipOp(), // NegateOp
1606  cellInterpolationMap().comm()
1607  );
1608 
1609  // All the donor cells will now have a valid cell centre. Construct a
1610  // stencil for these.
1611 
1612  DynamicList<label> donorCells(mesh_.nCells());
1613  forAll(samples, cellI)
1614  {
1615  if (samples[cellI] != greatPoint)
1616  {
1617  donorCells.append(cellI);
1618  }
1619  }
1620 
1621 
1622  // Get neighbours (global cell and centre) of donorCells.
1623  labelListList donorCellCells(mesh_.nCells());
1624  pointListList donorCellCentres(mesh_.nCells());
1625  globalCellCells
1626  (
1627  globalCells,
1628  mesh_,
1629  isValidDonor,
1630  donorCells,
1631  donorCellCells,
1632  donorCellCentres
1633  );
1634 
1635  // Determine the weights.
1636  scalarListList donorWeights(mesh_.nCells());
1637  forAll(donorCells, i)
1638  {
1639  label cellI = donorCells[i];
1640  const pointList& donorCentres = donorCellCentres[cellI];
1641  stencilWeights
1642  (
1643  samples[cellI],
1644  donorCentres,
1645  donorWeights[cellI]
1646  );
1647  }
1648 
1649  // Transfer the information back to the acceptor:
1650  // - donorCellCells : stencil (with first element the original donor)
1651  // - donorWeights : weights for donorCellCells
1652  cellInterpolationMap().distribute(donorCellCells);
1653  cellInterpolationMap().distribute(donorWeights);
1654  cellInterpolationMap().distribute(samples);
1655 
1656  // Check which acceptor has won and transfer
1657  forAll(interpolationCells_, i)
1658  {
1659  if (!doneAcceptor[i])
1660  {
1661  label cellI = interpolationCells_[i];
1662  const labelList& slots = cellStencil_[cellI];
1663 
1664  if (slots.size() > 1)
1665  {
1666  FatalErrorInFunction << "Problem:" << slots
1667  << abort(FatalError);
1668  }
1669  else if (slots.size() == 1)
1670  {
1671  const label sloti = slots[0];
1672 
1673  // Important: check if the stencil is actually for this cell
1674  if (samples[sloti] == mesh_.cellCentres()[cellI])
1675  {
1676  cellStencil_[cellI].transfer(donorCellCells[sloti]);
1677  cellInterpolationWeights_[cellI].transfer
1678  (
1679  donorWeights[sloti]
1680  );
1681  // Mark cell as being done so it does not get sent over
1682  // again.
1683  doneAcceptor.set(i);
1684  }
1685  }
1686  }
1687  }
1688  }
1689 
1690  // Re-do the mapDistribute
1691  List<Map<label>> compactMap;
1692  cellInterpolationMap_.reset
1693  (
1694  new mapDistribute
1695  (
1696  globalCells,
1697  cellStencil_,
1698  compactMap
1699  )
1700  );
1701 }
1702 
1704 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1705 
1706 Foam::cellCellStencils::inverseDistance::inverseDistance
1707 (
1708  const fvMesh& mesh,
1709  const dictionary& dict,
1710  const bool doUpdate
1711 )
1712 :
1714  dict_(dict),
1715  allowHoleDonors_(dict.getOrDefault("allowHoleDonors", false)),
1716  allowInterpolatedDonors_
1717  (
1718  dict.getOrDefault("allowInterpolatedDonors", true)
1719  ),
1720  smallVec_(Zero),
1721  cellTypes_(labelList(mesh.nCells(), CALCULATED)),
1722  interpolationCells_(0),
1723  cellInterpolationMap_(),
1724  cellStencil_(0),
1725  cellInterpolationWeights_(0),
1726  cellInterpolationWeight_
1727  (
1728  IOobject
1729  (
1730  "cellInterpolationWeight",
1731  mesh_.facesInstance(),
1732  mesh_,
1733  IOobject::NO_READ,
1734  IOobject::NO_WRITE,
1735  IOobject::NO_REGISTER
1736  ),
1737  mesh_,
1740  )
1741 {
1742  // Add motion-solver fields to non-interpolated fields
1744 
1745  // Read zoneID
1746  this->zoneID();
1747 
1748  // Read old-time cellTypes
1749  IOobject io
1750  (
1751  "cellTypes",
1752  mesh_.time().timeName(),
1753  mesh_,
1757  );
1758  if (io.typeHeaderOk<volScalarField>(true))
1759  {
1760  if (debug)
1761  {
1762  Pout<< "Reading cellTypes from time " << mesh_.time().timeName()
1763  << endl;
1764  }
1765 
1766  const volScalarField volCellTypes(io, mesh_);
1767  forAll(volCellTypes, celli)
1768  {
1769  // Round to integer
1770  cellTypes_[celli] = volCellTypes[celli];
1771  }
1772  }
1773 
1774  if (doUpdate)
1775  {
1776  update();
1777  }
1778 }
1780 
1781 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1782 
1784 {}
1786 
1787 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1788 
1790 {
1791  scalar layerRelax(dict_.getOrDefault("layerRelax", 1.0));
1792 
1793  scalar tol = dict_.getOrDefault("tolerance", 1e-10);
1794  smallVec_ = mesh_.bounds().span()*tol;
1795 
1796  const labelIOList& zoneID = this->zoneID();
1797 
1798  label nZones = gMax(zoneID)+1;
1799  labelList nCellsPerZone(nZones, Zero);
1800  forAll(zoneID, cellI)
1801  {
1802  nCellsPerZone[zoneID[cellI]]++;
1803  }
1804  Pstream::listCombineReduce(nCellsPerZone, plusEqOp<label>());
1805 
1806  const boundBox& allBb = mesh_.bounds();
1807 
1808  PtrList<fvMeshSubset> meshParts(nZones);
1810 
1811  // Determine zone meshes and bounding boxes
1812  {
1813  // Per processor, per zone the bounding box
1815  procBb[Pstream::myProcNo()].setSize(nZones);
1816 
1817  forAll(meshParts, zonei)
1818  {
1819  meshParts.set
1820  (
1821  zonei,
1822  new fvMeshSubset(mesh_, zonei, zoneID)
1823  );
1824  const fvMesh& subMesh = meshParts[zonei].subMesh();
1825 
1826  // Trigger early evaluation of mesh dimension (in case there are
1827  // zero cells in mesh)
1828  (void)subMesh.nGeometricD();
1829 
1830  if (subMesh.nPoints())
1831  {
1832  procBb[Pstream::myProcNo()][zonei] =
1833  treeBoundBox(subMesh.points());
1834  procBb[Pstream::myProcNo()][zonei].inflate(1e-6);
1835  }
1836  else
1837  {
1838  // No part of zone on this processor. Make up bb.
1839  procBb[Pstream::myProcNo()][zonei] = treeBoundBox
1840  (
1841  allBb.min() - 2*allBb.span(),
1842  allBb.min() - allBb.span()
1843  );
1844  procBb[Pstream::myProcNo()][zonei].inflate(1e-6);
1845  }
1846  }
1847 
1848  Pstream::allGatherList(procBb);
1849 
1850  // Move local bounding boxes to per-mesh indexing
1851  forAll(meshBb, zoneI)
1852  {
1853  treeBoundBoxList& bbs = meshBb[zoneI];
1854  bbs.setSize(Pstream::nProcs());
1855  forAll(procBb, procI)
1856  {
1857  bbs[procI] = procBb[procI][zoneI];
1858  }
1859  }
1860  }
1861 
1862 
1863  // Determine patch bounding boxes. These are either global and provided
1864  // by the user or processor-local as a copy of the mesh bounding box.
1865 
1866  List<treeBoundBoxList> patchBb(nZones);
1867  List<labelVector> patchDivisions(nZones);
1868  PtrList<PackedList<2>> patchParts(nZones);
1869  labelList allPatchTypes(mesh_.nCells(), OTHER);
1870 
1871  {
1872  treeBoundBox globalPatchBb;
1873  if (dict_.readIfPresent("searchBox", globalPatchBb))
1874  {
1875  // All processors, all zones have the same bounding box
1876  patchBb = treeBoundBoxList(Pstream::nProcs(), globalPatchBb);
1877  }
1878  else
1879  {
1880  // Use the meshBb (differing per zone, per processor)
1881  patchBb = meshBb;
1882  }
1883  }
1884 
1885  {
1886  labelVector globalDivs;
1887  if (dict_.readIfPresent("searchBoxDivisions", globalDivs))
1888  {
1889  patchDivisions = globalDivs;
1890  }
1891  else
1892  {
1893  const labelVector& dim = mesh_.geometricD();
1894  label nDivs = -1;
1895  if (mesh_.nGeometricD() == 1)
1896  {
1897  nDivs = mesh_.nCells();
1898  }
1899  else if (mesh_.nGeometricD() == 2)
1900  {
1901  nDivs = label(Foam::sqrt(scalar(mesh_.nCells())));
1902  }
1903  else
1904  {
1905  nDivs = label(Foam::cbrt(scalar(mesh_.nCells())));
1906  }
1907 
1908  labelVector v(nDivs, nDivs, nDivs);
1909  forAll(dim, i)
1910  {
1911  if (dim[i] == -1)
1912  {
1913  v[i] = 1;
1914  }
1915  }
1916  patchDivisions = v;
1917  }
1918  }
1919 
1920  forAll(patchParts, zoneI)
1921  {
1922  patchParts.set
1923  (
1924  zoneI,
1925  new PackedList<2>
1926  (
1927  patchDivisions[zoneI][0]
1928  *patchDivisions[zoneI][1]
1929  *patchDivisions[zoneI][2]
1930  )
1931  );
1932  markBoundaries
1933  (
1934  meshParts[zoneI].subMesh(),
1935  smallVec_,
1936 
1937  patchBb[zoneI][Pstream::myProcNo()],
1938  patchDivisions[zoneI],
1939  patchParts[zoneI],
1940 
1941  meshParts[zoneI].cellMap(),
1942  allPatchTypes
1943  );
1944  }
1945 
1946 
1947  // Print a bit
1948  {
1949  Info<< type() << " : detected " << nZones
1950  << " mesh regions" << endl;
1951  Info<< incrIndent;
1952  forAll(nCellsPerZone, zoneI)
1953  {
1954  Info<< indent<< "zone:" << zoneI
1955  << " nCells:" << nCellsPerZone[zoneI]
1956  << " voxels:" << patchDivisions[zoneI]
1957  << " bb:" << patchBb[zoneI][Pstream::myProcNo()]
1958  << endl;
1959  }
1960  Info<< decrIndent;
1961  }
1962 
1963 
1964  // Current best guess for cells. Includes best stencil. Weights should
1965  // add up to volume.
1966  labelList allCellTypes(mesh_.nCells(), CALCULATED);
1967  labelListList allStencil(mesh_.nCells());
1968  // zoneID of donor
1969  labelList allDonorID(mesh_.nCells(), -1);
1970 
1971  const globalIndex globalCells(mesh_.nCells());
1972 
1974 
1975  // Mark holes (in allCellTypes)
1976  for (label srcI = 0; srcI < meshParts.size()-1; srcI++)
1977  {
1978  for (label tgtI = srcI+1; tgtI < meshParts.size(); tgtI++)
1979  {
1980  markPatchesAsHoles
1981  (
1982  pBufs,
1983 
1984  meshParts,
1985 
1986  patchBb,
1987  patchDivisions,
1988  patchParts,
1989 
1990  srcI,
1991  tgtI,
1992  allCellTypes
1993  );
1994  markPatchesAsHoles
1995  (
1996  pBufs,
1997 
1998  meshParts,
1999 
2000  patchBb,
2001  patchDivisions,
2002  patchParts,
2003 
2004  tgtI,
2005  srcI,
2006  allCellTypes
2007  );
2008  }
2009  }
2010 
2011  // Find donors (which are not holes) in allStencil, allDonorID
2012  for (label srcI = 0; srcI < meshParts.size()-1; srcI++)
2013  {
2014  for (label tgtI = srcI+1; tgtI < meshParts.size(); tgtI++)
2015  {
2016  markDonors
2017  (
2018  globalCells,
2019  pBufs,
2020  meshParts,
2021  meshBb,
2022  allCellTypes,
2023 
2024  tgtI,
2025  srcI,
2026  allStencil,
2027  allDonorID
2028  );
2029  markDonors
2030  (
2031  globalCells,
2032  pBufs,
2033  meshParts,
2034  meshBb,
2035  allCellTypes,
2036 
2037  srcI,
2038  tgtI,
2039  allStencil,
2040  allDonorID
2041  );
2042  }
2043  }
2044 
2045  if ((debug & 2) && mesh_.time().writeTime())
2046  {
2047  tmp<volScalarField> tfld
2048  (
2049  createField(mesh_, "allCellTypes", allCellTypes)
2050  );
2051  tfld().write();
2052 
2053  tmp<volScalarField> tallDonorID
2054  (
2055  createField(mesh_, "allDonorID", allDonorID)
2056  );
2057  tallDonorID().write();
2058  }
2059 
2060  // Use the patch types and weights to decide what to do
2061  forAll(allPatchTypes, cellI)
2062  {
2063  if (allCellTypes[cellI] != HOLE)
2064  {
2065  switch (allPatchTypes[cellI])
2066  {
2067  case OVERSET:
2068  {
2069  // Require interpolation. See if possible.
2070  if (allStencil[cellI].size())
2071  {
2072  allCellTypes[cellI] = INTERPOLATED;
2073  }
2074  else
2075  {
2076  allCellTypes[cellI] = HOLE;
2077  }
2078  }
2079  }
2080  }
2081  }
2082 
2083  if ((debug & 2) && mesh_.time().writeTime())
2084  {
2085  tmp<volScalarField> tfld
2086  (
2087  createField(mesh_, "allCellTypes_patch", allCellTypes)
2088  );
2089  tfld().write();
2090 
2091  tmp<volScalarField> tfldOld
2092  (
2093  createField(mesh_, "allCellTypes_old", cellTypes_)
2094  );
2095  tfldOld().write();
2096  }
2097 
2098  // Mark unreachable bits
2099  findHoles(globalCells, mesh_, zoneID, allStencil, allCellTypes);
2100 
2101 
2102  if ((debug & 2) && mesh_.time().writeTime())
2103  {
2104  tmp<volScalarField> tfld
2105  (
2106  createField(mesh_, "allCellTypes_hole", allCellTypes)
2107  );
2108  tfld().write();
2109  }
2110  if ((debug & 2) && mesh_.time().writeTime())
2111  {
2112  labelList stencilSize(mesh_.nCells());
2113  forAll(allStencil, celli)
2114  {
2115  stencilSize[celli] = allStencil[celli].size();
2116  }
2117  tmp<volScalarField> tfld
2118  (
2119  createField(mesh_, "allStencil_hole", stencilSize)
2120  );
2121  tfld().write();
2122  }
2123 
2124  // Update allStencil with new fill HOLES
2125  forAll(allCellTypes, celli)
2126  {
2127  if (allCellTypes[celli] == HOLE && allStencil[celli].size())
2128  {
2129  allStencil[celli].clear();
2130  }
2131  }
2132 
2133 // // Check previous iteration cellTypes_ for any hole->calculated changes
2134 // // If so set the cell either to interpolated (if there are donors) or
2135 // // holes (if there are no donors). Note that any interpolated cell might
2136 // // still be overwritten by the flood filling
2137 // {
2138 // label nCalculated = 0;
2139 //
2140 // forAll(cellTypes_, celli)
2141 // {
2142 // if (allCellTypes[celli] == CALCULATED && cellTypes_[celli] == HOLE)
2143 // {
2144 // if (allStencil[celli].size() == 0)
2145 // {
2146 // // Reset to hole
2147 // allCellTypes[celli] = HOLE;
2148 // allStencil[celli].clear();
2149 // }
2150 // else
2151 // {
2152 // allCellTypes[celli] = INTERPOLATED;
2153 // nCalculated++;
2154 // }
2155 // }
2156 // }
2157 //
2158 // if (debug)
2159 // {
2160 // Pout<< "Detected " << nCalculated << " cells changing from hole"
2161 // << " to calculated. Changed to interpolated"
2162 // << endl;
2163 // }
2164 // }
2165 
2166  if ((debug & 2) && mesh_.time().writeTime())
2167  {
2168  tmp<volScalarField> tfld
2169  (
2170  createField(mesh_, "allCellTypes_pre_front", allCellTypes)
2171  );
2172  tfld().write();
2173  }
2174  // Add buffer interpolation layer(s) around holes
2175  scalarField allWeight(mesh_.nCells(), Zero);
2176 
2177  labelListList compactStencil(allStencil);
2178  List<Map<label>> compactStencilMap;
2179  mapDistribute map(globalCells, compactStencil, compactStencilMap);
2180 
2181  scalarList compactCellVol(mesh_.V());
2182  map.distribute(compactCellVol);
2183 
2184 
2185  label useLayer = dict_.getOrDefault("useLayer", -1);
2186  const dictionary& fvSchemes = mesh_.schemesDict();
2187  if (fvSchemes.found("oversetInterpolation"))
2188  {
2189  const dictionary& intDict = fvSchemes.subDict
2190  (
2191  "oversetInterpolation"
2192  );
2193  useLayer = intDict.getOrDefault("useLayer", -1);
2194  }
2195 
2196  walkFront
2197  (
2198  globalCells,
2199  layerRelax,
2200  allStencil,
2201  allCellTypes,
2202  allWeight,
2203  compactCellVol,
2204  compactStencil,
2205  zoneID,
2206  dict_.getOrDefault("holeLayers", 1),
2207  useLayer
2208  );
2209 
2210  if ((debug & 2) && mesh_.time().writeTime())
2211  {
2212  tmp<volScalarField> tfld
2213  (
2214  createField(mesh_, "allCellTypes_front", allCellTypes)
2215  );
2216  tfld().write();
2217  }
2218 
2219  // Check previous iteration cellTypes_ for any hole->calculated changes
2220  // If so set the cell either to interpolated (if there are donors) or
2221  // holes (if there are no donors). Note that any interpolated cell might
2222  // still be overwritten by the flood filling
2223  {
2224  label nCalculated = 0;
2225 
2226  forAll(cellTypes_, celli)
2227  {
2228  if (allCellTypes[celli] == CALCULATED && cellTypes_[celli] == HOLE)
2229  {
2230  if (allStencil[celli].size() == 0)
2231  {
2232  // Reset to hole
2233  allCellTypes[celli] = HOLE;
2234  allStencil[celli].clear();
2235  }
2236  else
2237  {
2238  allCellTypes[celli] = INTERPOLATED;
2239  nCalculated++;
2240  }
2241  }
2242  }
2243 
2244  if (debug)
2245  {
2246  Pout<< "Detected " << nCalculated << " cells changing from hole"
2247  << " to calculated. Changed to interpolated"
2248  << endl;
2249  }
2250  }
2251 
2252 
2253 
2254 
2255  // Convert cell-cell addressing to stencil in compact notation
2256  // Clean any potential INTERPOLATED with HOLE donors.
2257  // This is not eliminated in the front walk as the HOLE cells might
2258  // leak when inset region is overlapping backgroung mesh
2259  cellTypes_.transfer(allCellTypes);
2260  cellStencil_.setSize(mesh_.nCells());
2261  cellInterpolationWeights_.setSize(mesh_.nCells());
2262  DynamicList<label> interpolationCells;
2263 
2264 /*
2265  forAll(cellTypes_, cellI)
2266  {
2267  if (cellTypes_[cellI] == INTERPOLATED)
2268  {
2269  if (cellTypes_[allStencil[cellI][0]] != HOLE)
2270  {
2271  cellStencil_[cellI].transfer(allStencil[cellI]);
2272  cellInterpolationWeights_[cellI].setSize(1);
2273  cellInterpolationWeights_[cellI][0] = 1.0; interpolationCells.append(cellI);
2274  }
2275  else
2276  {
2277  cellTypes_[cellI] = POROUS;//CALCULATED;
2278  cellStencil_[cellI].clear();
2279  cellInterpolationWeights_[cellI].clear();
2280  }
2281  }
2282  else
2283  {
2284  cellStencil_[cellI].clear();
2285  cellInterpolationWeights_[cellI].clear();
2286  }
2287  }
2288 */
2289 
2290  //labelListList compactStencil(allStencil);
2291  //List<Map<label>> compactStencilMap;
2292  //mapDistribute map(globalCells, compactStencil, compactStencilMap);
2293 
2294  labelList compactCellTypes(cellTypes_);
2295  map.distribute(compactCellTypes);
2296 
2297  label nInterpolatedDonors = 0;
2298  forAll(cellTypes_, celli)
2299  {
2300  if (cellTypes_[celli] == INTERPOLATED)
2301  {
2302  const labelList& slots = compactStencil[celli];
2303  if (slots.size())
2304  {
2305  if
2306  (
2307  compactCellTypes[slots[0]] == HOLE
2308  ||
2309  (
2310  !allowInterpolatedDonors_
2311  && compactCellTypes[slots[0]] == INTERPOLATED
2312  )
2313  )
2314  {
2315  cellTypes_[celli] = POROUS;
2316  cellStencil_[celli].clear();
2317  cellInterpolationWeights_[celli].clear();
2318  nInterpolatedDonors++;
2319  }
2320  else
2321  {
2322  cellStencil_[celli].transfer(allStencil[celli]);
2323  cellInterpolationWeights_[celli].setSize(1);
2324  cellInterpolationWeights_[celli][0] = 1.0;
2325  interpolationCells.append(celli);
2326  }
2327  }
2328  }
2329  else
2330  {
2331  cellStencil_[celli].clear();
2332  cellInterpolationWeights_[celli].clear();
2333  }
2334  }
2335 
2336  reduce(nInterpolatedDonors, sumOp<label>());
2337 
2339  << "interpolate Donors : "
2340  << nInterpolatedDonors << endl;
2341 
2342  // Reset map with new cellStencil
2343  List<Map<label>> compactMap;
2344  cellInterpolationMap_.reset
2345  (
2346  new mapDistribute(globalCells, cellStencil_, compactMap)
2347  );
2348 
2349  interpolationCells_.transfer(interpolationCells);
2350 
2351  cellInterpolationWeight_.transfer(allWeight);
2353  <
2356  >(cellInterpolationWeight_.boundaryFieldRef(), false);
2357 
2358 
2359  if ((debug & 2) && mesh_.time().writeTime())
2360  {
2361  // Dump mesh
2362  mesh_.time().write();
2363 
2364  // Dump stencil
2365  mkDir(mesh_.time().timePath());
2366  OBJstream str(mesh_.time().timePath()/"injectionStencil.obj");
2367  Pout<< type() << " : dumping injectionStencil to "
2368  << str.name() << endl;
2369  pointField cc(mesh_.cellCentres());
2370  cellInterpolationMap().distribute(cc);
2371 
2372  forAll(cellStencil_, celli)
2373  {
2374  const labelList& slots = cellStencil_[celli];
2375  if (slots.size())
2376  {
2377  const point& accCc = mesh_.cellCentres()[celli];
2378  forAll(slots, i)
2379  {
2380  const point& donorCc = cc[slots[i]];
2381  str.writeLine(accCc, 0.1*accCc+0.9*donorCc);
2382  }
2383  }
2384  }
2385  }
2386 
2387 
2388  // Extend stencil to get inverse distance weighted neighbours
2389  createStencil(globalCells, allowHoleDonors_);
2390 
2391  // Optional: convert hole cells next to non-hole cells into
2392  // interpolate-from-neighbours (of cell type SPECIAL)
2393  if (allowHoleDonors_)
2394  {
2395  holeExtrapolationStencil(globalCells);
2396  }
2397 
2398  if ((debug & 2) && mesh_.time().writeTime())
2399  {
2400  // Dump weight
2401  cellInterpolationWeight_.instance() = mesh_.time().timeName();
2402  cellInterpolationWeight_.write();
2403 
2404  // Dump max weight
2405  {
2406  scalarField maxMagWeight(mesh_.nCells(), Zero);
2407  forAll(cellStencil_, celli)
2408  {
2409  const scalarList& wghts = cellInterpolationWeights_[celli];
2410  forAll(wghts, i)
2411  {
2412  if (mag(wghts[i]) > mag(maxMagWeight[celli]))
2413  {
2414  maxMagWeight[celli] = wghts[i];
2415  }
2416  }
2417  if (mag(maxMagWeight[celli]) > 1)
2418  {
2419  const pointField& cc = mesh_.cellCentres();
2420  Pout<< "cell:" << celli
2421  << " at:" << cc[celli]
2422  << " zone:" << zoneID[celli]
2423  << " donors:" << cellStencil_[celli]
2424  << " weights:" << wghts
2425  << " coords:"
2426  << UIndirectList<point>(cc, cellStencil_[celli])
2427  << " donorZone:"
2428  << UIndirectList<label>(zoneID, cellStencil_[celli])
2429  << endl;
2430  }
2431  }
2432  tmp<volScalarField> tfld
2433  (
2434  createField(mesh_, "maxMagWeight", maxMagWeight)
2435  );
2437  <
2440  >(tfld.ref().boundaryFieldRef(), false);
2441  tfld().write();
2442  }
2443 
2444  // Dump cell types
2445  {
2446  tmp<volScalarField> tfld
2447  (
2448  createField(mesh_, "cellTypes", cellTypes_)
2449  );
2450  //tfld.ref().correctBoundaryConditions();
2452  <
2455  >(tfld.ref().boundaryFieldRef(), false);
2456  tfld().write();
2457  }
2458 
2459 
2460  // Dump stencil, one per zone
2461  mkDir(mesh_.time().timePath());
2462  pointField cc(mesh_.cellCentres());
2463  cellInterpolationMap().distribute(cc);
2464  forAll(meshParts, zonei)
2465  {
2466  OBJstream str
2467  (
2468  mesh_.time().timePath()
2469  + "/stencil_" + name(zonei) + ".obj"
2470  );
2471  Pout<< type() << " : dumping to " << str.name() << endl;
2472 
2473  const labelList& subMeshCellMap = meshParts[zonei].cellMap();
2474 
2475  forAll(subMeshCellMap, subcelli)
2476  {
2477  const label celli = subMeshCellMap[subcelli];
2478  const labelList& slots = cellStencil_[celli];
2479  const point& accCc = mesh_.cellCentres()[celli];
2480  forAll(slots, i)
2481  {
2482  const point& donorCc = cc[slots[i]];
2483  str.writeLine(accCc, 0.1*accCc+0.9*donorCc);
2484  }
2485  }
2486  }
2487  }
2488 
2489  // Print some stats
2490  {
2491  labelList nCells(count(3, cellTypes_));
2492 
2493  label nLocal = 0;
2494  label nMixed = 0;
2495  label nRemote = 0;
2496  forAll(interpolationCells_, i)
2497  {
2498  label celli = interpolationCells_[i];
2499  const labelList& slots = cellStencil_[celli];
2500 
2501  bool hasLocal = false;
2502  bool hasRemote = false;
2503 
2504  forAll(slots, sloti)
2505  {
2506  if (slots[sloti] >= mesh_.nCells())
2507  {
2508  hasRemote = true;
2509  }
2510  else
2511  {
2512  hasLocal = true;
2513  }
2514  }
2515 
2516  if (hasRemote)
2517  {
2518  if (!hasLocal)
2519  {
2520  nRemote++;
2521  }
2522  else
2523  {
2524  nMixed++;
2525  }
2526  }
2527  else if (hasLocal)
2528  {
2529  nLocal++;
2530  }
2531  }
2532 
2533  Info<< "Overset analysis : nCells : "
2534  << returnReduce(cellTypes_.size(), sumOp<label>()) << nl
2535  << incrIndent
2536  << indent << "calculated : " << nCells[CALCULATED] << nl
2537  << indent << "interpolated : " << nCells[INTERPOLATED]
2538  << " (from local:" << returnReduce(nLocal, sumOp<label>())
2539  << " mixed local/remote:" << returnReduce(nMixed, sumOp<label>())
2540  << " remote:" << returnReduce(nRemote, sumOp<label>()) << ")" << nl
2541  << indent << "hole : " << nCells[HOLE] << nl
2542  << decrIndent << endl;
2543  }
2544 
2545  // Tbd: detect if anything changed. Most likely it did!
2546  return true;
2547 }
2548 
2549 
2550 // ************************************************************************* //
static point position(const boundBox &bb, const labelVector &nDivs, const label boxI)
Convert index back into coordinate.
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:136
Inverse-distance-weighted interpolation stencil.
const polyBoundaryMesh & pbm
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
uint8_t direction
Definition: direction.H:46
const word zeroGradientType
A zeroGradient patch field type.
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:493
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
static labelVector index3(const labelVector &nDivs, const label)
Convert single index into ijk.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition: List.C:326
label nPoints() const noexcept
Number of mesh points.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:517
void findHoles(const globalIndex &globalCells, const fvMesh &mesh, const labelList &zoneID, const labelListList &stencil, labelList &cellTypes) const
Do flood filling to detect unreachable (from patches) sections.
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute List data using default commsType, default flip/negate operator.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
static rangeType allProcs(const label communicator=worldComm)
Range of process indices for all processes.
Definition: UPstream.H:1176
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
wordList patchTypes(nPatches)
static void calculate(const polyMesh &src, const polyMesh &tgt, labelList &srcToTgtAddr)
Calculate addressing.
Definition: waveMethod.C:38
bool betterDonor(const label destMesh, const label currentDonorMesh, const label newDonorMesh) const
If multiple donors meshes: decide which is best.
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: IndirectList.H:65
scalarField samples(nIntervals, Zero)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static void markBoundaries(const fvMesh &mesh, const vector &smallVec, const boundBox &bb, const labelVector &nDivs, PackedList< 2 > &patchTypes, const labelList &cellMap, labelList &patchCellTypes)
Mark voxels of patchTypes with type of patch face.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
defineTypeNameAndDebug(cellVolumeWeight, 0)
static void fill(PackedList< 2 > &elems, const boundBox &bb, const labelVector &nDivs, const boundBox &subBb, const unsigned int val)
Fill all elements overlapping subBb with value val.
Template invariant parts for fvPatchField.
Definition: fvPatchField.H:77
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
labelList findIndices(const ListType &input, typename ListType::const_reference val, label start=0)
Linear search to find all occurrences of given element.
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:1229
List< treeBoundBox > treeBoundBoxList
A List of treeBoundBox.
Definition: treeBoundBox.H:83
constexpr label labelMin
Definition: label.H:54
label k
Boltzmann constant.
Ignore writing from objectRegistry::writeObject()
const point & min() const noexcept
Minimum describing the bounding box.
Definition: boundBoxI.H:162
const dimensionSet dimless
Dimensionless.
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1074
static const List< T > & null()
Return a null List.
Definition: ListI.H:130
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:441
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:876
Macros for easy insertion into run-time selection tables.
void inplaceToGlobal(const label proci, labelUList &labels) const
From local to global index on proci (inplace)
Definition: globalIndexI.H:351
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Definition: dictionaryI.H:104
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
void markDonors(const globalIndex &globalCells, PstreamBuffers &pBufs, const PtrList< fvMeshSubset > &meshParts, const List< treeBoundBoxList > &meshBb, const labelList &allCellTypes, const label srcI, const label tgtI, labelListList &allStencil, labelList &allDonor) const
Determine donors for all tgt cells.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
static void allGatherList(List< T > &values, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Gather data, but keep individual values separate. Uses linear/tree communication. ...
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
const point & max() const noexcept
Maximum describing the bounding box.
Definition: boundBoxI.H:168
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of &#39;true&#39; entries.
Definition: BitOps.H:73
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:61
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run. ...
Definition: UPstream.H:1065
PDRpatchDef PATCH
Definition: PDRpatchDef.H:120
void suppressMotionFields()
Helper: populate nonInterpolatedFields_ with motion solver.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:614
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
Calculation of interpolation stencils.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Definition: labelLists.C:44
virtual void stencilWeights(const point &sample, const pointList &donorCcs, scalarList &weights) const
Calculate inverse distance weights for a single acceptor.
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:137
List< scalarList > scalarListList
List of scalarList.
Definition: scalarList.H:35
void markPatchesAsHoles(PstreamBuffers &pBufs, const PtrList< fvMeshSubset > &meshParts, const List< treeBoundBoxList > &patchBb, const List< labelVector > &patchDivisions, const PtrList< PackedList< 2 >> &patchParts, const label srcI, const label tgtI, labelList &allCellTypes) const
Mark all cells overlapping (a voxel covered by) a src patch.
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:107
const label nSamples(pdfDictionary.get< label >("nSamples"))
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
dimensionedScalar cbrt(const dimensionedScalar &ds)
Reading is optional [identical to LAZY_READ].
const dictionary & schemesDict() const
The entire dictionary or the optional "select" sub-dictionary.
Vector< scalar > vector
Definition: vector.H:57
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicList.H:584
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
errorManip< error > abort(error &err)
Definition: errorManip.H:139
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: PtrList.H:159
const fvMesh & mesh_
Reference to the mesh.
const labelList & cellTypes
Definition: setCellMask.H:27
virtual void createStencil(const globalIndex &, const bool allowHoleDonors)
Create stencil starting from the donor containing the acceptor allowHoleDonors : allow donors of type...
#define DebugInfo
Report an information message using Foam::Info.
static void correctBoundaryConditions(typename GeoField::Boundary &bfld, const bool typeOnly)
Correct boundary conditions of certain type (typeOnly = true)
virtual bool update()
Update stencils. Return false if nothing changed.
addToRunTimeSelectionTable(cellCellStencil, cellVolumeWeight, mesh)
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files...
Definition: OBJstream.H:55
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
int debug
Static debugging option.
Type gMax(const FieldField< Field, Type > &f)
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
static label index(const labelVector &nDivs, const labelVector &)
Convert ijk indices into single index.
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:511
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
Definition: fvMeshSubset.H:75
Buffers for inter-processor communications streams (UOPstream, UIPstream).
static bool constraintType(const word &patchType)
Return true if the given type is a constraint type.
Definition: fvPatch.C:74
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData)
Swap to obtain neighbour cell values for all boundary faces.
const labelIOList & zoneID() const
Helper: get reference to registered zoneID. Loads volScalarField.
Class containing processor-to-processor mapping information.
vector point
Point is a vector.
Definition: point.H:37
label toGlobal(const label proci, const label i) const
From local to global on proci.
Definition: globalIndexI.H:308
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Foam::fvBoundaryMesh.
Selector class for finite volume differencing schemes. fvMesh is derived from fvSchemes so that all f...
Definition: fvSchemes.H:51
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A List with indirect addressing. Like IndirectList but does not store addressing. ...
Definition: faMatrix.H:56
label nRegions() const
Return total number of regions.
Definition: regionSplit.H:328
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:192
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:90
"nonBlocking" : (MPI_Isend, MPI_Irecv)
messageStream Info
Information stream (stdout output on master, null elsewhere)
Vector< label > labelVector
Vector of labels.
Definition: labelVector.H:47
virtual void write(Ostream &os) const
Write.
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
List< point > pointList
List of point.
Definition: pointList.H:32
List< label > labelList
A List of labels.
Definition: List.H:62
List< pointList > pointListList
List of pointList.
Definition: pointList.H:35
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
static bool overlaps(const boundBox &bb, const labelVector &nDivs, const PackedList< 2 > &voxels, const treeBoundBox &subBb, const unsigned int val)
Is any voxel inside subBb set to val.
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:502
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
labelList cellIDs
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Do not request registration (bool: false)
void holeExtrapolationStencil(const globalIndex &globalCells)
Make holes next to live ones type SPECIAL and include in interpolation stencil.
List< treeBoundBox > meshBb(1, treeBoundBox(coarseMesh.points()).extend(rndGen, 1e-3))
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
static tmp< volScalarField > createField(const fvMesh &mesh, const scalar val)
Boundary condition for use on overset patches. To be run in combination with special dynamicFvMesh ty...
static void listCombineReduce(List< T > &values, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
After completion all processors have the same data.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127