Loading...
Searching...
No Matches
meshRefinementProblemCells.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-2015 OpenFOAM Foundation
9 Copyright (C) 2015-2022,2024 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 "fvMesh.H"
31#include "syncTools.H"
32#include "Time.H"
33#include "refinementSurfaces.H"
34#include "pointSet.H"
35#include "faceSet.H"
37#include "cellSet.H"
38#include "searchableSurfaces.H"
39#include "polyMeshGeometry.H"
40#include "IOmanip.H"
41#include "unitConversion.H"
42#include "snappySnapDriver.H"
43
44#include "snapParameters.H"
45#include "motionSmoother.H"
46
47// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48
49void Foam::meshRefinement::markBoundaryFace
50(
51 const label facei,
52 boolList& isBoundaryFace,
53 boolList& isBoundaryEdge,
54 boolList& isBoundaryPoint
55) const
56{
57 isBoundaryFace[facei] = true;
58
59 const labelList& fEdges = mesh_.faceEdges(facei);
60
61 for (const label edgei : fEdges)
62 {
63 isBoundaryEdge[edgei] = true;
64 }
65
66 const face& f = mesh_.faces()[facei];
67
68 for (const label pointi : f)
69 {
70 isBoundaryPoint[pointi] = true;
71 }
72}
73
74
75void Foam::meshRefinement::findNearest
76(
77 const labelList& meshFaces,
78 List<pointIndexHit>& nearestInfo,
79 labelList& nearestSurface,
80 labelList& nearestRegion,
81 vectorField& nearestNormal
82) const
83{
84 pointField fc(meshFaces.size());
85 forAll(meshFaces, i)
86 {
87 fc[i] = mesh_.faceCentres()[meshFaces[i]];
88 }
89
90 const labelList allSurfaces(identity(surfaces_.surfaces().size()));
91
92 surfaces_.findNearest
93 (
94 allSurfaces,
95 fc,
96 scalarField(fc.size(), sqr(GREAT)), // sqr of attraction
97 nearestSurface,
98 nearestInfo
99 );
100
101 // Do normal testing per surface.
102 nearestNormal.setSize(nearestInfo.size());
103 nearestRegion.setSize(nearestInfo.size());
104
105 forAll(allSurfaces, surfI)
106 {
108
109 forAll(nearestSurface, i)
110 {
111 if (nearestSurface[i] == surfI)
112 {
113 localHits.append(nearestInfo[i]);
114 }
115 }
116
117 label geomI = surfaces_.surfaces()[surfI];
118
119 pointField localNormals;
120 surfaces_.geometry()[geomI].getNormal(localHits, localNormals);
121
122 labelList localRegion;
123 surfaces_.geometry()[geomI].getRegion(localHits, localRegion);
124
125 label localI = 0;
126 forAll(nearestSurface, i)
127 {
128 if (nearestSurface[i] == surfI)
129 {
130 nearestNormal[i] = localNormals[localI];
131 nearestRegion[i] = localRegion[localI];
132 localI++;
133 }
134 }
135 }
136}
137
138
139Foam::Map<Foam::label> Foam::meshRefinement::findEdgeConnectedProblemCells
140(
141 const scalarField& perpendicularAngle,
142 const labelList& globalToPatch
143) const
144{
145 // Construct addressing engine from all patches added for meshing.
147 (
149 (
150 mesh_,
151 meshedPatches()
152 )
153 );
154 const indirectPrimitivePatch& pp = ppPtr();
155
156
157 // 1. Collect faces to test
158 // ~~~~~~~~~~~~~~~~~~~~~~~~
159
160 DynamicList<label> candidateFaces(pp.size()/20);
161
162 const labelListList& edgeFaces = pp.edgeFaces();
163
164 const labelList& cellLevel = meshCutter_.cellLevel();
165
166 forAll(edgeFaces, edgeI)
167 {
168 const labelList& eFaces = edgeFaces[edgeI];
169
170 if (eFaces.size() == 2)
171 {
172 label face0 = pp.addressing()[eFaces[0]];
173 label face1 = pp.addressing()[eFaces[1]];
174
175 label cell0 = mesh_.faceOwner()[face0];
176 label cell1 = mesh_.faceOwner()[face1];
177
178 if (cellLevel[cell0] > cellLevel[cell1])
179 {
180 // cell0 smaller.
181 const vector& n0 = pp.faceNormals()[eFaces[0]];
182 const vector& n1 = pp.faceNormals()[eFaces[1]];
183
184 if (mag(n0 & n1) < 0.1)
185 {
186 candidateFaces.append(face0);
187 }
188 }
189 else if (cellLevel[cell1] > cellLevel[cell0])
190 {
191 // cell1 smaller.
192 const vector& n0 = pp.faceNormals()[eFaces[0]];
193 const vector& n1 = pp.faceNormals()[eFaces[1]];
194
195 if (mag(n0 & n1) < 0.1)
196 {
197 candidateFaces.append(face1);
198 }
199 }
200 }
201 }
202 candidateFaces.shrink();
203
204 Info<< "Testing " << returnReduce(candidateFaces.size(), sumOp<label>())
205 << " faces on edge-connected cells of differing level."
206 << endl;
207
209 {
210 faceSet fSet(mesh_, "edgeConnectedFaces", candidateFaces);
211 fSet.instance() = timeName();
212 Pout<< "Writing " << fSet.size()
213 << " with problematic topology to faceSet "
214 << fSet.objectPath() << endl;
215 fSet.write();
216 }
217
218
219 // 2. Find nearest surface on candidate faces
220 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
221
222 List<pointIndexHit> nearestInfo;
223 labelList nearestSurface;
224 labelList nearestRegion;
225 vectorField nearestNormal;
226 findNearest
227 (
228 candidateFaces,
229 nearestInfo,
230 nearestSurface,
231 nearestRegion,
232 nearestNormal
233 );
234
235
236 // 3. Test angle to surface
237 // ~~~~~~~~~~~~~~~~~~~~~~~~
238
239 Map<label> candidateCells(candidateFaces.size());
240
241 faceSet perpFaces(mesh_, "perpendicularFaces", pp.size()/100);
242
243 forAll(candidateFaces, i)
244 {
245 label facei = candidateFaces[i];
246
247 const vector n = normalised(mesh_.faceAreas()[facei]);
248
249 label region = surfaces_.globalRegion
250 (
251 nearestSurface[i],
252 nearestRegion[i]
253 );
254
255 scalar angle = degToRad(perpendicularAngle[region]);
256
257 if (angle >= 0)
258 {
259 if (mag(n & nearestNormal[i]) < Foam::sin(angle))
260 {
261 perpFaces.insert(facei);
262 candidateCells.insert
263 (
264 mesh_.faceOwner()[facei],
265 globalToPatch[region]
266 );
267 }
268 }
269 }
270
272 {
273 perpFaces.instance() = timeName();
274 Pout<< "Writing " << perpFaces.size()
275 << " faces that are perpendicular to the surface to set "
276 << perpFaces.objectPath() << endl;
277 perpFaces.write();
278 }
279 return candidateCells;
280}
281
282
283// Check if moving face to new points causes a 'collapsed' face.
284// Uses new point position only for the face, not the neighbouring
285// cell centres
286bool Foam::meshRefinement::isCollapsedFace
287(
288 const pointField& points,
289 const pointField& neiCc,
290 const scalar minFaceArea,
291 const scalar maxNonOrtho,
292 const label facei
293) const
294{
295 // Severe nonorthogonality threshold
296 const scalar severeNonorthogonalityThreshold =
297 ::cos(degToRad(maxNonOrtho));
298
299 vector s = mesh_.faces()[facei].areaNormal(points);
300 scalar magS = mag(s);
301
302 // Check face area
303 if (magS < minFaceArea)
304 {
305 return true;
306 }
307
308 // Check orthogonality
309 const point& ownCc = mesh_.cellCentres()[mesh_.faceOwner()[facei]];
310
311 if (mesh_.isInternalFace(facei))
312 {
313 label nei = mesh_.faceNeighbour()[facei];
314 vector d = mesh_.cellCentres()[nei] - ownCc;
315
316 scalar dDotS = (d & s)/(mag(d)*magS + VSMALL);
317
318 if (dDotS < severeNonorthogonalityThreshold)
319 {
320 return true;
321 }
322 else
323 {
324 return false;
325 }
326 }
327 else
328 {
329 label patchi = mesh_.boundaryMesh().whichPatch(facei);
330
331 if (mesh_.boundaryMesh()[patchi].coupled())
332 {
333 vector d = neiCc[facei-mesh_.nInternalFaces()] - ownCc;
334
335 scalar dDotS = (d & s)/(mag(d)*magS + VSMALL);
336
337 if (dDotS < severeNonorthogonalityThreshold)
338 {
339 return true;
340 }
341 else
342 {
343 return false;
344 }
345 }
346 else
347 {
348 // Collapsing normal boundary face does not cause problems with
349 // non-orthogonality
350 return false;
351 }
352 }
353}
354
355
356// Check if moving cell to new points causes it to collapse.
357bool Foam::meshRefinement::isCollapsedCell
358(
359 const pointField& points,
360 const scalar volFraction,
361 const label celli
362) const
363{
364 scalar vol = mesh_.cells()[celli].mag(points, mesh_.faces());
365
366 if (vol/mesh_.cellVolumes()[celli] < volFraction)
367 {
368 return true;
369 }
370 else
371 {
372 return false;
373 }
374}
375
376
377// Returns list with for every internal face -1 or the patch they should
378// be baffled into. Gets run after createBaffles so all the unzoned surface
379// intersections have already been turned into baffles. (Note: zoned surfaces
380// are not baffled at this stage)
381// Used to remove cells by baffling all their faces and have the
382// splitMeshRegions chuck away non used regions.
383void Foam::meshRefinement::markFacesOnProblemCells
384(
385 const dictionary& motionDict,
386 const bool removeEdgeConnectedCells,
387 const scalarField& perpendicularAngle,
388 const labelList& globalToPatch,
389
390 labelList& facePatch,
392) const
393{
394 const labelList& cellLevel = meshCutter_.cellLevel();
395 const labelList& pointLevel = meshCutter_.pointLevel();
396 const polyBoundaryMesh& patches = mesh_.boundaryMesh();
397
398
399 // Mark all points and edges on baffle patches (so not on any inlets,
400 // outlets etc.)
401 boolList isBoundaryPoint(mesh_.nPoints(), false);
402 boolList isBoundaryEdge(mesh_.nEdges(), false);
403 boolList isBoundaryFace(mesh_.nFaces(), false);
404
405 // Fill boundary data. All elements on meshed patches get marked.
406 // Get the labels of added patches.
407 const labelList adaptPatchIDs(meshedPatches());
408
409 forAll(adaptPatchIDs, i)
410 {
411 const polyPatch& pp = patches[adaptPatchIDs[i]];
412
413 label facei = pp.start();
414
415 forAll(pp, j)
416 {
417 markBoundaryFace
418 (
419 facei,
420 isBoundaryFace,
421 isBoundaryEdge,
422 isBoundaryPoint
423 );
424
425 facei++;
426 }
427 }
428
429
430 // Per mesh face the nearest adaptPatch and its faceZone (if any)
431 labelList nearestAdaptPatch;
432 labelList nearestAdaptZone;
433 nearestPatch(adaptPatchIDs, nearestAdaptPatch, nearestAdaptZone);
434
435
436 // Per internal face (boundary faces not used) the patch that the
437 // baffle should get (or -1)
438 facePatch.setSize(mesh_.nFaces());
439 facePatch = -1;
440 faceZone.setSize(mesh_.nFaces());
441 faceZone = -1;
442
443 // Swap neighbouring cell centres and cell level
444 labelList neiLevel(mesh_.nBoundaryFaces());
445 pointField neiCc(mesh_.nBoundaryFaces());
446 calcNeighbourData(neiLevel, neiCc);
447
448
449 // Count of faces marked for baffling
450 label nBaffleFaces = 0;
451 bitSet isMasterFace(syncTools::getMasterFaces(mesh_));
452
453 // Count of faces not baffled since would not cause a collapse
454 label nPrevented = 0;
455
456 if (removeEdgeConnectedCells && max(perpendicularAngle) >= 0)
457 {
458 Info<< "markFacesOnProblemCells :"
459 << " Checking for edge-connected cells of highly differing sizes."
460 << endl;
461
462 // Pick up the cells that need to be removed and (a guess for)
463 // the patch they should be patched with.
464 Map<label> problemCells
465 (
466 findEdgeConnectedProblemCells
467 (
468 perpendicularAngle,
469 globalToPatch
470 )
471 );
472
473 // Baffle all faces of cells that need to be removed
474 forAllConstIters(problemCells, iter)
475 {
476 for (const label facei : mesh_.cells()[iter.key()])
477 {
478 const label patchi = patches.whichPatch(facei);
479
480 if
481 (
482 facePatch[facei] == -1
483 // && mesh_.isInternalFace(facei)
484 && (patchi == -1 || patches[patchi].coupled())
485 )
486 {
487 facePatch[facei] = nearestAdaptPatch[facei];
488 faceZone[facei] = nearestAdaptZone[facei];
489 nBaffleFaces++;
490
491 // Mark face as a 'boundary'
492 markBoundaryFace
493 (
494 facei,
495 isBoundaryFace,
496 isBoundaryEdge,
497 isBoundaryPoint
498 );
499 }
500 }
501 }
502 Info<< "markFacesOnProblemCells : Marked "
503 << returnReduce(nBaffleFaces, sumOp<label>())
504 << " additional internal faces to be converted into baffles"
505 << " due to "
506 << returnReduce(problemCells.size(), sumOp<label>())
507 << " cells edge-connected to lower level cells." << endl;
508
510 {
511 cellSet problemCellSet(mesh_, "problemCells", problemCells.toc());
512 problemCellSet.instance() = timeName();
513 Pout<< "Writing " << problemCellSet.size()
514 << " cells that are edge connected to coarser cell to set "
515 << problemCellSet.objectPath() << endl;
516 problemCellSet.write();
517 }
518 }
519
521 (
522 mesh_,
523 isBoundaryPoint,
524 orEqOp<bool>(),
525 false // null value
526 );
527
529 (
530 mesh_,
531 isBoundaryEdge,
532 orEqOp<bool>(),
533 false // null value
534 );
535
537 (
538 mesh_,
539 isBoundaryFace,
541 );
542
543
544 // See if checking for collapse
545 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
546
547 // Collapse checking parameters
548 const scalar volFraction =
549 motionDict.getOrDefault<scalar>("minVolCollapseRatio", -1);
550
551 const bool checkCollapse = (volFraction > 0);
552 scalar minArea = -1;
553 scalar maxNonOrtho = -1;
554
555
556 // Find nearest (non-baffle) surface
557 pointField newPoints;
558
559 if (checkCollapse)
560 {
561 minArea = get<scalar>(motionDict, "minArea", dryRun_);
562 maxNonOrtho = get<scalar>(motionDict, "maxNonOrtho", dryRun_);
563
564 Info<< "markFacesOnProblemCells :"
565 << " Deleting all-anchor surface cells only if"
566 << " snapping them violates mesh quality constraints:" << nl
567 << " snapped/original cell volume < " << volFraction << nl
568 << " face area < " << minArea << nl
569 << " non-orthogonality > " << maxNonOrtho << nl
570 << endl;
571
572 // Construct addressing engine.
574 (
576 (
577 mesh_,
578 adaptPatchIDs
579 )
580 );
581 const indirectPrimitivePatch& pp = ppPtr();
582 const pointField& localPoints = pp.localPoints();
583 const labelList& meshPoints = pp.meshPoints();
584
585 List<pointIndexHit> hitInfo;
586 labelList hitSurface;
587 surfaces_.findNearest
588 (
589 surfaceZonesInfo::getUnnamedSurfaces(surfaces_.surfZones()),
590 localPoints,
591 scalarField(localPoints.size(), sqr(GREAT)), // sqr of attraction
592 hitSurface,
593 hitInfo
594 );
595
596 // Start off from current points
597 newPoints = mesh_.points();
598
599 forAll(hitInfo, i)
600 {
601 if (hitInfo[i].hit())
602 {
603 newPoints[meshPoints[i]] = hitInfo[i].point();
604 }
605 }
606
608 {
609 const_cast<Time&>(mesh_.time())++;
610 pointField oldPoints(mesh_.points());
611 mesh_.movePoints(newPoints);
612
613 // Unset any moving state
614 mesh_.moving(false);
615
616 Pout<< "Writing newPoints mesh to time " << timeName()
617 << endl;
618 write
619 (
620 debugType(debug),
621 writeType(writeLevel() | WRITEMESH),
622 mesh_.time().path()/"newPoints"
623 );
624 mesh_.movePoints(oldPoints);
625
626 // Unset any moving state
627 mesh_.moving(false);
628 }
629 }
630
631
632
633 // For each cell count the number of anchor points that are on
634 // the boundary:
635 // 8 : check the number of (baffle) boundary faces. If 3 or more block
636 // off the cell since the cell would get squeezed down to a diamond
637 // (probably; if the 3 or more faces are unrefined (only use the
638 // anchor points))
639 // 7 : store. Used to check later on whether there are points with
640 // 3 or more of these cells. (note that on a flat surface a boundary
641 // point will only have 4 cells connected to it)
642
643
644 // Does cell have exactly 7 of its 8 anchor points on the boundary?
645 bitSet hasSevenBoundaryAnchorPoints(mesh_.nCells());
646 // If so what is the remaining non-boundary anchor point?
647 labelHashSet nonBoundaryAnchors(mesh_.nCells()/10000);
648
649 // On-the-fly addressing storage.
650 DynamicList<label> dynFEdges;
651 DynamicList<label> dynCPoints;
652 labelHashSet pSet;
653
654 forAll(cellLevel, celli)
655 {
656 const labelList& cPoints = mesh_.cellPoints(celli, pSet, dynCPoints);
657
658 // Get number of anchor points (pointLevel <= cellLevel)
659
660 label nBoundaryAnchors = 0;
661 // label nNonAnchorBoundary = 0;
662 label nonBoundaryAnchor = -1;
663
664 for (const label pointi : cPoints)
665 {
666 if (pointLevel[pointi] <= cellLevel[celli])
667 {
668 // Anchor point
669 if (isBoundaryPoint[pointi])
670 {
671 nBoundaryAnchors++;
672 }
673 else
674 {
675 // Anchor point which is not on the surface
676 nonBoundaryAnchor = pointi;
677 }
678 }
679 // else if (isBoundaryPoint[pointi])
680 // {
681 // ++nNonAnchorBoundary;
682 // }
683 }
684
685 if (nBoundaryAnchors == 8)
686 {
687 const cell& cFaces = mesh_.cells()[celli];
688
689 // Count boundary faces.
690 // label nBfaces = 0;
691 //
692 // forAll(cFaces, cFacei)
693 // {
694 // if (isBoundaryFace[cFaces[cFacei]])
695 // {
696 // ++nBfaces;
697 // }
698 // }
699
700 // If nBfaces > 1 make all non-boundary non-baffle faces baffles.
701 // We assume that this situation is where there is a single
702 // cell sticking out which would get flattened.
703
704 // Eugene: delete cell no matter what.
705 //if (nBfaces > 1)
706 {
707 if
708 (
709 checkCollapse
710 && !isCollapsedCell(newPoints, volFraction, celli)
711 )
712 {
713 nPrevented++;
714 //Pout<< "Preventing baffling/removal of 8 anchor point"
715 // << " cell "
716 // << celli << " at " << mesh_.cellCentres()[celli]
717 // << " since new volume "
718 // << mesh_.cells()[celli].mag(newPoints, mesh_.faces())
719 // << " old volume " << mesh_.cellVolumes()[celli]
720 // << endl;
721 }
722 else
723 {
724 // Block all faces of cell
725 for (const label facei : cFaces)
726 {
727 const label patchi = patches.whichPatch(facei);
728
729 if
730 (
731 facePatch[facei] == -1
732 // && mesh_.isInternalFace(facei)
733 && (patchi == -1 || patches[patchi].coupled())
734 )
735 {
736 facePatch[facei] = nearestAdaptPatch[facei];
737 faceZone[facei] = nearestAdaptZone[facei];
738 nBaffleFaces++;
739
740 // Mark face as a 'boundary'
741 markBoundaryFace
742 (
743 facei,
744 isBoundaryFace,
745 isBoundaryEdge,
746 isBoundaryPoint
747 );
748 }
749 }
750 }
751 }
752 }
753 else if (nBoundaryAnchors == 7 && nonBoundaryAnchor != -1)
754 {
755 // Mark the cell. Store the (single!) non-boundary anchor point.
756 hasSevenBoundaryAnchorPoints.set(celli);
757 nonBoundaryAnchors.insert(nonBoundaryAnchor);
758 }
759 }
760
761
762 // Loop over all points. If a point is connected to 4 or more cells
763 // with 7 anchor points on the boundary set those cell's non-boundary faces
764 // to baffles
765
766 DynamicList<label> dynPCells;
767
768 for (const label pointi : nonBoundaryAnchors)
769 {
770 const labelList& pCells = mesh_.pointCells(pointi, dynPCells);
771
772 // Count number of 'hasSevenBoundaryAnchorPoints' cells.
773 label n = 0;
774
775 for (const label celli : pCells)
776 {
777 if (hasSevenBoundaryAnchorPoints.test(celli))
778 {
779 ++n;
780 }
781 }
782
783 if (n > 3)
784 {
785 // Point in danger of being what? Remove all 7-cells.
786 for (const label celli : pCells)
787 {
788 if (hasSevenBoundaryAnchorPoints.test(celli))
789 {
790 if
791 (
792 checkCollapse
793 && !isCollapsedCell(newPoints, volFraction, celli)
794 )
795 {
796 nPrevented++;
797 //Pout<< "Preventing baffling of 7 anchor cell "
798 // << celli
799 // << " at " << mesh_.cellCentres()[celli]
800 // << " since new volume "
801 // << mesh_.cells()[celli].mag
802 // (newPoints, mesh_.faces())
803 // << " old volume " << mesh_.cellVolumes()[celli]
804 // << endl;
805 }
806 else
807 {
808 const cell& cFaces = mesh_.cells()[celli];
809
810 for (const label facei : cFaces)
811 {
812 const label patchi = patches.whichPatch(facei);
813
814 if
815 (
816 facePatch[facei] == -1
817 // && mesh_.isInternalFace(facei)
818 && (patchi == -1 || patches[patchi].coupled())
819 )
820 {
821 facePatch[facei] = nearestAdaptPatch[facei];
822 faceZone[facei] = nearestAdaptZone[facei];
823 nBaffleFaces++;
824
825 // Mark face as a 'boundary'
826 markBoundaryFace
827 (
828 facei,
829 isBoundaryFace,
830 isBoundaryEdge,
831 isBoundaryPoint
832 );
833 }
834 }
835 }
836 }
837 }
838 }
839 }
840
841
842 // Sync all. (note that pointdata and facedata not used anymore but sync
843 // anyway)
844
846 (
847 mesh_,
848 isBoundaryPoint,
849 orEqOp<bool>(),
850 false // null value
851 );
852
854 (
855 mesh_,
856 isBoundaryEdge,
857 orEqOp<bool>(),
858 false // null value
859 );
860
862 (
863 mesh_,
864 isBoundaryFace,
866 );
867
868
869 // Find faces with all edges on the boundary and make them baffles
870 for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
871 {
872 if (facePatch[facei] == -1)
873 {
874 const labelList& fEdges = mesh_.faceEdges(facei, dynFEdges);
875 label nFaceBoundaryEdges = 0;
876
877 for (const label edgei : fEdges)
878 {
879 if (isBoundaryEdge[edgei])
880 {
881 ++nFaceBoundaryEdges;
882 }
883 }
884
885 if (nFaceBoundaryEdges == fEdges.size())
886 {
887 if
888 (
889 checkCollapse
890 && !isCollapsedFace
891 (
892 newPoints,
893 neiCc,
894 minArea,
895 maxNonOrtho,
896 facei
897 )
898 )
899 {
900 nPrevented++;
901 //Pout<< "Preventing baffling (to avoid collapse) of face "
902 // << facei
903 // << " with all boundary edges "
904 // << " at " << mesh_.faceCentres()[facei]
905 // << endl;
906 }
907 else
908 {
909 facePatch[facei] = nearestAdaptPatch[facei];
910 faceZone[facei] = nearestAdaptZone[facei];
911 nBaffleFaces++;
912
913 // Do NOT update boundary data since this would grow blocked
914 // faces across gaps.
915 }
916 }
917 }
918 }
919
920 for (const polyPatch& pp : patches)
921 {
922 if (pp.coupled())
923 {
924 label facei = pp.start();
925
926 forAll(pp, i)
927 {
928 if (facePatch[facei] == -1)
929 {
930 const labelList& fEdges = mesh_.faceEdges(facei, dynFEdges);
931 label nFaceBoundaryEdges = 0;
932
933 for (const label edgei : fEdges)
934 {
935 if (isBoundaryEdge[edgei])
936 {
937 ++nFaceBoundaryEdges;
938 }
939 }
940
941 if (nFaceBoundaryEdges == fEdges.size())
942 {
943 if
944 (
945 checkCollapse
946 && !isCollapsedFace
947 (
948 newPoints,
949 neiCc,
950 minArea,
951 maxNonOrtho,
952 facei
953 )
954 )
955 {
956 nPrevented++;
957 //Pout<< "Preventing baffling of coupled face "
958 // << facei
959 // << " with all boundary edges "
960 // << " at " << mesh_.faceCentres()[facei]
961 // << endl;
962 }
963 else
964 {
965 facePatch[facei] = nearestAdaptPatch[facei];
966 faceZone[facei] = nearestAdaptZone[facei];
967 if (isMasterFace.test(facei))
968 {
969 ++nBaffleFaces;
970 }
971
972 // Do NOT update boundary data since this would grow
973 // blocked faces across gaps.
974 }
975 }
976 }
977
978 facei++;
979 }
980 }
981 }
982
983
984 // Because of isCollapsedFace one side can decide not to baffle whereas
985 // the other side does so sync. Baffling is preferred over not baffling.
986 if (checkCollapse) // Or always?
987 {
989 (
990 mesh_,
991 facePatch,
993 );
995 (
996 mesh_,
997 faceZone,
999 );
1000 }
1001
1002 Info<< "markFacesOnProblemCells : marked "
1003 << returnReduce(nBaffleFaces, sumOp<label>())
1004 << " additional internal faces to be converted into baffles."
1005 << endl;
1006
1007 if (checkCollapse)
1008 {
1009 Info<< "markFacesOnProblemCells : prevented "
1010 << returnReduce(nPrevented, sumOp<label>())
1011 << " internal faces from getting converted into baffles."
1012 << endl;
1013 }
1014}
1015
1016
1017// Mark faces to be baffled to prevent snapping problems. Does
1018// test to find nearest surface and checks which faces would get squashed.
1019void Foam::meshRefinement::markFacesOnProblemCellsGeometric
1020(
1021 const snapParameters& snapParams,
1022 const dictionary& motionDict,
1023 const labelList& globalToMasterPatch,
1024 const labelList& globalToSlavePatch,
1025
1026 labelList& facePatch,
1028) const
1029{
1030 pointField oldPoints(mesh_.points());
1031
1032 // Repeat (most of) snappySnapDriver::doSnap
1033 {
1034 const labelList adaptPatchIDs(meshedPatches());
1035
1036 // Construct addressing engine.
1038 (
1040 (
1041 mesh_,
1042 adaptPatchIDs
1043 )
1044 );
1045 indirectPrimitivePatch& pp = ppPtr();
1046
1047 // Distance to attract to nearest feature on surface
1048 const scalarField snapDist
1049 (
1050 snappySnapDriver::calcSnapDistance(mesh_, snapParams, pp)
1051 );
1052
1053
1054 // Construct iterative mesh mover.
1055 Info<< "Constructing mesh displacer ..." << endl;
1056 Info<< "Using mesh parameters " << motionDict << nl << endl;
1057
1058 const pointMesh& pMesh = pointMesh::New(mesh_);
1059
1060 motionSmoother meshMover
1061 (
1062 mesh_,
1063 pp,
1064 adaptPatchIDs,
1065 meshRefinement::makeDisplacementField(pMesh, adaptPatchIDs)(),
1066 motionDict
1067 );
1068
1069
1070 // Check initial mesh
1071 Info<< "Checking initial mesh ..." << endl;
1072 labelHashSet wrongFaces(mesh_.nFaces()/100);
1074 (
1075 false,
1076 mesh_,
1077 motionDict,
1078 wrongFaces,
1079 dryRun_
1080 );
1081 const label nInitErrors = returnReduce
1082 (
1083 wrongFaces.size(),
1084 sumOp<label>()
1085 );
1086
1087 Info<< "Detected " << nInitErrors << " illegal faces"
1088 << " (concave, zero area or negative cell pyramid volume)"
1089 << endl;
1090
1091
1092 Info<< "Checked initial mesh in = "
1093 << mesh_.time().cpuTimeIncrement() << " s\n" << nl << endl;
1094
1095 // Pre-smooth patch vertices (so before determining nearest)
1097 (
1098 *this,
1099 snapParams,
1100 nInitErrors,
1101 List<labelPair>(0), //baffles
1102 meshMover
1103 );
1104
1105 pointField nearestPoint;
1106 vectorField nearestNormal;
1107 const vectorField disp
1108 (
1109 snappySnapDriver::calcNearestSurface
1110 (
1111 snapParams.strictRegionSnap(),
1112 *this,
1113 globalToMasterPatch,
1114 globalToSlavePatch,
1115 pp,
1116 pp.localPoints(),
1117 snapDist, // attraction
1118 nearestPoint,
1119 nearestNormal
1120 )
1121 );
1122
1123 const labelList& meshPoints = pp.meshPoints();
1124
1125 pointField newPoints(mesh_.points());
1126 forAll(meshPoints, i)
1127 {
1128 newPoints[meshPoints[i]] += disp[i];
1129 }
1130
1132 (
1133 mesh_,
1134 newPoints,
1135 minMagSqrEqOp<point>(), // combine op
1136 vector(GREAT, GREAT, GREAT) // null value (note: cannot use VGREAT)
1137 );
1138
1139 mesh_.movePoints(newPoints);
1140
1141 // Unset any moving state
1142 mesh_.moving(false);
1143 }
1144
1145
1146 // Per face the nearest adaptPatch
1147 //const labelList nearestAdaptPatch(nearestPatch(meshedPatches()));
1148 labelList nearestAdaptPatch;
1149 labelList nearestAdaptZone;
1150 nearestPatch(meshedPatches(), nearestAdaptPatch, nearestAdaptZone);
1151
1152 // Per face (internal or coupled!) the patch that the
1153 // baffle should get (or -1).
1154 facePatch.setSize(mesh_.nFaces());
1155 facePatch = -1;
1156 faceZone.setSize(mesh_.nFaces());
1157 faceZone = -1;
1158 // Count of baffled faces
1159 label nBaffleFaces = 0;
1160
1161 {
1162 faceSet wrongFaces(mesh_, "wrongFaces", 100);
1163 {
1164 //motionSmoother::checkMesh(false, mesh_, motionDict, wrongFaces);
1165
1166 // Just check the errors from squashing
1167 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1168
1169 const labelList allFaces(identity(mesh_.nFaces()));
1170 label nWrongFaces = 0;
1171
1172 //const scalar minV
1173 //(motionDict.get<scalar>("minVol", keyType::REGEX_RECURSIVE));
1174 //if (minV > -GREAT)
1175 //{
1176 // polyMeshGeometry::checkFacePyramids
1177 // (
1178 // false,
1179 // minV,
1180 // mesh_,
1181 // mesh_.cellCentres(),
1182 // mesh_.points(),
1183 // allFaces,
1184 // List<labelPair>(0),
1185 // &wrongFaces
1186 // );
1187 //
1188 // label nNewWrongFaces = returnReduce
1189 // (
1190 // wrongFaces.size(),
1191 // sumOp<label>()
1192 // );
1193 //
1194 // Info<< " faces with pyramid volume < "
1195 // << setw(5) << minV
1196 // << " m^3 : "
1197 // << nNewWrongFaces-nWrongFaces << endl;
1198 //
1199 // nWrongFaces = nNewWrongFaces;
1200 //}
1201
1202 scalar minArea(get<scalar>(motionDict, "minArea", dryRun_));
1203 if (minArea > -SMALL)
1204 {
1206 (
1207 false,
1208 minArea,
1209 mesh_,
1210 mesh_.faceAreas(),
1211 allFaces,
1212 &wrongFaces
1213 );
1214
1215 label nNewWrongFaces = returnReduce
1216 (
1217 wrongFaces.size(),
1218 sumOp<label>()
1219 );
1220
1221 Info<< " faces with area < "
1222 << setw(5) << minArea
1223 << " m^2 : "
1224 << nNewWrongFaces-nWrongFaces << endl;
1225
1226 nWrongFaces = nNewWrongFaces;
1227 }
1228
1229 scalar minDet(get<scalar>(motionDict, "minDeterminant", dryRun_));
1230 if (minDet > -1)
1231 {
1233 (
1234 false,
1235 minDet,
1236 mesh_,
1237 mesh_.faceAreas(),
1238 allFaces,
1239 polyMeshGeometry::affectedCells(mesh_, allFaces),
1240 &wrongFaces
1241 );
1242
1243 label nNewWrongFaces = returnReduce
1244 (
1245 wrongFaces.size(),
1246 sumOp<label>()
1247 );
1248
1249 Info<< " faces on cells with determinant < "
1250 << setw(5) << minDet << " : "
1251 << nNewWrongFaces-nWrongFaces << endl;
1252
1253 nWrongFaces = nNewWrongFaces;
1254 }
1255 }
1256
1257
1258 for (const label facei : wrongFaces)
1259 {
1260 const label patchi = mesh_.boundaryMesh().whichPatch(facei);
1261
1262 if (patchi == -1 || mesh_.boundaryMesh()[patchi].coupled())
1263 {
1264 facePatch[facei] = nearestAdaptPatch[facei];
1265 faceZone[facei] = nearestAdaptZone[facei];
1266 nBaffleFaces++;
1267
1268 //Pout<< " " << facei
1269 // //<< " on patch " << mesh_.boundaryMesh()[patchi].name()
1270 // << " is destined for patch " << facePatch[facei]
1271 // << endl;
1272 }
1273 }
1274 }
1275
1276
1277 // Restore points.
1278 mesh_.movePoints(oldPoints);
1279
1280 // Unset any moving state
1281 mesh_.moving(false);
1282
1283
1284 Info<< "markFacesOnProblemCellsGeometric : marked "
1285 << returnReduce(nBaffleFaces, sumOp<label>())
1286 << " additional internal and coupled faces"
1287 << " to be converted into baffles." << endl;
1288
1290 (
1291 mesh_,
1292 facePatch,
1294 );
1296 (
1297 mesh_,
1298 faceZone,
1300 );
1301}
1302
1303
1304// ************************************************************************* //
Istream and Ostream manipulators taking arguments.
label n
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
const Addr & addressing() const noexcept
The addressing used for the list.
label size() const noexcept
The number of elements in the list.
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
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
const labelList & meshPoints() const
Return labelList of mesh points in patch.
const Field< point_type > & localPoints() const
Return pointField of points in patch.
const Field< point_type > & faceNormals() const
Return face unit normals for patch.
const labelListList & edgeFaces() const
Return edge-face addressing.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
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
A collection of cell labels.
Definition cellSet.H:50
A cell is defined as a list of faces with extra functionality.
Definition cell.H:56
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A list of face labels.
Definition faceSet.H:50
A subset of mesh faces organised as a primitive patch.
Definition faceZone.H:63
virtual void movePoints(const pointField &pts)
Correct patch after moving points.
Definition faceZone.C:686
static tmp< pointVectorField > makeDisplacementField(const pointMesh &pMesh, const labelList &adaptPatchIDs)
Helper function to make a pointVectorField with correct.
static autoPtr< indirectPrimitivePatch > makePatch(const polyMesh &, const labelList &)
Create patch from set of patches.
static bool checkMesh(const bool report, const polyMesh &mesh, const dictionary &dict, labelHashSet &wrongFaces, const bool dryRun=false)
Check mesh with mesh settings in dict. Collects incorrect faces.
Given a displacement moves the mesh by scaling the displacement back until there are no more mesh err...
Mesh representing a set of points created from polyMesh.
Definition pointMesh.H:49
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
label whichPatch(const label meshFacei) const
Return patch index for a given mesh face index. Uses binary search.
static bool checkFaceArea(const bool report, const scalar minArea, const polyMesh &, const vectorField &faceAreas, const labelList &checkFaces, labelHashSet *setPtr)
Small faces.
static bool checkCellDeterminant(const bool report, const scalar minDet, const polyMesh &, const vectorField &faceAreas, const labelList &checkFaces, const labelList &affectedCells, labelHashSet *setPtr)
Area of internal faces v.s. boundary faces.
static labelList affectedCells(const polyMesh &, const labelList &changedFaces)
Helper function: get affected cells from faces.
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
Simple container to keep together snap specific information.
static void preSmoothPatch(const meshRefinement &meshRefiner, const snapParameters &snapParams, const label nInitErrors, const List< labelPair > &baffles, motionSmoother &)
Smooth the mesh (patch and internal) to increase visibility.
static scalarField calcSnapDistance(const fvMesh &mesh, const snapParameters &snapParams, const indirectPrimitivePatch &)
Calculate edge length per patch point.
static labelList getUnnamedSurfaces(const PtrList< surfaceZonesInfo > &surfList)
Get indices of unnamed surfaces (surfaces without faceZoneName).
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 syncPointPositions(const polyMesh &mesh, List< point > &positions, const CombineOp &cop, const point &nullValue)
Synchronize locations on all mesh points.
Definition syncTools.H:240
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
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
static void syncEdgeList(const polyMesh &mesh, List< T > &edgeValues, const CombineOp &cop, const T &nullValue, const TransformOp &top, const FlipOp &fop)
Synchronize values on all mesh edges.
const polyBoundaryMesh & patches
bool coupled
const pointField & points
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))
word timeName
Definition getTimeIndex.H:3
Namespace for handling debugging switches.
Definition debug.C:45
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
dimensionedScalar sin(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Omanip< int > setw(const int i)
Definition IOmanip.H:199
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
List< bool > boolList
A List of bools.
Definition List.H:60
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
dimensionedScalar cos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
runTime write()
labelList f(nPoints)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235
Unit conversion functions.