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-2020 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
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
44namespace Foam
45{
46namespace cellCellStencils
47{
50}
51}
52
53Foam::scalar
55
56
57// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
58
60(
61 const scalar layerRelax,
62 labelList& allCellTypes,
63 scalarField& allWeight
64) const
65{
66 // Current front
67 bitSet isFront(mesh_.nFaces());
68 // unused: bitSet doneCell(mesh_.nCells());
69
70 const fvBoundaryMesh& fvm = mesh_.boundary();
71
72
73 // 'overset' patches
74
75 forAll(fvm, patchI)
76 {
77 if (isA<oversetFvPatch>(fvm[patchI]))
78 {
79 //Pout<< "Storing faces on patch " << fvm[patchI].name() << endl;
80 forAll(fvm[patchI], i)
81 {
82 isFront.set(fvm[patchI].start()+i);
83 }
84 }
85 }
86
87
88 // Outside of 'hole' region
89 {
90 const labelList& own = mesh_.faceOwner();
91 const labelList& nei = mesh_.faceNeighbour();
92
93 for (label faceI = 0; faceI < mesh_.nInternalFaces(); faceI++)
94 {
95 label ownType = allCellTypes[own[faceI]];
96 label neiType = allCellTypes[nei[faceI]];
97 if
98 (
99 (ownType == HOLE && neiType != HOLE)
100 || (ownType != HOLE && neiType == HOLE)
101 )
102 {
103 //Pout<< "Front at face:" << faceI
104 // << " at:" << mesh_.faceCentres()[faceI] << endl;
105 isFront.set(faceI);
106 }
107 }
108
109 labelList nbrCellTypes;
110 syncTools::swapBoundaryCellList(mesh_, allCellTypes, nbrCellTypes);
111
112 for
113 (
114 label faceI = mesh_.nInternalFaces();
115 faceI < mesh_.nFaces();
116 faceI++
117 )
118 {
119 label ownType = allCellTypes[own[faceI]];
120 label neiType = nbrCellTypes[faceI-mesh_.nInternalFaces()];
121
122 if
123 (
124 (ownType == HOLE && neiType != HOLE)
125 || (ownType != HOLE && neiType == HOLE)
126 )
127 {
128 //Pout<< "Front at coupled face:" << faceI
129 // << " at:" << mesh_.faceCentres()[faceI] << endl;
130 isFront.set(faceI);
131 }
132 }
133 }
134
135
136
137 // Current interpolation fraction
138 scalar fraction = 1.0;
139
140 while (fraction > SMALL && returnReduce(isFront.count(), sumOp<label>()))
141 {
142 // Interpolate cells on front
143
144 Info<< "Front : fraction:" << fraction
145 << " size:" << returnReduce(isFront.count(), sumOp<label>())
146 << endl;
147
148 bitSet newIsFront(mesh_.nFaces());
149 forAll(isFront, faceI)
150 {
151 if (isFront.test(faceI))
152 {
153 label own = mesh_.faceOwner()[faceI];
154 if (allCellTypes[own] != HOLE)
155 {
156 if (allWeight[own] < fraction)
157 {
158 allWeight[own] = fraction;
159
160 //if (debug)
161 //{
162 // Pout<< " setting cell "
163 // << mesh_.cellCentres()[own]
164 // << " to " << fraction << endl;
165 //}
166 allCellTypes[own] = INTERPOLATED;
167 newIsFront.set(mesh_.cells()[own]);
168 }
169 }
170 if (mesh_.isInternalFace(faceI))
171 {
172 label nei = mesh_.faceNeighbour()[faceI];
173 if (allCellTypes[nei] != HOLE)
174 {
175 if (allWeight[nei] < fraction)
176 {
177 allWeight[nei] = fraction;
178
179 //if (debug)
180 //{
181 // Pout<< " setting cell "
182 // << mesh_.cellCentres()[nei]
183 // << " to " << fraction << endl;
184 //}
185
186 allCellTypes[nei] = INTERPOLATED;
187 newIsFront.set(mesh_.cells()[nei]);
188 }
189 }
190 }
191 }
192 }
193
195
196 isFront.transfer(newIsFront);
197
198 fraction -= layerRelax;
199 }
200}
201
202
204(
205 const globalIndex& globalCells,
206 const fvMesh& mesh,
207 const labelList& zoneID,
208 const labelListList& stencil,
210) const
211{
212 const fvBoundaryMesh& pbm = mesh.boundary();
213 const labelList& own = mesh.faceOwner();
214 const labelList& nei = mesh.faceNeighbour();
215
216
217 // The input cellTypes will be
218 // - HOLE : cell part covered by other-mesh patch
219 // - INTERPOLATED : cell fully covered by other-mesh patch
220 // or next to 'overset' patch
221 // - CALCULATED : otherwise
222 //
223 // so we start a walk from our patches and any cell we cannot reach
224 // (because we walk is stopped by other-mesh patch) is a hole.
225
226
227 boolList isBlockedFace(mesh.nFaces(), false);
228 for (label faceI = 0; faceI < mesh.nInternalFaces(); faceI++)
229 {
230 label ownType = cellTypes[own[faceI]];
231 label neiType = cellTypes[nei[faceI]];
232 if
233 (
234 (ownType == HOLE && neiType != HOLE)
235 || (ownType != HOLE && neiType == HOLE)
236 )
237 {
238 isBlockedFace[faceI] = true;
239 }
240 }
241
242 labelList nbrCellTypes;
244
245 for (label faceI = mesh.nInternalFaces(); faceI < mesh.nFaces(); faceI++)
246 {
247 label ownType = cellTypes[own[faceI]];
248 label neiType = nbrCellTypes[faceI-mesh.nInternalFaces()];
249
250 if
251 (
252 (ownType == HOLE && neiType != HOLE)
253 || (ownType != HOLE && neiType == HOLE)
254 )
255 {
256 isBlockedFace[faceI] = true;
257 }
258 }
259
260 regionSplit cellRegion(mesh, isBlockedFace);
261
262 Info<< typeName << " : detected " << cellRegion.nRegions()
263 << " mesh regions after overset" << nl << endl;
264
265
266
267 // Now we'll have a mesh split according to where there are cells
268 // covered by the other-side patches. See what we can reach from our
269 // real patches
270
271 // 0 : region not yet determined
272 // 1 : borders blockage so is not ok (but can be overridden by real
273 // patch)
274 // 2 : has real patch in it so is reachable
275 labelList regionType(cellRegion.nRegions(), Zero);
276
277
278 // See if any regions borders blockage. Note: isBlockedFace is already
279 // parallel synchronised.
280 {
281 for (label faceI = 0; faceI < mesh.nInternalFaces(); faceI++)
282 {
283 if (isBlockedFace[faceI])
284 {
285 label ownRegion = cellRegion[own[faceI]];
286
287 if (cellTypes[own[faceI]] != HOLE)
288 {
289 if (regionType[ownRegion] == 0)
290 {
291 //Pout<< "Mark region:" << ownRegion
292 // << " on zone:" << zoneID[own[faceI]]
293 // << " as next to blockage at:"
294 // << mesh.faceCentres()[faceI] << endl;
295 regionType[ownRegion] = 1;
296 }
297 }
298
299 label neiRegion = cellRegion[nei[faceI]];
300
301 if (cellTypes[nei[faceI]] != HOLE)
302 {
303 if (regionType[neiRegion] == 0)
304 {
305 //Pout<< "Mark region:" << neiRegion
306 // << " on zone:" << zoneID[nei[faceI]]
307 // << " as next to blockage at:"
308 // << mesh.faceCentres()[faceI] << endl;
309 regionType[neiRegion] = 1;
310 }
311 }
312 }
313 }
314 for
315 (
316 label faceI = mesh.nInternalFaces();
317 faceI < mesh.nFaces();
318 faceI++
319 )
320 {
321 if (isBlockedFace[faceI])
322 {
323 label ownRegion = cellRegion[own[faceI]];
324
325 if (regionType[ownRegion] == 0)
326 {
327 //Pout<< "Mark region:" << ownRegion
328 // << " on zone:" << zoneID[own[faceI]]
329 // << " as next to blockage at:"
330 // << mesh.faceCentres()[faceI] << endl;
331 regionType[ownRegion] = 1;
332 }
333 }
334 }
335 }
336
337
338 // Override with real patches
339 forAll(pbm, patchI)
340 {
341 const fvPatch& fvp = pbm[patchI];
342
343 if (isA<oversetFvPatch>(fvp))
344 {}
345 else if (!fvPatch::constraintType(fvp.type()))
346 {
347 //Pout<< "Proper patch " << fvp.name() << " of type " << fvp.type()
348 // << endl;
349
350 const labelList& fc = fvp.faceCells();
351 forAll(fc, i)
352 {
353 label regionI = cellRegion[fc[i]];
354
355 if (cellTypes[fc[i]] != HOLE && regionType[regionI] != 2)
356 {
357 //Pout<< "reachable region : " << regionI
358 // << " at cell " << mesh.cellCentres()[fc[i]]
359 // << " on zone " << zoneID[fc[i]] << endl;
360 regionType[regionI] = 2;
361 }
362 }
363 }
364 }
365
366 // Now we've handled
367 // - cells next to blocked cells
368 // - coupled boundaries
369 // Only thing to handle is the interpolation between regions
370
371
372 labelListList compactStencil(stencil);
373 List<Map<label>> compactMap;
374 mapDistribute map(globalCells, compactStencil, compactMap);
375
376 while (true)
377 {
378 // Synchronise region status on processors
379 // (could instead swap status through processor patches)
381
382 // Communicate region status through interpolative cells
383 labelList cellRegionType(labelUIndList(regionType, cellRegion));
384 map.distribute(cellRegionType);
385
386
387 label nChanged = 0;
388 forAll(pbm, patchI)
389 {
390 const fvPatch& fvp = pbm[patchI];
391
392 if (isA<oversetFvPatch>(fvp))
393 {
394 const labelUList& fc = fvp.faceCells();
395 forAll(fc, i)
396 {
397 label cellI = fc[i];
398 label regionI = cellRegion[cellI];
399
400 if (regionType[regionI] != 2)
401 {
402 const labelList& slots = compactStencil[cellI];
403 forAll(slots, i)
404 {
405 label otherType = cellRegionType[slots[i]];
406
407 if (otherType == 2)
408 {
409 //Pout<< "Reachable through interpolation : "
410 // << regionI << " at cell "
411 // << mesh.cellCentres()[cellI] << endl;
412 regionType[regionI] = 2;
413 nChanged++;
414 break;
415 }
416 }
417 }
418 }
419 }
420 }
421
422
423 reduce(nChanged, sumOp<label>());
424 if (nChanged == 0)
425 {
426 break;
427 }
428 }
429
430
431 // See which regions have not been visited (regionType == 1)
432 forAll(cellRegion, cellI)
433 {
434 label type = regionType[cellRegion[cellI]];
435 if (type == 1 && cellTypes[cellI] != HOLE)
436 {
437 cellTypes[cellI] = HOLE;
438 }
439 }
440}
441
442
444(
445 const fvMesh& mesh,
446 const labelList& cellMap,
447 labelList& patchCellTypes
448) const
449{
450 const fvBoundaryMesh& pbm = mesh.boundary();
451
452 forAll(pbm, patchI)
453 {
454 const fvPatch& fvp = pbm[patchI];
455 const labelList& fc = fvp.faceCells();
456
457 if (isA<oversetFvPatch>(fvp))
458 {
459 //Pout<< "Marking cells on overset patch " << fvp.name() << endl;
460 forAll(fc, i)
461 {
462 label cellI = fc[i];
463 patchCellTypes[cellMap[cellI]] = OVERSET;
464 }
465 }
466 else if (!fvPatch::constraintType(fvp.type()))
467 {
468 //Pout<< "Marking cells on proper patch " << fvp.name()
469 // << " with type " << fvp.type() << endl;
470 forAll(fc, i)
471 {
472 label cellI = fc[i];
473 if (patchCellTypes[cellMap[cellI]] != OVERSET)
474 {
475 patchCellTypes[cellMap[cellI]] = PATCH;
476 }
477 }
478 }
479 }
480}
481
482
484(
485 const labelListList& addressing,
486 const labelList& patchTypes,
487 labelList& result
488) const
489{
490 forAll(result, cellI)
491 {
492 const labelList& slots = addressing[cellI];
493 forAll(slots, i)
494 {
495 label type = patchTypes[slots[i]];
496
497 if (type == OVERSET)
498 {
499 // 'overset' overrides anything
500 result[cellI] = OVERSET;
501 break;
502 }
503 else if (type == PATCH)
504 {
505 // 'patch' overrides -1 and 'other'
506 result[cellI] = PATCH;
507 }
508 else if (result[cellI] == -1)
509 {
510 // 'other' overrides -1 only
511 result[cellI] = OTHER;
512 }
513 }
514 }
515}
516
517
519(
520 const autoPtr<mapDistribute>& mapPtr,
521 const labelListList& addressing,
522 const labelList& patchTypes,
523 labelList& result
524) const
525{
526 if (result.size() != addressing.size())
527 {
528 FatalErrorInFunction << "result:" << result.size()
529 << " addressing:" << addressing.size() << exit(FatalError);
530 }
531
532
533 // Initialise to not-mapped
534 result = -1;
535
536 if (mapPtr)
537 {
538 // Pull remote data into order of addressing
539 labelList work(patchTypes);
540 mapPtr().distribute(work);
541
542 interpolatePatchTypes(addressing, work, result);
543 }
544 else
545 {
546 interpolatePatchTypes(addressing, patchTypes, result);
547 }
548}
549
550
552(
553 const label subZoneID,
554 const fvMesh& subMesh,
555 const labelList& subCellMap,
556
557 const label donorZoneID,
558 const labelListList& addressing,
559 const List<scalarList>& weights,
560 const labelList& otherCells,
561 const labelList& interpolatedOtherPatchTypes,
562
563 labelListList& allStencil,
564 scalarListList& allWeights,
565 labelList& allCellTypes,
566 labelList& allDonorID
567) const
568{
569 forAll(subCellMap, subCellI)
570 {
571 label cellI = subCellMap[subCellI];
572
573 bool validDonors = true;
574 switch (interpolatedOtherPatchTypes[subCellI])
575 {
576 case -1:
577 {
578 validDonors = false;
579 }
580 break;
581
582 case OTHER:
583 {
584 // No patch interaction so keep valid
585 }
586 break;
587
588 case PATCH:
589 {
590 // Patch-patch interaction... For now disable always
591 allCellTypes[cellI] = HOLE;
592 validDonors = false;
593
594 // Alternative is to look at the amount of overlap but this
595 // is not very robust
596 //if (allCellTypes[cellI] != HOLE)
597 //{
598 // scalar overlapVol = sum(weights[subCellI]);
599 // scalar v = mesh_.V()[cellI];
600 // if (overlapVol < (1.0-overlapTolerance_)*v)
601 // {
602 // //Pout<< "** Patch overlap:" << cellI
603 // // << " at:" << mesh_.cellCentres()[cellI] << endl;
604 // allCellTypes[cellI] = HOLE;
605 // validDonors = false;
606 // }
607 //}
608 }
609 break;
610
611 case OVERSET:
612 {
613 validDonors = false;
614 }
615 break;
616 }
617
618
619 if (validDonors)
620 {
621 // There are a few possible choices how to choose between multiple
622 // donor candidates:
623 // 1 highest overlap volume. However this is generally already
624 // 99.9% so you're just measuring truncation error.
625 // 2 smallest donors cells or most donor cells. This is quite
626 // often done but can cause switching of donor zone from one
627 // time step to the other if the donor meshes are non-uniform
628 // and the acceptor cells just happens to be sweeping through
629 // some small donor cells.
630 // 3 nearest zoneID. So zone 0 preferentially interpolates from
631 // zone 1, zone 1 preferentially from zone 2 etc.
632
633 //- Option 1:
634 //scalar currentVol = sum(allWeights[cellI]);
635 //if (overlapVol[subCellI] > currentVol)
636
637 //- Option 3:
638 label currentDiff = mag(subZoneID-allDonorID[cellI]);
639 label thisDiff = mag(subZoneID-donorZoneID);
640
641 if
642 (
643 allDonorID[cellI] == -1
644 || (thisDiff < currentDiff)
645 || (thisDiff == currentDiff && donorZoneID > allDonorID[cellI])
646 )
647 {
648 allWeights[cellI] = weights[subCellI];
649 allStencil[cellI] =
650 labelUIndList(otherCells, addressing[subCellI]);
651 allDonorID[cellI] = donorZoneID;
652 }
653 }
654 }
655}
656
657
658// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
659
661(
662 const fvMesh& mesh,
663 const dictionary& dict,
664 const bool doUpdate
665)
666:
668 dict_(dict),
669 overlapTolerance_(defaultOverlapTolerance_),
670 cellTypes_(labelList(mesh.nCells(), CALCULATED)),
671 interpolationCells_(0),
672 cellInterpolationMap_(),
673 cellStencil_(0),
674 cellInterpolationWeights_(0),
675 cellInterpolationWeight_
676 (
678 (
679 "cellInterpolationWeight",
680 mesh_.facesInstance(),
681 mesh_,
682 IOobject::NO_READ,
683 IOobject::NO_WRITE,
684 false
685 ),
686 mesh_,
688 zeroGradientFvPatchScalarField::typeName
689 )
690{
691 // Protect local fields from interpolation
692 nonInterpolatedFields_.insert("cellTypes");
693 nonInterpolatedFields_.insert("cellInterpolationWeight");
694
695 // For convenience also suppress frequently used displacement field
696 nonInterpolatedFields_.insert("cellDisplacement");
697 nonInterpolatedFields_.insert("grad(cellDisplacement)");
698 const word w("snGradCorr(cellDisplacement)");
699 const word d("((viscosity*faceDiffusivity)*magSf)");
700 nonInterpolatedFields_.insert("surfaceIntegrate(("+d+"*"+w+"))");
701
702 // Read zoneID
703 this->zoneID();
704
705 // Read old-time cellTypes
707 (
708 "cellTypes",
709 mesh_.time().timeName(),
710 mesh_,
713 false
714 );
715 if (io.typeHeaderOk<volScalarField>(true))
716 {
717 if (debug)
718 {
719 Pout<< "Reading cellTypes from time " << mesh_.time().timeName()
720 << endl;
721 }
722
723 const volScalarField volCellTypes(io, mesh_);
724 forAll(volCellTypes, celli)
725 {
726 // Round to integer
727 cellTypes_[celli] = volCellTypes[celli];
728 }
729 }
730
731 if (doUpdate)
732 {
733 update();
734 }
735}
736
737
738// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
739
741{}
742
743
744// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
745
747{
748 scalar layerRelax(dict_.getOrDefault("layerRelax", 1.0));
749 const labelIOList& zoneID = this->zoneID();
750
751 label nZones = gMax(zoneID)+1;
753 forAll(zoneID, cellI)
754 {
755 nCellsPerZone[zoneID[cellI]]++;
756 }
758
759
760 Info<< typeName << " : detected " << nZones
761 << " mesh regions" << nl << endl;
762
763
765
767 forAll(meshParts, zonei)
768 {
769 Info<< indent<< "zone:" << zonei << " nCells:"
770 << nCellsPerZone[zonei] << nl;
771
772 meshParts.set
773 (
774 zonei,
775 new fvMeshSubset(mesh_, zonei, zoneID)
776 );
777 }
779
780
781
782 // Current best guess for cells. Includes best stencil. Weights should
783 // add up to volume.
784 labelList allCellTypes(mesh_.nCells(), CALCULATED);
785 labelList allPatchTypes(mesh_.nCells(), OTHER);
786 labelListList allStencil(mesh_.nCells());
787 scalarListList allWeights(mesh_.nCells());
788 // zoneID of donor
789 labelList allDonorID(mesh_.nCells(), -1);
790
791
792 // Marking patch cells
793 forAll(meshParts, partI)
794 {
795 const fvMesh& partMesh = meshParts[partI].subMesh();
796 const labelList& partCellMap = meshParts[partI].cellMap();
797
798 // Mark cells with
799 // - overset boundary
800 // - other, proper boundary
801 // - other cells
802 Info<< "Marking patch-cells on zone " << partI << endl;
803 markPatchCells(partMesh, partCellMap, allPatchTypes);
804 }
805
806
807 labelList nCells(count(3, allPatchTypes));
808 Info<< nl
809 << "After patch analysis : nCells : "
810 << returnReduce(allPatchTypes.size(), sumOp<label>()) << nl
811 << incrIndent
812 << indent << "other : " << nCells[OTHER] << nl
813 << indent << "patch : " << nCells[PATCH] << nl
814 << indent << "overset: " << nCells[OVERSET] << nl
815 << decrIndent << endl;
816
817 globalIndex globalCells(mesh_.nCells());
818
819
820 for (label srcI = 0; srcI < meshParts.size()-1; srcI++)
821 {
822 const fvMesh& srcMesh = meshParts[srcI].subMesh();
823 const labelList& srcCellMap = meshParts[srcI].cellMap();
824
825 for (label tgtI = srcI+1; tgtI < meshParts.size(); tgtI++)
826 {
827 const fvMesh& tgtMesh = meshParts[tgtI].subMesh();
828 const labelList& tgtCellMap = meshParts[tgtI].cellMap();
829
830 meshToMesh mapper
831 (
832 srcMesh,
833 tgtMesh,
835 HashTable<word>(0), // patchMap,
836 wordList(0), // cuttingPatches
838 false // do not normalise
839 );
840
841
842 {
843 // Get tgt patch types on src mesh
844 labelList interpolatedTgtPatchTypes(srcMesh.nCells(), -1);
845 interpolatePatchTypes
846 (
847 mapper.tgtMap(), // How to get remote data local
848 mapper.srcToTgtCellAddr(),
849 labelList(labelUIndList(allPatchTypes, tgtCellMap)),
850 interpolatedTgtPatchTypes
851 );
852
853 // Get target cell labels in global cell indexing (on overall
854 // mesh)
855 labelList tgtGlobalCells(tgtMesh.nCells());
856 {
857 forAll(tgtCellMap, tgtCellI)
858 {
859 label cellI = tgtCellMap[tgtCellI];
860 tgtGlobalCells[tgtCellI] = globalCells.toGlobal(cellI);
861 }
862 if (mapper.tgtMap())
863 {
864 mapper.tgtMap()->distribute(tgtGlobalCells);
865 }
866 }
867 combineCellTypes
868 (
869 srcI,
870 srcMesh,
871 srcCellMap,
872
873 tgtI,
874 mapper.srcToTgtCellAddr(),
875 mapper.srcToTgtCellWght(),
876 tgtGlobalCells,
877 interpolatedTgtPatchTypes,
878
879 // Overall mesh data
880 allStencil,
881 allWeights,
882 allCellTypes,
883 allDonorID
884 );
885 }
886
887 {
888 // Get src patch types on tgt mesh
889 labelList interpolatedSrcPatchTypes(tgtMesh.nCells(), -1);
890 interpolatePatchTypes
891 (
892 mapper.srcMap(), // How to get remote data local
893 mapper.tgtToSrcCellAddr(),
894 labelList(labelUIndList(allPatchTypes, srcCellMap)),
895 interpolatedSrcPatchTypes
896 );
897
898 labelList srcGlobalCells(srcMesh.nCells());
899 {
900 forAll(srcCellMap, srcCellI)
901 {
902 label cellI = srcCellMap[srcCellI];
903 srcGlobalCells[srcCellI] = globalCells.toGlobal(cellI);
904 }
905 if (mapper.srcMap())
906 {
907 mapper.srcMap()->distribute(srcGlobalCells);
908 }
909 }
910
911 combineCellTypes
912 (
913 tgtI,
914 tgtMesh,
915 tgtCellMap,
916
917 srcI,
918 mapper.tgtToSrcCellAddr(),
919 mapper.tgtToSrcCellWght(),
920 srcGlobalCells,
921 interpolatedSrcPatchTypes,
922
923 // Overall mesh data
924 allStencil,
925 allWeights,
926 allCellTypes,
927 allDonorID
928 );
929 }
930 }
931 }
932
933
934 if (debug)
935 {
937 (
938 createField(mesh_, "allCellTypes", allCellTypes)
939 );
940 tfld().write();
941 }
942
943
944 // Use the patch types and weights to decide what to do
945 forAll(allPatchTypes, cellI)
946 {
947 if (allCellTypes[cellI] != HOLE)
948 {
949 switch (allPatchTypes[cellI])
950 {
951 case OVERSET:
952 {
953 // Interpolate. Check if enough overlap
954 scalar v = mesh_.V()[cellI];
955 scalar overlapVol = sum(allWeights[cellI]);
956 if (overlapVol > overlapTolerance_*v)
957 {
958 allCellTypes[cellI] = INTERPOLATED;
959 }
960 else
961 {
962 //Pout<< "Holeing interpolated cell:" << cellI
963 // << " at:" << mesh_.cellCentres()[cellI] << endl;
964 allCellTypes[cellI] = HOLE;
965 allWeights[cellI].clear();
966 allStencil[cellI].clear();
967 }
968 break;
969 }
970 }
971 }
972 }
973
974/*
975 // Knock out cell with insufficient interpolation weights
976 forAll(allCellTypes, cellI)
977 {
978 if (allCellTypes[cellI] == INTERPOLATED)
979 {
980 scalar v = mesh_.V()[cellI];
981 scalar overlapVol = sum(allWeights[cellI]);
982 if (overlapVol < (1.0-overlapTolerance_)*v)
983 {
984 //Pout<< "Holeing cell:" << cellI
985 // << " at:" << mesh_.cellCentres()[cellI] << endl;
986 allCellTypes[cellI] = HOLE;
987 allWeights[cellI].clear();
988 allStencil[cellI].clear();
989 }
990 }
991 }
992*/
993 if (debug)
994 {
996 (
997 createField(mesh_, "allCellTypes_patch", allCellTypes)
998 );
999 //tfld.ref().correctBoundaryConditions();
1000 tfld().write();
1001 }
1002
1003
1004 // Check previous iteration cellTypes_ for any hole->calculated changes
1005 // If so set the cell either to interpolated (if there are donors) or
1006 // holes (if there are no donors). Note that any interpolated cell might
1007 // still be overwritten by the flood filling
1008 {
1009 label nCalculated = 0;
1010
1011 forAll(cellTypes_, celli)
1012 {
1013 if (allCellTypes[celli] == CALCULATED && cellTypes_[celli] == HOLE)
1014 {
1015 if (allStencil[celli].size() == 0)
1016 {
1017 // Reset to hole
1018 allCellTypes[celli] = HOLE;
1019 allWeights[celli].clear();
1020 allStencil[celli].clear();
1021 }
1022 else
1023 {
1024 allCellTypes[celli] = INTERPOLATED;
1025 nCalculated++;
1026 }
1027 }
1028 }
1029
1030 if (debug)
1031 {
1032 Pout<< "Detected " << nCalculated << " cells changing from hole"
1033 << " to calculated. Changed these to interpolated"
1034 << endl;
1035 }
1036 }
1037
1038
1039 // Mark unreachable bits
1040 findHoles(globalCells, mesh_, zoneID, allStencil, allCellTypes);
1041
1042 if (debug)
1043 {
1045 (
1046 createField(mesh_, "allCellTypes_hole", allCellTypes)
1047 );
1048 //tfld.ref().correctBoundaryConditions();
1049 tfld().write();
1050 }
1051
1052
1053 // Add buffer interpolation layer around holes
1054 scalarField allWeight(mesh_.nCells(), Zero);
1055 walkFront(layerRelax, allCellTypes, allWeight);
1056
1057 if (debug)
1058 {
1060 (
1061 createField(mesh_, "allCellTypes_front", allCellTypes)
1062 );
1063 //tfld.ref().correctBoundaryConditions();
1064 tfld().write();
1065 }
1066
1067 // Normalise weights, Clear storage
1068 forAll(allCellTypes, cellI)
1069 {
1070 if (allCellTypes[cellI] == INTERPOLATED)
1071 {
1072 if (allWeight[cellI] < SMALL || allStencil[cellI].size() == 0)
1073 {
1074 //Pout<< "Clearing cell:" << cellI
1075 // << " at:" << mesh_.cellCentres()[cellI] << endl;
1076 allWeights[cellI].clear();
1077 allStencil[cellI].clear();
1078 allWeight[cellI] = 0.0;
1079 }
1080 else
1081 {
1082 scalar s = sum(allWeights[cellI]);
1083 forAll(allWeights[cellI], i)
1084 {
1085 allWeights[cellI][i] /= s;
1086 }
1087 }
1088 }
1089 else
1090 {
1091 allWeights[cellI].clear();
1092 allStencil[cellI].clear();
1093 }
1094 }
1095
1096
1097 // Write to volField for debugging
1098 if (debug)
1099 {
1101 (
1102 IOobject
1103 (
1104 "patchTypes",
1105 mesh_.time().timeName(),
1106 mesh_,
1109 false
1110 ),
1111 mesh_,
1113 zeroGradientFvPatchScalarField::typeName
1114 );
1115
1116 forAll(patchTypes.internalField(), cellI)
1117 {
1118 patchTypes[cellI] = allPatchTypes[cellI];
1119 }
1120 //patchTypes.correctBoundaryConditions();
1122 <
1125 >(patchTypes.boundaryFieldRef(), false);
1126 patchTypes.write();
1127 }
1128 if (debug)
1129 {
1130 volScalarField volTypes
1131 (
1132 IOobject
1133 (
1134 "cellTypes",
1135 mesh_.time().timeName(),
1136 mesh_,
1139 false
1140 ),
1141 mesh_,
1143 zeroGradientFvPatchScalarField::typeName
1144 );
1145
1146 forAll(volTypes.internalField(), cellI)
1147 {
1148 volTypes[cellI] = allCellTypes[cellI];
1149 }
1150 //volTypes.correctBoundaryConditions();
1152 <
1155 >(volTypes.boundaryFieldRef(), false);
1156 volTypes.write();
1157 }
1158
1159
1160// // Check previous iteration cellTypes_ for any hole->calculated changes
1161// {
1162// label nCalculated = 0;
1163//
1164// forAll(cellTypes_, celli)
1165// {
1166// if (allCellTypes[celli] == CALCULATED && cellTypes_[celli] == HOLE)
1167// {
1168// if (allStencil[celli].size() == 0)
1169// {
1170// FatalErrorInFunction
1171// << "Cell:" << celli
1172// << " at:" << mesh_.cellCentres()[celli]
1173// << " zone:" << zoneID[celli]
1174// << " changed from hole to calculated"
1175// << " but there is no donor"
1176// << exit(FatalError);
1177// }
1178// else
1179// {
1180// allCellTypes[celli] = INTERPOLATED;
1181// nCalculated++;
1182// }
1183// }
1184// }
1185//
1186// if (debug)
1187// {
1188// Pout<< "Detected " << nCalculated << " cells changing from hole"
1189// << " to calculated. Changed these to interpolated"
1190// << endl;
1191// }
1192// }
1193
1194
1195 cellTypes_.transfer(allCellTypes);
1196 cellStencil_.transfer(allStencil);
1197 cellInterpolationWeights_.transfer(allWeights);
1198 cellInterpolationWeight_.transfer(allWeight);
1199 //cellInterpolationWeight_.correctBoundaryConditions();
1201 <
1204 >(cellInterpolationWeight_.boundaryFieldRef(), false);
1205
1206 DynamicList<label> interpolationCells;
1207 forAll(cellStencil_, cellI)
1208 {
1209 if (cellStencil_[cellI].size())
1210 {
1211 interpolationCells.append(cellI);
1212 }
1213 }
1214 interpolationCells_.transfer(interpolationCells);
1215
1216
1217 List<Map<label>> compactMap;
1218 cellInterpolationMap_.reset
1219 (
1220 new mapDistribute(globalCells, cellStencil_, compactMap)
1221 );
1222
1223 // Dump interpolation stencil
1224 if (debug)
1225 {
1226 // Dump weight
1227 cellInterpolationWeight_.instance() = mesh_.time().timeName();
1228 cellInterpolationWeight_.write();
1229
1230
1231 mkDir(mesh_.time().timePath());
1232 OBJstream str(mesh_.time().timePath()/"stencil2.obj");
1233 Info<< typeName << " : dumping to " << str.name() << endl;
1234 pointField cc(mesh_.cellCentres());
1235 cellInterpolationMap().distribute(cc);
1236
1237 forAll(interpolationCells_, compactI)
1238 {
1239 label cellI = interpolationCells_[compactI];
1240 const labelList& slots = cellStencil_[cellI];
1241
1242 Pout<< "cellI:" << cellI << " at:"
1243 << mesh_.cellCentres()[cellI]
1244 << " calculated from slots:" << slots
1245 << " cc:" << UIndirectList<point>(cc, slots)
1246 << " weights:" << cellInterpolationWeights_[cellI]
1247 << endl;
1248
1249 forAll(slots, i)
1250 {
1251 const point& donorCc = cc[slots[i]];
1252 const point& accCc = mesh_.cellCentres()[cellI];
1253
1254 str.write(linePointRef(accCc, 0.1*accCc+0.9*donorCc));
1255 }
1256 }
1257 }
1258
1259
1260 {
1261 labelList nCells(count(3, cellTypes_));
1262 Info<< "Overset analysis : nCells : "
1263 << returnReduce(cellTypes_.size(), sumOp<label>()) << nl
1264 << incrIndent
1265 << indent << "calculated : " << nCells[CALCULATED] << nl
1266 << indent << "interpolated : " << nCells[INTERPOLATED] << nl
1267 << indent << "hole : " << nCells[HOLE] << nl
1268 << decrIndent << endl;
1269 }
1270
1271 // Tbd: detect if anything changed. Most likely it did!
1272 return true;
1273}
1274
1275
1277(
1278 const point& sample,
1279 const pointList& donorCcs,
1280 scalarList& weights
1281) const
1282{
1283 // Inverse-distance weighting
1284
1285 weights.setSize(donorCcs.size());
1286 scalar sum = 0.0;
1287 forAll(donorCcs, i)
1288 {
1289 scalar d = mag(sample-donorCcs[i]);
1290
1291 if (d > ROOTVSMALL)
1292 {
1293 weights[i] = 1.0/d;
1294 sum += weights[i];
1295 }
1296 else
1297 {
1298 // Short circuit
1299 weights = 0.0;
1300 weights[i] = 1.0;
1301 return;
1302 }
1303 }
1304 forAll(weights, i)
1305 {
1306 weights[i] /= sum;
1307 }
1308}
1309
1310
1311// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
Minimal example by using system/controlDict.functions:
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.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
A HashTable similar to std::unordered_map.
Definition: HashTable.H:123
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 setSize(const label n)
Alias for resize()
Definition: List.H:218
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
Bookkeeping for patch definitions.
Definition: PDRpatchDef.H:54
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
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
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:66
unsigned int count(const bool on=true) const
Count number of bits set.
Definition: bitSetI.H:500
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
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.
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.
void walkFront(const scalar layerRelax, labelList &allCellTypes, scalarField &allWeight) const
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 bool update()
Update stencils. Return false if nothing changed.
void markPatchCells(const fvMesh &mesh, const labelList &cellMap, labelList &patchCellTypes) const
according to additionalDocumentation/MEJ_oversetMesh.txt
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
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
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.
Class to calculate the cell-addressing between two overlapping meshes.
Definition: meshToMesh.H:65
const labelListList & srcToTgtCellAddr() const
Return const access to the source to target cell addressing.
Definition: meshToMeshI.H:45
const autoPtr< mapDistribute > & srcMap() const
Source map pointer - valid if no singleMeshProc.
Definition: meshToMeshI.H:88
const scalarListList & srcToTgtCellWght() const
Return const access to the source to target cell weights.
Definition: meshToMeshI.H:57
const scalarListList & tgtToSrcCellWght() const
Return const access to the target to source cell weights.
Definition: meshToMeshI.H:63
const labelListList & tgtToSrcCellAddr() const
Return const access to the target to source cell addressing.
Definition: meshToMeshI.H:51
const autoPtr< mapDistribute > & tgtMap() const
Target map pointer - valid if no singleMeshProc.
Definition: meshToMeshI.H:95
Boundary condition for use on overset patches. To be run in combination with special dynamicFvMesh ty...
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1121
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1127
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
label nInternalFaces() const noexcept
Number of internal faces.
label nCells() const noexcept
Number of mesh cells.
label nFaces() const noexcept
Number of mesh faces.
const cellList & cells() const
virtual bool write(const bool valid=true) const
Write using setting from DB.
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
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 FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
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))
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)
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:63
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
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
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)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
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.
PDRpatchDef PATCH
Definition: PDRpatchDef.H:112
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)
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