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