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