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