Loading...
Searching...
No Matches
cellVolumeWeightCellCellStencil.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) 2014-2025 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify i
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
30#include "OBJstream.H"
31#include "Time.H"
32#include "meshToMesh.H"
34#include "fvMeshSubset.H"
35#include "regionSplit.H"
36#include "globalIndex.H"
37#include "oversetFvPatch.H"
39#include "syncTools.H"
41#include "oversetFvPatchField.H"
43// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44
45namespace Foam
46{
47namespace cellCellStencils
48{
51}
52}
54Foam::scalar
56
57
58// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
59
61(
62 const globalIndex& globalCells,
63 const fvMesh& mesh,
64 const labelList& zoneID,
65 const labelListList& stencil,
67) const
68{
69 const fvBoundaryMesh& pbm = mesh.boundary();
70 const labelList& own = mesh.faceOwner();
71 const labelList& nei = mesh.faceNeighbour();
72
73 // so we start a walk from our patches and any cell we cannot reach
74 // (because we walk is stopped by other-mesh patch) is a hole.
75
76
77 boolList isBlockedFace(mesh.nFaces(), false);
78 for (label faceI = 0; faceI < mesh.nInternalFaces(); faceI++)
79 {
80 label ownType = cellTypes[own[faceI]];
81 label neiType = cellTypes[nei[faceI]];
82 if
83 (
84 (ownType == HOLE && neiType != HOLE)
85 || (ownType != HOLE && neiType == HOLE)
86 )
87 {
88 isBlockedFace[faceI] = true;
89 }
90 }
91
92 labelList nbrCellTypes;
94
95 for (label faceI = mesh.nInternalFaces(); faceI < mesh.nFaces(); faceI++)
96 {
97 label ownType = cellTypes[own[faceI]];
98 label neiType = nbrCellTypes[faceI-mesh.nInternalFaces()];
99
100 if
101 (
102 (ownType == HOLE && neiType != HOLE)
103 || (ownType != HOLE && neiType == HOLE)
104 )
105 {
106 isBlockedFace[faceI] = true;
107 }
108 }
109
110 regionSplit cellRegion(mesh, isBlockedFace);
111
112 Info<< typeName << " : detected " << cellRegion.nRegions()
113 << " mesh regions after overset" << nl << endl;
114
115
116
117 // Now we'll have a mesh split according to where there are cells
118 // covered by the other-side patches. See what we can reach from our
119 // real patches
120
121 // 0 : region not yet determined
122 // 1 : borders blockage so is not ok (but can be overridden by real
123 // patch)
124 // 2 : has real patch in it so is reachable
125 labelList regionType(cellRegion.nRegions(), Zero);
126
127
128 // See if any regions borders blockage. Note: isBlockedFace is already
129 // parallel synchronised.
130 {
131 for (label faceI = 0; faceI < mesh.nInternalFaces(); faceI++)
132 {
133 if (isBlockedFace[faceI])
134 {
135 label ownRegion = cellRegion[own[faceI]];
136
137 if (cellTypes[own[faceI]] != HOLE)
138 {
139 if (regionType[ownRegion] == 0)
140 {
141 //Pout<< "Mark region:" << ownRegion
142 // << " on zone:" << zoneID[own[faceI]]
143 // << " as next to blockage at:"
144 // << mesh.faceCentres()[faceI] << endl;
145 regionType[ownRegion] = 1;
146 }
147 }
148
149 label neiRegion = cellRegion[nei[faceI]];
150
151 if (cellTypes[nei[faceI]] != HOLE)
152 {
153 if (regionType[neiRegion] == 0)
154 {
155 //Pout<< "Mark region:" << neiRegion
156 // << " on zone:" << zoneID[nei[faceI]]
157 // << " as next to blockage at:"
158 // << mesh.faceCentres()[faceI] << endl;
159 regionType[neiRegion] = 1;
160 }
161 }
162 }
163 }
164 for
165 (
166 label faceI = mesh.nInternalFaces();
167 faceI < mesh.nFaces();
168 faceI++
169 )
170 {
171 if (isBlockedFace[faceI])
172 {
173 label ownRegion = cellRegion[own[faceI]];
174
175 if (regionType[ownRegion] == 0)
176 {
177 //Pout<< "Mark region:" << ownRegion
178 // << " on zone:" << zoneID[own[faceI]]
179 // << " as next to blockage at:"
180 // << mesh.faceCentres()[faceI] << endl;
181 regionType[ownRegion] = 1;
182 }
183 }
184 }
185 }
186
187
188 // Override with real patches
189 forAll(pbm, patchI)
190 {
191 const fvPatch& fvp = pbm[patchI];
192
193 if (isA<oversetFvPatch>(fvp))
194 {}
195 else if (!fvPatch::constraintType(fvp.type()))
196 {
197 //Pout<< "Proper patch " << fvp.name() << " of type " << fvp.type()
198 // << endl;
199
200 for (const label celli : fvp.faceCells())
201 {
202 label regionI = cellRegion[celli];
203
204 if (cellTypes[celli] != HOLE && regionType[regionI] != 2)
205 {
206 //Pout<< "reachable region : " << regionI
207 // << " at cell " << mesh.cellCentres()[celli]
208 // << " on zone " << zoneID[celli] << endl;
209 regionType[regionI] = 2;
210 }
211 }
212 }
213 }
214
215 // Now we've handled
216 // - cells next to blocked cells
217 // - coupled boundaries
218 // Only thing to handle is the interpolation between regions
219
220
221 labelListList compactStencil(stencil);
222 List<Map<label>> compactMap;
223 mapDistribute map(globalCells, compactStencil, compactMap);
224
225 while (true)
226 {
227 // Synchronise region status on processors
228 // (could instead swap status through processor patches)
229 Pstream::listReduce(regionType, maxOp<label>());
230
231 // Communicate region status through interpolative cells
232 labelList cellRegionType(labelUIndList(regionType, cellRegion));
233 map.distribute(cellRegionType);
234
235
236 label nChanged = 0;
237 forAll(pbm, patchI)
238 {
239 const fvPatch& fvp = pbm[patchI];
240
241 if (isA<oversetFvPatch>(fvp))
242 {
243 const labelUList& fc = fvp.faceCells();
244 forAll(fc, i)
245 {
246 label cellI = fc[i];
247 label regionI = cellRegion[cellI];
248
249 if (regionType[regionI] != 2)
250 {
251 const labelList& slots = compactStencil[cellI];
252 forAll(slots, i)
253 {
254 label otherType = cellRegionType[slots[i]];
255
256 if (otherType == 2)
257 {
258 //Pout<< "Reachable through interpolation : "
259 // << regionI << " at cell "
260 // << mesh.cellCentres()[cellI] << endl;
261 regionType[regionI] = 2;
262 nChanged++;
263 break;
264 }
265 }
266 }
267 }
268 }
269 }
270
271
272 if (!returnReduceOr(nChanged))
273 {
274 break;
275 }
276 }
277
278
279 // See which regions have not been visited (regionType == 1)
280 forAll(cellRegion, cellI)
281 {
282 label type = regionType[cellRegion[cellI]];
283 if (type == 1 && cellTypes[cellI] != HOLE)
285 cellTypes[cellI] = HOLE;
286 }
287 }
288}
289
290
292(
293 const fvMesh& mesh,
294 const labelList& cellMap,
295 labelList& patchCellTypes
296) const
297{
298 const fvBoundaryMesh& pbm = mesh.boundary();
299
300 forAll(pbm, patchI)
301 {
302 const fvPatch& fvp = pbm[patchI];
303 const labelUList& fc = fvp.faceCells();
304
305 if (isA<oversetFvPatch>(fvp))
306 {
307 //Pout<< "Marking cells on overset patch " << fvp.name() << endl;
308 forAll(fc, i)
309 {
310 label cellI = fc[i];
311 patchCellTypes[cellMap[cellI]] = OVERSET;
312 }
313 }
314 else if (!fvPatch::constraintType(fvp.type()))
315 {
316 //Pout<< "Marking cells on proper patch " << fvp.name()
317 // << " with type " << fvp.type() << endl;
318 forAll(fc, i)
319 {
320 label cellI = fc[i];
321 if (patchCellTypes[cellMap[cellI]] != OVERSET)
322 {
323 patchCellTypes[cellMap[cellI]] = PATCH;
325 }
326 }
327 }
328}
329
330
332(
333 const labelListList& addressing,
334 const labelList& patchTypes,
335 labelList& result
336) const
337{
338 forAll(result, cellI)
339 {
340 const labelList& slots = addressing[cellI];
341 forAll(slots, i)
342 {
343 label type = patchTypes[slots[i]];
344
345 if (type == OVERSET)
346 {
347 // 'overset' overrides anything
348 result[cellI] = OVERSET;
349 break;
350 }
351 else if (type == PATCH)
352 {
353 // 'patch' overrides -1 and 'other'
354 result[cellI] = PATCH;
355 break;
356 }
357 else if (result[cellI] == -1)
358 {
359 // 'other' overrides -1 only
360 result[cellI] = OTHER;
361 }
362 }
363 }
364}
365
366
368(
369 const autoPtr<mapDistribute>& mapPtr,
370 const labelListList& addressing,
371 const labelList& patchTypes,
372 labelList& result
373) const
374{
375 if (result.size() != addressing.size())
376 {
377 FatalErrorInFunction << "result:" << result.size()
378 << " addressing:" << addressing.size() << exit(FatalError);
379 }
380
381
382 // Initialise to not-mapped
383 result = -1;
384
385 if (mapPtr)
386 {
387 // Pull remote data into order of addressing
388 labelList work(patchTypes);
389 mapPtr().distribute(work);
390
391 interpolatePatchTypes(addressing, work, result);
392 }
393 else
394 {
395 interpolatePatchTypes(addressing, patchTypes, result);
396 }
397}
398
399
401(
402 const label subZoneID,
403 const fvMesh& subMesh,
404 const labelList& subCellMap,
405
406 const label donorZoneID,
407 const labelListList& addressing,
408 const List<scalarList>& weights,
409 const labelList& otherCells,
410 const labelList& interpolatedOtherPatchTypes,
411
412 labelListList& allStencil,
413 scalarListList& allWeights,
414 labelList& allCellTypes,
415 labelList& allDonorID
416) const
417{
418 forAll(subCellMap, subCellI)
419 {
420 label cellI = subCellMap[subCellI];
421
422 bool validDonors = true;
423 switch (interpolatedOtherPatchTypes[subCellI])
424 {
425 case -1:
426 {
427 validDonors = false;
428 }
429 break;
430
431 case OTHER:
432 {
433 // No patch interaction so keep valid
434 }
435 break;
436
437 case PATCH:
438 {
439 // Patch-patch interaction... For now disable always
440 allCellTypes[cellI] = HOLE;
441 validDonors = false;
442 // Alternative is to look at the amount of overlap but this
443 // is not very robust
444// if (allCellTypes[cellI] != HOLE)
445// {
446// scalar overlapVol = sum(weights[subCellI]);
447// scalar v = mesh_.V()[cellI];
448// if (overlapVol < (1.0-overlapTolerance_)*v)
449// {
450// //Pout<< "** Patch overlap:" << cellI
451// // << " at:" << mesh_.cellCentres()[cellI] << endl;
452// allCellTypes[cellI] = HOLE;
453// validDonors = false;
454// }
455// }
456 }
457 break;
458
459 case OVERSET:
460 {
461 validDonors = true;
462 }
463 break;
464 }
465
466
467 if (validDonors)
468 {
469 // There are a few possible choices how to choose between multiple
470 // donor candidates:
471 // 1 highest overlap volume. However this is generally already
472 // 99.9% so you're just measuring truncation error.
473 // 2 smallest donors cells or most donor cells. This is quite
474 // often done but can cause switching of donor zone from one
475 // time step to the other if the donor meshes are non-uniform
476 // and the acceptor cells just happens to be sweeping through
477 // some small donor cells.
478 // 3 nearest zoneID. So zone 0 preferentially interpolates from
479 // zone 1, zone 1 preferentially from zone 2 etc.
480
481 //- Option 1:
482 //scalar currentVol = sum(allWeights[cellI]);
483 //if (overlapVol[subCellI] > currentVol)
484
485 //- Option 3:
486 label currentDiff = mag(subZoneID-allDonorID[cellI]);
487 label thisDiff = mag(subZoneID-donorZoneID);
488
489 if
490 (
491 allDonorID[cellI] == -1
492 || (thisDiff < currentDiff)
493 || (thisDiff == currentDiff && donorZoneID > allDonorID[cellI])
494 )
495 {
496 allWeights[cellI] = weights[subCellI];
497 allStencil[cellI] =
498 labelUIndList(otherCells, addressing[subCellI]);
499 allDonorID[cellI] = donorZoneID;
500 }
501 }
502 }
503}
504
505
506// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
507
508Foam::cellCellStencils::cellVolumeWeight::cellVolumeWeight
509(
510 const fvMesh& mesh,
511 const dictionary& dict,
512 const bool doUpdate
513)
514:
515 cellCellStencil(mesh),
516 dict_(dict),
521 cellStencil_(0),
524 (
526 (
527 "cellInterpolationWeight",
528 mesh_.facesInstance(),
529 mesh_,
530 IOobject::NO_READ,
531 IOobject::NO_WRITE,
532 IOobject::NO_REGISTER
533 ),
534 mesh_,
536 fvPatchFieldBase::zeroGradientType()
537 ),
539 (
540 dict.getOrDefault("allowInterpolatedDonors", true)
541 )
542{
543 // Add motion-solver fields to non-interpolated fields
545
546 // Read zoneID
547 this->zoneID();
548
550 dict_.getOrDefault("overlapTolerance", defaultOverlapTolerance_);
551
552 // Read old-time cellTypes
554 (
555 "cellTypes",
556 mesh_.time().timeName(),
557 mesh_,
561 );
562 if (io.typeHeaderOk<volScalarField>(true))
563 {
564 if (debug)
565 {
566 Pout<< "Reading cellTypes from time " << mesh_.time().timeName()
567 << endl;
568 }
569
570 const volScalarField volCellTypes(io, mesh_);
571 forAll(volCellTypes, celli)
572 {
573 // Round to integer
574 cellTypes_[celli] = volCellTypes[celli];
575 }
576 }
577
578 if (doUpdate)
579 {
580 update();
581 }
582}
583
584
585// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
586
588{}
589
590
591// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
592
594{
595 scalar layerRelax(dict_.getOrDefault("layerRelax", 1.0));
596 const labelIOList& zoneID = this->zoneID();
597
598 label nZones = gMax(zoneID)+1;
599 labelList nCellsPerZone(nZones, Zero);
600 forAll(zoneID, cellI)
601 {
602 nCellsPerZone[zoneID[cellI]]++;
603 }
604 Pstream::listReduce(nCellsPerZone, sumOp<label>());
605
606 Info<< typeName << " : detected " << nZones
607 << " mesh regions" << nl << endl;
608
609
610 PtrList<fvMeshSubset> meshParts(nZones);
611
613 forAll(meshParts, zonei)
614 {
615 Info<< indent<< "zone:" << zonei << " nCells:"
616 << nCellsPerZone[zonei] << nl;
617
618 meshParts.set
619 (
620 zonei,
621 new fvMeshSubset(mesh_, zonei, zoneID)
622 );
623 }
625
626
627 // Current best guess for cells. Includes best stencil. Weights should
628 // add up to volume.
629 labelList allCellTypes(mesh_.nCells(), CALCULATED);
630 labelList allPatchTypes(mesh_.nCells(), OTHER);
631 labelListList allStencil(mesh_.nCells());
632 scalarListList allWeights(mesh_.nCells());
633 // zoneID of donor
634 labelList allDonorID(mesh_.nCells(), -1);
635
636
637 // Marking patch cells
638 forAll(meshParts, partI)
639 {
640 const fvMesh& partMesh = meshParts[partI].subMesh();
641 const labelList& partCellMap = meshParts[partI].cellMap();
642
643 // Mark cells with
644 // - overset boundary
645 // - other, proper boundary
646 // - other cells
647 Info<< "Marking patch-cells on zone " << partI << endl;
648 markPatchCells(partMesh, partCellMap, allPatchTypes);
649 }
650
651 if ((debug & 2) && mesh_.time().writeTime())
652 {
654 (
655 createField(mesh_, "allPatchTypes", allPatchTypes)
656 );
657 tfld().write();
658 }
659
660
661 labelList nCells(count(3, allPatchTypes));
662 Info<< nl
663 << "After patch analysis : nCells : "
664 << returnReduce(allPatchTypes.size(), sumOp<label>()) << nl
665 << incrIndent
666 << indent << "other : " << nCells[OTHER] << nl
667 << indent << "patch : " << nCells[PATCH] << nl
668 << indent << "overset: " << nCells[OVERSET] << nl
669 << decrIndent << endl;
670
671 globalIndex globalCells(mesh_.nCells());
672
673 for (label srcI = 0; srcI < meshParts.size()-1; srcI++)
674 {
675 const fvMesh& srcMesh = meshParts[srcI].subMesh();
676 const labelList& srcCellMap = meshParts[srcI].cellMap();
677
678 for (label tgtI = srcI+1; tgtI < meshParts.size(); tgtI++)
679 {
680 const fvMesh& tgtMesh = meshParts[tgtI].subMesh();
681 const labelList& tgtCellMap = meshParts[tgtI].cellMap();
682
683 meshToMesh mapper
684 (
685 srcMesh,
686 tgtMesh,
688 HashTable<word>(0), // patchMap,
689 wordList(0), // cuttingPatches
691 false // do not normalise
692 );
693
694 {
695 // Get tgt patch types on src mesh
696 labelList interpolatedTgtPatchTypes(srcMesh.nCells(), -1);
698 (
699 mapper.tgtMap(), // How to get remote data local
700 mapper.srcToTgtCellAddr(),
701 labelList(labelUIndList(allPatchTypes, tgtCellMap)),
702 interpolatedTgtPatchTypes
703 );
704
705 // Get target cell labels in global cell indexing (on overall
706 // mesh)
707 labelList tgtGlobalCells(tgtMesh.nCells());
708 {
709 forAll(tgtCellMap, tgtCellI)
710 {
711 label cellI = tgtCellMap[tgtCellI];
712 tgtGlobalCells[tgtCellI] = globalCells.toGlobal(cellI);
713 }
714 if (mapper.tgtMap())
715 {
716 mapper.tgtMap()->distribute(tgtGlobalCells);
717 }
718 }
719
720
722 (
723 srcI,
724 srcMesh,
725 srcCellMap,
726
727 tgtI,
728 mapper.srcToTgtCellAddr(),
729 mapper.srcToTgtCellWght(),
730 tgtGlobalCells,
731 interpolatedTgtPatchTypes,
732
733 // Overall mesh data
734 allStencil,
735 allWeights,
736 allCellTypes,
737 allDonorID
738 );
739 }
740
741 {
742 // Get src patch types on tgt mesh
743 labelList interpolatedSrcPatchTypes(tgtMesh.nCells(), -1);
745 (
746 mapper.srcMap(), // How to get remote data local
747 mapper.tgtToSrcCellAddr(),
748 labelList(labelUIndList(allPatchTypes, srcCellMap)),
749 interpolatedSrcPatchTypes
750 );
751
752 labelList srcGlobalCells(srcMesh.nCells());
753 {
754 forAll(srcCellMap, srcCellI)
755 {
756 label cellI = srcCellMap[srcCellI];
757 srcGlobalCells[srcCellI] = globalCells.toGlobal(cellI);
758 }
759 if (mapper.srcMap())
760 {
761 mapper.srcMap()->distribute(srcGlobalCells);
762 }
763 }
764
766 (
767 tgtI,
768 tgtMesh,
769 tgtCellMap,
770
771 srcI,
772 mapper.tgtToSrcCellAddr(),
773 mapper.tgtToSrcCellWght(),
774 srcGlobalCells,
775 interpolatedSrcPatchTypes,
776
777 // Overall mesh data
778 allStencil,
779 allWeights,
780 allCellTypes,
781 allDonorID
782 );
783 }
784 }
785 }
786
787
788 if ((debug & 2) && mesh_.time().writeTime())
789 {
791 (
792 createField(mesh_, "allCellTypes", allCellTypes)
793 );
794 tfld().write();
795
796 tmp<volScalarField> tdonors
797 (
798 createField(mesh_, "allDonorID", allDonorID)
799 );
800 tdonors().write();
801
802 }
803
804
805 // Use the patch types and weights to decide what to do
806 forAll(allPatchTypes, cellI)
807 {
808 if (allCellTypes[cellI] != HOLE)
809 {
810 switch (allPatchTypes[cellI])
811 {
812 case OVERSET:
813 {
814 // Interpolate. Check if enough overlap
815 scalar v = mesh_.V()[cellI];
816 scalar overlapVol = sum(allWeights[cellI]);
817 if (overlapVol > overlapTolerance_*v)
818 {
819 allCellTypes[cellI] = INTERPOLATED;
820 }
821 else
822 {
823 allCellTypes[cellI] = HOLE;
824 allWeights[cellI].clear();
825 allStencil[cellI].clear();
826 }
827 break;
828 }
829 }
830 }
831 }
832
833
834 if ((debug & 2) && mesh_.time().writeTime())
835 {
837 (
838 createField(mesh_, "allCellTypes_patch", allCellTypes)
839 );
840 tfld().write();
841
842 labelList stencilSize(mesh_.nCells());
843 forAll(allStencil, celli)
844 {
845 stencilSize[celli] = allStencil[celli].size();
846 }
847 tmp<volScalarField> tfldStencil
848 (
849 createField(mesh_, "allStencil_patch", stencilSize)
850 );
851 tfldStencil().write();
852 }
853
854
855 // Mark unreachable bits
856 findHoles(globalCells, mesh_, zoneID, allStencil, allCellTypes);
857
858 if ((debug & 2) && mesh_.time().writeTime())
859 {
861 (
862 createField(mesh_, "allCellTypes_hole", allCellTypes)
863 );
864 tfld().write();
865
866 labelList stencilSize(mesh_.nCells());
867 forAll(allStencil, celli)
868 {
869 stencilSize[celli] = allStencil[celli].size();
870 }
871 tmp<volScalarField> tfldStencil
872 (
873 createField(mesh_, "allStencil_hole", stencilSize)
874 );
875 tfldStencil().write();
876 }
877
878
879 // Add buffer interpolation layer around holes
880 scalarField allWeight(mesh_.nCells(), Zero);
881
882 labelListList compactStencil(allStencil);
883 List<Map<label>> compactStencilMap;
884 mapDistribute map(globalCells, compactStencil, compactStencilMap);
885
886 scalarList compactCellVol(mesh_.V());
887 map.distribute(compactCellVol);
888
890 (
891 globalCells,
892 layerRelax,
893 allStencil,
894 allCellTypes,
895 allWeight,
896 compactCellVol,
897 compactStencil,
898 zoneID,
899 dict_.getOrDefault("holeLayers", 1),
900 dict_.getOrDefault("useLayer", -1)
901 );
902
903
904 if ((debug & 2) && mesh_.time().writeTime())
905 {
907 (
908 createField(mesh_, "allCellTypes_front", allCellTypes)
909 );
910 tfld().write();
911
912 labelList stencilSize(mesh_.nCells());
913 forAll(allStencil, celli)
914 {
915 stencilSize[celli] = allStencil[celli].size();
916 }
917 tmp<volScalarField> tfldStencil
918 (
919 createField(mesh_, "allStencil_front", stencilSize)
920 );
921 tfldStencil().write();
922 }
923
924 // Check previous iteration cellTypes_ for any hole->calculated changes
925 // If so set the cell either to interpolated (if there are donors) or
926 // holes (if there are no donors). Note that any interpolated cell might
927 // still be overwritten by the flood filling
928 {
929 label nCalculated = 0;
930
931 forAll(cellTypes_, celli)
932 {
933 if (allCellTypes[celli] == CALCULATED && cellTypes_[celli] == HOLE)
934 {
935 if (allStencil[celli].size() == 0)
936 {
937 // Reset to hole
938 allCellTypes[celli] = HOLE;
939 allWeights[celli].clear();
940 allStencil[celli].clear();
941 }
942 else
943 {
944 allCellTypes[celli] = INTERPOLATED;
945 nCalculated++;
946 }
947 }
948 }
949
950 if (debug)
951 {
952 Pout<< "Detected " << nCalculated << " cells changing from hole"
953 << " to calculated. Changed these to interpolated"
954 << endl;
955 }
956 }
957
958
959 labelList compactCellTypes(allCellTypes);
960 map.distribute(compactCellTypes);
961
962 label nHoleDonors = 0;
963 forAll(allCellTypes, cellI)
964 {
965 if (allCellTypes[cellI] == INTERPOLATED)
966 {
967 const labelList& slots = compactStencil[cellI];
968 forAll (slots, subCellI)
969 {
970 if
971 (
972 compactCellTypes[slots[0]] == HOLE
973 ||
974 (
976 && compactCellTypes[slots[0]] == INTERPOLATED
977 )
978 )
979 {
980 allWeights[cellI][subCellI] = 0;
981 nHoleDonors++;
982 }
983 }
984 }
985 else
986 {
987 allWeights[cellI].clear();
988 allStencil[cellI].clear();
989 }
990 }
991 reduce(nHoleDonors, sumOp<label>());
992
993
994 // Normalize weights
995 forAll(allCellTypes, cellI)
996 {
997 if (allCellTypes[cellI] == INTERPOLATED)
998 {
999 const scalar s = sum(allWeights[cellI]);
1000
1001 if (s < SMALL)
1002 {
1003 allCellTypes[cellI] = POROUS;
1004 allWeights[cellI].clear();
1005 allStencil[cellI].clear();
1006 }
1007 else
1008 {
1009 forAll(allWeights[cellI], i)
1010 {
1011 allWeights[cellI][i] /= s;
1012 }
1013 }
1014 }
1015 }
1016
1017 // Write to volField for debugging
1018 if ((debug & 2) && mesh_.time().writeTime())
1019 {
1021 (
1022 createField(mesh_, "allCellTypes_final", allCellTypes)
1023 );
1024 tfld().write();
1025 }
1026
1027
1028 cellTypes_.transfer(allCellTypes);
1029 cellStencil_.transfer(allStencil);
1030 cellInterpolationWeights_.transfer(allWeights);
1031 cellInterpolationWeight_.transfer(allWeight);
1032 //cellInterpolationWeight_.correctBoundaryConditions();
1034 <
1037 false
1038 >(cellInterpolationWeight_.boundaryFieldRef());
1039
1041 forAll(cellStencil_, cellI)
1042 {
1043 if (cellStencil_[cellI].size())
1044 {
1045 interpolationCells.append(cellI);
1046 }
1047 }
1049
1050
1051 List<Map<label>> compactMap;
1053 (
1054 new mapDistribute(globalCells, cellStencil_, compactMap)
1055 );
1056
1057 // Dump interpolation stencil
1058 if ((debug & 2) && mesh_.time().writeTime())
1059 {
1060 // Dump weight
1061 cellInterpolationWeight_.instance() = mesh_.time().timeName();
1063
1064
1065 mkDir(mesh_.time().timePath());
1066 OBJstream str(mesh_.time().timePath()/"stencil2.obj");
1067 Info<< typeName << " : dumping to " << str.name() << endl;
1068 pointField cc(mesh_.cellCentres());
1069 cellInterpolationMap().distribute(cc);
1070
1071 forAll(interpolationCells_, compactI)
1072 {
1073 label cellI = interpolationCells_[compactI];
1074 const labelList& slots = cellStencil_[cellI];
1075
1076 Pout<< "cellI:" << cellI << " at:"
1077 << mesh_.cellCentres()[cellI]
1078 << " calculated from slots:" << slots
1079 << " cc:" << UIndirectList<point>(cc, slots)
1080 << " weights:" << cellInterpolationWeights_[cellI]
1081 << endl;
1082
1083 forAll(slots, i)
1084 {
1085 if (cellInterpolationWeights_[cellI][slots[i]] > 0)
1086 {
1087 const point& donorCc = cc[slots[i]];
1088 const point& accCc = mesh_.cellCentres()[cellI];
1089 str.writeLine(accCc, 0.1*accCc+0.9*donorCc);
1090 }
1091 }
1092 }
1093 }
1094
1095 Info << this->info();
1096
1097 return true;
1098}
1100
1102(
1103 const point& sample,
1104 const pointList& donorCcs,
1105 scalarList& weights
1106) const
1107{
1108 // Inverse-distance weighting
1109
1110 weights.setSize(donorCcs.size());
1111 scalar sum = 0.0;
1112 forAll(donorCcs, i)
1113 {
1114 scalar d = mag(sample-donorCcs[i]);
1115
1116 if (d > ROOTVSMALL)
1117 {
1118 weights[i] = 1.0/d;
1119 sum += weights[i];
1120 }
1121 else
1122 {
1123 // Short circuit
1124 weights = 0.0;
1125 weights[i] = 1.0;
1126 return;
1127 }
1128 }
1129 forAll(weights, i)
1130 {
1131 weights[i] /= sum;
1132 }
1133}
1134
1135
1136// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const polyBoundaryMesh & pbm
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
Minimal example by using system/controlDict.functions:
A HashTable similar to std::unordered_map.
Definition HashTable.H:124
@ NO_REGISTER
Do not request registration (bool: false).
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void setSize(label n)
Alias for resize().
Definition List.H:536
void clear()
Clear the list, i.e. set size to zero.
Definition ListI.H:133
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files.
Definition OBJstream.H:58
Ostream & writeLine(const point &p0, const point &p1)
Write line joining two points.
Definition OBJstream.C:234
virtual const fileName & name() const override
Read/write access to the name of the stream.
Definition OSstream.H:134
static void listReduce(UList< T > &values, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce list elements (list must be equal size on all ranks), applying bop to each list element.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Calculation of interpolation stencils.
const fvMesh & mesh_
Reference to the mesh.
void walkFront(const globalIndex &globalCells, const scalar layerRelax, const labelListList &allStencil, labelList &allCellTypes, scalarField &allWeight, const scalarList &compactCellVol, const labelListList &compactStencil, const labelList &zoneID, const label holeLayers, const label useLayer) const
Surround holes with layer(s) of interpolated cells.
static labelList count(const label size, const labelUList &lst)
Count occurrences (in parallel).
const labelIOList & zoneID() const
Helper: get reference to registered zoneID. Loads volScalarField.
InfoProxy< cellCellStencil > info() const noexcept
Return info proxy, used to print stencil information to a stream.
static const labelIOList & zoneID(const fvMesh &)
Helper: get reference to registered zoneID. Loads volScalarField.
void suppressMotionFields()
Helper: populate nonInterpolatedFields_ with motion solver.
virtual void stencilWeights(const point &sample, const pointList &donorCcs, scalarList &weights) const
Calculate inverse distance weights for a single acceptor. Revert.
static scalar defaultOverlapTolerance_
Default overlap tolerance. Fraction of volume.
List< scalarList > cellInterpolationWeights_
Interpolation weights.
virtual const labelUList & interpolationCells() const
Indices of interpolated cells.
void combineCellTypes(const label subZoneID, const fvMesh &subMesh, const labelList &subCellMap, const label donorZoneID, const labelListList &toOtherCells, const List< scalarList > &weights, const labelList &otherCells, const labelList &interpolatedOtherPatchTypes, labelListList &allStencil, scalarListList &allWeights, labelList &allCellTypes, labelList &allDonorID) const
void interpolatePatchTypes(const labelListList &addressing, const labelList &patchTypes, labelList &result) const
interpolate (= combine) patch types
void findHoles(const globalIndex &globalCells, const fvMesh &mesh, const labelList &zoneID, const labelListList &stencil, labelList &cellTypes) const
Find cells next to cells of type PATCH.
virtual const mapDistribute & cellInterpolationMap() const
Return a communication schedule.
virtual const labelUList & cellTypes() const
Return the cell type list.
const bool allowInterpolatedDonors_
Allow interpolared as donors.
autoPtr< mapDistribute > cellInterpolationMap_
Fetch interpolated cells.
virtual bool update()
Update stencils. Return false if nothing changed.
const dictionary dict_
Dictionary of motion control parameters.
void markPatchCells(const fvMesh &mesh, const labelList &cellMap, labelList &patchCellTypes) const
according to additionalDocumentation/MEJ_oversetMesh.txt
volScalarField cellInterpolationWeight_
Amount of interpolation.
labelList interpolationCells_
Indices of interpolated cells.
scalar overlapTolerance_
Tolerance for volume overlap. Fraction of volume.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A fvBoundaryMesh is a fvPatch list with a reference to the associated fvMesh, with additional search ...
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Template invariant parts for fvPatchField.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
static bool constraintType(const word &patchType)
Return true if the given type is a constraint type.
Definition fvPatch.C:74
virtual const labelUList & faceCells() const
Return faceCells.
Definition fvPatch.C:107
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
label toGlobal(const label proci, const label i) const
From local to global on proci.
Class containing processor-to-processor mapping information.
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute List data using default commsType, default flip/negate operator.
Class to calculate the cell-addressing between two overlapping meshes.
Definition meshToMesh.H:61
const autoPtr< mapDistribute > & srcMap() const noexcept
Source map pointer - valid if no singleMeshProc.
const labelListList & srcToTgtCellAddr() const
Return const access to the source to target cell addressing.
Definition meshToMeshI.H:38
const scalarListList & srcToTgtCellWght() const
Return const access to the source to target cell weights.
Definition meshToMeshI.H:50
const scalarListList & tgtToSrcCellWght() const
Return const access to the target to source cell weights.
Definition meshToMeshI.H:56
const labelListList & tgtToSrcCellAddr() const
Return const access to the target to source cell addressing.
Definition meshToMeshI.H:44
const autoPtr< mapDistribute > & tgtMap() const noexcept
Target map pointer - valid if no singleMeshProc.
static void correctBoundaryConditions(typename GeoField::Boundary &bfld)
Correct boundary conditions of certain type (TypeOnly = true) or explicitly not of the type (TypeOnly...
Boundary condition for use on overset patches. To be run in combination with special dynamicFvMesh ty...
label nCells() const noexcept
Number of mesh cells.
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
label nRegions() const
Return total number of regions.
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData, const bool parRun=UPstream::parRun())
Extract and swap to obtain neighbour cell values for all boundary faces.
A class for managing temporary objects.
Definition tmp.H:75
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const auto & io
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const bitSet isBlockedFace(intersectedFaces())
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
List< scalarList > scalarListList
List of scalarList.
Definition scalarList.H:35
List< word > wordList
List of word.
Definition fileName.H:60
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
const dimensionSet dimless
Dimensionless.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition POSIX.C:616
GeometricField< scalar, fvPatchField, volMesh > volScalarField
UIndirectList< label > labelUIndList
UIndirectList of labels.
IOList< label > labelIOList
IO for a List of label.
Definition labelIOList.H:32
messageStream Info
Information stream (stdout output on master, null elsewhere).
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition Ostream.H:490
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Ostream & indent(Ostream &os)
Indent stream.
Definition Ostream.H:481
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
List< point > pointList
List of point.
Definition pointList.H:32
vector point
Point is a vector.
Definition point.H:37
List< bool > boolList
A List of bools.
Definition List.H:60
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
PDRpatchDef PATCH
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition Ostream.H:499
Type gMax(const FieldField< Field, Type > &f)
static tmp< volScalarField > createField(const fvMesh &mesh, const scalar val)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
wordList patchTypes(nPatches)
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299