Loading...
Searching...
No Matches
polyMeshFilter.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) 2012-2016 OpenFOAM Foundation
9 Copyright (C) 2022 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 "polyMeshFilter.H"
30#include "polyMesh.H"
31#include "fvMesh.H"
32#include "unitConversion.H"
33#include "edgeCollapser.H"
34#include "syncTools.H"
35#include "polyTopoChange.H"
36#include "globalIndex.H"
37#include "bitSet.H"
38#include "pointSet.H"
39#include "faceSet.H"
40#include "cellSet.H"
41
42// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
44namespace Foam
45{
46defineTypeNameAndDebug(polyMeshFilter, 0);
47}
48
49
50void Foam::polyMeshFilter::updateSets(const mapPolyMesh& map)
52 updateSets<pointSet>(map);
53 updateSets<faceSet>(map);
54 updateSets<cellSet>(map);
55}
56
57
58void Foam::polyMeshFilter::copySets
59(
60 const polyMesh& oldMesh,
61 const polyMesh& newMesh
62)
64 copySets<pointSet>(oldMesh, newMesh);
65 copySets<faceSet>(oldMesh, newMesh);
66 copySets<cellSet>(oldMesh, newMesh);
67}
68
69
71{
72 polyTopoChange originalMeshToNewMesh(mesh);
73
74 autoPtr<fvMesh> meshCopy;
75 autoPtr<mapPolyMesh> mapPtr = originalMeshToNewMesh.makeMesh
76 (
77 meshCopy,
79 (
80 mesh.name(),
81 mesh.polyMesh::instance(),
82 mesh.time(),
83 IOobject::READ_IF_PRESENT, // read fv* if present
86 ),
87 mesh,
88 true // parallel sync
89 );
90
91 const mapPolyMesh& map = mapPtr();
92
93 // Update fields
94 meshCopy().updateMesh(map);
95 if (map.hasMotionPoints())
96 {
97 meshCopy().movePoints(map.preMotionPoints());
98 }
99
100 copySets(mesh, meshCopy());
101
102 return meshCopy;
103}
104
105
106// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
107
108Foam::label Foam::polyMeshFilter::filterFacesLoop(const label nOriginalBadFaces)
109{
110 label nBadFaces = labelMax;
111 label nOuterIterations = 0;
112
113 // Maintain the number of times a point has been part of a bad face
114 labelList pointErrorCount(mesh_.nPoints(), Zero);
115
116 bitSet newErrorPoint(mesh_.nPoints());
118 (
119 mesh_,
120 meshQualityCoeffDict(),
121 newErrorPoint
122 );
123
124 bool newBadFaces = true;
125
126 // Main loop
127 // ~~~~~~~~~
128 // It tries and do some collapses, checks the resulting mesh and
129 // 'freezes' some edges (by marking in minEdgeLen) and tries again.
130 // This will iterate ultimately to the situation where every edge is
131 // frozen and nothing gets collapsed.
132 while
133 (
134 nOuterIterations < maxIterations()
135 //&& nBadFaces > nOriginalBadFaces
136 && newBadFaces
137 )
138 {
139 Info<< nl << "Outer Iteration = " << nOuterIterations++ << nl
140 << endl;
141
142 printScalarFieldStats("Edge Filter Factor", minEdgeLen_);
143 printScalarFieldStats("Face Filter Factor", faceFilterFactor_);
144
145 // Reset the new mesh to the old mesh
146 newMeshPtr_ = copyMesh(mesh_);
147 fvMesh& newMesh = newMeshPtr_();
148
149 scalarField newMeshFaceFilterFactor = faceFilterFactor_;
150 pointPriority_.reset(new labelList(originalPointPriority_));
151
152 labelList origToCurrentPointMap(identity(newMesh.nPoints()));
153
154 {
155 label nInnerIterations = 0;
156 label nPrevLocalCollapse = labelMax;
157
159
160 while (true)
161 {
162 Info<< nl << indent << "Inner iteration = "
163 << nInnerIterations++ << nl << incrIndent << endl;
164
165 label nLocalCollapse = filterFaces
166 (
167 newMesh,
168 newMeshFaceFilterFactor,
169 origToCurrentPointMap
170 );
172
173 if
174 (
175 nLocalCollapse >= nPrevLocalCollapse
176 || nLocalCollapse == 0
177 )
178 {
180 break;
181 }
182 else
183 {
184 nPrevLocalCollapse = nLocalCollapse;
185 }
186 }
187 }
188
189 scalarField newMeshMinEdgeLen = minEdgeLen_;
190
191 {
192 label nInnerIterations = 0;
193 label nPrevLocalCollapse = labelMax;
194
196
197 while (true)
198 {
199 Info<< nl << indent << "Inner iteration = "
200 << nInnerIterations++ << nl << incrIndent << endl;
201
202 label nLocalCollapse = filterEdges
203 (
204 newMesh,
205 newMeshMinEdgeLen,
206 origToCurrentPointMap
207 );
209
210 if
211 (
212 nLocalCollapse >= nPrevLocalCollapse
213 || nLocalCollapse == 0
214 )
215 {
217 break;
218 }
219 else
220 {
221 nPrevLocalCollapse = nLocalCollapse;
222 }
223 }
224 }
225
226
227 // Mesh check
228 // ~~~~~~~~~~~~~~~~~~
229 // Do not allow collapses in regions of error.
230 // Updates minEdgeLen, nRelaxedEdges
231
232 if (controlMeshQuality())
233 {
234 bitSet isErrorPoint(newMesh.nPoints());
236 (
237 newMesh,
238 meshQualityCoeffDict(),
239 isErrorPoint
240 );
241
242 Info<< nl << " Number of bad faces : " << nBadFaces << nl
243 << " Number of marked points : "
244 << returnReduce(isErrorPoint.count(), sumOp<unsigned int>())
245 << endl;
246
247 updatePointErrorCount
248 (
249 isErrorPoint,
250 origToCurrentPointMap,
251 pointErrorCount
252 );
253
254 checkMeshEdgesAndRelaxEdges
255 (
256 newMesh,
257 origToCurrentPointMap,
258 isErrorPoint,
259 pointErrorCount
260 );
261
262 checkMeshFacesAndRelaxEdges
263 (
264 newMesh,
265 origToCurrentPointMap,
266 isErrorPoint,
267 pointErrorCount
268 );
269
270 newBadFaces = false;
271 forAll(mesh_.points(), pI)
272 {
273 if
274 (
275 origToCurrentPointMap[pI] >= 0
276 && isErrorPoint[origToCurrentPointMap[pI]]
277 )
278 {
279 if (!newErrorPoint[pI])
280 {
281 newBadFaces = true;
282 break;
283 }
284 }
285 }
286
287 Pstream::reduceOr(newBadFaces);
288 }
289 else
290 {
291 return -1;
292 }
293 }
294
295 return nBadFaces;
296}
297
298
299Foam::label Foam::polyMeshFilter::filterFaces
300(
301 polyMesh& newMesh,
302 scalarField& newMeshFaceFilterFactor,
303 labelList& origToCurrentPointMap
304)
305{
306 // Per edge collapse status
307 bitSet collapseEdge(newMesh.nEdges());
308
309 Map<point> collapsePointToLocation(newMesh.nPoints());
310
311 edgeCollapser collapser(newMesh, collapseFacesCoeffDict());
312
313 {
314 // Collapse faces
315 labelPair nCollapsedPtEdge = collapser.markSmallSliverFaces
316 (
317 newMeshFaceFilterFactor,
318 pointPriority_(),
320 collapsePointToLocation
321 );
322
323 label nCollapsed = 0;
324 forAll(nCollapsedPtEdge, collapseTypeI)
325 {
326 nCollapsed += nCollapsedPtEdge[collapseTypeI];
327 }
328
329 reduce(nCollapsed, sumOp<label>());
330
331 label nToPoint = returnReduce(nCollapsedPtEdge.first(), sumOp<label>());
332 label nToEdge = returnReduce(nCollapsedPtEdge.second(), sumOp<label>());
333
334 Info<< indent
335 << "Collapsing " << nCollapsed << " faces "
336 << "(to point = " << nToPoint << ", to edge = " << nToEdge << ")"
337 << endl;
338
339 if (nCollapsed == 0)
340 {
341 return 0;
342 }
343 }
344
345 // Merge edge collapses into consistent collapse-network.
346 // Make sure no cells get collapsed.
347 List<pointEdgeCollapse> allPointInfo;
348 const globalIndex globalPoints(newMesh.nPoints());
349
350 collapser.consistentCollapse
351 (
353 pointPriority_(),
354 collapsePointToLocation,
356 allPointInfo
357 );
358
359 label nLocalCollapse = collapseEdge.count();
360
361 reduce(nLocalCollapse, sumOp<label>());
362 Info<< nl << indent << "Collapsing " << nLocalCollapse
363 << " edges after synchronisation and PointEdgeWave" << endl;
364
365 if (nLocalCollapse == 0)
366 {
367 return 0;
368 }
369
370 {
371 // Apply collapses to current mesh
372 polyTopoChange newMeshMod(newMesh);
373
374 // Insert mesh refinement into polyTopoChange.
375 collapser.setRefinement(allPointInfo, newMeshMod);
376
377 Info<< indent << "Apply changes to the current mesh" << endl;
378
379 // Apply changes to current mesh
380 autoPtr<mapPolyMesh> newMapPtr = newMeshMod.changeMesh
381 (
382 newMesh,
383 false
384 );
385 const mapPolyMesh& newMap = newMapPtr();
386
387 // Update fields
388 newMesh.updateMesh(newMap);
389 if (newMap.hasMotionPoints())
390 {
391 newMesh.movePoints(newMap.preMotionPoints());
392 }
393 updateSets(newMap);
394
395 updatePointPriorities(newMesh, newMap.pointMap());
396
397 mapOldMeshFaceFieldToNewMesh
398 (
399 newMesh,
400 newMap.faceMap(),
401 newMeshFaceFilterFactor
402 );
403
404 updateOldToNewPointMap
405 (
406 newMap.reversePointMap(),
407 origToCurrentPointMap
408 );
409 }
410
411 return nLocalCollapse;
412}
413
414
415Foam::label Foam::polyMeshFilter::filterEdges
416(
417 polyMesh& newMesh,
418 scalarField& newMeshMinEdgeLen,
419 labelList& origToCurrentPointMap
420)
421{
422 // Per edge collapse status
423 bitSet collapseEdge(newMesh.nEdges());
424
425 Map<point> collapsePointToLocation(newMesh.nPoints());
426
427 edgeCollapser collapser(newMesh, collapseFacesCoeffDict());
428
429 // Work out which edges to collapse
430 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
431 // This is by looking at minEdgeLen (to avoid frozen edges)
432 // and marking in collapseEdge.
433 label nSmallCollapsed = collapser.markSmallEdges
434 (
435 newMeshMinEdgeLen,
436 pointPriority_(),
438 collapsePointToLocation
439 );
440
441 reduce(nSmallCollapsed, sumOp<label>());
442 Info<< indent << "Collapsing " << nSmallCollapsed
443 << " small edges" << endl;
444
445 // Merge inline edges
446 label nMerged = collapser.markMergeEdges
447 (
448 maxCos(),
449 pointPriority_(),
451 collapsePointToLocation
452 );
453
454 reduce(nMerged, sumOp<label>());
455 Info<< indent << "Collapsing " << nMerged << " in line edges"
456 << endl;
457
458 if (nMerged + nSmallCollapsed == 0)
459 {
460 return 0;
461 }
462
463 // Merge edge collapses into consistent collapse-network.
464 // Make sure no cells get collapsed.
465 List<pointEdgeCollapse> allPointInfo;
466 const globalIndex globalPoints(newMesh.nPoints());
467
468 collapser.consistentCollapse
469 (
471 pointPriority_(),
472 collapsePointToLocation,
474 allPointInfo
475 );
476
477 label nLocalCollapse = collapseEdge.count();
478
479 reduce(nLocalCollapse, sumOp<label>());
480 Info<< nl << indent << "Collapsing " << nLocalCollapse
481 << " edges after synchronisation and PointEdgeWave" << endl;
482
483 if (nLocalCollapse == 0)
484 {
485 return 0;
486 }
487
488 // Apply collapses to current mesh
489 polyTopoChange newMeshMod(newMesh);
490
491 // Insert mesh refinement into polyTopoChange.
492 collapser.setRefinement(allPointInfo, newMeshMod);
493
494 Info<< indent << "Apply changes to the current mesh" << endl;
495
496 // Apply changes to current mesh
497 autoPtr<mapPolyMesh> newMapPtr = newMeshMod.changeMesh
498 (
499 newMesh,
500 false
501 );
502 const mapPolyMesh& newMap = newMapPtr();
503
504 // Update fields
505 newMesh.updateMesh(newMap);
506 if (newMap.hasMotionPoints())
507 {
508 newMesh.movePoints(newMap.preMotionPoints());
509 }
510 updateSets(newMap);
511
512 // Synchronise the factors
513 mapOldMeshEdgeFieldToNewMesh
514 (
515 newMesh,
516 newMap.pointMap(),
517 newMeshMinEdgeLen
518 );
519
520 updateOldToNewPointMap
521 (
522 newMap.reversePointMap(),
523 origToCurrentPointMap
524 );
525
526 updatePointPriorities(newMesh, newMap.pointMap());
527
528 return nLocalCollapse;
529}
530
531
532void Foam::polyMeshFilter::updatePointErrorCount
533(
534 const bitSet& isErrorPoint,
535 const labelList& oldToNewMesh,
536 labelList& pointErrorCount
537) const
538{
539 forAll(mesh_.points(), pI)
540 {
541 if (isErrorPoint[oldToNewMesh[pI]])
542 {
543 pointErrorCount[pI]++;
544 }
545 }
546}
547
548
549void Foam::polyMeshFilter::checkMeshEdgesAndRelaxEdges
550(
551 const polyMesh& newMesh,
552 const labelList& oldToNewMesh,
553 const bitSet& isErrorPoint,
554 const labelList& pointErrorCount
555)
556{
557 const edgeList& edges = mesh_.edges();
558
559 forAll(edges, edgeI)
560 {
561 const edge& e = edges[edgeI];
562 label newStart = oldToNewMesh[e[0]];
563 label newEnd = oldToNewMesh[e[1]];
564
565 if
566 (
567 pointErrorCount[e[0]] >= maxPointErrorCount()
568 || pointErrorCount[e[1]] >= maxPointErrorCount()
569 )
570 {
571 minEdgeLen_[edgeI] = -1;
572 }
573
574 if
575 (
576 (newStart >= 0 && isErrorPoint[newStart])
577 || (newEnd >= 0 && isErrorPoint[newEnd])
578 )
579 {
580 minEdgeLen_[edgeI] *= edgeReductionFactor();
581 }
582 }
583
584 syncTools::syncEdgeList(mesh_, minEdgeLen_, minEqOp<scalar>(), scalar(0));
585
586 for (label smoothIter = 0; smoothIter < maxSmoothIters(); ++smoothIter)
587 {
588 // Smooth minEdgeLen
589 forAll(mesh_.edges(), edgeI)
590 {
591 const edge& e = mesh_.edges()[edgeI];
592
593 scalar sumMinEdgeLen = 0;
594 label nEdges = 0;
595
596 forAll(e, pointi)
597 {
598 const labelList& pEdges = mesh_.pointEdges()[e[pointi]];
599
600 forAll(pEdges, pEdgeI)
601 {
602 const label pEdge = pEdges[pEdgeI];
603 sumMinEdgeLen += minEdgeLen_[pEdge];
604 nEdges++;
605 }
606 }
607
608 minEdgeLen_[edgeI] = min
609 (
610 minEdgeLen_[edgeI],
611 sumMinEdgeLen/nEdges
612 );
613 }
614
616 (
617 mesh_,
618 minEdgeLen_,
620 scalar(0)
621 );
622 }
623}
624
625
626void Foam::polyMeshFilter::checkMeshFacesAndRelaxEdges
627(
628 const polyMesh& newMesh,
629 const labelList& oldToNewMesh,
630 const bitSet& isErrorPoint,
631 const labelList& pointErrorCount
632)
633{
634 const faceList& faces = mesh_.faces();
635
636 forAll(faces, facei)
637 {
638 const face& f = faces[facei];
639
640 forAll(f, fpI)
641 {
642 const label ptIndex = oldToNewMesh[f[fpI]];
643
644 if (pointErrorCount[f[fpI]] >= maxPointErrorCount())
645 {
646 faceFilterFactor_[facei] = -1;
647 }
648
649 if (isErrorPoint[ptIndex])
650 {
651 faceFilterFactor_[facei] *= faceReductionFactor();
652
653 break;
654 }
655 }
656 }
657
658 syncTools::syncFaceList(mesh_, faceFilterFactor_, minEqOp<scalar>());
659
660 for (label smoothIter = 0; smoothIter < maxSmoothIters(); ++smoothIter)
661 {
662 // Smooth faceFilterFactor
663 forAll(faces, facei)
664 {
665 const labelList& fEdges = mesh_.faceEdges()[facei];
666
667 scalar sumFaceFilterFactors = 0;
668 label nFaces = 0;
669
670 // This is important: Only smooth around faces that share an
671 // edge with a bad face
672 bool skipFace = true;
673
674 forAll(fEdges, fEdgeI)
675 {
676 const labelList& eFaces = mesh_.edgeFaces()[fEdges[fEdgeI]];
677
678 forAll(eFaces, eFacei)
679 {
680 const label eFace = eFaces[eFacei];
681
682 const face& f = faces[eFace];
683
684 forAll(f, fpI)
685 {
686 const label ptIndex = oldToNewMesh[f[fpI]];
687
688 if (isErrorPoint[ptIndex])
689 {
690 skipFace = false;
691 break;
692 }
693 }
694
695 if (eFace != facei)
696 {
697 sumFaceFilterFactors += faceFilterFactor_[eFace];
698 nFaces++;
699 }
700 }
701 }
702
703 if (skipFace)
704 {
705 continue;
706 }
707
708 faceFilterFactor_[facei] = min
709 (
710 faceFilterFactor_[facei],
711 sumFaceFilterFactors/nFaces
712 );
713 }
714
715 // Face filter factor needs to be synchronised!
716 syncTools::syncFaceList(mesh_, faceFilterFactor_, minEqOp<scalar>());
717 }
718}
719
720
721void Foam::polyMeshFilter::updatePointPriorities
722(
723 const polyMesh& newMesh,
724 const labelList& pointMap
725)
726{
727 labelList newPointPriority(newMesh.nPoints(), labelMin);
728 const labelList& currPointPriority = pointPriority_();
729
730 forAll(newPointPriority, ptI)
731 {
732 const label newPointToOldPoint = pointMap[ptI];
733 const label origPointPriority = currPointPriority[newPointToOldPoint];
734
735 newPointPriority[ptI] = max(origPointPriority, newPointPriority[ptI]);
736 }
737
739 (
740 newMesh,
741 newPointPriority,
744 );
745
746 pointPriority_.reset(new labelList(newPointPriority));
747}
748
749
750void Foam::polyMeshFilter::printScalarFieldStats
751(
752 const string desc,
753 const scalarField& fld
754) const
755{
756 scalar sum = 0;
757 scalar validElements = 0;
758 scalar min = GREAT;
759 scalar max = -GREAT;
760
761 forAll(fld, i)
762 {
763 const scalar fldElement = fld[i];
764
765 if (fldElement >= 0)
766 {
767 sum += fldElement;
768
769 if (fldElement < min)
770 {
771 min = fldElement;
772 }
773
774 if (fldElement > max)
775 {
776 max = fldElement;
777 }
778
779 validElements++;
780 }
781 }
782
786 reduce(validElements, sumOp<label>());
787 const label totFieldSize = returnReduce(fld.size(), sumOp<label>());
788
789 Info<< incrIndent << indent << desc
790 << ": min = " << min
791 << " av = " << sum/(validElements + SMALL)
792 << " max = " << max << nl
793 << indent
794 << " " << validElements << " / " << totFieldSize << " elements used"
795 << decrIndent << endl;
796}
797
798
799void Foam::polyMeshFilter::mapOldMeshEdgeFieldToNewMesh
800(
801 const polyMesh& newMesh,
802 const labelList& pointMap,
803 scalarField& newMeshMinEdgeLen
804) const
805{
806 scalarField tmp(newMesh.nEdges());
807
808 const edgeList& newEdges = newMesh.edges();
809
810 forAll(newEdges, newEdgeI)
811 {
812 const edge& newEdge = newEdges[newEdgeI];
813 const label pStart = newEdge.start();
814 const label pEnd = newEdge.end();
815
816 tmp[newEdgeI] = min
817 (
818 newMeshMinEdgeLen[pointMap[pStart]],
819 newMeshMinEdgeLen[pointMap[pEnd]]
820 );
821 }
822
823 newMeshMinEdgeLen.transfer(tmp);
824
826 (
827 newMesh,
828 newMeshMinEdgeLen,
830 scalar(0)
831 );
832}
833
834
835void Foam::polyMeshFilter::mapOldMeshFaceFieldToNewMesh
836(
837 const polyMesh& newMesh,
838 const labelList& faceMap,
839 scalarField& newMeshFaceFilterFactor
840) const
841{
842 scalarField tmp(newMesh.nFaces());
843
844 forAll(faceMap, newFacei)
845 {
846 const label oldFacei = faceMap[newFacei];
847
848 tmp[newFacei] = newMeshFaceFilterFactor[oldFacei];
849 }
850
851 newMeshFaceFilterFactor.transfer(tmp);
852
854 (
855 newMesh,
856 newMeshFaceFilterFactor,
858 );
859}
860
861
862void Foam::polyMeshFilter::updateOldToNewPointMap
863(
864 const labelList& currToNew,
865 labelList& origToCurrentPointMap
866) const
867{
868 forAll(origToCurrentPointMap, origPointi)
869 {
870 label oldPointi = origToCurrentPointMap[origPointi];
871
872 if (oldPointi != -1)
873 {
874 label newPointi = currToNew[oldPointi];
875
876 if (newPointi >= 0)
877 {
878 origToCurrentPointMap[origPointi] = newPointi;
879 }
880 else if (newPointi == -1)
881 {
882 origToCurrentPointMap[origPointi] = -1;
883 }
884 else
885 {
886 origToCurrentPointMap[origPointi] = -newPointi-2;
887 }
889 }
890}
891
892
893// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
894
895Foam::polyMeshFilter::polyMeshFilter(const fvMesh& mesh)
896:
898 (
900 (
902 (
903 "collapseDict",
904 mesh.time().system(),
905 mesh.time(),
906 IOobject::MUST_READ,
907 IOobject::NO_WRITE
908 )
909 )
910 ),
911 mesh_(mesh),
912 newMeshPtr_(),
913 originalPointPriority_(mesh.nPoints(), labelMin),
914 pointPriority_(),
915 minEdgeLen_(),
916 faceFilterFactor_()
917{
919}
920
921
922Foam::polyMeshFilter::polyMeshFilter
923(
924 const fvMesh& mesh,
925 const labelList& pointPriority
926)
927:
929 (
931 (
933 (
934 "collapseDict",
935 mesh.time().system(),
936 mesh.time(),
937 IOobject::MUST_READ,
938 IOobject::NO_WRITE
939 )
940 )
941 ),
942 mesh_(mesh),
943 newMeshPtr_(),
944 originalPointPriority_(pointPriority),
945 pointPriority_(),
946 minEdgeLen_(),
947 faceFilterFactor_()
948{
950}
951
952Foam::polyMeshFilter::polyMeshFilter
953(
954 const fvMesh& mesh,
955 const labelList& pointPriority,
956 const dictionary& dict
957)
958:
960 mesh_(mesh),
961 newMeshPtr_(),
962 originalPointPriority_(pointPriority),
963 pointPriority_(),
964 minEdgeLen_(),
965 faceFilterFactor_()
968}
969
970
971// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
972
973Foam::label Foam::polyMeshFilter::filter(const label nOriginalBadFaces)
974{
975 minEdgeLen_.resize(mesh_.nEdges(), minLen());
976 faceFilterFactor_.resize(mesh_.nFaces(), initialFaceLengthFactor());
977
978 return filterFacesLoop(nOriginalBadFaces);
979}
980
981
982Foam::label Foam::polyMeshFilter::filter(const faceSet& fSet)
983{
984 minEdgeLen_.resize(mesh_.nEdges(), minLen());
985 faceFilterFactor_.resize(mesh_.nFaces(), initialFaceLengthFactor());
986
987 forAll(faceFilterFactor_, fI)
988 {
989 if (!fSet.found(fI))
990 {
991 faceFilterFactor_[fI] = -1;
993 }
994
995 return filterFacesLoop(0);
996}
997
998
999Foam::label Foam::polyMeshFilter::filterEdges
1000(
1001 const label nOriginalBadFaces
1002)
1003{
1004 label nBadFaces = labelMax/2;
1005 label nPreviousBadFaces = labelMax;
1006 label nOuterIterations = 0;
1007
1008 minEdgeLen_.resize(mesh_.nEdges(), minLen());
1009 faceFilterFactor_.clear();
1010
1011 labelList pointErrorCount(mesh_.nPoints(), Zero);
1012
1013 // Main loop
1014 // ~~~~~~~~~
1015 // It tries and do some collapses, checks the resulting mesh and
1016 // 'freezes' some edges (by marking in minEdgeLen) and tries again.
1017 // This will iterate ultimately to the situation where every edge is
1018 // frozen and nothing gets collapsed.
1019 while
1020 (
1021 nOuterIterations < maxIterations()
1022 && nBadFaces > nOriginalBadFaces
1023 && nBadFaces < nPreviousBadFaces
1024 )
1025 {
1026 Info<< nl << "Outer Iteration = " << nOuterIterations++ << nl
1027 << endl;
1028
1029 printScalarFieldStats("Edge Filter Factor", minEdgeLen_);
1030
1031 nPreviousBadFaces = nBadFaces;
1032
1033 // Reset the new mesh to the old mesh
1034 newMeshPtr_ = copyMesh(mesh_);
1035 fvMesh& newMesh = newMeshPtr_();
1036
1037 scalarField newMeshMinEdgeLen = minEdgeLen_;
1038 pointPriority_.reset(new labelList(originalPointPriority_));
1039
1040 labelList origToCurrentPointMap(identity(newMesh.nPoints()));
1041
1042 Info<< incrIndent;
1043
1044 label nInnerIterations = 0;
1045 label nPrevLocalCollapse = labelMax;
1046
1047 while (true)
1048 {
1049 Info<< nl
1050 << indent << "Inner iteration = " << nInnerIterations++ << nl
1051 << incrIndent << endl;
1052
1053 label nLocalCollapse = filterEdges
1054 (
1055 newMesh,
1056 newMeshMinEdgeLen,
1057 origToCurrentPointMap
1058 );
1059
1060 Info<< decrIndent;
1061
1062 if
1063 (
1064 nLocalCollapse >= nPrevLocalCollapse
1065 || nLocalCollapse == 0
1066 )
1067 {
1068 Info<< decrIndent;
1069 break;
1070 }
1071 else
1072 {
1073 nPrevLocalCollapse = nLocalCollapse;
1074 }
1075 }
1076
1077 // Mesh check
1078 // ~~~~~~~~~~~~~~~~~~
1079 // Do not allow collapses in regions of error.
1080 // Updates minEdgeLen, nRelaxedEdges
1081
1082 if (controlMeshQuality())
1083 {
1084 bitSet isErrorPoint(newMesh.nPoints());
1086 (
1087 newMesh,
1088 meshQualityCoeffDict(),
1089 isErrorPoint
1090 );
1091
1092 Info<< nl << " Number of bad faces : " << nBadFaces << nl
1093 << " Number of marked points : "
1094 << returnReduce(isErrorPoint.count(), sumOp<unsigned int>())
1095 << endl;
1096
1097 updatePointErrorCount
1098 (
1099 isErrorPoint,
1100 origToCurrentPointMap,
1101 pointErrorCount
1102 );
1103
1104 checkMeshEdgesAndRelaxEdges
1105 (
1106 newMesh,
1107 origToCurrentPointMap,
1108 isErrorPoint,
1109 pointErrorCount
1110 );
1111 }
1112 else
1113 {
1114 return -1;
1116 }
1117
1118 return nBadFaces;
1119}
1120
1121
1123{
1124 return newMeshPtr_;
1125}
1126
1127
1130{
1131 return pointPriority_;
1132}
1133
1134
1135// ************************************************************************* //
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_REGISTER
Do not request registration (bool: false).
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
static void reduceOr(bool &value, const int communicator=worldComm)
Logical (or) reduction (MPI_AllReduce).
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
unsigned int count(const bool on=true) const
Count number of bits set.
Definition bitSetI.H:420
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Does polyTopoChanges to remove edges. Can remove faces due to edge collapse but can not remove cells ...
static label checkMeshQuality(const polyMesh &mesh, const dictionary &meshQualityDict, bitSet &isErrorPoint)
Check mesh and mark points on faces in error.
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
label start() const noexcept
The start (first) vertex label.
Definition edge.H:155
A list of face labels.
Definition faceSet.H:50
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
Calculates points shared by more than two processor patches or cyclic patches.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const pointField & preMotionPoints() const noexcept
Pre-motion point positions.
bool hasMotionPoints() const noexcept
Has valid preMotionPoints?
Class to store the settings for the polyMeshFilter class.
polyMeshFilterSettings(const dictionary &dict)
Construct from dictionary.
void writeSettings(Ostream &os) const
Write the settings to a stream.
const scalar & initialFaceLengthFactor() const
const dictionary & meshQualityCoeffDict() const
Remove the edges and faces of a polyMesh whilst satisfying the given mesh quality criteria.
const autoPtr< labelList > & pointPriority() const
Return the new pointPriority list.
label filter(const label nOriginalBadFaces)
Filter edges and faces.
static autoPtr< fvMesh > copyMesh(const fvMesh &mesh)
Return a copy of an fvMesh.
const autoPtr< fvMesh > & filteredMesh() const
Return reference to the filtered mesh. Does not check if the.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
Direct mesh changes based on v1.3 polyTopoChange syntax.
autoPtr< mapPolyMesh > makeMesh(autoPtr< Type > &newMesh, const IOobject &io, const polyMesh &mesh, const labelUList &patchMap, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Create new mesh with old mesh patches. Additional dictionaries.
label nPoints() const noexcept
Number of mesh points.
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.
A class for managing temporary objects.
Definition tmp.H:75
virtual bool found(const label id) const
Has the given index?
Definition topoSet.C:539
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
label collapseEdge(triSurface &surf, const scalar minLen)
Keep collapsing all edges < minLen.
dynamicFvMesh & mesh
label nPoints
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Pair< label > labelPair
A pair of labels.
Definition Pair.H:54
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
int system(const std::string &command, const bool bg=false)
Execute the specified command via the shell.
Definition POSIX.C:1704
List< label > labelList
A List of labels.
Definition List.H:62
constexpr label labelMin
Definition label.H:54
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
constexpr label labelMax
Definition label.H:55
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition Ostream.H:490
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Ostream & indent(Ostream &os)
Indent stream.
Definition Ostream.H:481
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
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...
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition Ostream.H:499
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
label newPointi
labelList f(nPoints)
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Unit conversion functions.