Loading...
Searching...
No Matches
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-2025 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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"
43#include "topoDistanceData.H"
44#include "FaceCellWave.H"
45
46#include "OBJstream.H"
49// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
50
51namespace Foam
52{
53namespace cellCellStencils
54{
55 defineTypeNameAndDebug(inverseDistance, 0);
56 addToRunTimeSelectionTable(cellCellStencil, inverseDistance, mesh);
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 labelUList& 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 labelUList& 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);
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 }
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,
324 patchCellType::PATCH
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,
392 patchCellType::PATCH
393 )
394 )
395 {
396 allCellTypes[celli] = HOLE;
397 }
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
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 : "
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 for (const label celli : fvp.faceCells())
1078 {
1079 label regionI = cellRegion[celli];
1080
1081 if (cellTypes[celli] != HOLE && regionType[regionI] != 2)
1082 {
1083 regionType[regionI] = 2;
1084 }
1085 }
1086 }
1087 }
1088
1089 DebugInfo<< FUNCTION_NAME << " : Done local analysis" << endl;
1090
1091 // Now we've handled
1092 // - cells next to blocked cells
1093 // - coupled boundaries
1094 // Only thing to handle is the interpolation between regions
1095
1096
1097 labelListList compactStencil(stencil);
1098 List<Map<label>> compactMap;
1099 mapDistribute map(globalCells, compactStencil, compactMap);
1100
1101 DebugInfo<< FUNCTION_NAME << " : Converted stencil into compact form"
1102 << endl;
1103
1104
1105 while (true)
1106 {
1107 // Synchronise region status on processors
1108 // (could instead swap status through processor patches)
1109 Pstream::listReduce(regionType, maxOp<label>());
1110
1111 DebugInfo<< FUNCTION_NAME << " : Gathered region type" << endl;
1112
1113 // Communicate region status through interpolative cells
1114 labelList cellRegionType(labelUIndList(regionType, cellRegion));
1115 map.distribute(cellRegionType);
1116
1117 DebugInfo<< FUNCTION_NAME << " : Interpolated region type" << endl;
1118
1119
1120
1121 label nChanged = 0;
1122 forAll(pbm, patchI)
1123 {
1124 const fvPatch& fvp = pbm[patchI];
1125
1126 if (isA<oversetFvPatch>(fvp))
1127 {
1128 const labelUList& fc = fvp.faceCells();
1129 forAll(fc, i)
1130 {
1131 label cellI = fc[i];
1132 label regionI = cellRegion[cellI];
1133
1134 if (regionType[regionI] != 2)
1135 {
1136 const labelList& slots = compactStencil[cellI];
1137 forAll(slots, i)
1138 {
1139 label otherType = cellRegionType[slots[i]];
1140
1141 if (otherType == 2)
1142 {
1143 //Pout<< "Reachable through interpolation : "
1144 // << regionI << " at cell "
1145 // << mesh.cellCentres()[cellI] << endl;
1146 regionType[regionI] = 2;
1147 nChanged++;
1148 break;
1149 }
1150 }
1151 }
1152 }
1153 }
1154 }
1155
1156 reduce(nChanged, sumOp<label>());
1157 DebugInfo<< FUNCTION_NAME << " : Determined regions changed : "
1158 << nChanged << endl;
1159
1160 if (nChanged == 0)
1161 {
1162 break;
1163 }
1164 }
1165
1166
1167 // See which regions have not been visited (regionType == 1)
1168 // label count = 0;
1169 forAll(cellRegion, cellI)
1170 {
1171 label type = regionType[cellRegion[cellI]];
1172 if (type == 1 && cellTypes[cellI] != HOLE)
1173 {
1174 cellTypes[cellI] = HOLE;
1175 // ++count;
1176 }
1177 }
1178}
1179
1180
1182(
1183 const point& sample,
1184 const pointList& donorCcs,
1185 scalarList& weights
1186) const
1187{
1188 // Inverse-distance weighting
1189
1190 weights.setSize(donorCcs.size());
1191 scalar sum = 0;
1192 forAll(donorCcs, i)
1193 {
1194 const scalar d = mag(sample-donorCcs[i]);
1195
1196 if (d > ROOTVSMALL)
1197 {
1198 weights[i] = scalar(1)/d;
1199 sum += weights[i];
1200 }
1201 else
1202 {
1203 // Short circuit
1204 weights = scalar(0);
1205 weights[i] = scalar(1);
1206 return;
1207 }
1208 }
1209 forAll(weights, i)
1210 {
1211 weights[i] /= sum;
1212 }
1213}
1214
1215
1217(
1218 const globalIndex& globalCells
1219)
1220{
1221 // Detects holes that are used for interpolation. If so
1222 // - make type interpolated
1223 // - add interpolation to nearest non-hole cell(s)
1224
1225 const labelList& owner = mesh_.faceOwner();
1226 const labelList& neighbour = mesh_.faceNeighbour();
1227
1228
1229 // Get hole-cells that are used in an interpolation stencil
1230 boolList isDonor(cellInterpolationMap().constructSize(), false);
1231 label nHoleDonors = 0;
1232 {
1233 for (const label celli : interpolationCells_)
1234 {
1235 const labelList& slots = cellStencil_[celli];
1236 UIndirectList<bool>(isDonor, slots) = true;
1237 }
1238
1240 (
1242 List<labelPair>(),
1243 mesh_.nCells(),
1244 cellInterpolationMap().constructMap(),
1245 false,
1246 cellInterpolationMap().subMap(),
1247 false,
1248 isDonor,
1249 false,
1250 orEqOp<bool>(),
1251 flipOp() // negateOp
1252 );
1253
1254
1255 nHoleDonors = 0;
1256 forAll(cellTypes_, celli)
1257 {
1258 if (cellTypes_[celli] == cellCellStencil::HOLE && isDonor[celli])
1259 {
1260 nHoleDonors++;
1261 }
1262 }
1263 reduce(nHoleDonors, sumOp<label>());
1264 if (debug)
1265 {
1266 Pout<< "Detected " << nHoleDonors
1267 << " hole cells that are used as donors" << endl;
1268 }
1269 }
1270
1271
1272 if (nHoleDonors)
1273 {
1274 // Convert map and stencil back to global indices by 'interpolating'
1275 // global cells. Note: since we're only adding interpolations
1276 // (HOLE->SPECIAL) this could probably be done more efficiently.
1277
1278 labelList globalCellIDs(identity(mesh_.nCells()));
1279 globalCells.inplaceToGlobal(Pstream::myProcNo(), globalCellIDs);
1280 cellInterpolationMap().distribute(globalCellIDs);
1281
1282 forAll(interpolationCells_, i)
1283 {
1284 label cellI = interpolationCells_[i];
1285 const labelList& slots = cellStencil_[cellI];
1286 cellStencil_[cellI] = UIndirectList<label>(globalCellIDs, slots)();
1287 }
1288
1289
1290 labelList nbrTypes;
1291 syncTools::swapBoundaryCellList(mesh_, cellTypes_, nbrTypes);
1292
1293 label nSpecialNear = 0;
1294 label nSpecialFar = 0;
1295
1296
1297 // 1. Find hole cells next to live cells
1298 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1299
1300 {
1301 List<pointList> donorCcs(mesh_.nCells());
1302 labelListList donorCells(mesh_.nCells());
1303
1304 for (label facei = 0; facei < mesh_.nInternalFaces(); ++facei)
1305 {
1306 const label own = owner[facei];
1307 const bool ownHole = (cellTypes_[own] == cellCellStencil::HOLE);
1308 const label nbr = neighbour[facei];
1309 const bool nbrHole = (cellTypes_[nbr] == cellCellStencil::HOLE);
1310
1311 if (isDonor[own] && ownHole && !nbrHole)
1312 {
1313 donorCells[own].append(globalCells.toGlobal(nbr));
1314 donorCcs[own].append(mesh_.cellCentres()[nbr]);
1315 }
1316 else if (isDonor[nbr] && nbrHole && !ownHole)
1317 {
1318 donorCells[nbr].append(globalCells.toGlobal(own));
1319 donorCcs[nbr].append(mesh_.cellCentres()[own]);
1320 }
1321 }
1322 labelList globalCellIDs(identity(mesh_.nCells()));
1323 globalCells.inplaceToGlobal(Pstream::myProcNo(), globalCellIDs);
1324 labelList nbrCells;
1325 syncTools::swapBoundaryCellList(mesh_, globalCellIDs, nbrCells);
1326 pointField nbrCc;
1327 syncTools::swapBoundaryCellList(mesh_, mesh_.cellCentres(), nbrCc);
1328 forAll(nbrTypes, bFacei)
1329 {
1330 const label facei = bFacei+mesh_.nInternalFaces();
1331 const label own = owner[facei];
1332 const bool ownHole = (cellTypes_[own] == cellCellStencil::HOLE);
1333 const bool nbrHole = (nbrTypes[bFacei] == cellCellStencil::HOLE);
1334
1335 if (isDonor[own] && ownHole && !nbrHole)
1336 {
1337 donorCells[own].append(nbrCells[bFacei]);
1338 donorCcs[own].append(nbrCc[bFacei]);
1339 }
1340 }
1341
1342 forAll(donorCcs, celli)
1343 {
1344 if (donorCcs[celli].size())
1345 {
1346 cellStencil_[celli] = std::move(donorCells[celli]);
1347 stencilWeights
1348 (
1349 mesh_.cellCentres()[celli],
1350 donorCcs[celli],
1351 cellInterpolationWeights_[celli]
1352 );
1353 cellTypes_[celli] = SPECIAL;
1354 interpolationCells_.append(celli);
1355 cellInterpolationWeight_[celli] = scalar(1);
1356 nSpecialNear++;
1357 }
1358 }
1359 }
1360
1361
1362
1363
1364 // 2. Walk to find (topologically) nearest 'live' cell
1365 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1366 // This will do the remainder but will only find a single donor
1367
1368 // Field on cells and faces.
1369 List<topoDistanceData<int>> cellData(mesh_.nCells());
1370 List<topoDistanceData<int>> faceData(mesh_.nFaces());
1371
1372 int dummyTrackData(0);
1373
1374 {
1375 // Mark all non-hole cells to disable any wave across them
1376 forAll(cellTypes_, celli)
1377 {
1378 label cType = cellTypes_[celli];
1379 if
1380 (
1383 )
1384 {
1385 cellData[celli] = topoDistanceData<int>(labelMin, 0);
1386 }
1387 }
1388
1389 DynamicList<label> seedFaces(nHoleDonors);
1390 DynamicList<topoDistanceData<int>> seedData(nHoleDonors);
1391
1392 for (label facei = 0; facei < mesh_.nInternalFaces(); ++facei)
1393 {
1394 const label own = owner[facei];
1395 const bool ownLive =
1396 (
1397 cellTypes_[own] == cellCellStencil::CALCULATED
1398 || cellTypes_[own] == cellCellStencil::INTERPOLATED
1399 );
1400
1401 const label nbr = neighbour[facei];
1402 const bool nbrLive =
1403 (
1404 cellTypes_[nbr] == cellCellStencil::CALCULATED
1405 || cellTypes_[nbr] == cellCellStencil::INTERPOLATED
1406 );
1407
1408 if (ownLive && !nbrLive)
1409 {
1410 seedFaces.append(facei);
1411 const label globalOwn = globalCells.toGlobal(own);
1412 seedData.append(topoDistanceData<int>(globalOwn, 0));
1413 }
1414 else if (!ownLive && nbrLive)
1415 {
1416 seedFaces.append(facei);
1417 const label globalNbr = globalCells.toGlobal(nbr);
1418 seedData.append(topoDistanceData<int>(globalNbr, 0));
1419 }
1420 }
1421
1422 forAll(nbrTypes, bFacei)
1423 {
1424 const label facei = bFacei+mesh_.nInternalFaces();
1425 const label own = owner[facei];
1426 const bool ownLive =
1427 (
1428 cellTypes_[own] == cellCellStencil::CALCULATED
1429 || cellTypes_[own] == cellCellStencil::INTERPOLATED
1430 );
1431 const bool nbrLive =
1432 (
1433 nbrTypes[bFacei] == cellCellStencil::CALCULATED
1434 || nbrTypes[bFacei] == cellCellStencil::INTERPOLATED
1435 );
1436
1437 if (ownLive && !nbrLive)
1438 {
1439 seedFaces.append(facei);
1440 const label globalOwn = globalCells.toGlobal(own);
1441 seedData.append(topoDistanceData<int>(globalOwn, 0));
1442 }
1443 }
1444
1445 // Propagate information inwards
1446 FaceCellWave<topoDistanceData<int>> distanceCalc
1447 (
1448 mesh_,
1449 seedFaces,
1450 seedData,
1451 faceData,
1452 cellData,
1453 mesh_.globalData().nTotalCells()+1,
1454 dummyTrackData
1455 );
1456 }
1457
1458
1459 // Add the new donor ones
1460 forAll(cellData, celli)
1461 {
1462 if
1463 (
1464 cellTypes_[celli] == cellCellStencil::HOLE
1465 && isDonor[celli]
1466 && cellData[celli].valid(dummyTrackData)
1467 )
1468 {
1469 cellStencil_[celli].setSize(1);
1470 cellStencil_[celli] = cellData[celli].data(); // global cellID
1471 cellInterpolationWeights_[celli].setSize(1);
1472 cellInterpolationWeights_[celli] = 1.0;
1473 cellTypes_[celli] = SPECIAL; // handled later on
1474 interpolationCells_.append(celli);
1475 cellInterpolationWeight_[celli] = 1.0;
1476 nSpecialFar++;
1477 }
1478 }
1479
1480 if (debug & 2)
1481 {
1482 reduce(nSpecialNear, sumOp<label>());
1483 reduce(nSpecialFar, sumOp<label>());
1484
1485 Pout<< "Detected " << nHoleDonors
1486 << " hole cells that are used as donors of which" << nl
1487 << " next to live cells : " << nSpecialNear << nl
1488 << " other : " << nSpecialFar << nl
1489 << endl;
1490 }
1491
1492
1493 // Re-do the mapDistribute
1494 List<Map<label>> compactMap;
1495 cellInterpolationMap_.reset
1496 (
1497 new mapDistribute
1498 (
1499 globalCells,
1500 cellStencil_,
1501 compactMap
1502 )
1503 );
1504 }
1505}
1506
1507
1509(
1510 const globalIndex& globalCells,
1511 const bool allowHoleDonors
1512)
1513{
1514 // Send cell centre back to donor
1515 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1516 // The complication is that multiple acceptors need the same donor
1517 // (but with different weights obviously)
1518 // So we do multi-pass:
1519 // - send over cc of acceptor for which we want stencil.
1520 // Consistently choose the acceptor with smallest magSqr in case of
1521 // multiple acceptors for the containing cell/donor.
1522 // - find the cell-cells and weights for the donor
1523 // - send back together with the acceptor cc
1524 // - use the acceptor cc to see if it was 'me' that sent it. If so
1525 // mark me as complete so it doesn't get involved in the next loop.
1526 // - loop until all resolved.
1527
1528 // Special value for unused points
1529 const vector greatPoint(GREAT, GREAT, GREAT);
1530
1531 boolList isValidDonor(mesh_.nCells(), true);
1532 if (!allowHoleDonors)
1533 {
1534 forAll(cellTypes_, celli)
1535 {
1536 if (cellTypes_[celli] == HOLE)
1537 {
1538 isValidDonor[celli] = false;
1539 }
1540 }
1541 }
1542
1543 // Has acceptor been handled already?
1544 bitSet doneAcceptor(interpolationCells_.size());
1545
1546 while (true)
1547 {
1548 pointField samples(cellInterpolationMap().constructSize(), greatPoint);
1549
1550 // Fill remote slots (override old content). We'll find out later
1551 // on which one has won and mark this one in doneAcceptor.
1552 label nSamples = 0;
1553 forAll(interpolationCells_, i)
1554 {
1555 if (!doneAcceptor[i])
1556 {
1557 label cellI = interpolationCells_[i];
1558 const point& cc = mesh_.cellCentres()[cellI];
1559 const labelList& slots = cellStencil_[cellI];
1560
1561 //- Note: empty slots can happen for interpolated cells with
1562 // bad/insufficient donors. These are handled later on.
1563 if (slots.size() > 1)
1564 {
1565 FatalErrorInFunction<< "Problem:" << slots
1566 << abort(FatalError);
1567 }
1568 else if (slots.size())
1569 {
1570 forAll(slots, sloti)
1571 {
1572 const label elemi = slots[sloti];
1573 //Pout<< " acceptor:" << cellI
1574 // << " at:" << mesh_.cellCentres()[cellI]
1575 // << " global:" << globalCells.toGlobal(cellI)
1576 // << " found in donor:" << elemi << endl;
1577 minMagSqrEqOp<point>()(samples[elemi], cc);
1578 }
1579 nSamples++;
1580 }
1581 }
1582 }
1583
1584
1586 {
1587 break;
1588 }
1589
1590 // Send back to donor. Make sure valid point takes priority
1592 (
1595 mesh_.nCells(),
1596 cellInterpolationMap().constructMap(),
1597 false,
1598 cellInterpolationMap().subMap(),
1599 false,
1600 samples,
1601 greatPoint, // nullValue
1602 minMagSqrEqOp<point>(),
1603 flipOp(), // NegateOp
1605 cellInterpolationMap().comm()
1606 );
1607
1608 // All the donor cells will now have a valid cell centre. Construct a
1609 // stencil for these.
1610
1611 DynamicList<label> donorCells(mesh_.nCells());
1612 forAll(samples, cellI)
1613 {
1614 if (samples[cellI] != greatPoint)
1615 {
1616 donorCells.append(cellI);
1617 }
1618 }
1619
1620
1621 // Get neighbours (global cell and centre) of donorCells.
1622 labelListList donorCellCells(mesh_.nCells());
1623 pointListList donorCellCentres(mesh_.nCells());
1624 globalCellCells
1625 (
1626 globalCells,
1627 mesh_,
1628 isValidDonor,
1629 donorCells,
1630 donorCellCells,
1631 donorCellCentres
1632 );
1633
1634 // Determine the weights.
1635 scalarListList donorWeights(mesh_.nCells());
1636 forAll(donorCells, i)
1637 {
1638 label cellI = donorCells[i];
1639 const pointList& donorCentres = donorCellCentres[cellI];
1640 stencilWeights
1641 (
1642 samples[cellI],
1643 donorCentres,
1644 donorWeights[cellI]
1645 );
1646 }
1647
1648 // Transfer the information back to the acceptor:
1649 // - donorCellCells : stencil (with first element the original donor)
1650 // - donorWeights : weights for donorCellCells
1651 cellInterpolationMap().distribute(donorCellCells);
1652 cellInterpolationMap().distribute(donorWeights);
1653 cellInterpolationMap().distribute(samples);
1654
1655 // Check which acceptor has won and transfer
1656 forAll(interpolationCells_, i)
1657 {
1658 if (!doneAcceptor[i])
1659 {
1660 label cellI = interpolationCells_[i];
1661 const labelList& slots = cellStencil_[cellI];
1662
1663 if (slots.size() > 1)
1664 {
1665 FatalErrorInFunction << "Problem:" << slots
1666 << abort(FatalError);
1667 }
1668 else if (slots.size() == 1)
1669 {
1670 const label sloti = slots[0];
1671
1672 // Important: check if the stencil is actually for this cell
1673 if (samples[sloti] == mesh_.cellCentres()[cellI])
1674 {
1675 cellStencil_[cellI].transfer(donorCellCells[sloti]);
1676 cellInterpolationWeights_[cellI].transfer
1677 (
1678 donorWeights[sloti]
1679 );
1680 // Mark cell as being done so it does not get sent over
1681 // again.
1682 doneAcceptor.set(i);
1683 }
1684 }
1685 }
1686 }
1687 }
1688
1689 // Re-do the mapDistribute
1690 List<Map<label>> compactMap;
1691 cellInterpolationMap_.reset
1692 (
1693 new mapDistribute
1694 (
1695 globalCells,
1696 cellStencil_,
1697 compactMap
1698 )
1699 );
1700}
1702
1703// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1704
1705Foam::cellCellStencils::inverseDistance::inverseDistance
1706(
1707 const fvMesh& mesh,
1708 const dictionary& dict,
1709 const bool doUpdate
1710)
1711:
1712 cellCellStencil(mesh),
1713 dict_(dict),
1714 allowHoleDonors_(dict.getOrDefault("allowHoleDonors", false)),
1716 (
1717 dict.getOrDefault("allowInterpolatedDonors", true)
1718 ),
1719 smallVec_(Zero),
1720 cellTypes_(labelList(mesh.nCells(), CALCULATED)),
1723 cellStencil_(0),
1726 (
1727 IOobject
1728 (
1729 "cellInterpolationWeight",
1730 mesh_.facesInstance(),
1731 mesh_,
1732 IOobject::NO_READ,
1733 IOobject::NO_WRITE,
1734 IOobject::NO_REGISTER
1735 ),
1736 mesh_,
1738 fvPatchFieldBase::zeroGradientType()
1739 )
1740{
1741 // Add motion-solver fields to non-interpolated fields
1743
1744 // Read zoneID
1745 this->zoneID();
1746
1747 // Read old-time cellTypes
1748 IOobject io
1749 (
1750 "cellTypes",
1751 mesh_.time().timeName(),
1752 mesh_,
1756 );
1757 if (io.typeHeaderOk<volScalarField>(true))
1758 {
1759 if (debug)
1760 {
1761 Pout<< "Reading cellTypes from time " << mesh_.time().timeName()
1762 << endl;
1763 }
1764
1765 const volScalarField volCellTypes(io, mesh_);
1766 forAll(volCellTypes, celli)
1767 {
1768 // Round to integer
1769 cellTypes_[celli] = volCellTypes[celli];
1770 }
1771 }
1772
1773 if (doUpdate)
1774 {
1775 update();
1776 }
1777}
1779
1780// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1781
1783{}
1785
1786// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1787
1789{
1790 scalar layerRelax(dict_.getOrDefault("layerRelax", 1.0));
1791
1792 scalar tol = dict_.getOrDefault("tolerance", 1e-10);
1793 smallVec_ = mesh_.bounds().span()*tol;
1794
1795 const labelIOList& zoneID = this->zoneID();
1796
1797 label nZones = gMax(zoneID)+1;
1798 labelList nCellsPerZone(nZones, Zero);
1799 forAll(zoneID, cellI)
1800 {
1801 nCellsPerZone[zoneID[cellI]]++;
1802 }
1803 Pstream::listReduce(nCellsPerZone, sumOp<label>());
1804
1805 const boundBox& allBb = mesh_.bounds();
1806
1807 PtrList<fvMeshSubset> meshParts(nZones);
1809
1810 // Determine zone meshes and bounding boxes
1811 {
1812 // Per processor, per zone the bounding box
1814 procBb[Pstream::myProcNo()].setSize(nZones);
1815
1816 forAll(meshParts, zonei)
1817 {
1818 meshParts.set
1819 (
1820 zonei,
1821 new fvMeshSubset(mesh_, zonei, zoneID)
1822 );
1823 const fvMesh& subMesh = meshParts[zonei].subMesh();
1824
1825 // Trigger early evaluation of mesh dimension (in case there are
1826 // zero cells in mesh)
1827 (void)subMesh.nGeometricD();
1828
1829 if (subMesh.nPoints())
1830 {
1831 procBb[Pstream::myProcNo()][zonei] =
1832 treeBoundBox(subMesh.points());
1833 procBb[Pstream::myProcNo()][zonei].inflate(1e-6);
1834 }
1835 else
1836 {
1837 // No part of zone on this processor. Make up bb.
1838 procBb[Pstream::myProcNo()][zonei] = treeBoundBox
1839 (
1840 allBb.min() - 2*allBb.span(),
1841 allBb.min() - allBb.span()
1842 );
1843 procBb[Pstream::myProcNo()][zonei].inflate(1e-6);
1844 }
1845 }
1846
1847 Pstream::allGatherList(procBb);
1848
1849 // Move local bounding boxes to per-mesh indexing
1850 forAll(meshBb, zoneI)
1851 {
1852 treeBoundBoxList& bbs = meshBb[zoneI];
1853 bbs.setSize(Pstream::nProcs());
1854 forAll(procBb, procI)
1855 {
1856 bbs[procI] = procBb[procI][zoneI];
1857 }
1858 }
1859 }
1860
1861
1862 // Determine patch bounding boxes. These are either global and provided
1863 // by the user or processor-local as a copy of the mesh bounding box.
1864
1865 List<treeBoundBoxList> patchBb(nZones);
1866 List<labelVector> patchDivisions(nZones);
1867 PtrList<PackedList<2>> patchParts(nZones);
1868 labelList allPatchTypes(mesh_.nCells(), OTHER);
1869
1870 {
1871 treeBoundBox globalPatchBb;
1872 if (dict_.readIfPresent("searchBox", globalPatchBb))
1873 {
1874 // All processors, all zones have the same bounding box
1875 patchBb = treeBoundBoxList(Pstream::nProcs(), globalPatchBb);
1876 }
1877 else
1878 {
1879 // Use the meshBb (differing per zone, per processor)
1880 patchBb = meshBb;
1881 }
1882 }
1883
1884 {
1885 labelVector globalDivs;
1886 if (dict_.readIfPresent("searchBoxDivisions", globalDivs))
1887 {
1888 patchDivisions = globalDivs;
1889 }
1890 else
1891 {
1892 const labelVector& dim = mesh_.geometricD();
1893 label nDivs = -1;
1894 if (mesh_.nGeometricD() == 1)
1895 {
1896 nDivs = mesh_.nCells();
1897 }
1898 else if (mesh_.nGeometricD() == 2)
1899 {
1900 nDivs = label(Foam::sqrt(scalar(mesh_.nCells())));
1901 }
1902 else
1903 {
1904 nDivs = label(Foam::cbrt(scalar(mesh_.nCells())));
1905 }
1906
1907 labelVector v(nDivs, nDivs, nDivs);
1908 forAll(dim, i)
1909 {
1910 if (dim[i] == -1)
1911 {
1912 v[i] = 1;
1913 }
1914 }
1915 patchDivisions = v;
1916 }
1917 }
1918
1919 forAll(patchParts, zoneI)
1920 {
1921 patchParts.set
1922 (
1923 zoneI,
1924 new PackedList<2>
1925 (
1926 patchDivisions[zoneI][0]
1927 *patchDivisions[zoneI][1]
1928 *patchDivisions[zoneI][2]
1929 )
1930 );
1932 (
1933 meshParts[zoneI].subMesh(),
1934 smallVec_,
1935
1936 patchBb[zoneI][Pstream::myProcNo()],
1937 patchDivisions[zoneI],
1938 patchParts[zoneI],
1939
1940 meshParts[zoneI].cellMap(),
1941 allPatchTypes
1942 );
1943 }
1944
1945
1946 // Print a bit
1947 {
1948 Info<< type() << " : detected " << nZones
1949 << " mesh regions" << endl;
1950 Info<< incrIndent;
1951 forAll(nCellsPerZone, zoneI)
1952 {
1953 Info<< indent<< "zone:" << zoneI
1954 << " nCells:" << nCellsPerZone[zoneI]
1955 << " voxels:" << patchDivisions[zoneI]
1956 << " bb:" << patchBb[zoneI][Pstream::myProcNo()]
1957 << endl;
1958 }
1959 Info<< decrIndent;
1960 }
1961
1962
1963 // Current best guess for cells. Includes best stencil. Weights should
1964 // add up to volume.
1965 labelList allCellTypes(mesh_.nCells(), CALCULATED);
1966 labelListList allStencil(mesh_.nCells());
1967 // zoneID of donor
1968 labelList allDonorID(mesh_.nCells(), -1);
1969
1970 const globalIndex globalCells(mesh_.nCells());
1971
1972 PstreamBuffers pBufs;
1973
1974 // Mark holes (in allCellTypes)
1975 for (label srcI = 0; srcI < meshParts.size()-1; srcI++)
1976 {
1977 for (label tgtI = srcI+1; tgtI < meshParts.size(); tgtI++)
1978 {
1980 (
1981 pBufs,
1982
1983 meshParts,
1984
1985 patchBb,
1986 patchDivisions,
1987 patchParts,
1988
1989 srcI,
1990 tgtI,
1991 allCellTypes
1992 );
1994 (
1995 pBufs,
1996
1997 meshParts,
1998
1999 patchBb,
2000 patchDivisions,
2001 patchParts,
2002
2003 tgtI,
2004 srcI,
2005 allCellTypes
2006 );
2007 }
2008 }
2009
2010 // Find donors (which are not holes) in allStencil, allDonorID
2011 for (label srcI = 0; srcI < meshParts.size()-1; srcI++)
2012 {
2013 for (label tgtI = srcI+1; tgtI < meshParts.size(); tgtI++)
2014 {
2016 (
2017 globalCells,
2018 pBufs,
2019 meshParts,
2020 meshBb,
2021 allCellTypes,
2022
2023 tgtI,
2024 srcI,
2025 allStencil,
2026 allDonorID
2027 );
2029 (
2030 globalCells,
2031 pBufs,
2032 meshParts,
2033 meshBb,
2034 allCellTypes,
2035
2036 srcI,
2037 tgtI,
2038 allStencil,
2039 allDonorID
2040 );
2041 }
2042 }
2043
2044 if ((debug & 2) && mesh_.time().writeTime())
2045 {
2047 (
2048 createField(mesh_, "allCellTypes", allCellTypes)
2049 );
2050 tfld().write();
2051
2052 tmp<volScalarField> tallDonorID
2053 (
2054 createField(mesh_, "allDonorID", allDonorID)
2055 );
2056 tallDonorID().write();
2057 }
2058
2059 // Use the patch types and weights to decide what to do
2060 forAll(allPatchTypes, cellI)
2061 {
2062 if (allCellTypes[cellI] != HOLE)
2063 {
2064 switch (allPatchTypes[cellI])
2065 {
2066 case OVERSET:
2067 {
2068 // Require interpolation. See if possible.
2069 if (allStencil[cellI].size())
2070 {
2071 allCellTypes[cellI] = INTERPOLATED;
2072 }
2073 else
2074 {
2075 allCellTypes[cellI] = HOLE;
2076 }
2077 }
2078 }
2079 }
2080 }
2081
2082 if ((debug & 2) && mesh_.time().writeTime())
2083 {
2085 (
2086 createField(mesh_, "allCellTypes_patch", allCellTypes)
2087 );
2088 tfld().write();
2089
2090 tmp<volScalarField> tfldOld
2091 (
2092 createField(mesh_, "allCellTypes_old", cellTypes_)
2093 );
2094 tfldOld().write();
2095 }
2096
2097 // Mark unreachable bits
2098 findHoles(globalCells, mesh_, zoneID, allStencil, allCellTypes);
2099
2100
2101 if ((debug & 2) && mesh_.time().writeTime())
2102 {
2104 (
2105 createField(mesh_, "allCellTypes_hole", allCellTypes)
2106 );
2107 tfld().write();
2108 }
2109 if ((debug & 2) && mesh_.time().writeTime())
2110 {
2111 labelList stencilSize(mesh_.nCells());
2112 forAll(allStencil, celli)
2113 {
2114 stencilSize[celli] = allStencil[celli].size();
2115 }
2117 (
2118 createField(mesh_, "allStencil_hole", stencilSize)
2119 );
2120 tfld().write();
2121 }
2122
2123 // Update allStencil with new fill HOLES
2124 forAll(allCellTypes, celli)
2125 {
2126 if (allCellTypes[celli] == HOLE && allStencil[celli].size())
2127 {
2128 allStencil[celli].clear();
2129 }
2130 }
2131
2132// // Check previous iteration cellTypes_ for any hole->calculated changes
2133// // If so set the cell either to interpolated (if there are donors) or
2134// // holes (if there are no donors). Note that any interpolated cell might
2135// // still be overwritten by the flood filling
2136// {
2137// label nCalculated = 0;
2138//
2139// forAll(cellTypes_, celli)
2140// {
2141// if (allCellTypes[celli] == CALCULATED && cellTypes_[celli] == HOLE)
2142// {
2143// if (allStencil[celli].size() == 0)
2144// {
2145// // Reset to hole
2146// allCellTypes[celli] = HOLE;
2147// allStencil[celli].clear();
2148// }
2149// else
2150// {
2151// allCellTypes[celli] = INTERPOLATED;
2152// nCalculated++;
2153// }
2154// }
2155// }
2156//
2157// if (debug)
2158// {
2159// Pout<< "Detected " << nCalculated << " cells changing from hole"
2160// << " to calculated. Changed to interpolated"
2161// << endl;
2162// }
2163// }
2164
2165 if ((debug & 2) && mesh_.time().writeTime())
2166 {
2168 (
2169 createField(mesh_, "allCellTypes_pre_front", allCellTypes)
2170 );
2171 tfld().write();
2172 }
2173 // Add buffer interpolation layer(s) around holes
2174 scalarField allWeight(mesh_.nCells(), Zero);
2175
2176 labelListList compactStencil(allStencil);
2177 List<Map<label>> compactStencilMap;
2178 mapDistribute map(globalCells, compactStencil, compactStencilMap);
2179
2180 scalarList compactCellVol(mesh_.V());
2181 map.distribute(compactCellVol);
2182
2183
2184 label useLayer = dict_.getOrDefault("useLayer", -1);
2185 const dictionary& fvSchemes = mesh_.schemesDict();
2186 if (fvSchemes.found("oversetInterpolation"))
2187 {
2188 const dictionary& intDict = fvSchemes.subDict
2189 (
2190 "oversetInterpolation"
2191 );
2192 useLayer = intDict.getOrDefault("useLayer", -1);
2193 }
2194
2195 walkFront
2196 (
2197 globalCells,
2198 layerRelax,
2199 allStencil,
2200 allCellTypes,
2201 allWeight,
2202 compactCellVol,
2203 compactStencil,
2204 zoneID,
2205 dict_.getOrDefault("holeLayers", 1),
2206 useLayer
2207 );
2208
2209 if ((debug & 2) && mesh_.time().writeTime())
2210 {
2212 (
2213 createField(mesh_, "allCellTypes_front", allCellTypes)
2214 );
2215 tfld().write();
2216 }
2217
2218 // Check previous iteration cellTypes_ for any hole->calculated changes
2219 // If so set the cell either to interpolated (if there are donors) or
2220 // holes (if there are no donors). Note that any interpolated cell might
2221 // still be overwritten by the flood filling
2222 {
2223 label nCalculated = 0;
2224
2225 forAll(cellTypes_, celli)
2226 {
2227 if (allCellTypes[celli] == CALCULATED && cellTypes_[celli] == HOLE)
2228 {
2229 if (allStencil[celli].size() == 0)
2230 {
2231 // Reset to hole
2232 allCellTypes[celli] = HOLE;
2233 allStencil[celli].clear();
2234 }
2235 else
2236 {
2237 allCellTypes[celli] = INTERPOLATED;
2238 nCalculated++;
2239 }
2240 }
2241 }
2242
2243 if (debug)
2244 {
2245 Pout<< "Detected " << nCalculated << " cells changing from hole"
2246 << " to calculated. Changed to interpolated"
2247 << endl;
2248 }
2249 }
2250
2251
2252
2253
2254 // Convert cell-cell addressing to stencil in compact notation
2255 // Clean any potential INTERPOLATED with HOLE donors.
2256 // This is not eliminated in the front walk as the HOLE cells might
2257 // leak when inset region is overlapping backgroung mesh
2258 cellTypes_.transfer(allCellTypes);
2259 cellStencil_.setSize(mesh_.nCells());
2260 cellInterpolationWeights_.setSize(mesh_.nCells());
2262
2263/*
2264 forAll(cellTypes_, cellI)
2265 {
2266 if (cellTypes_[cellI] == INTERPOLATED)
2267 {
2268 if (cellTypes_[allStencil[cellI][0]] != HOLE)
2269 {
2270 cellStencil_[cellI].transfer(allStencil[cellI]);
2271 cellInterpolationWeights_[cellI].setSize(1);
2272 cellInterpolationWeights_[cellI][0] = 1.0; interpolationCells.append(cellI);
2273 }
2274 else
2275 {
2276 cellTypes_[cellI] = POROUS;//CALCULATED;
2277 cellStencil_[cellI].clear();
2278 cellInterpolationWeights_[cellI].clear();
2279 }
2280 }
2281 else
2282 {
2283 cellStencil_[cellI].clear();
2284 cellInterpolationWeights_[cellI].clear();
2285 }
2286 }
2287*/
2288
2289 //labelListList compactStencil(allStencil);
2290 //List<Map<label>> compactStencilMap;
2291 //mapDistribute map(globalCells, compactStencil, compactStencilMap);
2292
2293 labelList compactCellTypes(cellTypes_);
2294 map.distribute(compactCellTypes);
2295
2296 label nInterpolatedDonors = 0;
2297 forAll(cellTypes_, celli)
2298 {
2299 if (cellTypes_[celli] == INTERPOLATED)
2300 {
2301 const labelList& slots = compactStencil[celli];
2302 if (slots.size())
2303 {
2304 if
2305 (
2306 compactCellTypes[slots[0]] == HOLE
2307 ||
2308 (
2310 && compactCellTypes[slots[0]] == INTERPOLATED
2311 )
2312 )
2313 {
2314 cellTypes_[celli] = POROUS;
2315 cellStencil_[celli].clear();
2316 cellInterpolationWeights_[celli].clear();
2317 nInterpolatedDonors++;
2318 }
2319 else
2320 {
2321 cellStencil_[celli].transfer(allStencil[celli]);
2322 cellInterpolationWeights_[celli].setSize(1);
2323 cellInterpolationWeights_[celli][0] = 1.0;
2324 interpolationCells.append(celli);
2325 }
2326 }
2327 }
2328 else
2329 {
2330 cellStencil_[celli].clear();
2331 cellInterpolationWeights_[celli].clear();
2332 }
2333 }
2334
2335 reduce(nInterpolatedDonors, sumOp<label>());
2336
2338 << "interpolate Donors : "
2339 << nInterpolatedDonors << endl;
2340
2341 // Reset map with new cellStencil
2342 List<Map<label>> compactMap;
2344 (
2345 new mapDistribute(globalCells, cellStencil_, compactMap)
2346 );
2347
2349
2350 cellInterpolationWeight_.transfer(allWeight);
2352 <
2355 false
2356 >(cellInterpolationWeight_.boundaryFieldRef());
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();
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:"
2429 << endl;
2430 }
2431 }
2433 (
2434 createField(mesh_, "maxMagWeight", maxMagWeight)
2435 );
2437 <
2440 false
2441 >(tfld.ref().boundaryFieldRef());
2442 tfld().write();
2443 }
2444
2445 // Dump cell types
2446 {
2448 (
2449 createField(mesh_, "cellTypes", cellTypes_)
2450 );
2451 //tfld.ref().correctBoundaryConditions();
2453 <
2456 false
2457 >(tfld.ref().boundaryFieldRef());
2458 tfld().write();
2459 }
2460
2461
2462 // Dump stencil, one per zone
2463 mkDir(mesh_.time().timePath());
2464 pointField cc(mesh_.cellCentres());
2465 cellInterpolationMap().distribute(cc);
2466 forAll(meshParts, zonei)
2467 {
2468 OBJstream str
2469 (
2470 mesh_.time().timePath()
2471 + "/stencil_" + name(zonei) + ".obj"
2472 );
2473 Pout<< type() << " : dumping to " << str.name() << endl;
2474
2475 const labelList& subMeshCellMap = meshParts[zonei].cellMap();
2476
2477 forAll(subMeshCellMap, subcelli)
2478 {
2479 const label celli = subMeshCellMap[subcelli];
2480 const labelList& slots = cellStencil_[celli];
2481 const point& accCc = mesh_.cellCentres()[celli];
2482 forAll(slots, i)
2483 {
2484 const point& donorCc = cc[slots[i]];
2485 str.writeLine(accCc, 0.1*accCc+0.9*donorCc);
2486 }
2487 }
2488 }
2489 }
2490
2491 // Print some stats
2492 {
2493 labelList nCells(count(3, cellTypes_));
2494
2495 label nLocal = 0;
2496 label nMixed = 0;
2497 label nRemote = 0;
2499 {
2500 label celli = interpolationCells_[i];
2501 const labelList& slots = cellStencil_[celli];
2502
2503 bool hasLocal = false;
2504 bool hasRemote = false;
2505
2506 forAll(slots, sloti)
2507 {
2508 if (slots[sloti] >= mesh_.nCells())
2509 {
2510 hasRemote = true;
2511 }
2512 else
2513 {
2514 hasLocal = true;
2515 }
2516 }
2517
2518 if (hasRemote)
2519 {
2520 if (!hasLocal)
2521 {
2522 nRemote++;
2523 }
2524 else
2525 {
2526 nMixed++;
2527 }
2528 }
2529 else if (hasLocal)
2530 {
2531 nLocal++;
2532 }
2533 }
2534
2535 Info<< "Overset analysis : nCells : "
2536 << returnReduce(cellTypes_.size(), sumOp<label>()) << nl
2537 << incrIndent
2538 << indent << "calculated : " << nCells[CALCULATED] << nl
2539 << indent << "interpolated : " << nCells[INTERPOLATED]
2540 << " (from local:" << returnReduce(nLocal, sumOp<label>())
2541 << " mixed local/remote:" << returnReduce(nMixed, sumOp<label>())
2542 << " remote:" << returnReduce(nRemote, sumOp<label>()) << ")" << nl
2543 << indent << "hole : " << nCells[HOLE] << nl
2544 << decrIndent << endl;
2545 }
2546
2547 // Tbd: detect if anything changed. Most likely it did!
2548 return true;
2549}
2550
2551
2552// ************************************************************************* //
label k
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void append(const T &val)
Copy append an element to the end of this list.
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Minimal example by using system/controlDict.functions:
@ NO_REGISTER
Do not request registration (bool: false).
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void append(const T &val)
Append an element at the end of the list.
Definition List.H:497
void setSize(label n)
Alias for resize().
Definition List.H:536
void clear()
Clear the list, i.e. set size to zero.
Definition ListI.H:133
static const List< T > & null() noexcept
Return a null List (reference to a nullObject). Behaves like an empty List.
Definition List.H:138
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files.
Definition OBJstream.H:58
Ostream & writeLine(const point &p0, const point &p1)
Write line joining two points.
Definition OBJstream.C:234
virtual const fileName & name() const override
Read/write access to the name of the stream.
Definition OSstream.H:134
A dynamic list of packed unsigned integers, with the number of bits per item specified by the <Width>...
Definition PackedList.H:146
Buffers for inter-processor communications streams (UOPstream, UIPstream).
void finishedSends(const bool wait=true)
Mark the send phase as being finished.
void clear()
Clear all send/recv buffers and reset states.
static void allGatherList(UList< T > &values, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Gather data, but keep individual values separate. Uses MPI_Allgather or manual communication.
static void listReduce(UList< T > &values, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce list elements (list must be equal size on all ranks), applying bop to each list element.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
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:171
Input inter-processor communications stream using MPI send/recv etc. - operating on external buffer.
Definition UIPstream.H:313
A List with indirect addressing. Like IndirectList but does not store addressing.
T * data() noexcept
Return pointer to the underlying array serving as data storage.
Definition UListI.H:274
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Output inter-processor communications stream using MPI send/recv etc. - operating on external buffer.
Definition UOPstream.H:408
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
@ nonBlocking
"nonBlocking" (immediate) : (MPI_Isend, MPI_Irecv)
Definition UPstream.H:84
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
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:1697
static rangeType allProcs(const label communicator=worldComm)
Range of process indices for all processes.
Definition UPstream.H:1857
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
static const Form zero
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition bitSetI.H:502
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
const point & max() const noexcept
Maximum describing the bounding box.
Definition boundBoxI.H:168
void grow(const scalar delta)
Expand box by adjusting min/max by specified amount in each dimension.
Definition boundBoxI.H:367
bool overlaps(const boundBox &bb) const
Overlaps/touches boundingBox?
Definition boundBoxI.H:439
const point & min() const noexcept
Minimum describing the bounding box.
Definition boundBoxI.H:162
vector span() const
The bounding box span (from minimum to maximum).
Definition boundBoxI.H:192
Calculation of interpolation stencils.
const fvMesh & mesh_
Reference to the mesh.
void walkFront(const globalIndex &globalCells, const scalar layerRelax, const labelListList &allStencil, labelList &allCellTypes, scalarField &allWeight, const scalarList &compactCellVol, const labelListList &compactStencil, const labelList &zoneID, const label holeLayers, const label useLayer) const
Surround holes with layer(s) of interpolated cells.
static labelList count(const label size, const labelUList &lst)
Count occurrences (in parallel).
const labelIOList & zoneID() const
Helper: get reference to registered zoneID. Loads volScalarField.
static const labelIOList & zoneID(const fvMesh &)
Helper: get reference to registered zoneID. Loads volScalarField.
void suppressMotionFields()
Helper: populate nonInterpolatedFields_ with motion solver.
static void globalCellCells(const globalIndex &gi, const polyMesh &mesh, const boolList &isValidDonor, const labelList &selectedCells, labelListList &cellCells, pointListList &cellCellCentres)
Helper: create cell-cell addressing in global numbering.
Inverse-distance-weighted interpolation stencil.
vector smallVec_
Small perturbation vector for geometric tests.
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.
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.
virtual void stencilWeights(const point &sample, const pointList &donorCcs, scalarList &weights) const
Calculate inverse distance weights for a single acceptor.
static point position(const boundBox &bb, const labelVector &nDivs, const label boxI)
Convert index back into coordinate.
static label index(const labelVector &nDivs, const labelVector &)
Convert ijk indices into single index.
virtual const labelUList & interpolationCells() const
Indices of interpolated cells.
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.
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.
static labelVector index3(const labelVector &nDivs, const label)
Convert single index into ijk.
bool betterDonor(const label destMesh, const label currentDonorMesh, const label newDonorMesh) const
If multiple donors meshes: decide which is best.
virtual const mapDistribute & cellInterpolationMap() const
Return a communication schedule.
virtual const labelUList & cellTypes() const
Return the cell type list.
const bool allowInterpolatedDonors_
Allow interpolared as donors.
autoPtr< mapDistribute > cellInterpolationMap_
Fetch interpolated cells.
virtual bool update()
Update stencils. Return false if nothing changed.
void holeExtrapolationStencil(const globalIndex &globalCells)
Make holes next to live ones type SPECIAL and include in interpolation stencil.
const dictionary dict_
Dictionary of motion control parameters.
volScalarField cellInterpolationWeight_
Amount of interpolation.
scalarListList cellInterpolationWeights_
Interpolation weights.
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.
labelList interpolationCells_
Indices of interpolated cells.
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 void createStencil(const globalIndex &, const bool allowHoleDonors)
Create stencil starting from the donor containing the acceptor allowHoleDonors : allow donors of type...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
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...
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
A fvBoundaryMesh is a fvPatch list with a reference to the associated fvMesh, with additional search ...
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Template invariant parts for fvPatchField.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
const polyPatch & patch() const noexcept
Return the polyPatch.
Definition fvPatch.H:202
static bool constraintType(const word &patchType)
Return true if the given type is a constraint type.
Definition fvPatch.C:74
virtual const labelUList & faceCells() const
Return faceCells.
Definition fvPatch.C:107
Selector class for finite volume differencing schemes. fvMesh is derived from fvSchemes so that all f...
Definition fvSchemes.H:54
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
void inplaceToGlobal(const label proci, labelUList &labels) const
From local to global index on proci (inplace).
label toGlobal(const label proci, const label i) const
From local to global on proci.
static void distribute(const UPstream::commsTypes commsType, const UList< labelPair > &schedule, const label constructSize, const labelListList &subMap, const bool subHasFlip, const labelListList &constructMap, const bool constructHasFlip, List< T > &field, const T &nullValue, const CombineOp &cop, const NegateOp &negOp, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Distribute combine data with specified combine operation and negate operator (for flips).
Class containing processor-to-processor mapping information.
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.
static void correctBoundaryConditions(typename GeoField::Boundary &bfld)
Correct boundary conditions of certain type (TypeOnly = true) or explicitly not of the type (TypeOnly...
Boundary condition for use on overset patches. To be run in combination with special dynamicFvMesh ty...
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition polyMesh.C:861
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition polyMesh.C:1496
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition polyMesh.C:884
virtual const pointField & points() const
Return raw points.
Definition polyMesh.C:1063
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
const vectorField & cellCentres() const
label nPoints() const noexcept
Number of mesh points.
label nCells() const noexcept
Number of mesh cells.
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
label nRegions() const
Return total number of regions.
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData, const bool parRun=UPstream::parRun())
Extract and swap to obtain neighbour cell values for all boundary faces.
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
For use with FaceCellWave. Determines topological distance to starting faces. Templated on passive tr...
Standard boundBox with extra functionality for use in octree.
bool overlaps(const boundBox &bb) const
Overlaps with other bounding box, sphere etc?
Definition boundBoxI.H:439
static void calculate(const polyMesh &src, const polyMesh &tgt, labelList &srcToTgtAddr)
Calculate addressing.
Definition waveMethod.C:38
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const auto & io
const bitSet isBlockedFace(intersectedFaces())
#define DebugInfo
Report an information message using Foam::Info.
#define FUNCTION_NAME
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
List< scalarList > scalarListList
List of scalarList.
Definition scalarList.H:35
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
const dimensionSet dimless
Dimensionless.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
Vector< label > labelVector
Vector of labels.
Definition labelVector.H:47
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:616
constexpr label labelMin
Definition label.H:54
GeometricField< scalar, fvPatchField, volMesh > volScalarField
UIndirectList< label > labelUIndList
UIndirectList of labels.
IOList< label > labelIOList
IO for a List of label.
Definition labelIOList.H:32
List< treeBoundBox > treeBoundBoxList
A List of treeBoundBox.
messageStream Info
Information stream (stdout output on master, null elsewhere).
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition Ostream.H:490
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Ostream & indent(Ostream &os)
Indent stream.
Definition Ostream.H:481
dimensionedScalar sqrt(const dimensionedScalar &ds)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
List< point > pointList
List of point.
Definition pointList.H:32
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...
errorManip< error > abort(error &err)
Definition errorManip.H:139
uint8_t direction
Definition direction.H:49
vector point
Point is a vector.
Definition point.H:37
List< bool > boolList
A List of bools.
Definition List.H:60
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensionedScalar cbrt(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
labelList findIndices(const ListType &input, typename ListType::const_reference val, label start=0)
Linear search to find all occurrences of given element.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition Ostream.H:499
Type gMax(const FieldField< Field, Type > &f)
List< pointList > pointListList
List of pointList.
Definition pointList.H:35
static tmp< volScalarField > createField(const fvMesh &mesh, const scalar val)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
wordList patchTypes(nPatches)
dictionary dict
List< treeBoundBox > meshBb(1, treeBoundBox(coarseMesh.points()).extend(rndGen, 1e-3))
volScalarField & e
Us boundaryFieldRef().evaluateCoupled< coupledFaPatch >()
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Functor to negate primitives. Dummy for most other types.
Definition flipOp.H:67
const label nSamples(pdfDictionary.get< label >("nSamples"))
scalarField samples(nIntervals, Zero)