trackingInverseDistanceCellCellStencil.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#include "globalIndex.H"
34#include "oversetFvPatch.H"
36#include "syncTools.H"
37#include "treeBoundBoxList.H"
38#include "voxelMeshSearch.H"
40
41// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42
43namespace Foam
44{
45namespace cellCellStencils
46{
49}
50}
51
52// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53
55(
56 const fvMesh& mesh,
57 const vector& smallVec,
58
59 const boundBox& bb,
60 labelVector& nDivs,
62
63 const labelList& cellMap,
64 labelList& patchCellTypes
65)
66{
67 // Mark all voxels that overlap the bounding box of any patch
68
69 static bool hasWarned = false;
70
71 const fvBoundaryMesh& pbm = mesh.boundary();
72
74
75 // Mark wall boundaries
76 forAll(pbm, patchi)
77 {
78 const fvPatch& fvp = pbm[patchi];
79 const labelList& fc = fvp.faceCells();
80
81 if (!fvPatch::constraintType(fvp.type()))
82 {
83 //Info<< "Marking cells on proper patch " << fvp.name()
84 // << " with type " << fvp.type() << endl;
85 const polyPatch& pp = fvp.patch();
86 forAll(pp, i)
87 {
88 // Mark in overall patch types
89 patchCellTypes[cellMap[fc[i]]] = patchCellType::PATCH;
90
91 // Mark in voxel mesh
92 boundBox faceBb(pp.points(), pp[i], false);
93 faceBb.min() -= smallVec;
94 faceBb.max() += smallVec;
95
96 if (bb.overlaps(faceBb))
97 {
99 (
101 bb,
102 nDivs,
103 faceBb,
105 );
106 }
107 }
108 }
109 }
110
111 // Override with overset boundaries
112 forAll(pbm, patchi)
113 {
114 const fvPatch& fvp = pbm[patchi];
115 const labelList& fc = fvp.faceCells();
116
117 if (isA<oversetFvPatch>(fvp))
118 {
119 //Info<< "Marking cells on overset patch " << fvp.name() << endl;
120 const polyPatch& pp = fvp.patch();
121 const vectorField::subField faceCentres(pp.faceCentres());
122 forAll(pp, i)
123 {
124 // Mark in overall patch types
125 patchCellTypes[cellMap[fc[i]]] = patchCellType::OVERSET;
126
127 // Mark in voxel mesh
128 boundBox faceBb(pp.points(), pp[i], false);
129 faceBb.min() -= smallVec;
130 faceBb.max() += smallVec;
131 if (!bb.contains(faceCentres[i]))
132 {
133 if (!hasWarned)
134 {
135 hasWarned = true;
136 WarningInFunction << "Patch " << pp.name()
137 << " : face at " << faceCentres[i]
138 << " is not inside search bounding box of"
139 << " voxel mesh " << bb << endl
140 << " Is your 'searchBox' specification correct?"
141 << endl;
142 }
143 }
144 else
145 {
146 // Test for having voxels already marked as patch
147 // -> voxel mesh is too coarse
148 if
149 (
151 (
152 bb,
153 nDivs,
154 faceBb,
156 static_cast<unsigned int>(patchCellType::PATCH)
157 )
158 )
159 {
160 // Determine number of voxels from number of cells
161 // in mesh
162 const labelVector& dim = mesh.geometricD();
163 forAll(dim, i)
164 {
165 if (dim[i] != -1)
166 {
167 nDivs[i] *= 2;
168 }
169 }
170
171 Pout<< "Voxel mesh too coarse. Bounding box "
172 << faceBb
173 << " contains both non-overset and overset patches"
174 << ". Refining voxel mesh to " << nDivs << endl;
175
176 return false;
177 }
178
180 (
182 bb,
183 nDivs,
184 faceBb,
186 );
187 }
188 }
189 }
190 }
191
192 return true;
193}
194
195
197(
198 PstreamBuffers& pBufs,
199 const List<treeBoundBoxList>& patchBb,
200 const List<labelVector>& patchDivisions,
201 const PtrList<PackedList<2>>& patchParts,
202
203 const label srcI,
204 const label tgtI,
205 labelList& allCellTypes
206) const
207{
208 const pointField& allPoints = mesh_.points();
209 const labelListList& allCellPoints = mesh_.cellPoints();
210
211 const treeBoundBoxList& srcPatchBbs = patchBb[srcI];
212 const treeBoundBoxList& tgtPatchBbs = patchBb[tgtI];
213 const labelList& tgtCellMap = meshParts_[tgtI].cellMap();
214
215 // 1. do processor-local src-tgt patch overlap
216 {
217 const treeBoundBox& srcPatchBb = srcPatchBbs[Pstream::myProcNo()];
218 const treeBoundBox& tgtPatchBb = tgtPatchBbs[Pstream::myProcNo()];
219
220 if (srcPatchBb.overlaps(tgtPatchBb))
221 {
222 const PackedList<2>& srcPatchTypes = patchParts[srcI];
223 const labelVector& srcDivs = patchDivisions[srcI];
224
225 forAll(tgtCellMap, tgtCelli)
226 {
227 label celli = tgtCellMap[tgtCelli];
228 boundBox cBb(allPoints, allCellPoints[celli], false);
229 cBb.min() -= smallVec_;
230 cBb.max() += smallVec_;
231
232 if
233 (
235 (
236 srcPatchBb,
237 srcDivs,
238 cBb,
239 srcPatchTypes,
240 static_cast<unsigned int>(patchCellType::PATCH)
241 )
242 )
243 {
244 allCellTypes[celli] = HOLE;
245 }
246 }
247 }
248 }
249
250
251 // 2. Send over srcMesh bits that overlap tgt and do calculation
252 pBufs.clear();
253 for (const int procI : Pstream::allProcs())
254 {
255 if (procI != Pstream::myProcNo())
256 {
257 const treeBoundBox& srcPatchBb = srcPatchBbs[Pstream::myProcNo()];
258 const treeBoundBox& tgtPatchBb = tgtPatchBbs[procI];
259
260 if (srcPatchBb.overlaps(tgtPatchBb))
261 {
262 // Send over complete patch voxel map. Tbd: could
263 // subset
264 UOPstream os(procI, pBufs);
265 os << srcPatchBb << patchDivisions[srcI] << patchParts[srcI];
266 }
267 }
268 }
269 pBufs.finishedSends();
270 for (const int procI : Pstream::allProcs())
271 {
272 if (procI != Pstream::myProcNo())
273 {
274 const treeBoundBox& srcPatchBb = srcPatchBbs[procI];
275 const treeBoundBox& tgtPatchBb = tgtPatchBbs[Pstream::myProcNo()];
276
277 if (srcPatchBb.overlaps(tgtPatchBb))
278 {
279 UIPstream is(procI, pBufs);
280 {
281 treeBoundBox receivedBb(is);
282 if (srcPatchBb != receivedBb)
283 {
285 << "proc:" << procI
286 << " srcPatchBb:" << srcPatchBb
287 << " receivedBb:" << receivedBb
288 << exit(FatalError);
289 }
290 }
291 const labelVector srcDivs(is);
292 const PackedList<2> srcPatchTypes(is);
293
294 forAll(tgtCellMap, tgtCelli)
295 {
296 label celli = tgtCellMap[tgtCelli];
297 boundBox cBb(allPoints, allCellPoints[celli], false);
298 cBb.min() -= smallVec_;
299 cBb.max() += smallVec_;
300
301 if
302 (
304 (
305 srcPatchBb,
306 srcDivs,
307 cBb,
308 srcPatchTypes,
309 static_cast<unsigned int>(patchCellType::PATCH)
310 )
311 )
312 {
313 allCellTypes[celli] = HOLE;
314 }
315 }
316 }
317 }
318 }
319}
320
321
323(
324 PstreamBuffers& pBufs,
326 const PtrList<voxelMeshSearch>& meshSearches,
327 const labelList& allCellTypes,
328
329 const label srcI,
330 const label tgtI,
331 labelListList& allStencil,
332 labelList& allDonor
333) const
334{
335 const treeBoundBoxList& srcBbs = meshBb[srcI];
336 const treeBoundBoxList& tgtBbs = meshBb[tgtI];
337
338 const fvMesh& srcMesh = meshParts_[srcI].subMesh();
339 const labelList& srcCellMap = meshParts_[srcI].cellMap();
340 const voxelMeshSearch& meshSearch = meshSearches[srcI];
341 const fvMesh& tgtMesh = meshParts_[tgtI].subMesh();
342 const pointField& tgtCc = tgtMesh.cellCentres();
343 const labelList& tgtCellMap = meshParts_[tgtI].cellMap();
344
345 // 1. do processor-local src/tgt overlap
346 {
347 forAll(tgtCellMap, tgtCelli)
348 {
349 label srcCelli = meshSearch.findCell(tgtCc[tgtCelli]);
350 if (srcCelli != -1 && allCellTypes[srcCellMap[srcCelli]] != HOLE)
351 {
352 label celli = tgtCellMap[tgtCelli];
353
354 // TBD: check for multiple donors. Maybe better one? For
355 // now check 'nearer' mesh
356 if (betterDonor(tgtI, allDonor[celli], srcI))
357 {
358 allStencil[celli].setSize(1);
359 allStencil[celli][0] =
360 globalCells_.toGlobal(srcCellMap[srcCelli]);
361 allDonor[celli] = srcI;
362 }
363 }
364 }
365 }
366
367
368 // 2. Send over tgtMesh bits that overlap src and do calculation on
369 // srcMesh.
370
371
372 // (remote) processors where the tgt overlaps my src
373 DynamicList<label> tgtOverlapProcs(Pstream::nProcs());
374 // (remote) processors where the src overlaps my tgt
375 DynamicList<label> srcOverlapProcs(Pstream::nProcs());
376 for (const int procI : Pstream::allProcs())
377 {
378 if (procI != Pstream::myProcNo())
379 {
380 if (tgtBbs[procI].overlaps(srcBbs[Pstream::myProcNo()]))
381 {
382 tgtOverlapProcs.append(procI);
383 }
384 if (srcBbs[procI].overlaps(tgtBbs[Pstream::myProcNo()]))
385 {
386 srcOverlapProcs.append(procI);
387 }
388 }
389 }
390
391
392
393 // Indices of tgtcells to send over to each processor
395 forAll(srcOverlapProcs, i)
396 {
397 label procI = srcOverlapProcs[i];
398 tgtSendCells[procI].reserve(tgtMesh.nCells()/srcOverlapProcs.size());
399 }
400
401
402 forAll(tgtCellMap, tgtCelli)
403 {
404 label celli = tgtCellMap[tgtCelli];
405 if (srcOverlapProcs.size())
406 {
407 treeBoundBox subBb(cellBb(mesh_, celli));
408 subBb.min() -= smallVec_;
409 subBb.max() += smallVec_;
410
411 forAll(srcOverlapProcs, i)
412 {
413 label procI = srcOverlapProcs[i];
414 if (subBb.overlaps(srcBbs[procI]))
415 {
416 tgtSendCells[procI].append(tgtCelli);
417 }
418 }
419 }
420 }
421
422 // Send target cell centres to overlapping processors
423 pBufs.clear();
424
425 forAll(srcOverlapProcs, i)
426 {
427 label procI = srcOverlapProcs[i];
428 const labelList& cellIDs = tgtSendCells[procI];
429
430 UOPstream os(procI, pBufs);
431 os << UIndirectList<point>(tgtCc, cellIDs);
432 }
433 pBufs.finishedSends();
434
435 // Receive bits of target processors; find; send back
436 (void)srcMesh.tetBasePtIs();
437 forAll(tgtOverlapProcs, i)
438 {
439 label procI = tgtOverlapProcs[i];
440
441 UIPstream is(procI, pBufs);
442 pointList samples(is);
443
444 labelList donors(samples.size(), -1);
445 forAll(samples, sampleI)
446 {
447 label srcCelli = meshSearch.findCell(samples[sampleI]);
448 if (srcCelli != -1 && allCellTypes[srcCellMap[srcCelli]] != HOLE)
449 {
450 donors[sampleI] = globalCells_.toGlobal(srcCellMap[srcCelli]);
451 }
452 }
453
454 // Use same pStreamBuffers to send back.
455 UOPstream os(procI, pBufs);
456 os << donors;
457 }
458 pBufs.finishedSends();
459
460 forAll(srcOverlapProcs, i)
461 {
462 label procI = srcOverlapProcs[i];
463 const labelList& cellIDs = tgtSendCells[procI];
464
465 UIPstream is(procI, pBufs);
466 labelList donors(is);
467
468 if (donors.size() != cellIDs.size())
469 {
470 FatalErrorInFunction<< "problem : cellIDs:" << cellIDs.size()
471 << " donors:" << donors.size() << abort(FatalError);
472 }
473
474 forAll(donors, donorI)
475 {
476 label globalDonor = donors[donorI];
477
478 if (globalDonor != -1)
479 {
480 label celli = tgtCellMap[cellIDs[donorI]];
481
482 // TBD: check for multiple donors. Maybe better one?
483 if (betterDonor(tgtI, allDonor[celli], srcI))
484 {
485 allStencil[celli].setSize(1);
486 allStencil[celli][0] = globalDonor;
487 allDonor[celli] = srcI;
488 }
489 }
490 }
491 }
492}
493
494
495// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
496
498(
499 const fvMesh& mesh,
500 const dictionary& dict,
501 const bool doUpdate
502)
503:
504 inverseDistance(mesh, dict, false),
505 globalCells_(mesh_.nCells())
506{
507 // Initialise donor cell
509 globalDonor_ = -1;
510
511 // Initialise mesh partitions
512 const labelIOList& zoneID = this->zoneID();
513 label nZones = gMax(zoneID)+1;
514
516 forAll(zoneID, celli)
517 {
518 nCellsPerZone[zoneID[celli]]++;
519 }
521
522 meshParts_.setSize(nZones);
523 forAll(meshParts_, zonei)
524 {
525 meshParts_.set
526 (
527 zonei,
528 new fvMeshSubset(mesh_, zonei, zoneID)
529 );
530
531 // Trigger early evaluation of mesh dimension
532 // (in case there are locally zero cells in mesh)
533 (void)meshParts_[zonei].subMesh().nGeometricD();
534 }
535
536
537 // Print a bit
538 {
539 Info<< typeName << " : detected " << nZones
540 << " mesh regions" << endl;
542 forAll(nCellsPerZone, zonei)
543 {
544 Info<< indent<< "zone:" << zonei
545 << " nCells:" << nCellsPerZone[zonei]
546 << endl;
547 }
549 }
550
551
552 // Do geometry update
553 if (doUpdate)
554 {
555 update();
556 }
557}
558
559
560// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
561
563{}
564
565
566// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
567
569{
570 DebugInfo<< FUNCTION_NAME << " : Start of analysis" << endl;
571
572 scalar layerRelax(dict_.getOrDefault("layerRelax", 1.0));
573 const labelIOList& zoneID = this->zoneID();
574 label nZones = meshParts_.size();
575
576 // Update stored mesh partitions for geometry changes
577 forAll(meshParts_, zonei)
578 {
579 pointField subPoints(mesh_.points(), meshParts_[zonei].pointMap());
580
581 fvMesh& subMesh = meshParts_[zonei].subMesh();
582 subMesh.movePoints(subPoints);
583 }
584
585 DebugInfo<< FUNCTION_NAME << " : Moved zone sub-meshes" << endl;
586
587 // Calculate fast search structure for each zone
588 List<labelVector> searchBoxDivisions(dict_.lookup("searchBoxDivisions"));
589 if (searchBoxDivisions.size() != nZones)
590 {
591 FatalIOErrorInFunction(dict_) << "Number of searchBoxDivisions "
592 << searchBoxDivisions.size()
593 << " should equal the number of zones " << nZones
594 << exit(FatalIOError);
595 }
596
597 const boundBox& allBb(mesh_.bounds());
598
600
601 // Determine zone meshes and bounding boxes
602 {
603 // Per processor, per zone the bounding box
606
607 forAll(meshParts_, zonei)
608 {
609 const fvMesh& subMesh = meshParts_[zonei].subMesh();
610
611 if (subMesh.nPoints())
612 {
613 procBb[Pstream::myProcNo()][zonei] =
614 treeBoundBox(subMesh.points());
615 procBb[Pstream::myProcNo()][zonei].inflate(1e-6);
616 }
617 else
618 {
619 // No part of zone on this processor. Make up bb.
620 procBb[Pstream::myProcNo()][zonei] = treeBoundBox
621 (
622 allBb.min() - 2*allBb.span(),
623 allBb.min() - allBb.span()
624 );
625 procBb[Pstream::myProcNo()][zonei].inflate(1e-6);
626 }
627 }
628
630
631 // Move local bounding boxes to per-mesh indexing
632 forAll(meshBb, zonei)
633 {
634 treeBoundBoxList& bbs = meshBb[zonei];
636 forAll(procBb, proci)
637 {
638 bbs[proci] = procBb[proci][zonei];
639 }
640 }
641 }
642
643 DebugInfo<< FUNCTION_NAME << " : Calculated bounding boxes" << endl;
644
645
646 // Determine patch bounding boxes. These are either global and provided
647 // by the user or processor-local as a copy of the mesh bounding box.
648
650 List<labelVector> patchDivisions(searchBoxDivisions);
651 PtrList<PackedList<2>> patchParts(nZones);
652 labelList allPatchTypes(mesh_.nCells(), OTHER);
653
654 {
655 treeBoundBox globalPatchBb;
656 if (dict_.readIfPresent("searchBox", globalPatchBb))
657 {
658 // All processors, all zones have the same bounding box
659 patchBb = treeBoundBoxList(Pstream::nProcs(), globalPatchBb);
660 }
661 else
662 {
663 // Use the meshBb (differing per zone, per processor)
664 patchBb = meshBb;
665 }
666 }
667
668 forAll(patchParts, zonei)
669 {
670 while (true)
671 {
672 patchParts.set
673 (
674 zonei,
675 new PackedList<2>(cmptProduct(patchDivisions[zonei]))
676 );
677
678 bool ok = markBoundaries
679 (
680 meshParts_[zonei].subMesh(),
681 smallVec_,
682
683 patchBb[zonei][Pstream::myProcNo()],
684 patchDivisions[zonei],
685 patchParts[zonei],
686
687 meshParts_[zonei].cellMap(),
688 allPatchTypes
689 );
690
691 //if (returnReduce(ok, andOp<bool>()))
692 if (ok)
693 {
694 break;
695 }
696 }
697 }
698 DebugInfo<< FUNCTION_NAME << " : Calculated boundary voxel meshes" << endl;
699
700
701 PtrList<voxelMeshSearch> meshSearches(meshParts_.size());
702 forAll(meshParts_, zonei)
703 {
704 meshSearches.set
705 (
706 zonei,
708 (
709 meshParts_[zonei].subMesh(),
710 patchBb[zonei][Pstream::myProcNo()],
711 patchDivisions[zonei],
712 true
713 )
714 );
715 }
716
717 DebugInfo<< FUNCTION_NAME << " : Constructed cell voxel meshes" << endl;
718
719 PtrList<fvMesh> voxelMeshes;
720 if (debug)
721 {
722 voxelMeshes.setSize(meshSearches.size());
723 forAll(meshSearches, zonei)
724 {
725 IOobject meshIO
726 (
727 word("voxelMesh" + Foam::name(zonei)),
728 mesh_.time().timeName(),
729 mesh_,
731 );
732
733 Pout<< "Constructing voxel mesh " << meshIO.objectPath() << endl;
734 voxelMeshes.set(zonei, meshSearches[zonei].makeMesh(meshIO));
735 }
736 }
737
738
739 if (debug)
740 {
741 forAll(patchParts, zonei)
742 {
743 const labelList voxelTypes(patchParts[zonei].values());
744 const fvMesh& vm = voxelMeshes[zonei];
746 (
747 createField<label>
748 (
749 vm,
750 "patchTypes",
751 voxelTypes
752 )
753 );
754 tfld().write();
755 }
756 }
757
758 // Current best guess for cells
759 labelList allCellTypes(mesh_.nCells(), CALCULATED);
760 labelListList allStencil(mesh_.nCells());
761 // zoneID of donor
762 labelList allDonorID(mesh_.nCells(), -1);
763
764 const globalIndex globalCells(mesh_.nCells());
765
767
768 DebugInfo<< FUNCTION_NAME << " : Allocated donor-cell structures" << endl;
769
770 for (label srci = 0; srci < meshParts_.size()-1; srci++)
771 {
772 for (label tgti = srci+1; tgti < meshParts_.size(); tgti++)
773 {
774 markPatchesAsHoles
775 (
776 pBufs,
777
778 patchBb,
779 patchDivisions,
780 patchParts,
781
782 srci,
783 tgti,
784 allCellTypes
785 );
786 markPatchesAsHoles
787 (
788 pBufs,
789
790 patchBb,
791 patchDivisions,
792 patchParts,
793
794 tgti,
795 srci,
796 allCellTypes
797 );
798 }
799 }
800
801 for (label srci = 0; srci < meshParts_.size()-1; srci++)
802 {
803 for (label tgti = srci+1; tgti < meshParts_.size(); tgti++)
804 {
805 markDonors
806 (
807 pBufs,
808 meshBb,
809 meshSearches,
810 allCellTypes, // to exclude hole donors
811
812 tgti,
813 srci,
814 allStencil,
815 allDonorID
816 );
817 markDonors
818 (
819 pBufs,
820 meshBb,
821 meshSearches,
822 allCellTypes, // to exclude hole donors
823
824 srci,
825 tgti,
826 allStencil,
827 allDonorID
828 );
829 }
830 }
831
832 DebugInfo<< FUNCTION_NAME << " : Determined holes and donor-acceptors"
833 << endl;
834
835 if (debug)
836 {
838 (
839 createField(mesh_, "allCellTypes", allCellTypes)
840 );
841 tfld().write();
842 }
843
844
845 // Use the patch types and weights to decide what to do
846 forAll(allPatchTypes, celli)
847 {
848 if (allCellTypes[celli] != HOLE)
849 {
850 switch (allPatchTypes[celli])
851 {
852 case OVERSET:
853 {
854 // Require interpolation. See if possible.
855
856 if (allStencil[celli].size())
857 {
858 allCellTypes[celli] = INTERPOLATED;
859 }
860 else
861 {
862 allCellTypes[celli] = HOLE;
863 }
864 }
865 }
866 }
867 }
868 DebugInfo<< FUNCTION_NAME << " : Removed bad donors" << endl;
869
870 if (debug)
871 {
873 (
874 createField(mesh_, "allCellTypes_patch", allCellTypes)
875 );
876 tfld().write();
877 }
878
879
880 // Check previous iteration cellTypes_ for any hole->calculated changes
881 // If so set the cell either to interpolated (if there are donors) or
882 // holes (if there are no donors). Note that any interpolated cell might
883 // still be overwritten by the flood filling
884 {
885 label nCalculated = 0;
886
887 forAll(cellTypes_, celli)
888 {
889 if (allCellTypes[celli] == CALCULATED && cellTypes_[celli] == HOLE)
890 {
891 if (allStencil[celli].size() == 0)
892 {
893 // Reset to hole
894 allCellTypes[celli] = HOLE;
895 allStencil[celli].clear();
896 }
897 else
898 {
899 allCellTypes[celli] = INTERPOLATED;
900 nCalculated++;
901 }
902 }
903 }
904
905 if (debug)
906 {
907 Pout<< "Detected " << nCalculated << " cells changing from hole"
908 << " to calculated. Changed to interpolated"
909 << endl;
910 }
911 }
912
913
914 // Mark unreachable bits
915 findHoles(globalCells_, mesh_, zoneID, allStencil, allCellTypes);
916 DebugInfo<< FUNCTION_NAME << " : Flood-filled holes" << endl;
917
918 if (debug)
919 {
921 (
922 createField(mesh_, "allCellTypes_hole", allCellTypes)
923 );
924 tfld().write();
925 }
926
927 // Add buffer interpolation layer(s) around holes
928 scalarField allWeight(mesh_.nCells(), Zero);
929 walkFront(layerRelax, allStencil, allCellTypes, allWeight);
930 DebugInfo<< FUNCTION_NAME << " : Implemented layer relaxation" << endl;
931
932 if (debug)
933 {
935 (
936 createField(mesh_, "allCellTypes_front", allCellTypes)
937 );
938 tfld().write();
939 }
940
941
942 // Convert cell-cell addressing to stencil in compact notation
943
944 cellTypes_.transfer(allCellTypes);
945 cellStencil_.setSize(mesh_.nCells());
946 cellInterpolationWeights_.setSize(mesh_.nCells());
947 DynamicList<label> interpolationCells;
948 forAll(cellTypes_, celli)
949 {
950 if (cellTypes_[celli] == INTERPOLATED)
951 {
952 cellStencil_[celli].transfer(allStencil[celli]);
953 cellInterpolationWeights_[celli].setSize(1);
954 cellInterpolationWeights_[celli][0] = 1.0;
955 interpolationCells.append(celli);
956 }
957 else
958 {
959 cellStencil_[celli].clear();
960 cellInterpolationWeights_[celli].clear();
961 }
962 }
963 interpolationCells_.transfer(interpolationCells);
964
965 List<Map<label>> compactMap;
966 cellInterpolationMap_.reset
967 (
968 new mapDistribute
969 (
970 globalCells,
971 cellStencil_,
972 compactMap
973 )
974 );
975 cellInterpolationWeight_.transfer(allWeight);
976 //cellInterpolationWeight_.correctBoundaryConditions();
978 <
981 >(cellInterpolationWeight_.boundaryFieldRef(), false);
982
983
984 if (debug & 2)
985 {
986 // Dump stencil
987 mkDir(mesh_.time().timePath());
988 OBJstream str(mesh_.time().timePath()/"injectionStencil.obj");
989 Pout<< typeName << " : dumping injectionStencil to "
990 << str.name() << endl;
991 pointField cc(mesh_.cellCentres());
992 cellInterpolationMap().distribute(cc);
993
994 forAll(cellStencil_, celli)
995 {
996 const labelList& slots = cellStencil_[celli];
997 if (slots.size())
998 {
999 const point& accCc = mesh_.cellCentres()[celli];
1000 forAll(slots, i)
1001 {
1002 const point& donorCc = cc[slots[i]];
1003 str.write(linePointRef(accCc, 0.1*accCc+0.9*donorCc));
1004 }
1005 }
1006 }
1007 }
1008
1009 DebugInfo<< FUNCTION_NAME << " : Transferred donor to stencil" << endl;
1010
1011
1012 // Extend stencil to get inverse distance weighted neighbours
1013 createStencil(globalCells);
1014 DebugInfo<< FUNCTION_NAME << " : Extended stencil" << endl;
1015
1016
1017 if (debug & 2)
1018 {
1019 // Dump weight
1020 cellInterpolationWeight_.instance() = mesh_.time().timeName();
1021 cellInterpolationWeight_.write();
1022
1023 // Dump cell types
1024 volScalarField volTypes
1025 (
1026 IOobject
1027 (
1028 "cellTypes",
1029 mesh_.time().timeName(),
1030 mesh_,
1033 false
1034 ),
1035 mesh_,
1037 zeroGradientFvPatchScalarField::typeName
1038 );
1039
1040 forAll(volTypes.internalField(), celli)
1041 {
1042 volTypes[celli] = cellTypes_[celli];
1043 }
1044 //volTypes.correctBoundaryConditions();
1046 <
1049 >(volTypes.boundaryFieldRef(), false);
1050 volTypes.write();
1051
1052 // Dump stencil
1053 mkDir(mesh_.time().timePath());
1054 OBJstream str(mesh_.time().timePath()/"stencil.obj");
1055 Pout<< typeName << " : dumping to " << str.name() << endl;
1056 pointField cc(mesh_.cellCentres());
1057 cellInterpolationMap().distribute(cc);
1058
1059 forAll(cellStencil_, celli)
1060 {
1061 const labelList& slots = cellStencil_[celli];
1062 if (slots.size())
1063 {
1064 const point& accCc = mesh_.cellCentres()[celli];
1065 forAll(slots, i)
1066 {
1067 const point& donorCc = cc[slots[i]];
1068 str.write(linePointRef(accCc, 0.1*accCc+0.9*donorCc));
1069 }
1070 }
1071 }
1072 }
1073
1074
1075 // Print some stats
1076 {
1077 labelList nCells(count(3, cellTypes_));
1078
1079 label nLocal = 0;
1080 label nMixed = 0;
1081 label nRemote = 0;
1082 forAll(interpolationCells_, i)
1083 {
1084 label celli = interpolationCells_[i];
1085 const labelList& slots = cellStencil_[celli];
1086
1087 bool hasLocal = false;
1088 bool hasRemote = false;
1089
1090 forAll(slots, sloti)
1091 {
1092 if (slots[sloti] >= mesh_.nCells())
1093 {
1094 hasRemote = true;
1095 }
1096 else
1097 {
1098 hasLocal = true;
1099 }
1100 }
1101
1102 if (hasRemote)
1103 {
1104 if (!hasLocal)
1105 {
1106 nRemote++;
1107 }
1108 else
1109 {
1110 nMixed++;
1111 }
1112 }
1113 else if (hasLocal)
1114 {
1115 nLocal++;
1116 }
1117 }
1118 reduce(nLocal, sumOp<label>());
1119 reduce(nMixed, sumOp<label>());
1120 reduce(nRemote, sumOp<label>());
1121
1122 Info<< "Overset analysis : nCells : "
1123 << returnReduce(cellTypes_.size(), sumOp<label>()) << nl
1124 << incrIndent
1125 << indent << "calculated : " << nCells[CALCULATED] << nl
1126 << indent << "interpolated : " << nCells[INTERPOLATED]
1127 << " (interpolated from local:" << nLocal
1128 << " mixed local/remote:" << nMixed
1129 << " remote:" << nRemote << ")" << nl
1130 << indent << "hole : " << nCells[HOLE] << nl
1131 << decrIndent << endl;
1132 }
1133 DebugInfo<< FUNCTION_NAME << " : Finished analysis" << endl;
1134
1135 // Tbd: detect if anything changed. Most likely it did!
1136 return true;
1137}
1138
1139
1140// ************************************************************************* //
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
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
const Internal & internalField() const
Return a const-reference to the dimensioned internal field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
fileName objectPath() const
The complete path + object name.
Definition: IOobjectI.H:214
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
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.
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
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
void setSize(const label newLen)
Same as resize()
Definition: PtrList.H:151
SubField is a Field obtained as a section of another Field, without its own allocation....
Definition: SubField.H:62
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
@ nonBlocking
"nonBlocking"
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
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
bool contains(const point &pt) const
Contains point? (inside or on edge)
Definition: boundBoxI.H:271
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.
Inverse-distance-weighted interpolation stencil.
Inverse-distance-weighted interpolation stencil.
PtrList< fvMeshSubset > meshParts_
Subset according to zone.
static bool markBoundaries(const fvMesh &mesh, const vector &smallVec, const boundBox &bb, labelVector &nDivs, PackedList< 2 > &patchTypes, const labelList &cellMap, labelList &patchCellTypes)
Mark voxels of patchTypes with type of patch face.
void markPatchesAsHoles(PstreamBuffers &pBufs, 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 bool update()
Update stencils. Return false if nothing changed.
void markDonors(PstreamBuffers &pBufs, const List< treeBoundBoxList > &meshBb, const PtrList< voxelMeshSearch > &meshSearches, const labelList &allCellTypes, const label srcI, const label tgtI, labelListList &allStencil, labelList &allDonor) const
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
Foam::fvBoundaryMesh.
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 fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:712
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:895
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
Class containing processor-to-processor mapping information.
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:61
label findCell(const point &location, const label seedCelli=-1, const bool useTreeSearch=true) const
Find cell containing location.
Definition: meshSearch.C:780
Boundary condition for use on overset patches. To be run in combination with special dynamicFvMesh ty...
const word & name() const noexcept
The patch name.
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:906
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:872
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:321
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.
virtual bool write(const bool valid=true) const
Write using setting from DB.
A class for managing temporary objects.
Definition: tmp.H:65
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
Fast, non-parallel searching in mesh without use of octree.
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
dynamicFvMesh & mesh
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
label nZones
labelList nCellsPerZone(nZones, Zero)
const labelIOList & zoneID
#define DebugInfo
Report an information message using Foam::Info.
#define WarningInFunction
Report a warning using Foam::Warning.
#define FUNCTION_NAME
Namespace for OpenFOAM.
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
messageStream Info
Information stream (stdout output on master, null elsewhere)
Cmpt cmptProduct(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H: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
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
List< treeBoundBox > treeBoundBoxList
List of bounding boxes.
errorManip< error > abort(error &err)
Definition: errorManip.H:144
IOerror FatalIOError
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
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.
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)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
wordList patchTypes(nPatches)
List< treeBoundBox > meshBb(1, treeBoundBox(boundBox(coarseMesh.points(), false)).extend(rndGen, 1e-3))
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
scalarField samples(nIntervals, Zero)