polyMeshFilter.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2012-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "polyMeshFilter.H"
27 #include "polyMesh.H"
28 #include "fvMesh.H"
29 #include "edgeCollapser.H"
30 #include "syncTools.H"
31 #include "polyTopoChange.H"
32 #include "globalIndex.H"
33 #include "PackedBoolList.H"
34 #include "pointSet.H"
35 #include "faceSet.H"
36 #include "cellSet.H"
37 
38 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
43 }
44 
45 
46 void Foam::polyMeshFilter::updateSets(const polyTopoChangeMap& map)
47 {
48  updateSets<pointSet>(map);
49  updateSets<faceSet>(map);
50  updateSets<cellSet>(map);
51 }
52 
53 
54 void Foam::polyMeshFilter::copySets
55 (
56  const polyMesh& oldMesh,
57  const polyMesh& newMesh
58 )
59 {
60  copySets<pointSet>(oldMesh, newMesh);
61  copySets<faceSet>(oldMesh, newMesh);
62  copySets<cellSet>(oldMesh, newMesh);
63 }
64 
65 
67 {
68  polyTopoChange originalMeshToNewMesh(mesh);
69 
70  autoPtr<fvMesh> meshCopy;
71  autoPtr<polyTopoChangeMap> mapPtr = originalMeshToNewMesh.makeMesh
72  (
73  meshCopy,
74  IOobject
75  (
76  mesh.name(),
77  mesh.polyMesh::instance(),
78  mesh.time(),
81  false
82  ),
83  mesh,
84  true // parallel sync
85  );
86 
87  const polyTopoChangeMap& map = mapPtr();
88 
89  // Update fields
90  meshCopy().topoChange(map);
91 
92  copySets(mesh, meshCopy());
93 
94  return meshCopy;
95 }
96 
97 
98 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
99 
100 Foam::label Foam::polyMeshFilter::filterFacesLoop(const label nOriginalBadFaces)
101 {
102  label nBadFaces = labelMax;
103  label nOuterIterations = 0;
104 
105  // Maintain the number of times a point has been part of a bad face
106  labelList pointErrorCount(mesh_.nPoints(), 0);
107 
108  PackedBoolList newErrorPoint(mesh_.nPoints());
110  (
111  mesh_,
112  meshQualityCoeffDict(),
113  newErrorPoint
114  );
115 
116  bool newBadFaces = true;
117 
118  // Main loop
119  // ~~~~~~~~~
120  // It tries and do some collapses, checks the resulting mesh and
121  // 'freezes' some edges (by marking in minEdgeLen) and tries again.
122  // This will iterate ultimately to the situation where every edge is
123  // frozen and nothing gets collapsed.
124  while
125  (
126  nOuterIterations < maxIterations()
127  //&& nBadFaces > nOriginalBadFaces
128  && newBadFaces
129  )
130  {
131  Info<< nl << "Outer Iteration = " << nOuterIterations++ << nl
132  << endl;
133 
134  printScalarFieldStats("Edge Filter Factor", minEdgeLen_);
135  printScalarFieldStats("Face Filter Factor", faceFilterFactor_);
136 
137  // Reset the new mesh to the old mesh
138  newMeshPtr_ = copyMesh(mesh_);
139  fvMesh& newMesh = newMeshPtr_();
140 
141  scalarField newMeshFaceFilterFactor = faceFilterFactor_;
142  pointPriority_.reset(new labelList(originalPointPriority_));
143 
144  labelList origToCurrentPointMap(identityMap(newMesh.nPoints()));
145 
146  {
147  label nInnerIterations = 0;
148  label nPrevLocalCollapse = labelMax;
149 
150  Info<< incrIndent;
151 
152  while (true)
153  {
154  Info<< nl << indent << "Inner iteration = "
155  << nInnerIterations++ << nl << incrIndent << endl;
156 
157  label nLocalCollapse = filterFaces
158  (
159  newMesh,
160  newMeshFaceFilterFactor,
161  origToCurrentPointMap
162  );
163  Info<< decrIndent;
164 
165  if
166  (
167  nLocalCollapse >= nPrevLocalCollapse
168  || nLocalCollapse == 0
169  )
170  {
171  Info<< decrIndent;
172  break;
173  }
174  else
175  {
176  nPrevLocalCollapse = nLocalCollapse;
177  }
178  }
179  }
180 
181  scalarField newMeshMinEdgeLen = minEdgeLen_;
182 
183  {
184  label nInnerIterations = 0;
185  label nPrevLocalCollapse = labelMax;
186 
187  Info<< incrIndent;
188 
189  while (true)
190  {
191  Info<< nl << indent << "Inner iteration = "
192  << nInnerIterations++ << nl << incrIndent << endl;
193 
194  label nLocalCollapse = filterEdges
195  (
196  newMesh,
197  newMeshMinEdgeLen,
198  origToCurrentPointMap
199  );
200  Info<< decrIndent;
201 
202  if
203  (
204  nLocalCollapse >= nPrevLocalCollapse
205  || nLocalCollapse == 0
206  )
207  {
208  Info<< decrIndent;
209  break;
210  }
211  else
212  {
213  nPrevLocalCollapse = nLocalCollapse;
214  }
215  }
216  }
217 
218 
219  // Mesh check
220  // ~~~~~~~~~~~~~~~~~~
221  // Do not allow collapses in regions of error.
222  // Updates minEdgeLen, nRelaxedEdges
223 
224  if (controlMeshQuality())
225  {
226  PackedBoolList isErrorPoint(newMesh.nPoints());
228  (
229  newMesh,
230  meshQualityCoeffDict(),
231  isErrorPoint
232  );
233 
234  Info<< nl << " Number of bad faces : " << nBadFaces << nl
235  << " Number of marked points : "
236  << returnReduce(isErrorPoint.count(), sumOp<unsigned int>())
237  << endl;
238 
239  updatePointErrorCount
240  (
241  isErrorPoint,
242  origToCurrentPointMap,
243  pointErrorCount
244  );
245 
246  checkMeshEdgesAndRelaxEdges
247  (
248  newMesh,
249  origToCurrentPointMap,
250  isErrorPoint,
251  pointErrorCount
252  );
253 
254  checkMeshFacesAndRelaxEdges
255  (
256  newMesh,
257  origToCurrentPointMap,
258  isErrorPoint,
259  pointErrorCount
260  );
261 
262  newBadFaces = false;
263  forAll(mesh_.points(), pI)
264  {
265  if
266  (
267  origToCurrentPointMap[pI] >= 0
268  && isErrorPoint[origToCurrentPointMap[pI]]
269  )
270  {
271  if (!newErrorPoint[pI])
272  {
273  newBadFaces = true;
274  break;
275  }
276  }
277  }
278 
279  reduce(newBadFaces, orOp<bool>());
280  }
281  else
282  {
283  return -1;
284  }
285  }
286 
287  return nBadFaces;
288 }
289 
290 
291 Foam::label Foam::polyMeshFilter::filterFaces
292 (
293  polyMesh& newMesh,
294  scalarField& newMeshFaceFilterFactor,
295  labelList& origToCurrentPointMap
296 )
297 {
298  // Per edge collapse status
299  PackedBoolList collapseEdge(newMesh.nEdges());
300 
301  Map<point> collapsePointToLocation(newMesh.nPoints());
302 
303  edgeCollapser collapser(newMesh, collapseFacesCoeffDict());
304 
305  {
306  // Collapse faces
307  labelPair nCollapsedPtEdge = collapser.markSmallSliverFaces
308  (
309  newMeshFaceFilterFactor,
310  pointPriority_(),
311  collapseEdge,
312  collapsePointToLocation
313  );
314 
315  label nCollapsed = 0;
316  forAll(nCollapsedPtEdge, collapseTypeI)
317  {
318  nCollapsed += nCollapsedPtEdge[collapseTypeI];
319  }
320 
321  reduce(nCollapsed, sumOp<label>());
322 
323  label nToPoint = returnReduce(nCollapsedPtEdge.first(), sumOp<label>());
324  label nToEdge = returnReduce(nCollapsedPtEdge.second(), sumOp<label>());
325 
326  Info<< indent
327  << "Collapsing " << nCollapsed << " faces "
328  << "(to point = " << nToPoint << ", to edge = " << nToEdge << ")"
329  << endl;
330 
331  if (nCollapsed == 0)
332  {
333  return 0;
334  }
335  }
336 
337  // Merge edge collapses into consistent collapse-network.
338  // Make sure no cells get collapsed.
339  List<pointEdgeCollapse> allPointInfo;
340  const globalIndex globalPoints(newMesh.nPoints());
341 
342  collapser.consistentCollapse
343  (
344  globalPoints,
345  pointPriority_(),
346  collapsePointToLocation,
347  collapseEdge,
348  allPointInfo
349  );
350 
351  label nLocalCollapse = collapseEdge.count();
352 
353  reduce(nLocalCollapse, sumOp<label>());
354  Info<< nl << indent << "Collapsing " << nLocalCollapse
355  << " edges after synchronisation and PointEdgeWave" << endl;
356 
357  if (nLocalCollapse == 0)
358  {
359  return 0;
360  }
361 
362  {
363  // Apply collapses to current mesh
364  polyTopoChange newMeshMod(newMesh);
365 
366  // Insert mesh refinement into polyTopoChange.
367  collapser.setRefinement(allPointInfo, newMeshMod);
368 
369  Info<< indent << "Apply changes to the current mesh" << endl;
370 
371  // Apply changes to current mesh
372  autoPtr<polyTopoChangeMap> newMapPtr = newMeshMod.changeMesh(newMesh);
373  const polyTopoChangeMap& newMap = newMapPtr();
374 
375  // Update fields
376  newMesh.topoChange(newMap);
377  updateSets(newMap);
378 
379  updatePointPriorities(newMesh, newMap.pointMap());
380 
381  mapOldMeshFaceFieldToNewMesh
382  (
383  newMesh,
384  newMap.faceMap(),
385  newMeshFaceFilterFactor
386  );
387 
388  updateOldToNewPointMap
389  (
390  newMap.reversePointMap(),
391  origToCurrentPointMap
392  );
393  }
394 
395  return nLocalCollapse;
396 }
397 
398 
399 Foam::label Foam::polyMeshFilter::filterEdges
400 (
401  polyMesh& newMesh,
402  scalarField& newMeshMinEdgeLen,
403  labelList& origToCurrentPointMap
404 )
405 {
406  // Per edge collapse status
407  PackedBoolList collapseEdge(newMesh.nEdges());
408 
409  Map<point> collapsePointToLocation(newMesh.nPoints());
410 
411  edgeCollapser collapser(newMesh, collapseFacesCoeffDict());
412 
413  // Work out which edges to collapse
414  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
415  // This is by looking at minEdgeLen (to avoid frozen edges)
416  // and marking in collapseEdge.
417  label nSmallCollapsed = collapser.markSmallEdges
418  (
419  newMeshMinEdgeLen,
420  pointPriority_(),
421  collapseEdge,
422  collapsePointToLocation
423  );
424 
425  reduce(nSmallCollapsed, sumOp<label>());
426  Info<< indent << "Collapsing " << nSmallCollapsed
427  << " small edges" << endl;
428 
429  // Merge inline edges
430  label nMerged = collapser.markMergeEdges
431  (
432  maxCos(),
433  pointPriority_(),
434  collapseEdge,
435  collapsePointToLocation
436  );
437 
438  reduce(nMerged, sumOp<label>());
439  Info<< indent << "Collapsing " << nMerged << " in line edges"
440  << endl;
441 
442  if (nMerged + nSmallCollapsed == 0)
443  {
444  return 0;
445  }
446 
447  // Merge edge collapses into consistent collapse-network.
448  // Make sure no cells get collapsed.
449  List<pointEdgeCollapse> allPointInfo;
450  const globalIndex globalPoints(newMesh.nPoints());
451 
452  collapser.consistentCollapse
453  (
454  globalPoints,
455  pointPriority_(),
456  collapsePointToLocation,
457  collapseEdge,
458  allPointInfo
459  );
460 
461  label nLocalCollapse = collapseEdge.count();
462 
463  reduce(nLocalCollapse, sumOp<label>());
464  Info<< nl << indent << "Collapsing " << nLocalCollapse
465  << " edges after synchronisation and PointEdgeWave" << endl;
466 
467  if (nLocalCollapse == 0)
468  {
469  return 0;
470  }
471 
472  // Apply collapses to current mesh
473  polyTopoChange newMeshMod(newMesh);
474 
475  // Insert mesh refinement into polyTopoChange.
476  collapser.setRefinement(allPointInfo, newMeshMod);
477 
478  Info<< indent << "Apply changes to the current mesh" << endl;
479 
480  // Apply changes to current mesh
481  autoPtr<polyTopoChangeMap> newMapPtr = newMeshMod.changeMesh(newMesh);
482  const polyTopoChangeMap& newMap = newMapPtr();
483 
484  // Update fields
485  newMesh.topoChange(newMap);
486  updateSets(newMap);
487 
488  // Synchronise the factors
489  mapOldMeshEdgeFieldToNewMesh
490  (
491  newMesh,
492  newMap.pointMap(),
493  newMeshMinEdgeLen
494  );
495 
496  updateOldToNewPointMap
497  (
498  newMap.reversePointMap(),
499  origToCurrentPointMap
500  );
501 
502  updatePointPriorities(newMesh, newMap.pointMap());
503 
504  return nLocalCollapse;
505 }
506 
507 
508 void Foam::polyMeshFilter::updatePointErrorCount
509 (
510  const PackedBoolList& isErrorPoint,
511  const labelList& oldToNewMesh,
512  labelList& pointErrorCount
513 ) const
514 {
515  forAll(mesh_.points(), pI)
516  {
517  if (isErrorPoint[oldToNewMesh[pI]])
518  {
519  pointErrorCount[pI]++;
520  }
521  }
522 }
523 
524 
525 void Foam::polyMeshFilter::checkMeshEdgesAndRelaxEdges
526 (
527  const polyMesh& newMesh,
528  const labelList& oldToNewMesh,
529  const PackedBoolList& isErrorPoint,
530  const labelList& pointErrorCount
531 )
532 {
533  const edgeList& edges = mesh_.edges();
534 
535  forAll(edges, edgeI)
536  {
537  const edge& e = edges[edgeI];
538  label newStart = oldToNewMesh[e[0]];
539  label newEnd = oldToNewMesh[e[1]];
540 
541  if
542  (
543  pointErrorCount[e[0]] >= maxPointErrorCount()
544  || pointErrorCount[e[1]] >= maxPointErrorCount()
545  )
546  {
547  minEdgeLen_[edgeI] = -1;
548  }
549 
550  if
551  (
552  (newStart >= 0 && isErrorPoint[newStart])
553  || (newEnd >= 0 && isErrorPoint[newEnd])
554  )
555  {
556  minEdgeLen_[edgeI] *= edgeReductionFactor();
557  }
558  }
559 
560  syncTools::syncEdgeList(mesh_, minEdgeLen_, minEqOp<scalar>(), scalar(0));
561 
562  for (label smoothIter = 0; smoothIter < maxSmoothIters(); ++smoothIter)
563  {
564  // Smooth minEdgeLen
565  forAll(mesh_.edges(), edgeI)
566  {
567  const edge& e = mesh_.edges()[edgeI];
568 
569  scalar sumMinEdgeLen = 0;
570  label nEdges = 0;
571 
572  forAll(e, pointi)
573  {
574  const labelList& pEdges = mesh_.pointEdges()[e[pointi]];
575 
576  forAll(pEdges, pEdgeI)
577  {
578  const label pEdge = pEdges[pEdgeI];
579  sumMinEdgeLen += minEdgeLen_[pEdge];
580  nEdges++;
581  }
582  }
583 
584  minEdgeLen_[edgeI] = min
585  (
586  minEdgeLen_[edgeI],
587  sumMinEdgeLen/nEdges
588  );
589  }
590 
592  (
593  mesh_,
594  minEdgeLen_,
595  minEqOp<scalar>(),
596  scalar(0)
597  );
598  }
599 }
600 
601 
602 void Foam::polyMeshFilter::checkMeshFacesAndRelaxEdges
603 (
604  const polyMesh& newMesh,
605  const labelList& oldToNewMesh,
606  const PackedBoolList& isErrorPoint,
607  const labelList& pointErrorCount
608 )
609 {
610  const faceList& faces = mesh_.faces();
611 
612  forAll(faces, facei)
613  {
614  const face& f = faces[facei];
615 
616  forAll(f, fpI)
617  {
618  const label ptIndex = oldToNewMesh[f[fpI]];
619 
620  if (pointErrorCount[f[fpI]] >= maxPointErrorCount())
621  {
622  faceFilterFactor_[facei] = -1;
623  }
624 
625  if (isErrorPoint[ptIndex])
626  {
627  faceFilterFactor_[facei] *= faceReductionFactor();
628 
629  break;
630  }
631  }
632  }
633 
634  syncTools::syncFaceList(mesh_, faceFilterFactor_, minEqOp<scalar>());
635 
636  for (label smoothIter = 0; smoothIter < maxSmoothIters(); ++smoothIter)
637  {
638  // Smooth faceFilterFactor
639  forAll(faces, facei)
640  {
641  const labelList& fEdges = mesh_.faceEdges()[facei];
642 
643  scalar sumFaceFilterFactors = 0;
644  label nFaces = 0;
645 
646  // This is important: Only smooth around faces that share an
647  // edge with a bad face
648  bool skipFace = true;
649 
650  forAll(fEdges, fEdgeI)
651  {
652  const labelList& eFaces = mesh_.edgeFaces()[fEdges[fEdgeI]];
653 
654  forAll(eFaces, eFacei)
655  {
656  const label eFace = eFaces[eFacei];
657 
658  const face& f = faces[eFace];
659 
660  forAll(f, fpI)
661  {
662  const label ptIndex = oldToNewMesh[f[fpI]];
663 
664  if (isErrorPoint[ptIndex])
665  {
666  skipFace = false;
667  break;
668  }
669  }
670 
671  if (eFace != facei)
672  {
673  sumFaceFilterFactors += faceFilterFactor_[eFace];
674  nFaces++;
675  }
676  }
677  }
678 
679  if (skipFace)
680  {
681  continue;
682  }
683 
684  faceFilterFactor_[facei] = min
685  (
686  faceFilterFactor_[facei],
687  sumFaceFilterFactors/nFaces
688  );
689  }
690 
691  // Face filter factor needs to be synchronised!
692  syncTools::syncFaceList(mesh_, faceFilterFactor_, minEqOp<scalar>());
693  }
694 }
695 
696 
697 void Foam::polyMeshFilter::updatePointPriorities
698 (
699  const polyMesh& newMesh,
700  const labelList& pointMap
701 )
702 {
703  labelList newPointPriority(newMesh.nPoints(), labelMin);
704  const labelList& currPointPriority = pointPriority_();
705 
706  forAll(newPointPriority, ptI)
707  {
708  const label newPointToOldPoint = pointMap[ptI];
709  const label origPointPriority = currPointPriority[newPointToOldPoint];
710 
711  newPointPriority[ptI] = max(origPointPriority, newPointPriority[ptI]);
712  }
713 
715  (
716  newMesh,
717  newPointPriority,
718  maxEqOp<label>(),
719  labelMin
720  );
721 
722  pointPriority_.reset(new labelList(newPointPriority));
723 }
724 
725 
726 void Foam::polyMeshFilter::printScalarFieldStats
727 (
728  const string desc,
729  const scalarField& fld
730 ) const
731 {
732  scalar sum = 0;
733  scalar validElements = 0;
734  scalar min = great;
735  scalar max = -great;
736 
737  forAll(fld, i)
738  {
739  const scalar fldElement = fld[i];
740 
741  if (fldElement >= 0)
742  {
743  sum += fldElement;
744 
745  if (fldElement < min)
746  {
747  min = fldElement;
748  }
749 
750  if (fldElement > max)
751  {
752  max = fldElement;
753  }
754 
755  validElements++;
756  }
757  }
758 
759  reduce(sum, sumOp<scalar>());
760  reduce(min, minOp<scalar>());
761  reduce(max, maxOp<scalar>());
762  reduce(validElements, sumOp<label>());
763  const label totFieldSize = returnReduce(fld.size(), sumOp<label>());
764 
765  Info<< incrIndent << indent << desc
766  << ": min = " << min
767  << " av = " << sum/(validElements + small)
768  << " max = " << max << nl
769  << indent
770  << " " << validElements << " / " << totFieldSize << " elements used"
771  << decrIndent << endl;
772 }
773 
774 
775 void Foam::polyMeshFilter::mapOldMeshEdgeFieldToNewMesh
776 (
777  const polyMesh& newMesh,
778  const labelList& pointMap,
779  scalarField& newMeshMinEdgeLen
780 ) const
781 {
782  scalarField tmp(newMesh.nEdges());
783 
784  const edgeList& newEdges = newMesh.edges();
785 
786  forAll(newEdges, newEdgeI)
787  {
788  const edge& newEdge = newEdges[newEdgeI];
789  const label pStart = newEdge.start();
790  const label pEnd = newEdge.end();
791 
792  tmp[newEdgeI] = min
793  (
794  newMeshMinEdgeLen[pointMap[pStart]],
795  newMeshMinEdgeLen[pointMap[pEnd]]
796  );
797  }
798 
799  newMeshMinEdgeLen.transfer(tmp);
800 
802  (
803  newMesh,
804  newMeshMinEdgeLen,
805  maxEqOp<scalar>(),
806  scalar(0)
807  );
808 }
809 
810 
811 void Foam::polyMeshFilter::mapOldMeshFaceFieldToNewMesh
812 (
813  const polyMesh& newMesh,
814  const labelList& faceMap,
815  scalarField& newMeshFaceFilterFactor
816 ) const
817 {
818  scalarField tmp(newMesh.nFaces());
819 
820  forAll(faceMap, newFacei)
821  {
822  const label oldFacei = faceMap[newFacei];
823 
824  tmp[newFacei] = newMeshFaceFilterFactor[oldFacei];
825  }
826 
827  newMeshFaceFilterFactor.transfer(tmp);
828 
830  (
831  newMesh,
832  newMeshFaceFilterFactor,
833  maxEqOp<scalar>()
834  );
835 }
836 
837 
838 void Foam::polyMeshFilter::updateOldToNewPointMap
839 (
840  const labelList& currToNew,
841  labelList& origToCurrentPointMap
842 ) const
843 {
844  forAll(origToCurrentPointMap, origPointi)
845  {
846  label oldPointi = origToCurrentPointMap[origPointi];
847 
848  if (oldPointi != -1)
849  {
850  label newPointi = currToNew[oldPointi];
851 
852  if (newPointi >= 0)
853  {
854  origToCurrentPointMap[origPointi] = newPointi;
855  }
856  else if (newPointi == -1)
857  {
858  origToCurrentPointMap[origPointi] = -1;
859  }
860  else
861  {
862  origToCurrentPointMap[origPointi] = -newPointi-2;
863  }
864  }
865  }
866 }
867 
868 
869 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
870 
872 :
874  (
876  (
877  IOobject
878  (
879  "collapseDict",
880  mesh.time().system(),
881  mesh.time(),
882  IOobject::MUST_READ,
883  IOobject::NO_WRITE
884  )
885  )
886  ),
887  mesh_(mesh),
888  newMeshPtr_(),
889  originalPointPriority_(mesh.nPoints(), labelMin),
890  pointPriority_(),
891  minEdgeLen_(),
892  faceFilterFactor_()
893 {
895 }
896 
897 
899 (
900  const fvMesh& mesh,
901  const labelList& pointPriority
902 )
903 :
905  (
907  (
908  IOobject
909  (
910  "collapseDict",
911  mesh.time().system(),
912  mesh.time(),
913  IOobject::MUST_READ,
914  IOobject::NO_WRITE
915  )
916  )
917  ),
918  mesh_(mesh),
919  newMeshPtr_(),
920  originalPointPriority_(pointPriority),
921  pointPriority_(),
922  minEdgeLen_(),
923  faceFilterFactor_()
924 {
926 }
927 
929 (
930  const fvMesh& mesh,
931  const labelList& pointPriority,
932  const dictionary& dict
933 )
934 :
936  mesh_(mesh),
937  newMeshPtr_(),
938  originalPointPriority_(pointPriority),
939  pointPriority_(),
940  minEdgeLen_(),
941  faceFilterFactor_()
942 {
944 }
945 
946 
947 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
948 
950 {}
951 
952 
953 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
954 
956 {
957  minEdgeLen_.resize(mesh_.nEdges(), minLen());
958  faceFilterFactor_.resize(mesh_.nFaces(), initialFaceLengthFactor());
959 
960  return filterFacesLoop(nOriginalBadFaces);
961 }
962 
963 
965 {
966  minEdgeLen_.resize(mesh_.nEdges(), minLen());
967  faceFilterFactor_.resize(mesh_.nFaces(), initialFaceLengthFactor());
968 
969  forAll(faceFilterFactor_, fI)
970  {
971  if (!fSet.found(fI))
972  {
973  faceFilterFactor_[fI] = -1;
974  }
975  }
976 
977  return filterFacesLoop(0);
978 }
979 
980 
981 Foam::label Foam::polyMeshFilter::filterEdges
982 (
983  const label nOriginalBadFaces
984 )
985 {
986  label nBadFaces = labelMax/2;
987  label nPreviousBadFaces = labelMax;
988  label nOuterIterations = 0;
989 
990  minEdgeLen_.resize(mesh_.nEdges(), minLen());
991  faceFilterFactor_.resize(0);
992 
993  labelList pointErrorCount(mesh_.nPoints(), 0);
994 
995  // Main loop
996  // ~~~~~~~~~
997  // It tries and do some collapses, checks the resulting mesh and
998  // 'freezes' some edges (by marking in minEdgeLen) and tries again.
999  // This will iterate ultimately to the situation where every edge is
1000  // frozen and nothing gets collapsed.
1001  while
1002  (
1003  nOuterIterations < maxIterations()
1004  && nBadFaces > nOriginalBadFaces
1005  && nBadFaces < nPreviousBadFaces
1006  )
1007  {
1008  Info<< nl << "Outer Iteration = " << nOuterIterations++ << nl
1009  << endl;
1010 
1011  printScalarFieldStats("Edge Filter Factor", minEdgeLen_);
1012 
1013  nPreviousBadFaces = nBadFaces;
1014 
1015  // Reset the new mesh to the old mesh
1016  newMeshPtr_ = copyMesh(mesh_);
1017  fvMesh& newMesh = newMeshPtr_();
1018 
1019  scalarField newMeshMinEdgeLen = minEdgeLen_;
1020  pointPriority_.reset(new labelList(originalPointPriority_));
1021 
1022  labelList origToCurrentPointMap(identityMap(newMesh.nPoints()));
1023 
1024  Info<< incrIndent;
1025 
1026  label nInnerIterations = 0;
1027  label nPrevLocalCollapse = labelMax;
1028 
1029  while (true)
1030  {
1031  Info<< nl
1032  << indent << "Inner iteration = " << nInnerIterations++ << nl
1033  << incrIndent << endl;
1034 
1035  label nLocalCollapse = filterEdges
1036  (
1037  newMesh,
1038  newMeshMinEdgeLen,
1039  origToCurrentPointMap
1040  );
1041 
1042  Info<< decrIndent;
1043 
1044  if
1045  (
1046  nLocalCollapse >= nPrevLocalCollapse
1047  || nLocalCollapse == 0
1048  )
1049  {
1050  Info<< decrIndent;
1051  break;
1052  }
1053  else
1054  {
1055  nPrevLocalCollapse = nLocalCollapse;
1056  }
1057  }
1058 
1059  // Mesh check
1060  // ~~~~~~~~~~~~~~~~~~
1061  // Do not allow collapses in regions of error.
1062  // Updates minEdgeLen, nRelaxedEdges
1063 
1064  if (controlMeshQuality())
1065  {
1066  PackedBoolList isErrorPoint(newMesh.nPoints());
1068  (
1069  newMesh,
1070  meshQualityCoeffDict(),
1071  isErrorPoint
1072  );
1073 
1074  Info<< nl << " Number of bad faces : " << nBadFaces << nl
1075  << " Number of marked points : "
1076  << returnReduce(isErrorPoint.count(), sumOp<unsigned int>())
1077  << endl;
1078 
1079  updatePointErrorCount
1080  (
1081  isErrorPoint,
1082  origToCurrentPointMap,
1083  pointErrorCount
1084  );
1085 
1086  checkMeshEdgesAndRelaxEdges
1087  (
1088  newMesh,
1089  origToCurrentPointMap,
1090  isErrorPoint,
1091  pointErrorCount
1092  );
1093  }
1094  else
1095  {
1096  return -1;
1097  }
1098  }
1099 
1100  return nBadFaces;
1101 }
1102 
1103 
1105 {
1106  return newMeshPtr_;
1107 }
1108 
1109 
1112 {
1113  return pointPriority_;
1114 }
1115 
1116 
1117 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
void reset(const Field< Type > &)
Reset the field values to the given field.
Definition: Field.C:456
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:138
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
A bit-packed bool list.
unsigned int count() const
Count number of bits set, O(log(n))
Definition: PackedList.C:55
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
static label checkMeshQuality(const polyMesh &mesh, const dictionary &meshQualityDict, PackedBoolList &isErrorPoint)
Check mesh and mark points on faces in error.
Definition: edgeCollapser.C:82
A list of face labels.
Definition: faceSet.H:51
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
const word & name() const
Return reference to name.
Definition: fvMesh.H:432
Class to store the settings for the polyMeshFilter class.
void writeSettings(Ostream &os) const
Write the settings to a stream.
Remove the edges and faces of a polyMesh whilst satisfying the given mesh quality criteria.
polyMeshFilter(const fvMesh &mesh)
Construct from fvMesh.
const autoPtr< labelList > & pointPriority() const
Return the new pointPriority list.
~polyMeshFilter()
Destructor.
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:80
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Direct mesh changes based on v1.3 polyTopoChange syntax.
autoPtr< polyTopoChangeMap > makeMesh(autoPtr< fvMesh > &newMesh, const IOobject &io, const polyMesh &mesh, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Create new mesh with old mesh patches.
label nPoints() const
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronise values on all mesh points.
static void syncEdgeList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronise values on all mesh edges.
static void syncFaceList(const polyMesh &mesh, UList< T > &l, const CombineOp &cop)
Synchronise values on all mesh faces.
Definition: syncTools.H:387
label collapseEdge(triSurface &surf, const scalar minLen)
Keep collapsing all edges < minLen.
label nPoints
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Namespace for OpenFOAM.
int system(const std::string &command)
Execute the specified command.
Definition: POSIX.C:1230
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
const doubleScalar e
Definition: doubleScalar.H:106
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:241
List< label > labelList
A List of labels.
Definition: labelList.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
messageStream Info
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:234
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
List< edge > edgeList
Definition: edgeList.H:38
static const label labelMax
Definition: label.H:62
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
List< face > faceList
Definition: faceListFwd.H:41
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:227
static const char nl
Definition: Ostream.H:266
static const label labelMin
Definition: label.H:61
label newPointi
Definition: readKivaGrid.H:495
labelList f(nPoints)
dictionary dict