meshRefinementMerge.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-2014 OpenFOAM Foundation
9  Copyright (C) 2016-2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "combineFaces.H"
31 #include "polyTopoChange.H"
32 #include "removePoints.H"
33 #include "faceSet.H"
34 #include "Time.H"
35 #include "motionSmoother.H"
36 #include "syncTools.H"
37 
38 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
39 
40 // Merge faces that are in-line.
42 (
43  const scalar minCos,
44  const scalar concaveCos,
45  const label mergeSize,
46  const labelList& patchIDs,
47  const meshRefinement::FaceMergeType mergeType
48 )
49 {
50  // Patch face merging engine
51  combineFaces faceCombiner(mesh_, false);
52 
53  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
54 
55  // Pick up all candidate cells on boundary
56  labelHashSet boundaryCells(mesh_.nBoundaryFaces());
57 
58  for (const label patchi : patchIDs)
59  {
60  const polyPatch& patch = patches[patchi];
61 
62  if (!patch.coupled())
63  {
64  forAll(patch, i)
65  {
66  boundaryCells.insert(mesh_.faceOwner()[patch.start()+i]);
67  }
68  }
69  }
70 
71  // Get all sets of faces that can be merged
72  labelListList mergeSets
73  (
74  faceCombiner.getMergeSets
75  (
76  minCos,
77  concaveCos,
78  boundaryCells,
79  (mergeType == FaceMergeType::IGNOREPATCH) // merge across patches?
80  )
81  );
82 
83  if (mergeSize != -1)
84  {
85  // Keep only those that are mergeSize faces
86  label compactI = 0;
87  forAll(mergeSets, setI)
88  {
89  if (mergeSets[setI].size() == mergeSize && compactI != setI)
90  {
91  mergeSets[compactI++] = mergeSets[setI];
92  }
93  }
94  mergeSets.setSize(compactI);
95  }
96 
97 
98  label nFaceSets = returnReduce(mergeSets.size(), sumOp<label>());
99 
100  Info<< "Merging " << nFaceSets << " sets of faces." << nl << endl;
101 
102  if (nFaceSets > 0)
103  {
104  // Topology changes container
105  polyTopoChange meshMod(mesh_);
106 
107  // Merge all faces of a set into the first face of the set. Remove
108  // unused points.
109  faceCombiner.setRefinement(mergeSets, meshMod);
110 
111  // Change the mesh (no inflation)
112  autoPtr<mapPolyMesh> mapPtr = meshMod.changeMesh(mesh_, false, true);
113  mapPolyMesh& map = *mapPtr;
114 
115  // Update fields
116  mesh_.updateMesh(map);
117 
118  // Move mesh (since morphing does not do this)
119  if (map.hasMotionPoints())
120  {
121  mesh_.movePoints(map.preMotionPoints());
122  }
123  else
124  {
125  // Delete mesh volumes. No other way to do this?
126  mesh_.clearOut();
127  }
128 
129  // Reset the instance for if in overwrite mode
130  mesh_.setInstance(timeName());
131 
132  faceCombiner.updateMesh(map);
133 
134  // Get the kept faces that need to be recalculated.
135  // Merging two boundary faces might shift the cell centre
136  // (unless the faces are absolutely planar)
137  labelHashSet retestFaces(2*mergeSets.size());
138 
139  forAll(mergeSets, setI)
140  {
141  label oldMasterI = mergeSets[setI][0];
142  retestFaces.insert(map.reverseFaceMap()[oldMasterI]);
143  }
144  updateMesh(map, growFaceCellFace(retestFaces));
145  }
146 
147  return nFaceSets;
148 }
149 
150 
153 //Foam::autoPtr<Foam::mapPolyMesh> Foam::meshRefinement::mergeEdges
154 //(
155 // const scalar minCos
156 //)
157 //{
158 // // Point removal analysis engine
159 // removePoints pointRemover(mesh_);
160 //
161 // // Count usage of points
162 // boolList pointCanBeDeleted;
163 // label nRemove = pointRemover.countPointUsage(minCos, pointCanBeDeleted);
164 //
165 // Info<< "Removing " << nRemove
166 // << " straight edge points." << endl;
167 //
168 // autoPtr<mapPolyMesh> map;
169 //
170 // if (nRemove > 0)
171 // {
172 // // Save my local faces that will change. These changed faces might
173 // // cause a shift in the cell centre which needs to be retested.
174 // // Have to do this before changing mesh since point will be removed.
175 // labelHashSet retestOldFaces(nRemove / Pstream::nProcs());
176 //
177 // {
178 // const faceList& faces = mesh_.faces();
179 //
180 // forAll(faces, faceI)
181 // {
182 // const face& f = faces[faceI];
183 //
184 // forAll(f, fp)
185 // {
186 // if (pointCanBeDeleted[f[fp]])
187 // {
188 // retestOldFaces.insert(faceI);
189 // break;
190 // }
191 // }
192 // }
193 // }
194 //
195 // // Topology changes container
196 // polyTopoChange meshMod(mesh_);
197 //
198 // pointRemover.setRefinement(pointCanBeDeleted, meshMod);
199 //
200 // // Change the mesh (no inflation)
201 // map = meshMod.changeMesh(mesh_, false, true);
202 //
203 // // Update fields
204 // mesh_.updateMesh(map);
205 //
206 // // Move mesh (since morphing does not do this)
207 // if (map().hasMotionPoints())
208 // {
209 // mesh_.movePoints(map().preMotionPoints());
210 // }
211 // else
212 // {
213 // // Delete mesh volumes. No other way to do this?
214 // mesh_.clearOut();
215 // }
216 //
217 // // Reset the instance for if in overwrite mode
218 // mesh_.setInstance(timeName());
219 //
220 // pointRemover.updateMesh(map);
221 //
222 // // Get the kept faces that need to be recalculated.
223 // labelHashSet retestFaces(6*retestOldFaces.size());
224 //
225 // const cellList& cells = mesh_.cells();
226 //
227 // for (const label oldFacei : retestOldFaces)
228 // {
229 // const label facei = map().reverseFaceMap()[oldFacei];
230 //
231 // const cell& ownFaces = cells[mesh_.faceOwner()[facei]];
232 //
233 // retestFaces.insert(ownFaces);
234 //
235 // if (mesh_.isInternalFace(facei))
236 // {
237 // const cell& neiFaces = cells[mesh_.faceNeighbour()[facei]];
238 //
239 // retestFaces.insert(neiFaces);
240 // }
241 // }
242 // updateMesh(map, retestFaces.toc());
243 // }
244 //
245 // return map;
246 //}
247 
248 
250 (
251  const scalar minCos,
252  const scalar concaveCos,
253  const labelList& patchIDs,
254  const dictionary& motionDict,
255  const labelList& preserveFaces,
256  const meshRefinement::FaceMergeType mergeType
257 )
258 {
259  // Patch face merging engine
260  combineFaces faceCombiner(mesh_, true);
261 
262  // Pick up all candidate cells on boundary
263  labelHashSet boundaryCells(mesh_.nBoundaryFaces());
264 
265  {
266  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
267 
268  for (const label patchi : patchIDs)
269  {
270  const polyPatch& patch = patches[patchi];
271 
272  if (!patch.coupled())
273  {
274  forAll(patch, i)
275  {
276  boundaryCells.insert(mesh_.faceOwner()[patch.start()+i]);
277  }
278  }
279  }
280  }
281 
282  // Get all sets of faces that can be merged. Since only faces on the same
283  // patch get merged there is no risk of e.g. patchID faces getting merged
284  // with original patches (or even processor patches). There is a risk
285  // though of original patch faces getting merged if they make a small
286  // angle. Would be pretty weird starting mesh though.
287  labelListList allFaceSets
288  (
289  faceCombiner.getMergeSets
290  (
291  minCos,
292  concaveCos,
293  boundaryCells,
294  (mergeType == FaceMergeType::IGNOREPATCH) // merge across patches?
295  )
296  );
297 
298  // Filter out any set that contains any preserveFace
299  label compactI = 0;
300  forAll(allFaceSets, i)
301  {
302  const labelList& set = allFaceSets[i];
303 
304  bool keep = true;
305  forAll(set, j)
306  {
307  if (preserveFaces[set[j]] != -1)
308  {
309  keep = false;
310  break;
311  }
312  }
313 
314  if (keep)
315  {
316  if (compactI != i)
317  {
318  allFaceSets[compactI] = set;
319  }
320 
321  compactI++;
322  }
323  }
324  allFaceSets.setSize(compactI);
325 
326  label nFaceSets = returnReduce(allFaceSets.size(), sumOp<label>());
327 
328  Info<< "Merging " << nFaceSets << " sets of faces." << nl << endl;
329 
330  if (nFaceSets > 0)
331  {
332  if (debug&meshRefinement::MESH)
333  {
334  faceSet allSets(mesh_, "allFaceSets", allFaceSets.size());
335  forAll(allFaceSets, setI)
336  {
337  allSets.insert(allFaceSets[setI]);
338  }
339  Pout<< "Writing all faces to be merged to set "
340  << allSets.objectPath() << endl;
341  allSets.instance() = timeName();
342  allSets.write();
343 
344  const_cast<Time&>(mesh_.time())++;
345  }
346 
347 
348  // Topology changes container
349  polyTopoChange meshMod(mesh_);
350 
351  // Merge all faces of a set into the first face of the set.
352  faceCombiner.setRefinement(allFaceSets, meshMod);
353 
354  // Experimental: store data for all the points that have been deleted
355  storeData
356  (
357  faceCombiner.savedPointLabels(), // points to store
358  labelList(0), // faces to store
359  labelList(0) // cells to store
360  );
361 
362  // Change the mesh (no inflation)
363  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh_, false, true);
364 
365  // Update fields
366  mesh_.updateMesh(map());
367 
368  // Move mesh (since morphing does not do this)
369  if (map().hasMotionPoints())
370  {
371  mesh_.movePoints(map().preMotionPoints());
372  }
373  else
374  {
375  // Delete mesh volumes.
376  mesh_.clearOut();
377  }
378 
379  // Reset the instance for if in overwrite mode
380  mesh_.setInstance(timeName());
381 
382  faceCombiner.updateMesh(map());
383 
384  // Get the kept faces that need to be recalculated.
385  // Merging two boundary faces might shift the cell centre
386  // (unless the faces are absolutely planar)
387  labelHashSet retestFaces(2*allFaceSets.size());
388 
389  forAll(allFaceSets, setI)
390  {
391  label oldMasterI = allFaceSets[setI][0];
392  retestFaces.insert(map().reverseFaceMap()[oldMasterI]);
393  }
394  updateMesh(map(), growFaceCellFace(retestFaces));
395 
396  if (debug&meshRefinement::MESH)
397  {
398  // Check sync
399  Pout<< "Checking sync after initial merging " << nFaceSets
400  << " faces." << endl;
401  checkData();
402 
403  Pout<< "Writing initial merged-faces mesh to time "
404  << timeName() << nl << endl;
405  write();
406  }
407 
408  for (label iteration = 0; iteration < 100; iteration++)
409  {
410  Info<< nl
411  << "Undo iteration " << iteration << nl
412  << "----------------" << endl;
413 
414 
415  // Check mesh for errors
416  // ~~~~~~~~~~~~~~~~~~~~~
417 
418  faceSet errorFaces(mesh_, "errorFaces", mesh_.nBoundaryFaces());
419  bool hasErrors = motionSmoother::checkMesh
420  (
421  false, // report
422  mesh_,
423  motionDict,
424  errorFaces,
425  dryRun_
426  );
427 
428  //if (checkEdgeConnectivity)
429  //{
430  // Info<< "Checking edge-face connectivity (duplicate faces"
431  // << " or non-consecutive shared vertices)" << endl;
432  //
433  // label nOldSize = errorFaces.size();
434  //
435  // hasErrors =
436  // mesh_.checkFaceFaces
437  // (
438  // false,
439  // &errorFaces
440  // )
441  // || hasErrors;
442  //
443  // Info<< "Detected additional "
444  // << returnReduce
445  // (
446  // errorFaces.size() - nOldSize,
447  // sumOp<label>()
448  // )
449  // << " faces with illegal face-face connectivity"
450  // << endl;
451  //}
452 
453  if (!hasErrors)
454  {
455  break;
456  }
457 
458 
459  if (debug&meshRefinement::MESH)
460  {
461  errorFaces.instance() = timeName();
462  Pout<< "Writing all faces in error to faceSet "
463  << errorFaces.objectPath() << nl << endl;
464  errorFaces.write();
465  }
466 
467 
468  // Check any master cells for using any of the error faces
469  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
470 
471  DynamicList<label> mastersToRestore(allFaceSets.size());
472 
473  forAll(allFaceSets, setI)
474  {
475  label masterFaceI = faceCombiner.masterFace()[setI];
476 
477  if (masterFaceI != -1)
478  {
479  label masterCellII = mesh_.faceOwner()[masterFaceI];
480 
481  const cell& cFaces = mesh_.cells()[masterCellII];
482 
483  forAll(cFaces, i)
484  {
485  if (errorFaces.found(cFaces[i]))
486  {
487  mastersToRestore.append(masterFaceI);
488  break;
489  }
490  }
491  }
492  }
493  mastersToRestore.shrink();
494 
495  label nRestore = returnReduce
496  (
497  mastersToRestore.size(),
498  sumOp<label>()
499  );
500 
501  Info<< "Masters that need to be restored:"
502  << nRestore << endl;
503 
504  if (debug&meshRefinement::MESH)
505  {
506  faceSet restoreSet(mesh_, "mastersToRestore", mastersToRestore);
507  restoreSet.instance() = timeName();
508  Pout<< "Writing all " << mastersToRestore.size()
509  << " masterfaces to be restored to set "
510  << restoreSet.objectPath() << endl;
511  restoreSet.write();
512  }
513 
514 
515  if (nRestore == 0)
516  {
517  break;
518  }
519 
520 
521  // Undo
522  // ~~~~
523 
524  if (debug)
525  {
526  const_cast<Time&>(mesh_.time())++;
527  }
528 
529  // Topology changes container
530  polyTopoChange meshMod(mesh_);
531 
532  // Merge all faces of a set into the first face of the set.
533  // Experimental:mark all points/faces/cells that have been restored.
534  Map<label> restoredPoints(0);
535  Map<label> restoredFaces(0);
536  Map<label> restoredCells(0);
537 
538  faceCombiner.setUnrefinement
539  (
540  mastersToRestore,
541  meshMod,
542  restoredPoints,
543  restoredFaces,
544  restoredCells
545  );
546 
547  // Change the mesh (no inflation)
548  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh_, false, true);
549 
550  // Update fields
551  mesh_.updateMesh(map());
552 
553  // Move mesh (since morphing does not do this)
554  if (map().hasMotionPoints())
555  {
556  mesh_.movePoints(map().preMotionPoints());
557  }
558  else
559  {
560  // Delete mesh volumes.
561  mesh_.clearOut();
562  }
563 
564 
565  // Reset the instance for if in overwrite mode
566  mesh_.setInstance(timeName());
567 
568  faceCombiner.updateMesh(map());
569 
570  // Renumber restore maps
571  inplaceMapKey(map().reversePointMap(), restoredPoints);
572  inplaceMapKey(map().reverseFaceMap(), restoredFaces);
573  inplaceMapKey(map().reverseCellMap(), restoredCells);
574 
575 
576  // Get the kept faces that need to be recalculated.
577  // Merging two boundary faces might shift the cell centre
578  // (unless the faces are absolutely planar)
579  labelHashSet retestFaces(2*restoredFaces.size());
580 
581  forAllConstIters(restoredFaces, iter)
582  {
583  retestFaces.insert(iter.key());
584  }
585 
586  // Experimental:restore all points/face/cells in maps
587  updateMesh
588  (
589  map(),
590  growFaceCellFace(retestFaces),
591  restoredPoints,
592  restoredFaces,
593  restoredCells
594  );
595 
596  if (debug&meshRefinement::MESH)
597  {
598  // Check sync
599  Pout<< "Checking sync after restoring " << retestFaces.size()
600  << " faces." << endl;
601  checkData();
602 
603  Pout<< "Writing merged-faces mesh to time "
604  << timeName() << nl << endl;
605  write();
606  }
607 
608  Info<< endl;
609  }
610  }
611  else
612  {
613  Info<< "No faces merged ..." << endl;
614  }
615 
616  return nFaceSets;
617 }
618 
619 
620 // Remove points. pointCanBeDeleted is parallel synchronised.
622 (
623  removePoints& pointRemover,
624  const boolList& pointCanBeDeleted
625 )
626 {
627  // Topology changes container
628  polyTopoChange meshMod(mesh_);
629 
630  pointRemover.setRefinement(pointCanBeDeleted, meshMod);
631 
632  // Change the mesh (no inflation)
633  autoPtr<mapPolyMesh> mapPtr = meshMod.changeMesh(mesh_, false, true);
634  mapPolyMesh& map = *mapPtr;
635 
636  // Update fields
637  mesh_.updateMesh(map);
638 
639  // Move mesh (since morphing does not do this)
640  if (map.hasMotionPoints())
641  {
642  mesh_.movePoints(map.preMotionPoints());
643  }
644  else
645  {
646  // Delete mesh volumes.
647  mesh_.clearOut();
648  }
649 
650  // Reset the instance for if in overwrite mode
651  mesh_.setInstance(timeName());
652 
653  pointRemover.updateMesh(map);
654 
655 
656  // Retest all affected faces and all the cells using them
657  labelHashSet retestFaces(pointRemover.savedFaceLabels().size());
658  for (const label facei : pointRemover.savedFaceLabels())
659  {
660  if (facei >= 0)
661  {
662  retestFaces.insert(facei);
663  }
664  }
665  updateMesh(map, growFaceCellFace(retestFaces));
666 
667  if (debug)
668  {
669  // Check sync
670  Pout<< "Checking sync after removing points." << endl;
671  checkData();
672  }
673 
674  return mapPtr;
675 }
676 
677 
678 // Restore faces (which contain removed points)
680 (
681  removePoints& pointRemover,
682  const labelList& facesToRestore
683 )
684 {
685  // Topology changes container
686  polyTopoChange meshMod(mesh_);
687 
688  // Determine sets of points and faces to restore
689  labelList localFaces, localPoints;
690  pointRemover.getUnrefimentSet
691  (
692  facesToRestore,
693  localFaces,
694  localPoints
695  );
696 
697  // Undo the changes on the faces that are in error.
698  pointRemover.setUnrefinement
699  (
700  localFaces,
701  localPoints,
702  meshMod
703  );
704 
705  // Change the mesh (no inflation)
706  autoPtr<mapPolyMesh> mapPtr = meshMod.changeMesh(mesh_, false, true);
707  mapPolyMesh& map = *mapPtr;
708 
709  // Update fields
710  mesh_.updateMesh(map);
711 
712  // Move mesh (since morphing does not do this)
713  if (map.hasMotionPoints())
714  {
715  mesh_.movePoints(map.preMotionPoints());
716  }
717  else
718  {
719  // Delete mesh volumes.
720  mesh_.clearOut();
721  }
722 
723  // Reset the instance for if in overwrite mode
724  mesh_.setInstance(timeName());
725 
726  pointRemover.updateMesh(map);
727 
728  labelHashSet retestFaces(2*facesToRestore.size());
729  forAll(facesToRestore, i)
730  {
731  label faceI = map.reverseFaceMap()[facesToRestore[i]];
732  if (faceI >= 0)
733  {
734  retestFaces.insert(faceI);
735  }
736  }
737  updateMesh(map, growFaceCellFace(retestFaces));
738 
739  if (debug)
740  {
741  // Check sync
742  Pout<< "Checking sync after restoring points on "
743  << facesToRestore.size() << " faces." << endl;
744  checkData();
745  }
746 
747  return mapPtr;
748 }
749 
750 
751 // Collect all faces that are both in candidateFaces and in the set.
752 // If coupled face also collects the coupled face.
754 (
755  const labelList& candidateFaces,
756  const labelHashSet& set
757 ) const
758 {
759  // Has face been selected?
760  boolList selected(mesh_.nFaces(), false);
761 
762  for (const label facei : candidateFaces)
763  {
764  if (set.found(facei))
765  {
766  selected[facei] = true;
767  }
768  }
769  syncTools::syncFaceList
770  (
771  mesh_,
772  selected,
773  orEqOp<bool>() // combine operator
774  );
775 
776  labelList selectedFaces(findIndices(selected, true));
777 
778  return selectedFaces;
779 }
780 
781 
782 namespace Foam
783 {
784  // Pick up faces of cells of faces in set.
785  // file-scope
786  static inline void markGrowFaceCellFace
787  (
788  const polyMesh& pMesh,
789  const label faceI,
790  boolList& selected
791  )
792  {
793  const label own = pMesh.faceOwner()[faceI];
794 
795  const cell& ownFaces = pMesh.cells()[own];
796  for (const label facei : ownFaces)
797  {
798  selected[facei] = true;
799  }
800 
801  if (pMesh.isInternalFace(faceI))
802  {
803  const label nbr = pMesh.faceNeighbour()[faceI];
804 
805  const cell& nbrFaces = pMesh.cells()[nbr];
806  forAll(nbrFaces, nbrFaceI)
807  {
808  selected[nbrFaces[nbrFaceI]] = true;
809  }
810  }
811  }
812 }
813 
814 
815 // Pick up faces of cells of faces in set.
817 (
818  const labelUList& set
819 ) const
820 {
821  boolList selected(mesh_.nFaces(), false);
822 
823  for (const label facei : set)
824  {
825  markGrowFaceCellFace(mesh_, facei, selected);
826  }
827 
828  syncTools::syncFaceList
829  (
830  mesh_,
831  selected,
832  orEqOp<bool>() // combine operator
833  );
834  return findIndices(selected, true);
835 }
836 
837 
838 // Pick up faces of cells of faces in set.
840 (
841  const labelHashSet& set
842 ) const
843 {
844  boolList selected(mesh_.nFaces(), false);
845 
846  for (const label facei : set)
847  {
848  markGrowFaceCellFace(mesh_, facei, selected);
849  }
850 
851  syncTools::syncFaceList
852  (
853  mesh_,
854  selected,
855  orEqOp<bool>() // combine operator
856  );
857  return findIndices(selected, true);
858 }
859 
860 
861 // Remove points not used by any face or points used by only two faces where
862 // the edges are in line
864 (
865  const scalar minCos,
866  const dictionary& motionDict
867 )
868 {
869  Info<< nl
870  << "Merging all points on surface that" << nl
871  << "- are used by only two boundary faces and" << nl
872  << "- make an angle with a cosine of more than " << minCos
873  << "." << nl << endl;
874 
875  // Point removal analysis engine with undo
876  removePoints pointRemover(mesh_, true);
877 
878  // Count usage of points
879  boolList pointCanBeDeleted;
880  label nRemove = pointRemover.countPointUsage(minCos, pointCanBeDeleted);
881 
882  if (nRemove > 0)
883  {
884  Info<< "Removing " << nRemove
885  << " straight edge points ..." << nl << endl;
886 
887  // Remove points
888  // ~~~~~~~~~~~~~
889 
890  doRemovePoints(pointRemover, pointCanBeDeleted);
891 
892 
893  for (label iteration = 0; iteration < 100; iteration++)
894  {
895  Info<< nl
896  << "Undo iteration " << iteration << nl
897  << "----------------" << endl;
898 
899 
900  // Check mesh for errors
901  // ~~~~~~~~~~~~~~~~~~~~~
902 
903  faceSet errorFaces(mesh_, "errorFaces", mesh_.nBoundaryFaces());
904  bool hasErrors = motionSmoother::checkMesh
905  (
906  false, // report
907  mesh_,
908  motionDict,
909  errorFaces,
910  dryRun_
911  );
912  //if (checkEdgeConnectivity)
913  //{
914  // Info<< "Checking edge-face connectivity (duplicate faces"
915  // << " or non-consecutive shared vertices)" << endl;
916  //
917  // label nOldSize = errorFaces.size();
918  //
919  // hasErrors =
920  // mesh_.checkFaceFaces
921  // (
922  // false,
923  // &errorFaces
924  // )
925  // || hasErrors;
926  //
927  // Info<< "Detected additional "
928  // << returnReduce(errorFaces.size()-nOldSize,sumOp<label>())
929  // << " faces with illegal face-face connectivity"
930  // << endl;
931  //}
932 
933  if (!hasErrors)
934  {
935  break;
936  }
937 
938  if (debug&meshRefinement::MESH)
939  {
940  errorFaces.instance() = timeName();
941  Pout<< "**Writing all faces in error to faceSet "
942  << errorFaces.objectPath() << nl << endl;
943  errorFaces.write();
944  }
945 
946  labelList masterErrorFaces
947  (
948  collectFaces
949  (
950  pointRemover.savedFaceLabels(),
951  errorFaces
952  )
953  );
954 
955  label n = returnReduce(masterErrorFaces.size(), sumOp<label>());
956 
957  Info<< "Detected " << n
958  << " error faces on boundaries that have been merged."
959  << " These will be restored to their original faces." << nl
960  << endl;
961 
962  if (n == 0)
963  {
964  if (hasErrors)
965  {
966  Info<< "Detected "
967  << returnReduce(errorFaces.size(), sumOp<label>())
968  << " error faces in mesh."
969  << " Restoring neighbours of faces in error." << nl
970  << endl;
971 
972  labelList expandedErrorFaces
973  (
974  growFaceCellFace
975  (
976  errorFaces
977  )
978  );
979 
980  doRestorePoints(pointRemover, expandedErrorFaces);
981  }
982 
983  break;
984  }
985 
986  doRestorePoints(pointRemover, masterErrorFaces);
987  }
988 
989  if (debug&meshRefinement::MESH)
990  {
991  const_cast<Time&>(mesh_.time())++;
992  Pout<< "Writing merged-edges mesh to time "
993  << timeName() << nl << endl;
994  write();
995  }
996  }
997  else
998  {
999  Info<< "No straight edges simplified and no points removed ..." << endl;
1000  }
1001 
1002  return nRemove;
1003 }
1004 
1005 
1006 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::meshRefinement::mergePatchFacesUndo
label mergePatchFacesUndo(const scalar minCos, const scalar concaveCos, const labelList &patchIDs, const dictionary &motionDict, const labelList &preserveFaces, const meshRefinement::FaceMergeType mergeType)
Merge coplanar faces. preserveFaces is != -1 for faces.
Definition: meshRefinementMerge.C:250
Foam::meshRefinement::growFaceCellFace
labelList growFaceCellFace(const labelUList &set) const
Definition: meshRefinementMerge.C:817
Foam::BitOps::set
void set(List< bool > &bools, const labelRange &range)
Set the specified range 'on' in a boolList.
Definition: BitOps.C:37
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::combineFaces
Combines boundary faces into single face. The faces get the patch of the first face ('the master')
Definition: combineFaces.H:58
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:63
Foam::combineFaces::savedPointLabels
const labelList & savedPointLabels() const
If undoable: set of original point labels of stored points.
Definition: combineFaces.H:148
Foam::IOobject::instance
const fileName & instance() const noexcept
Definition: IOobjectI.H:196
Foam::removePoints
Removes selected points from mesh and updates faces using these points.
Definition: removePoints.H:60
Foam::DynamicList< label >
Foam::combineFaces::getMergeSets
labelListList getMergeSets(const scalar featureCos, const scalar minConcaveCos, const labelHashSet &boundaryCells, const bool mergeAcrossPatches=false) const
Extract lists of all (non-coupled) boundary faces on selected.
Definition: combineFaces.C:305
Foam::checkData
label checkData(const fvMesh &mesh, const instantList &timeDirs, wordList &objectNames)
Check if fields are good to use (available at all times)
Foam::meshRefinement::collectFaces
labelList collectFaces(const labelList &candidateFaces, const labelHashSet &set) const
Definition: meshRefinementMerge.C:754
motionSmoother.H
polyTopoChange.H
combineFaces.H
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:138
Foam::polyTopoChange
Direct mesh changes based on v1.3 polyTopoChange syntax.
Definition: polyTopoChange.H:99
Foam::removePoints::getUnrefimentSet
void getUnrefimentSet(const labelList &undoFaces, labelList &localFaces, labelList &localPoints) const
Given set of faces to restore calculates a consistent set of.
Definition: removePoints.C:553
Foam::Map< label >
Foam::faceSet
A list of face labels.
Definition: faceSet.H:51
Foam::orEqOp
Definition: ops.H:86
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::mapPolyMesh::preMotionPoints
const pointField & preMotionPoints() const
Pre-motion point positions.
Definition: mapPolyMesh.H:613
Foam::meshRefinement::mergePatchFaces
label mergePatchFaces(const scalar minCos, const scalar concaveCos, const label mergeSize, const labelList &patchIDs, const meshRefinement::FaceMergeType mergeType)
Merge coplanar faces if sets are of size mergeSize.
Definition: meshRefinementMerge.C:42
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::HashSet< label, Hash< label > >
syncTools.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::meshRefinement::FaceMergeType
FaceMergeType
Enumeration for what to do with co-planar patch faces on a single.
Definition: meshRefinement.H:133
n
label n
Definition: TABSMDCalcMethod2.H:31
removePoints.H
Foam::inplaceMapKey
void inplaceMapKey(const labelUList &oldToNew, Container &input)
Rewrite with mapped keys. Ignore elements with negative key.
Definition: ListOpsTemplates.C:240
Foam::regIOobject::write
virtual bool write(const bool valid=true) const
Write using setting from DB.
Definition: regIOobjectWrite.C:132
Foam::markGrowFaceCellFace
static void markGrowFaceCellFace(const polyMesh &pMesh, const label faceI, boolList &selected)
Definition: meshRefinementMerge.C:787
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::combineFaces::updateMesh
void updateMesh(const mapPolyMesh &)
Force recalculation of locally stored data on topological change.
Definition: combineFaces.C:813
Foam::combineFaces::setUnrefinement
void setUnrefinement(const labelList &masterFaces, polyTopoChange &meshMod, Map< label > &restoredPoints, Map< label > &restoredFaces, Map< label > &restoredCells)
Play commands into polyTopoChange to reinsert original faces.
Definition: combineFaces.C:862
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1107
faceSet.H
Foam::meshRefinement::doRemovePoints
autoPtr< mapPolyMesh > doRemovePoints(removePoints &pointRemover, const boolList &pointCanBeDeleted)
Definition: meshRefinementMerge.C:622
Foam::removePoints::savedFaceLabels
const labelList & savedFaceLabels() const
Definition: removePoints.H:115
Foam::meshRefinement::mergeEdgesUndo
label mergeEdgesUndo(const scalar minCos, const dictionary &motionDict)
Merge edges, maintain mesh quality. Return global number.
Definition: meshRefinementMerge.C:864
Foam::polyTopoChange::changeMesh
autoPtr< mapPolyMesh > changeMesh(polyMesh &mesh, const labelUList &patchMap, const bool inflate, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Inplace changes mesh without change of patches.
Definition: polyTopoChange.C:2980
Foam::removePoints::setUnrefinement
void setUnrefinement(const labelList &localFaces, const labelList &localPoints, polyTopoChange &)
Restore selected faces and vertices.
Definition: removePoints.C:765
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
meshRefinement.H
Foam::Ostream::write
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::combineFaces::masterFace
const labelList & masterFace() const
If undoable: masterface for every set.
Definition: combineFaces.H:142
Foam::mapPolyMesh::reverseFaceMap
const labelList & reverseFaceMap() const
Reverse face map.
Definition: mapPolyMesh.H:501
Time.H
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::primitiveMesh::isInternalFace
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
Definition: primitiveMeshI.H:103
Foam::nl
constexpr char nl
Definition: Ostream.H:404
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< label >
Foam::UList< label >
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:36
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::mapPolyMesh::hasMotionPoints
bool hasMotionPoints() const
Has valid preMotionPoints?
Definition: mapPolyMesh.H:619
Foam::findIndices
labelList findIndices(const ListType &input, typename ListType::const_reference val, label start=0)
Linear search to find all occurrences of given element.
Foam::combineFaces::setRefinement
void setRefinement(const labelListList &, polyTopoChange &)
Play commands into polyTopoChange to combine faces. Gets.
Definition: combineFaces.C:588
Foam::meshRefinement::doRestorePoints
autoPtr< mapPolyMesh > doRestorePoints(removePoints &pointRemover, const labelList &facesToRestore)
Definition: meshRefinementMerge.C:680
Foam::IOobject::objectPath
fileName objectPath() const
The complete path + object name.
Definition: IOobjectI.H:214
Foam::removePoints::updateMesh
void updateMesh(const mapPolyMesh &)
Force recalculation of locally stored data on topological change.
Definition: removePoints.C:451
Foam::removePoints::setRefinement
void setRefinement(const boolList &, polyTopoChange &)
Play commands into polyTopoChange to remove points. Gets.
Definition: removePoints.C:296
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:54
Foam::removePoints::countPointUsage
label countPointUsage(const scalar minCos, boolList &pointCanBeDeleted) const
Mark in pointCanBeDeleted the points that can be deleted.
Definition: removePoints.C:148
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1113