Loading...
Searching...
No Matches
meshRefinementRefine.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2015-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "meshRefinement.H"
30#include "trackedParticle.H"
31#include "syncTools.H"
32#include "Time.H"
33#include "refinementSurfaces.H"
34#include "refinementFeatures.H"
35#include "shellSurfaces.H"
36#include "faceSet.H"
37#include "decompositionMethod.H"
38#include "fvMeshDistribute.H"
39#include "polyTopoChange.H"
41#include "Cloud.H"
42//#include "globalIndex.H"
43#include "cellSet.H"
44#include "treeDataCell.H"
45
46#include "cellCuts.H"
47#include "refineCell.H"
48#include "hexCellLooper.H"
49#include "meshCutter.H"
50
51// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
52
53namespace Foam
54{
55 //- To compare normals
56 static bool less(const vector& x, const vector& y)
57 {
58 for (direction i = 0; i < vector::nComponents; i++)
59 {
60 if (x[i] < y[i])
61 {
62 return true;
63 }
64 else if (x[i] > y[i])
65 {
66 return false;
67 }
68 }
69 // All components the same
70 return false;
71 }
73
74 //- To compare normals
75 class normalLess
76 {
77 const vectorList& values_;
79 public:
80
81 normalLess(const vectorList& values)
82 :
83 values_(values)
84 {}
85
86 bool operator()(const label a, const label b) const
87 {
88 return less(values_[a], values_[b]);
89 }
90 };
91
92
93 //- Template specialization for pTraits<labelList> so we can have fields
94 template<>
95 class pTraits<labelList>
96 {
97 public:
98
99 //- Component type
100 typedef labelList cmptType;
101 };
102
103 //- Template specialization for pTraits<vectorList> so we can have fields
104 template<>
105 class pTraits<vectorList>
106 {
107 public:
109 //- Component type
110 typedef vectorList cmptType;
111 };
112}
113
114
115// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
116
117// Get faces (on the new mesh) that have in some way been affected by the
118// mesh change. Picks up all faces but those that are between two
119// unrefined faces. (Note that of an unchanged face the edge still might be
120// split but that does not change any face centre or cell centre.
121Foam::labelList Foam::meshRefinement::getChangedFaces
122(
123 const mapPolyMesh& map,
124 const labelList& oldCellsToRefine
125)
126{
127 const polyMesh& mesh = map.mesh();
128
129 labelList changedFaces;
130
131 // For reporting: number of masterFaces changed
132 label nMasterChanged = 0;
134
135 {
136 // Mark any face on a cell which has been added or changed
137 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
138 // Note that refining a face changes the face centre (for a warped face)
139 // which changes the cell centre. This again changes the cellcentre-
140 // cellcentre edge across all faces of the cell.
141 // Note: this does not happen for unwarped faces but unfortunately
142 // we don't have this information.
143
144 const labelList& faceOwner = mesh.faceOwner();
145 const labelList& faceNeighbour = mesh.faceNeighbour();
146 const cellList& cells = mesh.cells();
147 const label nInternalFaces = mesh.nInternalFaces();
148
149 // Mark refined cells on old mesh
150 bitSet oldRefineCell(map.nOldCells(), oldCellsToRefine);
151
152 // Mark refined faces
153 bitSet refinedInternalFace(nInternalFaces);
154
155 // 1. Internal faces
156
157 for (label faceI = 0; faceI < nInternalFaces; faceI++)
158 {
159 label oldOwn = map.cellMap()[faceOwner[faceI]];
160 label oldNei = map.cellMap()[faceNeighbour[faceI]];
161
162 if
163 (
164 oldOwn >= 0 && !oldRefineCell.test(oldOwn)
165 && oldNei >= 0 && !oldRefineCell.test(oldNei)
166 )
167 {
168 // Unaffected face since both neighbours were not refined.
169 }
170 else
171 {
172 refinedInternalFace.set(faceI);
173 }
174 }
175
176
177 // 2. Boundary faces
178
179 boolList refinedBoundaryFace(mesh.nBoundaryFaces(), false);
180
181 forAll(mesh.boundaryMesh(), patchI)
182 {
183 const polyPatch& pp = mesh.boundaryMesh()[patchI];
184
185 label faceI = pp.start();
186
187 forAll(pp, i)
188 {
189 label oldOwn = map.cellMap()[faceOwner[faceI]];
190
191 if (oldOwn >= 0 && !oldRefineCell.test(oldOwn))
192 {
193 // owner did exist and wasn't refined.
194 }
195 else
196 {
197 refinedBoundaryFace[faceI-nInternalFaces] = true;
198 }
199 faceI++;
200 }
201 }
202
203 // Synchronise refined face status
205 (
206 mesh,
207 refinedBoundaryFace,
209 );
210
211
212 // 3. Mark all faces affected by refinement. Refined faces are in
213 // - refinedInternalFace
214 // - refinedBoundaryFace
215 boolList changedFace(mesh.nFaces(), false);
216
217 forAll(refinedInternalFace, faceI)
218 {
219 if (refinedInternalFace.test(faceI))
220 {
221 const cell& ownFaces = cells[faceOwner[faceI]];
222 forAll(ownFaces, ownI)
223 {
224 changedFace[ownFaces[ownI]] = true;
225 }
226 const cell& neiFaces = cells[faceNeighbour[faceI]];
227 forAll(neiFaces, neiI)
228 {
229 changedFace[neiFaces[neiI]] = true;
230 }
231 }
232 }
233
234 forAll(refinedBoundaryFace, i)
235 {
236 if (refinedBoundaryFace[i])
237 {
238 const cell& ownFaces = cells[faceOwner[i+nInternalFaces]];
239 forAll(ownFaces, ownI)
240 {
241 changedFace[ownFaces[ownI]] = true;
242 }
243 }
244 }
245
247 (
248 mesh,
249 changedFace,
251 );
252
253
254 // Now we have in changedFace marked all affected faces. Pack.
255 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
256
257 changedFaces = findIndices(changedFace, true);
258
259 // Count changed master faces.
260 nMasterChanged = 0;
261
262 forAll(changedFace, faceI)
263 {
264 if (changedFace[faceI] && isMasterFace.test(faceI))
265 {
266 nMasterChanged++;
267 }
268 }
269
270 }
271
273 {
274 Pout<< "getChangedFaces : Detected "
275 << " local:" << changedFaces.size()
276 << " global:" << returnReduce(nMasterChanged, sumOp<label>())
277 << " changed faces out of " << mesh.globalData().nTotalFaces()
278 << endl;
279
280 faceSet changedFacesSet(mesh, "changedFaces", changedFaces);
281 Pout<< "getChangedFaces : Writing " << changedFaces.size()
282 << " changed faces to faceSet " << changedFacesSet.name()
283 << endl;
284 changedFacesSet.write();
285 }
286
287 return changedFaces;
288}
289
290
291// Mark cell for refinement (if not already marked). Return false if
292// refinelimit hit. Keeps running count (in nRefine) of cells marked for
293// refinement
294bool Foam::meshRefinement::markForRefine
295(
296 const label markValue,
297 const label nAllowRefine,
298
299 label& cellValue,
300 label& nRefine
301)
302{
303 if (cellValue == -1)
304 {
305 cellValue = markValue;
306 nRefine++;
307 }
308
309 return nRefine <= nAllowRefine;
310}
311
312
313void Foam::meshRefinement::markFeatureCellLevel
314(
315 const pointField& keepPoints,
316 labelList& maxFeatureLevel
317) const
318{
319 // We want to refine all cells containing a feature edge.
320 // - don't want to search over whole mesh
321 // - don't want to build octree for whole mesh
322 // - so use tracking to follow the feature edges
323 //
324 // 1. find non-manifold points on feature edges (i.e. start of feature edge
325 // or 'knot')
326 // 2. seed particle starting at keepPoint going to this non-manifold point
327 // 3. track particles to their non-manifold point
328 //
329 // 4. track particles across their connected feature edges, marking all
330 // visited cells with their level (through trackingData)
331 // 5. do 4 until all edges have been visited.
332
333
334 // Find all start cells of features. Is done by tracking from keepPoint.
335 Cloud<trackedParticle> startPointCloud
336 (
337 mesh_,
338 Foam::zero{},
339 "startPointCloud"
340 );
341
342
343 // Features are identical on all processors. Number them so we know
344 // what to seed. Do this on only the processor that
345 // holds the keepPoint.
346
347 for (const point& keepPoint : keepPoints)
348 {
349 const label celli = mesh_.cellTree().findInside(keepPoint);
350
351 if (celli != -1)
352 {
353 // I am the processor that holds the keepPoint
354
355 forAll(features_, feati)
356 {
357 const edgeMesh& featureMesh = features_[feati];
358 const label featureLevel = features_.levels()[feati][0];
359 const labelListList& pointEdges = featureMesh.pointEdges();
360
361 // Find regions on edgeMesh
362 labelList edgeRegion;
363 label nRegions = featureMesh.regions(edgeRegion);
364
365
366 bitSet regionVisited(nRegions);
367
368
369 // 1. Seed all 'knots' in edgeMesh
370
371
372 forAll(pointEdges, pointi)
373 {
374 if (pointEdges[pointi].size() != 2)
375 {
377 {
378 Pout<< "Adding particle from point:" << pointi
379 << " coord:" << featureMesh.points()[pointi]
380 << " since number of emanating edges:"
381 << pointEdges[pointi].size()
382 << endl;
383 }
384
385 // Non-manifold point. Create particle.
386 startPointCloud.addParticle
387 (
389 (
390 mesh_,
391 keepPoint,
392 celli,
393 featureMesh.points()[pointi], // endpos
394 featureLevel, // level
395 feati, // featureMesh
396 pointi, // end point
397 -1 // feature edge
398 )
399 );
400
401 // Mark
402 if (pointEdges[pointi].size() > 0)
403 {
404 label e0 = pointEdges[pointi][0];
405 label regioni = edgeRegion[e0];
406 regionVisited.set(regioni);
407 }
408 }
409 }
410
411
412 // 2. Any regions that have not been visited at all? These can
413 // only be circular regions!
414 forAll(featureMesh.edges(), edgei)
415 {
416 if (regionVisited.set(edgeRegion[edgei]))
417 {
418 const edge& e = featureMesh.edges()[edgei];
419 label pointi = e.start();
421 {
422 Pout<< "Adding particle from point:" << pointi
423 << " coord:" << featureMesh.points()[pointi]
424 << " on circular region:" << edgeRegion[edgei]
425 << endl;
426 }
427
428 // Non-manifold point. Create particle.
429 startPointCloud.addParticle
430 (
432 (
433 mesh_,
434 keepPoint,
435 celli,
436 featureMesh.points()[pointi], // endpos
437 featureLevel, // level
438 feati, // featureMesh
439 pointi, // end point
440 -1 // feature edge
441 )
442 );
443 }
444 }
445 }
446 }
447 }
448
449
450 // Largest refinement level of any feature passed through
451 maxFeatureLevel = labelList(mesh_.nCells(), -1);
452
453 // Whether edge has been visited.
454 List<bitSet> featureEdgeVisited(features_.size());
455
456 forAll(features_, featI)
457 {
458 featureEdgeVisited[featI].setSize(features_[featI].edges().size());
459 featureEdgeVisited[featI] = false;
460 }
461
462 // Database to pass into trackedParticle::move
464 (
465 startPointCloud,
466 maxFeatureLevel,
467 featureEdgeVisited
468 );
469
470
471 // Track all particles to their end position (= starting feature point)
472 // Note that the particle might have started on a different processor
473 // so this will transfer across nicely until we can start tracking proper.
474 scalar maxTrackLen = 2.0*mesh_.bounds().mag();
475
477 {
478 Pout<< "Tracking " << startPointCloud.size()
479 << " particles over distance " << maxTrackLen
480 << " to find the starting cell" << endl;
481 }
482 startPointCloud.move(startPointCloud, td, maxTrackLen);
483
484
485 // Reset levels
486 maxFeatureLevel = -1;
487 forAll(features_, featI)
488 {
489 featureEdgeVisited[featI] = false;
490 }
491
492
493 // Start with empty cloud
494 Cloud<trackedParticle> cloud(mesh_, Foam::zero{}, "featureCloud");
495
497 {
498 Pout<< "Constructing cloud for cell marking" << endl;
499 }
500
501 for (const trackedParticle& startTp : startPointCloud)
502 {
503 const label featI = startTp.i();
504 const label pointI = startTp.j();
505
506 const edgeMesh& featureMesh = features_[featI];
507 const labelList& pEdges = featureMesh.pointEdges()[pointI];
508
509 // Now shoot particles down all pEdges.
510 forAll(pEdges, pEdgeI)
511 {
512 label edgeI = pEdges[pEdgeI];
513
514 if (featureEdgeVisited[featI].set(edgeI))
515 {
516 // Unvisited edge. Make the particle go to the other point
517 // on the edge.
518
519 const edge& e = featureMesh.edges()[edgeI];
520 label otherPointi = e.otherVertex(pointI);
521
522 trackedParticle* tp(new trackedParticle(startTp));
523 tp->start() = tp->position();
524 tp->end() = featureMesh.points()[otherPointi];
525 tp->j() = otherPointi;
526 tp->k() = edgeI;
527
529 {
530 Pout<< "Adding particle for point:" << pointI
531 << " coord:" << tp->position()
532 << " feature:" << featI
533 << " to track to:" << tp->end()
534 << endl;
535 }
536
537 cloud.addParticle(tp);
538 }
539 }
540 }
541
542 startPointCloud.clear();
543
544
545 while (true)
546 {
547 // Track all particles to their end position.
549 {
550 Pout<< "Tracking " << cloud.size()
551 << " particles over distance " << maxTrackLen
552 << " to mark cells" << endl;
553 }
554 cloud.move(cloud, td, maxTrackLen);
555
556 // Make particle follow edge.
557 for (trackedParticle& tp : cloud)
558 {
559 const label featI = tp.i();
560 const label pointI = tp.j();
561
562 const edgeMesh& featureMesh = features_[featI];
563 const labelList& pEdges = featureMesh.pointEdges()[pointI];
564
565 // Particle now at pointI. Check connected edges to see which one
566 // we have to visit now.
567
568 bool keepParticle = false;
569
570 forAll(pEdges, i)
571 {
572 label edgeI = pEdges[i];
573
574 if (featureEdgeVisited[featI].set(edgeI))
575 {
576 // Unvisited edge. Make the particle go to the other point
577 // on the edge.
578
579 const edge& e = featureMesh.edges()[edgeI];
580 label otherPointi = e.otherVertex(pointI);
581
582 tp.start() = tp.position();
583 tp.end() = featureMesh.points()[otherPointi];
584 tp.j() = otherPointi;
585 tp.k() = edgeI;
586 keepParticle = true;
587 break;
588 }
589 }
590
591 if (!keepParticle)
592 {
593 // Particle at 'knot' where another particle already has been
594 // seeded. Delete particle.
595 cloud.deleteParticle(tp);
596 }
597 }
598
599
601 {
602 Pout<< "Remaining particles " << cloud.size() << endl;
603 }
604
605 if (!returnReduceOr(cloud.size()))
606 {
607 break;
608 }
609 }
610
611
612
613 //if (debug&meshRefinement::FEATURESEEDS)
614 //{
615 // forAll(maxFeatureLevel, cellI)
616 // {
617 // if (maxFeatureLevel[cellI] != -1)
618 // {
619 // Pout<< "Feature went through cell:" << cellI
620 // << " coord:" << mesh_.cellCentres()[cellI]
621 // << " level:" << maxFeatureLevel[cellI]
622 // << endl;
623 // }
624 // }
625 //}
626}
627
628
629// Calculates list of cells to refine based on intersection with feature edge.
630Foam::label Foam::meshRefinement::markFeatureRefinement
631(
632 const pointField& keepPoints,
633 const label nAllowRefine,
634
636 label& nRefine
637) const
638{
639 // Largest refinement level of any feature passed through
640 labelList maxFeatureLevel;
641 markFeatureCellLevel(keepPoints, maxFeatureLevel);
642
643 // See which cells to refine. maxFeatureLevel will hold highest level
644 // of any feature edge that passed through.
645
646 const labelList& cellLevel = meshCutter_.cellLevel();
647
648 label oldNRefine = nRefine;
649
650 forAll(maxFeatureLevel, cellI)
651 {
652 if (maxFeatureLevel[cellI] > cellLevel[cellI])
653 {
654 // Mark
655 if
656 (
657 !markForRefine
658 (
659 0, // surface (n/a)
660 nAllowRefine,
661 refineCell[cellI],
662 nRefine
663 )
664 )
665 {
666 // Reached limit
667 break;
668 }
669 }
670 }
671
672 if
673 (
674 returnReduce(nRefine, sumOp<label>())
675 > returnReduce(nAllowRefine, sumOp<label>())
676 )
677 {
678 Info<< "Reached refinement limit." << endl;
679 }
680
681 return returnReduce(nRefine-oldNRefine, sumOp<label>());
682}
683
684
685// Mark cells for distance-to-feature based refinement.
686Foam::label Foam::meshRefinement::markInternalDistanceToFeatureRefinement
687(
688 const label nAllowRefine,
689
691 label& nRefine
692) const
693{
694 const labelList& cellLevel = meshCutter_.cellLevel();
695 const pointField& cellCentres = mesh_.cellCentres();
696
697 // Detect if there are any distance shells
698 if (features_.maxDistance() <= 0.0)
699 {
700 return 0;
701 }
702
703 label oldNRefine = nRefine;
704
705 // Collect cells to test
706 pointField testCc(cellLevel.size()-nRefine);
707 labelList testLevels(cellLevel.size()-nRefine);
708 label testI = 0;
709
710 forAll(cellLevel, cellI)
711 {
712 if (refineCell[cellI] == -1)
713 {
714 testCc[testI] = cellCentres[cellI];
715 testLevels[testI] = cellLevel[cellI];
716 testI++;
717 }
718 }
719
720 // Do test to see whether cells are inside/outside shell with higher level
721 labelList maxLevel;
722 features_.findHigherLevel(testCc, testLevels, maxLevel);
723
724 // Mark for refinement. Note that we didn't store the original cellID so
725 // now just reloop in same order.
726 testI = 0;
727 forAll(cellLevel, cellI)
728 {
729 if (refineCell[cellI] == -1)
730 {
731 if (maxLevel[testI] > testLevels[testI])
732 {
733 bool reachedLimit = !markForRefine
734 (
735 maxLevel[testI], // mark with any positive value
736 nAllowRefine,
737 refineCell[cellI],
738 nRefine
739 );
740
741 if (reachedLimit)
742 {
743 if (debug)
744 {
745 Pout<< "Stopped refining internal cells"
746 << " since reaching my cell limit of "
747 << mesh_.nCells()+7*nRefine << endl;
748 }
749 break;
750 }
751 }
752 testI++;
753 }
754 }
755
756 if
757 (
758 returnReduce(nRefine, sumOp<label>())
759 > returnReduce(nAllowRefine, sumOp<label>())
760 )
761 {
762 Info<< "Reached refinement limit." << endl;
763 }
764
765 return returnReduce(nRefine-oldNRefine, sumOp<label>());
766}
767
768
769// Mark cells for non-surface intersection based refinement.
770Foam::label Foam::meshRefinement::markInternalRefinement
771(
772 const label nAllowRefine,
773
775 label& nRefine
776) const
777{
778 const labelList& cellLevel = meshCutter_.cellLevel();
779 const pointField& cellCentres = mesh_.cellCentres();
780
781 label oldNRefine = nRefine;
782
783 // Collect cells to test
784 pointField testCc(cellLevel.size()-nRefine);
785 labelList testLevels(cellLevel.size()-nRefine);
786 label testI = 0;
787
788 forAll(cellLevel, cellI)
789 {
790 if (refineCell[cellI] == -1)
791 {
792 testCc[testI] = cellCentres[cellI];
793 testLevels[testI] = cellLevel[cellI];
794 testI++;
795 }
796 }
797
798 // Do test to see whether cells are inside/outside shell with higher level
799 labelList maxLevel;
800 shells_.findHigherLevel(testCc, testLevels, maxLevel);
801
802 // Mark for refinement. Note that we didn't store the original cellID so
803 // now just reloop in same order.
804 testI = 0;
805 forAll(cellLevel, cellI)
806 {
807 if (refineCell[cellI] == -1)
808 {
809 if (maxLevel[testI] > testLevels[testI])
810 {
811 bool reachedLimit = !markForRefine
812 (
813 maxLevel[testI], // mark with any positive value
814 nAllowRefine,
815 refineCell[cellI],
816 nRefine
817 );
818
819 if (reachedLimit)
820 {
821 if (debug)
822 {
823 Pout<< "Stopped refining internal cells"
824 << " since reaching my cell limit of "
825 << mesh_.nCells()+7*nRefine << endl;
826 }
827 break;
828 }
829 }
830 testI++;
831 }
832 }
833
834 if
835 (
836 returnReduce(nRefine, sumOp<label>())
837 > returnReduce(nAllowRefine, sumOp<label>())
838 )
839 {
840 Info<< "Reached refinement limit." << endl;
841 }
842
843 return returnReduce(nRefine-oldNRefine, sumOp<label>());
844}
845
846
847Foam::label Foam::meshRefinement::unmarkInternalRefinement
848(
850 label& nRefine
851) const
852{
853 const labelList& cellLevel = meshCutter_.cellLevel();
854 const pointField& cellCentres = mesh_.cellCentres();
855
856 label oldNRefine = nRefine;
857
858 // Collect cells to test
859 pointField testCc(nRefine);
860 labelList testLevels(nRefine);
861 label testI = 0;
862
863 forAll(cellLevel, cellI)
864 {
865 if (refineCell[cellI] >= 0)
866 {
867 testCc[testI] = cellCentres[cellI];
868 testLevels[testI] = cellLevel[cellI];
869 testI++;
870 }
871 }
872
873 // Do test to see whether cells are inside/outside shell with higher level
874 labelList levelShell;
875 limitShells_.findLevel(testCc, testLevels, levelShell);
876
877 // Mark for refinement. Note that we didn't store the original cellID so
878 // now just reloop in same order.
879 testI = 0;
880 forAll(cellLevel, cellI)
881 {
882 if (refineCell[cellI] >= 0)
883 {
884 if (levelShell[testI] != -1)
885 {
886 refineCell[cellI] = -1;
887 nRefine--;
888 }
889 testI++;
890 }
891 }
892
893 return returnReduce(oldNRefine-nRefine, sumOp<label>());
894}
895
896
897// Collect faces that are intersected and whose neighbours aren't yet marked
898// for refinement.
899Foam::labelList Foam::meshRefinement::getRefineCandidateFaces
900(
901 const labelList& refineCell
902) const
903{
904 labelList testFaces(mesh_.nFaces());
905
906 label nTest = 0;
907
908 const labelList& surfIndex = surfaceIndex();
909
910 labelList boundaryRefineCell;
911 syncTools::swapBoundaryCellList(mesh_, refineCell, boundaryRefineCell);
912
913 forAll(surfIndex, faceI)
914 {
915 if (surfIndex[faceI] != -1)
916 {
917 label own = mesh_.faceOwner()[faceI];
918
919 if (mesh_.isInternalFace(faceI))
920 {
921 label nei = mesh_.faceNeighbour()[faceI];
922
923 if (refineCell[own] == -1 || refineCell[nei] == -1)
924 {
925 testFaces[nTest++] = faceI;
926 }
927 }
928 else
929 {
930 const label bFacei = faceI - mesh_.nInternalFaces();
931 if (refineCell[own] == -1 || boundaryRefineCell[bFacei] == -1)
932 {
933 testFaces[nTest++] = faceI;
934 }
935 }
936 }
937 }
938 testFaces.setSize(nTest);
939
940 return testFaces;
941}
942
943
944// Mark cells for surface intersection based refinement.
945Foam::label Foam::meshRefinement::markSurfaceRefinement
946(
947 const label nAllowRefine,
948 const labelList& neiLevel,
949 const pointField& neiCc,
950
952 label& nRefine
953) const
954{
955 const labelList& cellLevel = meshCutter_.cellLevel();
956
957 label oldNRefine = nRefine;
958
959 // Use cached surfaceIndex_ to detect if any intersection. If so
960 // re-intersect to determine level wanted.
961
962 // Collect candidate faces
963 // ~~~~~~~~~~~~~~~~~~~~~~~
964
965 labelList testFaces(getRefineCandidateFaces(refineCell));
966
967 // Collect segments
968 // ~~~~~~~~~~~~~~~~
969
970 pointField start(testFaces.size());
971 pointField end(testFaces.size());
972 labelList minLevel(testFaces.size());
973
974 calcCellCellRays
975 (
976 neiCc,
977 neiLevel,
978 testFaces,
979 start,
980 end,
981 minLevel
982 );
983
984
985 // Do test for higher intersections
986 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
987
988 labelList surfaceHit;
989 labelList surfaceMinLevel;
990 surfaces_.findHigherIntersection
991 (
992 shells_,
993 start,
994 end,
995 minLevel,
996
997 surfaceHit,
998 surfaceMinLevel
999 );
1000
1001
1002 // Mark cells for refinement
1003 // ~~~~~~~~~~~~~~~~~~~~~~~~~
1004
1005 forAll(testFaces, i)
1006 {
1007 label faceI = testFaces[i];
1008
1009 label surfI = surfaceHit[i];
1010
1011 if (surfI != -1)
1012 {
1013 // Found intersection with surface with higher wanted
1014 // refinement. Check if the level field on that surface
1015 // specifies an even higher level. Note:this is weird. Should
1016 // do the check with the surfaceMinLevel whilst intersecting the
1017 // surfaces?
1018
1019 label own = mesh_.faceOwner()[faceI];
1020
1021 if (surfaceMinLevel[i] > cellLevel[own])
1022 {
1023 // Owner needs refining
1024 if
1025 (
1026 !markForRefine
1027 (
1028 surfI,
1029 nAllowRefine,
1030 refineCell[own],
1031 nRefine
1032 )
1033 )
1034 {
1035 break;
1036 }
1037 }
1038
1039 if (mesh_.isInternalFace(faceI))
1040 {
1041 label nei = mesh_.faceNeighbour()[faceI];
1042 if (surfaceMinLevel[i] > cellLevel[nei])
1043 {
1044 // Neighbour needs refining
1045 if
1046 (
1047 !markForRefine
1048 (
1049 surfI,
1050 nAllowRefine,
1051 refineCell[nei],
1052 nRefine
1053 )
1054 )
1055 {
1056 break;
1057 }
1058 }
1059 }
1060 }
1061 }
1062
1063 if
1064 (
1065 returnReduce(nRefine, sumOp<label>())
1066 > returnReduce(nAllowRefine, sumOp<label>())
1067 )
1068 {
1069 Info<< "Reached refinement limit." << endl;
1070 }
1071
1072 return returnReduce(nRefine-oldNRefine, sumOp<label>());
1073}
1074
1075
1076// Count number of matches of first argument in second argument
1077Foam::label Foam::meshRefinement::countMatches
1078(
1079 const List<point>& normals1,
1080 const List<point>& normals2,
1081 const scalar tol
1082) const
1083{
1084 label nMatches = 0;
1085
1086 forAll(normals1, i)
1087 {
1088 const vector& n1 = normals1[i];
1089
1090 forAll(normals2, j)
1091 {
1092 const vector& n2 = normals2[j];
1093
1094 if (magSqr(n1-n2) < tol)
1095 {
1096 nMatches++;
1097 break;
1098 }
1099 }
1100 }
1101
1102 return nMatches;
1103}
1104
1105
1106//bool Foam::meshRefinement::highCurvature
1107//(
1108// const scalar minCosAngle,
1109// const point& p0,
1110// const vector& n0,
1111// const point& p1,
1112// const vector& n1
1113//) const
1114//{
1115// return ((n0&n1) < minCosAngle);
1116//}
1117bool Foam::meshRefinement::highCurvature
1118(
1119 const scalar minCosAngle,
1120 const scalar lengthScale,
1121 const point& p0,
1122 const vector& n0,
1123 const point& p1,
1124 const vector& n1
1125) const
1126{
1127 const scalar cosAngle = (n0&n1);
1128
1129 if (cosAngle < minCosAngle)
1130 {
1131 // Sharp feature, independent of intersection points
1132 return true;
1133 }
1134 else if (cosAngle > 1-1e-6)
1135 {
1136 // Co-planar
1137 return false;
1138 }
1139 else if (lengthScale > SMALL)
1140 {
1141 // Calculate radius of curvature
1142
1143 vector axis(n0 ^ n1);
1144 const plane pl0(p0, (n0 ^ axis));
1145 const scalar r1 = pl0.normalIntersect(p1, n1);
1146
1147 //const point radiusPoint(p1+r1*n1);
1148 //DebugVar(radiusPoint);
1149 const plane pl1(p1, (n1 ^ axis));
1150 const scalar r0 = pl1.normalIntersect(p0, n0);
1151
1152 //const point radiusPoint(p0+r0*n0);
1153 //DebugVar(radiusPoint);
1154 //- Take average radius. Bit ad hoc
1155 const scalar radius = 0.5*(mag(r1)+mag(r0));
1156
1157 if (radius < lengthScale)
1158 {
1159 return true;
1160 }
1161 else
1162 {
1163 return false;
1164 }
1165 }
1166 else
1167 {
1168 return false;
1169 }
1170}
1171//XXXXX
1172
1173
1174// Mark cells for surface curvature based refinement.
1175Foam::label Foam::meshRefinement::markSurfaceCurvatureRefinement
1176(
1177 const scalar curvature,
1178 const label nAllowRefine,
1179 const labelList& neiLevel,
1180 const pointField& neiCc,
1181
1183 label& nRefine
1184) const
1185{
1186 const labelList& cellLevel = meshCutter_.cellLevel();
1187 const pointField& cellCentres = mesh_.cellCentres();
1188
1189 label oldNRefine = nRefine;
1190
1191 // 1. local test: any cell on more than one surface gets refined
1192 // (if its current level is < max of the surface max level)
1193
1194 // 2. 'global' test: any cell on only one surface with a neighbour
1195 // on a different surface gets refined (if its current level etc.)
1196
1197
1198 const bitSet isMasterFace(syncTools::getMasterFaces(mesh_));
1199
1200
1201 // Collect candidate faces (i.e. intersecting any surface and
1202 // owner/neighbour not yet refined.
1203 labelList testFaces(getRefineCandidateFaces(refineCell));
1204
1205 // Collect segments
1206 pointField start(testFaces.size());
1207 pointField end(testFaces.size());
1208 labelList minLevel(testFaces.size());
1209
1210 // Note: uses isMasterFace otherwise could be call to calcCellCellRays
1211 forAll(testFaces, i)
1212 {
1213 label faceI = testFaces[i];
1214
1215 label own = mesh_.faceOwner()[faceI];
1216
1217 if (mesh_.isInternalFace(faceI))
1218 {
1219 label nei = mesh_.faceNeighbour()[faceI];
1220
1221 start[i] = cellCentres[own];
1222 end[i] = cellCentres[nei];
1223 minLevel[i] = min(cellLevel[own], cellLevel[nei]);
1224 }
1225 else
1226 {
1227 label bFaceI = faceI - mesh_.nInternalFaces();
1228
1229 start[i] = cellCentres[own];
1230 end[i] = neiCc[bFaceI];
1231
1232 if (!isMasterFace[faceI])
1233 {
1234 std::swap(start[i], end[i]);
1235 }
1236
1237 minLevel[i] = min(cellLevel[own], neiLevel[bFaceI]);
1238 }
1239 }
1240
1241 // Extend segments a bit
1242 {
1243 const vectorField smallVec(ROOTSMALL*(end-start));
1244 start -= smallVec;
1245 end += smallVec;
1246 }
1247
1248
1249 // If no curvature supplied behave as before
1250 const bool hasCurvatureLevels = (max(surfaces_.maxCurvatureLevel()) > 0);
1251
1252
1253 // Test for all intersections (with surfaces of higher max level than
1254 // minLevel) and cache per cell the interesting inter
1255 labelListList cellSurfLevels(mesh_.nCells());
1256 List<pointList> cellSurfLocations(mesh_.nCells());
1257 List<vectorList> cellSurfNormals(mesh_.nCells());
1258
1259 {
1260 // Per segment the intersection point of the surfaces
1262 // Per segment the normals of the surfaces
1263 List<vectorList> surfaceNormal;
1264 // Per segment the list of levels of the surfaces
1265 labelListList surfaceLevel;
1266
1267 surfaces_.findAllIntersections
1268 (
1269 start,
1270 end,
1271 minLevel, // max level of surface has to be bigger
1272 // than min level of neighbouring cells
1273
1274 labelList(surfaces_.maxLevel().size(), 0), // min level
1275 surfaces_.maxLevel(),
1276
1278 surfaceNormal,
1279 surfaceLevel
1280 );
1281
1282
1283 // Sort the data according to surface location. This will guarantee
1284 // that on coupled faces both sides visit the intersections in
1285 // the same order so will decide the same
1286 labelList visitOrder;
1287 forAll(surfaceNormal, pointI)
1288 {
1289 pointList& pLocations = surfaceLocation[pointI];
1290 vectorList& pNormals = surfaceNormal[pointI];
1291 labelList& pLevel = surfaceLevel[pointI];
1292
1293 sortedOrder(pNormals, visitOrder, normalLess(pLocations));
1294
1295 pLocations = List<point>(pLocations, visitOrder);
1296 pNormals = List<vector>(pNormals, visitOrder);
1297 pLevel = labelUIndList(pLevel, visitOrder);
1298 }
1299
1300
1301 //- At some point could just return the intersected surface+region
1302 // and derive all the surface information (maxLevel, curvatureLevel)
1303 // from that - we're now doing it inside refinementSurfaces itself.
1305 //List<labelList> hitSurface;
1306 //List<pointList> hitLocation;
1307 //List<labelList> hitRegion;
1308 //List<vectorField> hitNormal;
1309 //surfaces_.findAllIntersections
1310 //(
1311 // identity(surfaces_.surfaces()), // all refinement geometries
1312 // start,
1313 // end,
1314 //
1315 // hitSurface,
1316 // hitLocation,
1317 // hitRegion,
1318 // hitNormal
1319 //);
1320 //
1323 //const auto& maxLevel = surfaces_.maxLevel();
1324 //labelList visitOrder;
1325 //DynamicList<label> valid;
1326 //forAll(hitSurface, segmenti)
1327 //{
1328 // const label meshLevel = minLevel[segmenti];
1329 //
1330 // auto& fSurface = hitSurface[segmenti];
1331 // auto& fLocation = hitLocation[segmenti];
1332 // auto& fRegion = hitRegion[segmenti];
1333 // auto& fNormal = hitNormal[segmenti];
1334 //
1335 // // Sort the data according to intersection location. This will
1336 // // guarantee
1337 // // that on coupled faces both sides visit the intersections in
1338 // // the same order so will decide the same
1339 // sortedOrder(fLocation, visitOrder, normalLess(hfLocation));
1340 // fLocation = List<point>(fLocation, visitOrder);
1341 // fSurface = labelUIndList(fSurface, visitOrder);
1342 // fRegion = labelUIndList(fRegion, visitOrder);
1343 // fNormal = List<vector>(fNormal, visitOrder);
1344 //
1345 // // Filter out any intersections with surfaces outside cell level.
1346 // // Note that min refinement level of surfaces is ignored.
1347 // valid.clear();
1348 // forAll(fSurface, hiti)
1349 // {
1350 // const label regioni =
1351 // surfaces_.globalRegion(fSurface[hiti], fRegion[hiti]);
1352 // if (meshLevel < maxLevel[regioni]) //&& >= minLevel(regioni)
1353 // {
1354 // valid.append(hiti);
1355 // }
1356 // }
1357 // fLocation = List<point>(fLocation, valid);
1358 // fSurface = labelUIndList(fSurface, valid);
1359 // fRegion = labelUIndList(fRegion, valid);
1360 // fNormal = List<vector>(fNormal, valid);
1361 //}
1362
1363 // Clear out unnecessary data
1364 start.clear();
1365 end.clear();
1366 minLevel.clear();
1367
1368 // Convert face-wise data to cell.
1369 forAll(testFaces, i)
1370 {
1371 label faceI = testFaces[i];
1372 label own = mesh_.faceOwner()[faceI];
1373
1374 const pointList& fPoints = surfaceLocation[i];
1375 const vectorList& fNormals = surfaceNormal[i];
1376 const labelList& fLevels = surfaceLevel[i];
1377
1378 forAll(fNormals, hitI)
1379 {
1380 if (fLevels[hitI] > cellLevel[own])
1381 {
1382 cellSurfLevels[own].append(fLevels[hitI]);
1383 cellSurfLocations[own].append(fPoints[hitI]);
1384 cellSurfNormals[own].append(fNormals[hitI]);
1385 }
1386
1387 if (mesh_.isInternalFace(faceI))
1388 {
1389 label nei = mesh_.faceNeighbour()[faceI];
1390 if (fLevels[hitI] > cellLevel[nei])
1391 {
1392 cellSurfLevels[nei].append(fLevels[hitI]);
1393 cellSurfLocations[nei].append(fPoints[hitI]);
1394 cellSurfNormals[nei].append(fNormals[hitI]);
1395 }
1396 }
1397 }
1398 }
1399 }
1400
1401
1402
1403 // Bit of statistics
1404 if (debug)
1405 {
1406 label nSet = 0;
1407 label nNormals = 0;
1408 forAll(cellSurfNormals, cellI)
1409 {
1410 const vectorList& normals = cellSurfNormals[cellI];
1411 if (normals.size())
1412 {
1413 nSet++;
1414 nNormals += normals.size();
1415 }
1416 }
1417 reduce(nSet, sumOp<label>());
1418 reduce(nNormals, sumOp<label>());
1419 Info<< "markSurfaceCurvatureRefinement :"
1420 << " cells:" << mesh_.globalData().nTotalCells()
1421 << " of which with normals:" << nSet
1422 << " ; total normals stored:" << nNormals
1423 << endl;
1424 }
1425
1426
1427
1428 bool reachedLimit = false;
1429
1430
1431 // 1. Check normals per cell
1432 // ~~~~~~~~~~~~~~~~~~~~~~~~~
1433
1434 for
1435 (
1436 label cellI = 0;
1437 !reachedLimit && cellI < cellSurfLocations.size();
1438 cellI++
1439 )
1440 {
1441 const pointList& points = cellSurfLocations[cellI];
1442 const vectorList& normals = cellSurfNormals[cellI];
1443 const labelList& levels = cellSurfLevels[cellI];
1444
1445 // Current cell size
1446 const scalar cellSize =
1447 meshCutter_.level0EdgeLength()/pow(2.0, cellLevel[cellI]);
1448
1449 // n^2 comparison of all normals in a cell
1450 for (label i = 0; !reachedLimit && i < normals.size(); i++)
1451 {
1452 for (label j = i+1; !reachedLimit && j < normals.size(); j++)
1453 {
1454 // TBD: calculate curvature size (if curvatureLevel specified)
1455 // and pass in instead of cellSize
1456
1457 //if ((normals[i] & normals[j]) < curvature)
1458 if
1459 (
1460 highCurvature
1461 (
1462 curvature,
1463 (hasCurvatureLevels ? cellSize : scalar(0)),
1464 points[i],
1465 normals[i],
1466 points[j],
1467 normals[j]
1468 )
1469 )
1470 {
1471 const label maxLevel = max(levels[i], levels[j]);
1472
1473 if (cellLevel[cellI] < maxLevel)
1474 {
1475 if
1476 (
1477 !markForRefine
1478 (
1479 maxLevel,
1480 nAllowRefine,
1481 refineCell[cellI],
1482 nRefine
1483 )
1484 )
1485 {
1486 if (debug)
1487 {
1488 Pout<< "Stopped refining since reaching my cell"
1489 << " limit of " << mesh_.nCells()+7*nRefine
1490 << endl;
1491 }
1492 reachedLimit = true;
1493 break;
1494 }
1495 }
1496 }
1497 }
1498 }
1499 }
1500
1501
1502
1503 // 2. Find out a measure of surface curvature
1504 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1505 // Look at normals between neighbouring surfaces
1506 // Loop over all faces. Could only be checkFaces, except if they're coupled
1507
1508 // Internal faces
1509 for
1510 (
1511 label faceI = 0;
1512 !reachedLimit && faceI < mesh_.nInternalFaces();
1513 faceI++
1514 )
1515 {
1516 label own = mesh_.faceOwner()[faceI];
1517 label nei = mesh_.faceNeighbour()[faceI];
1518
1519 const pointList& ownPoints = cellSurfLocations[own];
1520 const vectorList& ownNormals = cellSurfNormals[own];
1521 const labelList& ownLevels = cellSurfLevels[own];
1522 const pointList& neiPoints = cellSurfLocations[nei];
1523 const vectorList& neiNormals = cellSurfNormals[nei];
1524 const labelList& neiLevels = cellSurfLevels[nei];
1525
1526
1527 // Special case: owner normals set is a subset of the neighbour
1528 // normals. Do not do curvature refinement since other cell's normals
1529 // do not add any information. Avoids weird corner extensions of
1530 // refinement regions.
1531 bool ownIsSubset =
1532 countMatches(ownNormals, neiNormals)
1533 == ownNormals.size();
1534
1535 bool neiIsSubset =
1536 countMatches(neiNormals, ownNormals)
1537 == neiNormals.size();
1538
1539
1540 if (!ownIsSubset && !neiIsSubset)
1541 {
1542 // Current average cell size. Is min? max? average?
1543 const scalar cellSize =
1544 meshCutter_.level0EdgeLength()
1545 / pow(2.0, min(cellLevel[own], cellLevel[nei]));
1546
1547 // n^2 comparison of between ownNormals and neiNormals
1548 for (label i = 0; !reachedLimit && i < ownNormals.size(); i++)
1549 {
1550 for (label j = 0; !reachedLimit && j < neiNormals.size(); j++)
1551 {
1552 // Have valid data on both sides. Check curvature.
1553 //if ((ownNormals[i] & neiNormals[j]) < curvature)
1554 if
1555 (
1556 highCurvature
1557 (
1558 curvature,
1559 (hasCurvatureLevels ? cellSize : scalar(0)),
1560 ownPoints[i],
1561 ownNormals[i],
1562 neiPoints[j],
1563 neiNormals[j]
1564 )
1565 )
1566 {
1567 // See which side to refine.
1568 if (cellLevel[own] < ownLevels[i])
1569 {
1570 if
1571 (
1572 !markForRefine
1573 (
1574 ownLevels[i],
1575 nAllowRefine,
1576 refineCell[own],
1577 nRefine
1578 )
1579 )
1580 {
1581 if (debug)
1582 {
1583 Pout<< "Stopped refining since reaching"
1584 << " my cell limit of "
1585 << mesh_.nCells()+7*nRefine << endl;
1586 }
1587 reachedLimit = true;
1588 break;
1589 }
1590 }
1591 if (cellLevel[nei] < neiLevels[j])
1592 {
1593 if
1594 (
1595 !markForRefine
1596 (
1597 neiLevels[j],
1598 nAllowRefine,
1599 refineCell[nei],
1600 nRefine
1601 )
1602 )
1603 {
1604 if (debug)
1605 {
1606 Pout<< "Stopped refining since reaching"
1607 << " my cell limit of "
1608 << mesh_.nCells()+7*nRefine << endl;
1609 }
1610 reachedLimit = true;
1611 break;
1612 }
1613 }
1614 }
1615 }
1616 }
1617 }
1618 }
1619
1620
1621 // Send over surface point/normal to neighbour cell.
1622// labelListList neiSurfaceLevel;
1623// syncTools::swapBoundaryCellList(mesh_, cellSurfLevels, neiSurfaceLevel);
1624 List<vectorList> neiSurfacePoints;
1625 syncTools::swapBoundaryCellList(mesh_, cellSurfLocations, neiSurfacePoints);
1626 List<vectorList> neiSurfaceNormals;
1627 syncTools::swapBoundaryCellList(mesh_, cellSurfNormals, neiSurfaceNormals);
1628
1629 // Boundary faces
1630 for
1631 (
1632 label faceI = mesh_.nInternalFaces();
1633 !reachedLimit && faceI < mesh_.nFaces();
1634 faceI++
1635 )
1636 {
1637 label own = mesh_.faceOwner()[faceI];
1638 label bFaceI = faceI - mesh_.nInternalFaces();
1639
1640 const pointList& ownPoints = cellSurfLocations[own];
1641 const vectorList& ownNormals = cellSurfNormals[own];
1642 const labelList& ownLevels = cellSurfLevels[own];
1643
1644 const pointList& neiPoints = neiSurfacePoints[bFaceI];
1645 const vectorList& neiNormals = neiSurfaceNormals[bFaceI];
1646 //const labelList& neiLevels = neiSurfaceLevel[bFacei];
1647
1648 // Special case: owner normals set is a subset of the neighbour
1649 // normals. Do not do curvature refinement since other cell's normals
1650 // do not add any information. Avoids weird corner extensions of
1651 // refinement regions.
1652 bool ownIsSubset =
1653 countMatches(ownNormals, neiNormals)
1654 == ownNormals.size();
1655
1656 bool neiIsSubset =
1657 countMatches(neiNormals, ownNormals)
1658 == neiNormals.size();
1659
1660
1661 if (!ownIsSubset && !neiIsSubset)
1662 {
1663 // Current average cell size. Is min? max? average?
1664 const scalar cellSize =
1665 meshCutter_.level0EdgeLength()
1666 / pow(2.0, min(cellLevel[own], neiLevel[bFaceI]));
1667
1668 // n^2 comparison of between ownNormals and neiNormals
1669 for (label i = 0; !reachedLimit && i < ownNormals.size(); i++)
1670 {
1671 for (label j = 0; !reachedLimit && j < neiNormals.size(); j++)
1672 {
1673 // Have valid data on both sides. Check curvature.
1674 //if ((ownNormals[i] & neiNormals[j]) < curvature)
1675 if
1676 (
1677 highCurvature
1678 (
1679 curvature,
1680 (hasCurvatureLevels ? cellSize : scalar(0)),
1681 ownPoints[i],
1682 ownNormals[i],
1683 neiPoints[j],
1684 neiNormals[j]
1685 )
1686 )
1687 {
1688 if (cellLevel[own] < ownLevels[i])
1689 {
1690 if
1691 (
1692 !markForRefine
1693 (
1694 ownLevels[i],
1695 nAllowRefine,
1696 refineCell[own],
1697 nRefine
1698 )
1699 )
1700 {
1701 if (debug)
1702 {
1703 Pout<< "Stopped refining since reaching"
1704 << " my cell limit of "
1705 << mesh_.nCells()+7*nRefine
1706 << endl;
1707 }
1708 reachedLimit = true;
1709 break;
1710 }
1711 }
1712 }
1713 }
1714 }
1715 }
1716 }
1717
1718
1719 if
1720 (
1721 returnReduce(nRefine, sumOp<label>())
1722 > returnReduce(nAllowRefine, sumOp<label>())
1723 )
1724 {
1725 Info<< "Reached refinement limit." << endl;
1726 }
1727
1728 return returnReduce(nRefine-oldNRefine, sumOp<label>());
1729}
1730
1731
1733(
1734 const scalar planarCos,
1735 const vector& point0,
1736 const vector& normal0,
1737
1738 const vector& point1,
1739 const vector& normal1
1740) const
1741{
1742 //- Hits differ and angles are oppositeish and
1743 // hits have a normal distance
1744 vector d = point1-point0;
1745 scalar magD = mag(d);
1746
1747 if (magD > mergeDistance())
1748 {
1749 scalar cosAngle = (normal0 & normal1);
1750
1751 vector avg = Zero;
1752 if (cosAngle < (-1+planarCos))
1753 {
1754 // Opposite normals
1755 avg = 0.5*(normal0-normal1);
1756 }
1757 else if (cosAngle > (1-planarCos))
1758 {
1759 avg = 0.5*(normal0+normal1);
1760 }
1761
1762 if (avg != vector::zero)
1763 {
1764 avg /= mag(avg);
1765
1766 // Check normal distance of intersection locations
1767 if (mag(avg&d) > mergeDistance())
1768 {
1769 return true;
1770 }
1771 }
1772 }
1773
1774 return false;
1775}
1776
1777
1778// Mark small gaps
1780(
1781 const scalar planarCos,
1782 const label level0,
1783 const vector& point0,
1784 const vector& normal0,
1785
1786 const label level1,
1787 const vector& point1,
1788 const vector& normal1
1789) const
1790{
1791 //const scalar edge0Len = meshCutter_.level0EdgeLength();
1793 //- Hits differ and angles are oppositeish and
1794 // hits have a normal distance
1795 vector d = point1-point0;
1796 scalar magD = mag(d);
1797
1798 if (magD > mergeDistance())
1799 {
1800 scalar cosAngle = (normal0 & normal1);
1801
1802 vector avgNormal = Zero;
1803 if (cosAngle < (-1+planarCos))
1804 {
1805 // Opposite normals
1806 avgNormal = 0.5*(normal0-normal1);
1807 }
1808 else if (cosAngle > (1-planarCos))
1809 {
1810 avgNormal = 0.5*(normal0+normal1);
1811 }
1812
1813 if (avgNormal != vector::zero)
1814 {
1815 avgNormal /= mag(avgNormal);
1816
1817 //scalar normalDist = mag(d&avgNormal);
1818 //const scalar avgCellSize = edge0Len/pow(2.0, 0.5*(level0+level1));
1819 //if (normalDist > 1*avgCellSize)
1820 //{
1821 // Pout<< "** skipping large distance " << endl;
1822 // return false;
1823 //}
1824
1825 // Check average normal with respect to intersection locations
1826 if (mag(avgNormal&d/magD) > Foam::cos(degToRad(45.0)))
1827 {
1828 return true;
1829 }
1830 }
1831 }
1832
1833 return false;
1834}
1835
1836
1837bool Foam::meshRefinement::checkProximity
1838(
1839 const scalar planarCos,
1840 const label nAllowRefine,
1841
1842 const label surfaceLevel, // current intersection max level
1843 const vector& surfaceLocation, // current intersection location
1844 const vector& surfaceNormal, // current intersection normal
1845
1846 const label cellI,
1847
1848 label& cellMaxLevel, // cached max surface level for this cell
1849 vector& cellMaxLocation, // cached surface location for this cell
1850 vector& cellMaxNormal, // cached surface normal for this cell
1851
1853 label& nRefine
1854) const
1855{
1856 const labelList& cellLevel = meshCutter_.cellLevel();
1857
1858 // Test if surface applicable
1859 if (surfaceLevel > cellLevel[cellI])
1860 {
1861 if (cellMaxLevel == -1)
1862 {
1863 // First visit of cell. Store
1864 cellMaxLevel = surfaceLevel;
1865 cellMaxLocation = surfaceLocation;
1866 cellMaxNormal = surfaceNormal;
1867 }
1868 else
1869 {
1870 // Second or more visit.
1871 // Check if
1872 // - different location
1873 // - opposite surface
1874
1875 bool closeSurfaces = isNormalGap
1876 (
1877 planarCos,
1878 cellLevel[cellI], //cellMaxLevel,
1879 cellMaxLocation,
1880 cellMaxNormal,
1881 cellLevel[cellI], //surfaceLevel,
1883 surfaceNormal
1884 );
1885
1886 // Set normal to that of highest surface. Not really necessary
1887 // over here but we reuse cellMax info when doing coupled faces.
1888 if (surfaceLevel > cellMaxLevel)
1889 {
1890 cellMaxLevel = surfaceLevel;
1891 cellMaxLocation = surfaceLocation;
1892 cellMaxNormal = surfaceNormal;
1893 }
1894
1895
1896 if (closeSurfaces)
1897 {
1898 //Pout<< "Found gap:" << nl
1899 // << " location:" << surfaceLocation
1900 // << "\tnormal :" << surfaceNormal << nl
1902 // << "\tnormal :" << cellMaxNormal << nl
1903 // << "\tcos :" << (surfaceNormal&cellMaxNormal) << nl
1904 // << endl;
1905
1906 return markForRefine
1907 (
1908 surfaceLevel, // mark with any non-neg number.
1909 nAllowRefine,
1910 refineCell[cellI],
1911 nRefine
1912 );
1913 }
1914 }
1915 }
1916
1917 // Did not reach refinement limit.
1918 return true;
1919}
1920
1921
1922Foam::label Foam::meshRefinement::markProximityRefinement
1923(
1924 const scalar planarCos,
1925
1926 const labelList& surfaceMinLevel,
1927 const labelList& surfaceMaxLevel,
1928
1929 const label nAllowRefine,
1930 const labelList& neiLevel,
1931 const pointField& neiCc,
1932
1934 label& nRefine
1935) const
1936{
1937 const labelList& cellLevel = meshCutter_.cellLevel();
1938
1939 label oldNRefine = nRefine;
1940
1941 // 1. local test: any cell on more than one surface gets refined
1942 // (if its current level is < max of the surface max level)
1943
1944 // 2. 'global' test: any cell on only one surface with a neighbour
1945 // on a different surface gets refined (if its current level etc.)
1946
1947
1948 // Collect candidate faces (i.e. intersecting any surface and
1949 // owner/neighbour not yet refined.
1950 const labelList testFaces(getRefineCandidateFaces(refineCell));
1951
1952 // Collect segments
1953 pointField start(testFaces.size());
1954 pointField end(testFaces.size());
1955 labelList minLevel(testFaces.size());
1956
1957 calcCellCellRays
1958 (
1959 neiCc,
1960 neiLevel,
1961 testFaces,
1962 start,
1963 end,
1964 minLevel
1965 );
1966
1967
1968 // Test for all intersections (with surfaces of higher gap level than
1969 // minLevel) and cache per cell the max surface level and the local normal
1970 // on that surface.
1971 labelList cellMaxLevel(mesh_.nCells(), -1);
1972 vectorField cellMaxNormal(mesh_.nCells(), Zero);
1973 pointField cellMaxLocation(mesh_.nCells(), Zero);
1974
1975 {
1976 // Per segment the normals of the surfaces
1978 List<vectorList> surfaceNormal;
1979 // Per segment the list of levels of the surfaces
1980 labelListList surfaceLevel;
1981
1982 surfaces_.findAllIntersections
1983 (
1984 start,
1985 end,
1986 minLevel, // gap level of surface has to be bigger
1987 // than min level of neighbouring cells
1988
1989 surfaceMinLevel,
1990 surfaceMaxLevel,
1991
1993 surfaceNormal,
1994 surfaceLevel
1995 );
1996 // Clear out unnecessary data
1997 start.clear();
1998 end.clear();
1999 minLevel.clear();
2000
2001
2003
2004 forAll(testFaces, i)
2005 {
2006 label faceI = testFaces[i];
2007 label own = mesh_.faceOwner()[faceI];
2008
2009 const labelList& fLevels = surfaceLevel[i];
2010 const vectorList& fPoints = surfaceLocation[i];
2011 const vectorList& fNormals = surfaceNormal[i];
2012
2013 forAll(fLevels, hitI)
2014 {
2015 checkProximity
2016 (
2017 planarCos,
2018 nAllowRefine,
2019
2020 fLevels[hitI],
2021 fPoints[hitI],
2022 fNormals[hitI],
2023
2024 own,
2025 cellMaxLevel[own],
2026 cellMaxLocation[own],
2027 cellMaxNormal[own],
2028
2029 refineCell,
2030 nRefine
2031 );
2032 }
2033
2034 if (mesh_.isInternalFace(faceI))
2035 {
2036 label nei = mesh_.faceNeighbour()[faceI];
2037
2038 forAll(fLevels, hitI)
2039 {
2040 checkProximity
2041 (
2042 planarCos,
2043 nAllowRefine,
2044
2045 fLevels[hitI],
2046 fPoints[hitI],
2047 fNormals[hitI],
2048
2049 nei,
2050 cellMaxLevel[nei],
2051 cellMaxLocation[nei],
2052 cellMaxNormal[nei],
2053
2054 refineCell,
2055 nRefine
2056 );
2057 }
2058 }
2059 }
2060 }
2061
2062 // 2. Find out a measure of surface curvature and region edges.
2063 // Send over surface region and surface normal to neighbour cell.
2064
2065 labelList neiBndMaxLevel(mesh_.nBoundaryFaces());
2066 pointField neiBndMaxLocation(mesh_.nBoundaryFaces());
2067 vectorField neiBndMaxNormal(mesh_.nBoundaryFaces());
2068
2069 for (label faceI = mesh_.nInternalFaces(); faceI < mesh_.nFaces(); faceI++)
2070 {
2071 label bFaceI = faceI-mesh_.nInternalFaces();
2072 label own = mesh_.faceOwner()[faceI];
2073
2074 neiBndMaxLevel[bFaceI] = cellMaxLevel[own];
2075 neiBndMaxLocation[bFaceI] = cellMaxLocation[own];
2076 neiBndMaxNormal[bFaceI] = cellMaxNormal[own];
2077 }
2078 syncTools::swapBoundaryFaceList(mesh_, neiBndMaxLevel);
2079 syncTools::swapBoundaryFaceList(mesh_, neiBndMaxLocation);
2080 syncTools::swapBoundaryFaceList(mesh_, neiBndMaxNormal);
2081
2082 // Loop over all faces. Could only be checkFaces.. except if they're coupled
2083
2084 // Internal faces
2085 for (label faceI = 0; faceI < mesh_.nInternalFaces(); faceI++)
2086 {
2087 label own = mesh_.faceOwner()[faceI];
2088 label nei = mesh_.faceNeighbour()[faceI];
2089
2090 if (cellMaxLevel[own] != -1 && cellMaxLevel[nei] != -1)
2091 {
2092 // Have valid data on both sides. Check planarCos.
2093 if
2094 (
2095 isNormalGap
2096 (
2097 planarCos,
2098 cellLevel[own], //cellMaxLevel[own],
2099 cellMaxLocation[own],
2100 cellMaxNormal[own],
2101 cellLevel[nei], //cellMaxLevel[nei],
2102 cellMaxLocation[nei],
2103 cellMaxNormal[nei]
2104 )
2105 )
2106 {
2107 // See which side to refine
2108 if (cellLevel[own] < cellMaxLevel[own])
2109 {
2110 if
2111 (
2112 !markForRefine
2113 (
2114 cellMaxLevel[own],
2115 nAllowRefine,
2116 refineCell[own],
2117 nRefine
2118 )
2119 )
2120 {
2121 if (debug)
2122 {
2123 Pout<< "Stopped refining since reaching my cell"
2124 << " limit of " << mesh_.nCells()+7*nRefine
2125 << endl;
2126 }
2127 break;
2128 }
2129 }
2130
2131 if (cellLevel[nei] < cellMaxLevel[nei])
2132 {
2133 if
2134 (
2135 !markForRefine
2136 (
2137 cellMaxLevel[nei],
2138 nAllowRefine,
2139 refineCell[nei],
2140 nRefine
2141 )
2142 )
2143 {
2144 if (debug)
2145 {
2146 Pout<< "Stopped refining since reaching my cell"
2147 << " limit of " << mesh_.nCells()+7*nRefine
2148 << endl;
2149 }
2150 break;
2151 }
2152 }
2153 }
2154 }
2155 }
2156 // Boundary faces
2157 for (label faceI = mesh_.nInternalFaces(); faceI < mesh_.nFaces(); faceI++)
2158 {
2159 label own = mesh_.faceOwner()[faceI];
2160 label bFaceI = faceI - mesh_.nInternalFaces();
2161
2162 if (cellLevel[own] < cellMaxLevel[own] && neiBndMaxLevel[bFaceI] != -1)
2163 {
2164 // Have valid data on both sides. Check planarCos.
2165 if
2166 (
2167 isNormalGap
2168 (
2169 planarCos,
2170 cellLevel[own], //cellMaxLevel[own],
2171 cellMaxLocation[own],
2172 cellMaxNormal[own],
2173 neiLevel[bFaceI], //neiBndMaxLevel[bFaceI],
2174 neiBndMaxLocation[bFaceI],
2175 neiBndMaxNormal[bFaceI]
2176 )
2177 )
2178 {
2179 if
2180 (
2181 !markForRefine
2182 (
2183 cellMaxLevel[own],
2184 nAllowRefine,
2185 refineCell[own],
2186 nRefine
2187 )
2188 )
2189 {
2190 if (debug)
2191 {
2192 Pout<< "Stopped refining since reaching my cell"
2193 << " limit of " << mesh_.nCells()+7*nRefine
2194 << endl;
2195 }
2196 break;
2197 }
2198 }
2199 }
2200 }
2201
2202 if
2203 (
2204 returnReduce(nRefine, sumOp<label>())
2205 > returnReduce(nAllowRefine, sumOp<label>())
2206 )
2207 {
2208 Info<< "Reached refinement limit." << endl;
2209 }
2210
2211 return returnReduce(nRefine-oldNRefine, sumOp<label>());
2212}
2213
2214
2215// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
2216
2217// Calculate list of cells to refine. Gets for any edge (start - end)
2218// whether it intersects the surface. Does more accurate test and checks
2219// the wanted level on the surface intersected.
2220// Does approximate precalculation of how many cells can be refined before
2221// hitting overall limit maxGlobalCells.
2223(
2224 const pointField& keepPoints,
2225 const scalar curvature,
2226 const scalar planarAngle,
2227
2228 const bool featureRefinement,
2229 const bool featureDistanceRefinement,
2230 const bool internalRefinement,
2231 const bool surfaceRefinement,
2232 const bool curvatureRefinement,
2233 const bool smallFeatureRefinement,
2234 const bool gapRefinement,
2235 const bool bigGapRefinement,
2236 const bool spreadGapSize,
2237 const label maxGlobalCells,
2238 const label maxLocalCells
2239) const
2240{
2241 label totNCells = mesh_.globalData().nTotalCells();
2242
2243 labelList cellsToRefine;
2244
2245 if (totNCells >= maxGlobalCells)
2246 {
2247 Info<< "No cells marked for refinement since reached limit "
2248 << maxGlobalCells << '.' << endl;
2249 }
2250 else
2251 {
2252 // Every cell I refine adds 7 cells. Estimate the number of cells
2253 // I am allowed to refine.
2254 // Assume perfect distribution so can only refine as much the fraction
2255 // of the mesh I hold. This prediction step prevents us having to do
2256 // lots of reduces to keep count of the total number of cells selected
2257 // for refinement.
2258
2259 //scalar fraction = scalar(mesh_.nCells())/totNCells;
2260 //label myMaxCells = label(maxGlobalCells*fraction);
2261 //label nAllowRefine = (myMaxCells - mesh_.nCells())/7;
2263 //
2264 //Pout<< "refineCandidates:" << nl
2265 // << " total cells:" << totNCells << nl
2266 // << " local cells:" << mesh_.nCells() << nl
2267 // << " local fraction:" << fraction << nl
2268 // << " local allowable cells:" << myMaxCells << nl
2269 // << " local allowable refinement:" << nAllowRefine << nl
2270 // << endl;
2271
2272 //- Disable refinement shortcut. nAllowRefine is per processor limit.
2273 label nAllowRefine = labelMax / Pstream::nProcs();
2274
2275 // Marked for refinement (>= 0) or not (-1). Actual value is the
2276 // index of the surface it intersects / shell it is inside.
2277 labelList refineCell(mesh_.nCells(), -1);
2278 label nRefine = 0;
2279
2280
2281 // Swap neighbouring cell centres and cell level
2282 labelList neiLevel(mesh_.nBoundaryFaces());
2283 pointField neiCc(mesh_.nBoundaryFaces());
2284 calcNeighbourData(neiLevel, neiCc);
2285
2286
2287 const scalar planarCos = Foam::cos(degToRad(planarAngle));
2288
2289
2290 // Cells pierced by feature lines
2291 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2292
2293 if (featureRefinement)
2294 {
2295 label nFeatures = markFeatureRefinement
2296 (
2297 keepPoints,
2298 nAllowRefine,
2299
2300 refineCell,
2301 nRefine
2302 );
2303
2304 Info<< "Marked for refinement due to explicit features "
2305 << ": " << nFeatures << " cells." << endl;
2306 }
2307
2308 // Inside distance-to-feature shells
2309 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2310
2311 if (featureDistanceRefinement)
2312 {
2313 label nShell = markInternalDistanceToFeatureRefinement
2314 (
2315 nAllowRefine,
2316
2317 refineCell,
2318 nRefine
2319 );
2320 Info<< "Marked for refinement due to distance to explicit features "
2321 ": " << nShell << " cells." << endl;
2322 }
2323
2324 // Inside refinement shells
2325 // ~~~~~~~~~~~~~~~~~~~~~~~~
2326
2327 if (internalRefinement)
2328 {
2329 label nShell = markInternalRefinement
2330 (
2331 nAllowRefine,
2332
2333 refineCell,
2334 nRefine
2335 );
2336 Info<< "Marked for refinement due to refinement shells "
2337 << ": " << nShell << " cells." << endl;
2338 }
2339
2340 // Refinement based on intersection of surface
2341 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2342
2343 if (surfaceRefinement)
2344 {
2345 label nSurf = markSurfaceRefinement
2346 (
2347 nAllowRefine,
2348 neiLevel,
2349 neiCc,
2350
2351 refineCell,
2352 nRefine
2353 );
2354 Info<< "Marked for refinement due to surface intersection "
2355 << ": " << nSurf << " cells." << endl;
2356
2357 // Refine intersected-cells only inside gaps. See
2358 // markInternalGapRefinement to refine all cells inside gaps.
2359 if
2360 (
2361 planarCos >= -1
2362 && planarCos <= 1
2363 && max(shells_.maxGapLevel()) > 0
2364 )
2365 {
2366 label nGapSurf = markSurfaceGapRefinement
2367 (
2368 planarCos,
2369 nAllowRefine,
2370 neiLevel,
2371 neiCc,
2372
2373 refineCell,
2374 nRefine
2375 );
2376 Info<< "Marked for refinement due to surface intersection"
2377 << " (at gaps)"
2378 << ": " << nGapSurf << " cells." << endl;
2379 }
2380 }
2381
2382 // Refinement based on curvature of surface
2383 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2384
2385 if
2386 (
2387 curvatureRefinement
2388 && (curvature >= -1 && curvature <= 1)
2389 && (surfaces_.minLevel() != surfaces_.maxLevel())
2390 )
2391 {
2392 label nCurv = markSurfaceCurvatureRefinement
2393 (
2394 curvature,
2395 nAllowRefine,
2396 neiLevel,
2397 neiCc,
2398
2399 refineCell,
2400 nRefine
2401 );
2402 Info<< "Marked for refinement due to curvature/regions "
2403 << ": " << nCurv << " cells." << endl;
2404 }
2405
2406
2407 // Refinement based on features smaller than cell size
2408 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2409
2410 if (smallFeatureRefinement)
2411 {
2412 label nGap = 0;
2413 if
2414 (
2415 planarCos >= -1
2416 && planarCos <= 1
2417 && max(shells_.maxGapLevel()) > 0
2418 )
2419 {
2420 nGap = markSmallFeatureRefinement
2421 (
2422 planarCos,
2423 nAllowRefine,
2424 neiLevel,
2425 neiCc,
2426
2427 refineCell,
2428 nRefine
2429 );
2430 }
2431 Info<< "Marked for refinement due to close opposite surfaces "
2432 << ": " << nGap << " cells." << endl;
2433
2434 label nCurv = 0;
2435 if (max(surfaces_.maxCurvatureLevel()) > 0)
2436 {
2437 nCurv = markSurfaceFieldRefinement
2438 (
2439 nAllowRefine,
2440 neiLevel,
2441 neiCc,
2442
2443 refineCell,
2444 nRefine
2445 );
2446 Info<< "Marked for refinement"
2447 << " due to curvature "
2448 << ": " << nCurv << " cells." << endl;
2449 }
2450 }
2451
2452
2453 // Refinement based on gap (only neighbouring cells)
2454 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2455
2456 const labelList& surfaceGapLevel = surfaces_.gapLevel();
2457
2458 if
2459 (
2460 gapRefinement
2461 && (planarCos >= -1 && planarCos <= 1)
2462 && (max(surfaceGapLevel) > -1)
2463 )
2464 {
2465 Info<< "Specified gap level : " << max(surfaceGapLevel)
2466 << ", planar angle " << planarAngle << endl;
2467
2468 label nGap = markProximityRefinement
2469 (
2470 planarCos,
2471
2472 labelList(surfaceGapLevel.size(), -1), // surface min level
2473 surfaceGapLevel, // surface max level
2474
2475 nAllowRefine,
2476 neiLevel,
2477 neiCc,
2478
2479 refineCell,
2480 nRefine
2481 );
2482 Info<< "Marked for refinement due to close opposite surfaces "
2483 << ": " << nGap << " cells." << endl;
2484 }
2485
2486
2487 // Refinement based on gaps larger than cell size
2488 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2489
2490 if
2491 (
2492 bigGapRefinement
2493 && (planarCos >= -1 && planarCos <= 1)
2494 && max(shells_.maxGapLevel()) > 0
2495 )
2496 {
2497 // Refine based on gap information provided by shell and nearest
2498 // surface
2499 labelList numGapCells;
2500 scalarField gapSize;
2501 label nGap = markInternalGapRefinement
2502 (
2503 planarCos,
2504 spreadGapSize,
2505 nAllowRefine,
2506
2507 refineCell,
2508 nRefine,
2509 numGapCells,
2510 gapSize
2511 );
2512 Info<< "Marked for refinement due to opposite surfaces"
2513 << " "
2514 << ": " << nGap << " cells." << endl;
2515 }
2516
2517
2518 // Limit refinement
2519 // ~~~~~~~~~~~~~~~~
2520
2521 if (limitShells_.shells().size())
2522 {
2523 label nUnmarked = unmarkInternalRefinement(refineCell, nRefine);
2524 if (nUnmarked > 0)
2525 {
2526 Info<< "Unmarked for refinement due to limit shells"
2527 << " : " << nUnmarked << " cells." << endl;
2528 }
2529 }
2530
2531
2532
2533 // Pack cells-to-refine
2534 // ~~~~~~~~~~~~~~~~~~~~
2535
2536 cellsToRefine.setSize(nRefine);
2537 nRefine = 0;
2538
2539 forAll(refineCell, cellI)
2540 {
2541 if (refineCell[cellI] != -1)
2542 {
2543 cellsToRefine[nRefine++] = cellI;
2544 }
2545 }
2546 }
2547
2548 return cellsToRefine;
2549}
2550
2551
2553(
2554 const labelList& cellsToRefine
2555)
2556{
2557 // Mesh changing engine.
2558 polyTopoChange meshMod(mesh_);
2559
2560 // Play refinement commands into mesh changer.
2561 meshCutter_.setRefinement(cellsToRefine, meshMod);
2562
2563 // Remove any unnecessary fields
2564 mesh_.clearOut();
2565 mesh_.moving(false);
2566
2567 // Create mesh (no inflation), return map from old to new mesh.
2568 autoPtr<mapPolyMesh> mapPtr = meshMod.changeMesh(mesh_, false);
2569 mapPolyMesh& map = *mapPtr;
2571 // Update fields
2572 mesh_.updateMesh(map);
2573
2574 // Optionally inflate mesh
2575 if (map.hasMotionPoints())
2576 {
2577 mesh_.movePoints(map.preMotionPoints());
2578 }
2579 else
2580 {
2581 // Delete mesh volumes.
2582 mesh_.clearOut();
2583 }
2584
2585 // Reset the instance for if in overwrite mode
2586 mesh_.setInstance(timeName());
2587
2588 // Update intersection info
2589 updateMesh(map, getChangedFaces(map, cellsToRefine));
2590
2591 return mapPtr;
2592}
2593
2594
2595// Load balancing
2596Foam::autoPtr<Foam::mapDistributePolyMesh> Foam::meshRefinement::balance
2597(
2598 const string& msg,
2599 decompositionMethod& decomposer,
2600 fvMeshDistribute& distributor,
2601 labelList& cellsToRefine,
2602 const scalar maxLoadUnbalance,
2603 const label maxCellUnbalance
2604)
2605{
2607
2608 if (Pstream::nProcs() > 1)
2609 {
2610 // First check if we need to balance at all. Precalculate number of
2611 // cells after refinement and see what maximum difference is.
2612 const scalar nNewCells =
2613 scalar(mesh_.nCells() + 7*cellsToRefine.size());
2614 const scalar nNewCellsAll = returnReduce(nNewCells, sumOp<scalar>());
2615 const scalar nIdealNewCells = nNewCellsAll / Pstream::nProcs();
2616 const scalar unbalance = returnReduce
2617 (
2618 mag(1.0-nNewCells/nIdealNewCells),
2620 );
2621
2622 // Trigger the balancing to avoid too early balancing for better
2623 // scaling performance.
2624 const scalar nNewCellsOnly = scalar(7*cellsToRefine.size());
2625
2626 const label maxNewCells =
2627 label(returnReduce(nNewCellsOnly, maxOp<scalar>()));
2628
2629 const label maxDeltaCells =
2630 label(mag(returnReduce(nNewCells, maxOp<scalar>())-nIdealNewCells));
2631
2632 // New trigger to avoid too early balancing
2633 // 1. Check if globally one proc exceeds the maxCellUnbalance value
2634 // related to the new added cells at the refinement loop
2635 // 2. Check if globally one proc exceeds the maxCellUnbalance based on
2636 // the average cell count a proc should have
2637 if
2638 (
2639 (maxNewCells <= maxCellUnbalance)
2640 && (maxDeltaCells <= maxCellUnbalance)
2641 )
2642 {
2643 Info<< "Skipping balancing since trigger value not reached:" << "\n"
2644 << " Trigger cell count: " << maxCellUnbalance << nl
2645 << " Max new cell count in proc: " << maxNewCells << nl
2646 << " Max difference between new cells and balanced: "
2647 << maxDeltaCells << nl
2648 << " Max load unbalance " << maxLoadUnbalance
2649 << nl <<endl;
2650 }
2651 else
2652 {
2653 if (unbalance <= maxLoadUnbalance)
2654 {
2655 Info<< "Skipping balancing since max unbalance " << unbalance
2656 << " is less than allowable " << maxLoadUnbalance
2657 << endl;
2658 }
2659 else
2660 {
2661 Info<< "Balancing since max unbalance " << unbalance
2662 << " is larger than allowable " << maxLoadUnbalance
2663 << endl;
2664
2665 scalarField cellWeights(mesh_.nCells(), 1);
2666 forAll(cellsToRefine, i)
2667 {
2668 cellWeights[cellsToRefine[i]] += 7;
2669 }
2670
2671 distMap = balance
2672 (
2673 false, //keepZoneFaces
2674 false, //keepBaffles
2675 labelList::null(), //singleProcPoints
2676 cellWeights,
2677 decomposer,
2678 distributor
2679 );
2680
2681 // Update cells to refine
2682 distMap().distributeCellIndices(cellsToRefine);
2683
2684 Info<< "Balanced mesh in = "
2685 << mesh_.time().cpuTimeIncrement() << " s" << endl;
2686
2687 printMeshInfo(debug, "After balancing " + msg, true);
2688
2689
2691 {
2692 Pout<< "Writing balanced " << msg
2693 << " mesh to time " << timeName() << endl;
2694 write
2695 (
2696 debugType(debug),
2697 writeType(writeLevel() | WRITEMESH),
2698 mesh_.time().path()/timeName()
2699 );
2700 Pout<< "Dumped debug data in = "
2701 << mesh_.time().cpuTimeIncrement() << " s" << endl;
2702
2703 // test all is still synced across proc patches
2704 checkData();
2705 }
2706 }
2707 }
2708 }
2709
2710 return distMap;
2711}
2712
2713
2714// Do refinement of consistent set of cells followed by truncation and
2715// load balancing.
2716Foam::autoPtr<Foam::mapDistributePolyMesh>
2718(
2719 const string& msg,
2720 decompositionMethod& decomposer,
2721 fvMeshDistribute& distributor,
2722 const labelList& cellsToRefine,
2723 const scalar maxLoadUnbalance,
2724 const label maxCellUnbalance
2725)
2726{
2727 // Refinement
2728 // ~~~~~~~~~~
2729
2730 refine(cellsToRefine);
2731
2733 {
2734 Pout<< "Writing refined but unbalanced " << msg
2735 << " mesh to time " << timeName() << endl;
2736 write
2737 (
2740 mesh_.time().path()/timeName()
2741 );
2742 Pout<< "Dumped debug data in = "
2743 << mesh_.time().cpuTimeIncrement() << " s" << endl;
2744
2745 // test all is still synced across proc patches
2746 checkData();
2747 }
2748
2749 Info<< "Refined mesh in = "
2750 << mesh_.time().cpuTimeIncrement() << " s" << endl;
2751 printMeshInfo(debug, "After refinement " + msg, true);
2752
2753
2754 // Load balancing
2755 // ~~~~~~~~~~~~~~
2756
2757 labelList noCellsToRefine;
2758
2759 auto distMap = balance
2760 (
2761 msg,
2762 decomposer,
2763 distributor,
2764 noCellsToRefine, // mesh is already refined; no need to predict
2765 maxLoadUnbalance,
2766 maxCellUnbalance
2767 );
2768
2769 return distMap;
2770}
2771
2772
2773// Do load balancing followed by refinement of consistent set of cells.
2776(
2777 const string& msg,
2778 decompositionMethod& decomposer,
2779 fvMeshDistribute& distributor,
2780 const labelList& initCellsToRefine,
2781 const scalar maxLoadUnbalance,
2782 const label maxCellUnbalance
2783)
2784{
2785 labelList cellsToRefine(initCellsToRefine);
2786
2787 //{
2788 // globalIndex globalCells(mesh_.nCells());
2789 //
2790 // Info<< "** Distribution before balancing/refining:" << endl;
2791 // for (const int procI : Pstream::allProcs())
2792 // {
2793 // Info<< " " << procI << '\t'
2794 // << globalCells.localSize(procI) << endl;
2795 // }
2796 // Info<< endl;
2797 //}
2798 //{
2799 // globalIndex globalCells(cellsToRefine.size());
2800 //
2801 // Info<< "** Cells to be refined:" << endl;
2802 // for (const int procI : Pstream::allProcs())
2803 // {
2804 // Info<< " " << procI << '\t'
2805 // << globalCells.localSize(procI) << endl;
2806 // }
2807 // Info<< endl;
2808 //}
2809
2810
2811
2812 // Load balancing
2813 // ~~~~~~~~~~~~~~
2814
2815 auto distMap = balance
2816 (
2817 msg,
2818 decomposer,
2819 distributor,
2820 cellsToRefine,
2821 maxLoadUnbalance,
2822 maxCellUnbalance
2823 );
2824
2825
2826 // Refinement
2827 // ~~~~~~~~~~
2828 // Note: uses updated cellsToRefine
2829
2830 refine(cellsToRefine);
2831
2833 {
2834 Pout<< "Writing refined " << msg
2835 << " mesh to time " << timeName() << endl;
2836 write
2837 (
2840 mesh_.time().path()/timeName()
2841 );
2842 Pout<< "Dumped debug data in = "
2843 << mesh_.time().cpuTimeIncrement() << " s" << endl;
2844
2845 // test all is still synced across proc patches
2846 checkData();
2847 }
2848
2849 Info<< "Refined mesh in = "
2850 << mesh_.time().cpuTimeIncrement() << " s" << endl;
2851
2852 //{
2853 // globalIndex globalCells(mesh_.nCells());
2854 //
2855 // Info<< "** After refinement distribution:" << endl;
2856 // for (const int procI : Pstream::allProcs())
2857 // {
2858 // Info<< " " << procI << '\t'
2859 // << globalCells.localSize(procI) << endl;
2860 // }
2861 // Info<< endl;
2862 //}
2863
2864 printMeshInfo(debug, "After refinement " + msg, true);
2865
2866 return distMap;
2867}
2868
2869
2871(
2872 const label maxGlobalCells,
2873 const label maxLocalCells,
2874 const labelList& currentLevel,
2875 const direction dir
2876) const
2877{
2878 const labelList& cellLevel = meshCutter_.cellLevel();
2879 const pointField& cellCentres = mesh_.cellCentres();
2880
2881 label totNCells = mesh_.globalData().nTotalCells();
2882
2883 labelList cellsToRefine;
2884
2885 if (totNCells >= maxGlobalCells)
2886 {
2887 Info<< "No cells marked for refinement since reached limit "
2888 << maxGlobalCells << '.' << endl;
2889 }
2890 else
2891 {
2892 // Disable refinement shortcut. nAllowRefine is per processor limit.
2893 label nAllowRefine = labelMax / Pstream::nProcs();
2894
2895 // Marked for refinement (>= 0) or not (-1). Actual value is the
2896 // index of the surface it intersects / shell it is inside
2897 labelList refineCell(mesh_.nCells(), -1);
2898 label nRefine = 0;
2899
2900 // Find cells inside the shells with directional levels
2901 labelList insideShell;
2902 shells_.findDirectionalLevel
2903 (
2904 cellCentres,
2905 cellLevel,
2906 currentLevel, // current directional level
2907 dir,
2908 insideShell
2909 );
2910
2911 // Mark for refinement
2912 forAll(insideShell, celli)
2913 {
2914 if (insideShell[celli] >= 0)
2915 {
2916 bool reachedLimit = !markForRefine
2917 (
2918 insideShell[celli], // mark with any positive value
2919 nAllowRefine,
2920 refineCell[celli],
2921 nRefine
2922 );
2923
2924 if (reachedLimit)
2925 {
2926 if (debug)
2927 {
2928 Pout<< "Stopped refining cells"
2929 << " since reaching my cell limit of "
2930 << mesh_.nCells()+nAllowRefine << endl;
2931 }
2932 break;
2933 }
2934 }
2935 }
2936
2937 // Limit refinement
2938 // ~~~~~~~~~~~~~~~~
2939
2940 {
2941 label nUnmarked = unmarkInternalRefinement(refineCell, nRefine);
2942 if (nUnmarked > 0)
2943 {
2944 Info<< "Unmarked for refinement due to limit shells"
2945 << " : " << nUnmarked << " cells." << endl;
2946 }
2947 }
2948
2949
2950
2951 // Pack cells-to-refine
2952 // ~~~~~~~~~~~~~~~~~~~~
2953
2954 cellsToRefine.setSize(nRefine);
2955 nRefine = 0;
2956
2957 forAll(refineCell, cellI)
2958 {
2959 if (refineCell[cellI] != -1)
2960 {
2961 cellsToRefine[nRefine++] = cellI;
2962 }
2963 }
2964 }
2965
2966 return cellsToRefine;
2967}
2968
2969
2970Foam::autoPtr<Foam::mapPolyMesh> Foam::meshRefinement::directionalRefine
2971(
2972 const string& msg,
2973 const direction cmpt,
2974 const labelList& cellsToRefine
2975)
2976{
2977 // Set splitting direction
2978 vector refDir(Zero);
2979 refDir[cmpt] = 1;
2980 List<refineCell> refCells(cellsToRefine.size());
2981 forAll(cellsToRefine, i)
2982 {
2983 refCells[i] = refineCell(cellsToRefine[i], refDir);
2984 }
2985
2986 // How to walk circumference of cells
2987 hexCellLooper cellWalker(mesh_);
2989 // Analyse cuts
2990 cellCuts cuts(mesh_, cellWalker, refCells);
2991
2992 // Cell cutter
2993 Foam::meshCutter meshRefiner(mesh_);
2994
2995 polyTopoChange meshMod(mesh_);
2996
2997 // Insert mesh refinement into polyTopoChange.
2998 meshRefiner.setRefinement(cuts, meshMod);
2999
3000 // Remove any unnecessary fields
3001 mesh_.clearOut();
3002 mesh_.moving(false);
3003
3004 autoPtr<mapPolyMesh> morphMap = meshMod.changeMesh(mesh_, false);
3005
3006 // Update fields
3007 mesh_.updateMesh(*morphMap);
3008
3009 // Move mesh (since morphing does not do this)
3010 if (morphMap().hasMotionPoints())
3011 {
3012 mesh_.movePoints(morphMap().preMotionPoints());
3013 }
3014 else
3015 {
3016 // Delete mesh volumes.
3017 mesh_.clearOut();
3018 }
3019
3020 // Reset the instance for if in overwrite mode
3021 mesh_.setInstance(timeName());
3022
3023 // Update stored refinement pattern
3024 meshRefiner.updateMesh(*morphMap);
3025
3026 // Update intersection info
3027 updateMesh(*morphMap, getChangedFaces(*morphMap, cellsToRefine));
3028
3029 return morphMap;
3030}
3031
3032
3033// ************************************************************************* //
scalar y
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Base cloud calls templated on particle type.
Definition Cloud.H:64
label size() const noexcept
The number of elements in table.
Definition HashTable.H:358
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
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
static const List< label > & null() noexcept
Definition List.H:138
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name).
Definition OSstream.H:134
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
iterator end() noexcept
Return an iterator to end traversing the UList.
Definition UListI.H:454
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run.
Definition UPstream.H:1697
static const Form zero
static constexpr direction nComponents
Number of components in this vector space.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition bitSetI.H:502
bool test(const label pos) const
Test for True value at specified position, never auto-vivify entries.
Definition bitSet.H:334
Description of cuts across cells.
Definition cellCuts.H:109
A cell is defined as a list of faces with extra functionality.
Definition cell.H:56
A cloud is a registry collection of lagrangian particles.
Definition cloud.H:56
Abstract base class for domain decomposition.
Mesh data needed to do the Finite Area discretisation.
Definition edgeFaMesh.H:50
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
A list of face labels.
Definition faceSet.H:50
Sends/receives parts of mesh+fvfields to neighbouring processors. Used in load balancing.
Implementation of cellLooper.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const pointField & preMotionPoints() const noexcept
Pre-motion point positions.
label nOldCells() const noexcept
Number of old cells.
bool hasMotionPoints() const noexcept
Has valid preMotionPoints?
const polyMesh & mesh() const noexcept
Return polyMesh.
const labelList & cellMap() const noexcept
Old cell map.
Cuts (splits) cells.
Definition meshCutter.H:137
void updateMesh(const mapPolyMesh &)
Force recalculation of locally stored data on topological change.
Definition meshCutter.C:988
void setRefinement(const cellCuts &cuts, polyTopoChange &meshMod)
Do actual cutting with cut description. Inserts mesh changes.
Definition meshCutter.C:515
void printMeshInfo(const bool debug, const string &msg, const bool printCellLevel) const
Print some mesh stats.
void updateMesh(const mapPolyMesh &, const labelList &changedFaces)
Update for external change to mesh. changedFaces are in new mesh.
writeType
Enumeration for what to write. Used as a bit-pattern.
labelList refineCandidates(const pointField &keepPoints, const scalar curvature, const scalar planarAngle, const bool featureRefinement, const bool featureDistanceRefinement, const bool internalRefinement, const bool surfaceRefinement, const bool curvatureRefinement, const bool smallFeatureRefinement, const bool gapRefinement, const bool bigGapRefinement, const bool spreadGapSize, const label maxGlobalCells, const label maxLocalCells) const
Calculate list of cells to refine.
word timeName() const
Replacement for Time::timeName() that returns oldInstance (if overwrite_).
debugType
Enumeration for what to debug. Used as a bit-pattern.
bool isNormalGap(const scalar planarCos, const label level0, const vector &point0, const vector &normal0, const label level1, const vector &point1, const vector &normal1) const
Is local topology a small gap normal to the test vector.
autoPtr< mapDistributePolyMesh > refineAndBalance(const string &msg, decompositionMethod &decomposer, fvMeshDistribute &distributor, const labelList &cellsToRefine, const scalar maxLoadUnbalance, const label maxCellUnbalance)
Refine some cells and rebalance.
const fvMesh & mesh() const
Reference to mesh.
bool isGap(const scalar, const vector &, const vector &, const vector &, const vector &) const
Is local topology a small gap?
scalar mergeDistance() const
autoPtr< mapPolyMesh > directionalRefine(const string &msg, const direction cmpt, const labelList &cellsToRefine)
Directionally refine in direction cmpt.
autoPtr< mapPolyMesh > refine(const labelList &cellsToRefine)
Refine some cells.
labelList directionalRefineCandidates(const label maxGlobalCells, const label maxLocalCells, const labelList &currentLevel, const direction dir) const
Calculate list of cells to directionally refine.
autoPtr< mapDistributePolyMesh > balance(const bool keepZoneFaces, const bool keepBaffles, const labelList &singleProcPoints, const scalarField &cellWeights, decompositionMethod &decomposer, fvMeshDistribute &distributor)
Redecompose according to cell count.
autoPtr< mapDistributePolyMesh > balanceAndRefine(const string &msg, decompositionMethod &decomposer, fvMeshDistribute &distributor, const labelList &cellsToRefine, const scalar maxLoadUnbalance, const label maxCellUnbalance)
Balance before refining some cells.
static writeType writeLevel()
Get/set write level.
To compare normals.
normalLess(const vectorList &values)
bool operator()(const label a, const label b) const
void clear()
Clear all entries from the registry.
labelList cmptType
Component type.
vectorList cmptType
Component type.
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
pTraits(const labelList &obj)
Definition pTraits.H:72
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition plane.H:91
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
Direct mesh changes based on v1.3 polyTopoChange syntax.
autoPtr< mapPolyMesh > changeMesh(polyMesh &mesh, const labelUList &patchMap, const bool inflate, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Inplace changes mesh without change of patches.
Container with cells to refine. Refinement given as single direction.
Definition refineCell.H:53
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
Contains information about location on a triSurface.
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.
static bitSet getMasterFaces(const polyMesh &mesh)
Get per face whether it is uncoupled or a master of a coupled set of faces.
Definition syncTools.C:119
static void syncBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop, const TransformOp &top, const bool parRun=UPstream::parRun())
Synchronize values on boundary faces only.
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues, const bool parRun=UPstream::parRun())
Swap coupled boundary face values. Uses eqOp.
Definition syncTools.H:524
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop, const bool parRun=UPstream::parRun())
Synchronize values on all mesh faces.
Definition syncTools.H:465
Class used to pass tracking data to the trackToFace function.
Particle class that marks cells it passes through. Used to mark cells visited by feature edges.
const volScalarField & p0
Definition EEqn.H:36
dynamicFvMesh & mesh
const pointField & points
const cellShapeList & cells
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
const label oldNRefine
word timeName
Definition getTimeIndex.H:3
void set(List< bool > &bools, const labelUList &locations)
Set the listed locations (assign 'true').
Definition BitOps.C:35
const char * end
Definition SVGTools.H:223
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
UIndirectList< label > labelUIndList
UIndirectList of labels.
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< vector > vectorList
List of vector.
Definition vectorList.H:32
label checkData(const objectRegistry &obr, const instantList &timeDirs, wordList &objectNames, const fileName &local=fileName::null)
Check if fields are good to use (available at all times).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
constexpr label labelMax
Definition label.H:55
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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 & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
List< cell > cellList
List of cell.
Definition cellListFwd.H:41
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static bool less(const vector &x, const vector &y)
To compare normals.
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
List< point > pointList
List of point.
Definition pointList.H:32
Field< vector > vectorField
Specialisation of Field<T> for vector.
uint8_t direction
Definition direction.H:49
vector point
Point is a vector.
Definition point.H:37
List< bool > boolList
A List of bools.
Definition List.H:60
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
labelList findIndices(const ListType &input, typename ListType::const_reference val, label start=0)
Linear search to find all occurrences of given element.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensionedScalar cos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
runTime write()
volScalarField & b
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299