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