removePoints.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-2013 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 "BiIndirectList.H"
27 #include "removePoints.H"
28 #include "PstreamReduceOps.H"
29 #include "polyMesh.H"
30 #include "polyTopoChange.H"
31 #include "polyRemovePoint.H"
32 #include "polyAddPoint.H"
33 #include "polyModifyFace.H"
34 #include "syncTools.H"
35 #include "faceSet.H"
36 #include "dummyTransform.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 
43 defineTypeNameAndDebug(removePoints, 0);
44 
45 //- Combine-reduce operator to combine data on faces. Takes care
46 // of reverse orientation on coupled face.
47 template<class T, template<class> class CombineOp>
48 class faceEqOp
49 {
50 
51 public:
52 
53  void operator()(List<T>& x, const List<T>& y) const
54  {
55  if (y.size() > 0)
56  {
57  if (x.size() == 0)
58  {
59  x = y;
60  }
61  else
62  {
63  label j = 0;
64  forAll(x, i)
65  {
66  CombineOp<T>()(x[i], y[j]);
67  j = y.rcIndex(j);
68  }
69  }
70  }
71  }
72 };
73 
74 
76 //template<class T>
77 //class dummyTransformList
78 //{
79 //public:
80 // void operator()(const coupledPolyPatch&, Field<List<T> >&) const
81 // {}
82 //};
84 //template<>
85 //class pTraits<boolList>
86 //{
87 //public:
88 //
89 // //- Component type
90 // typedef label cmptType;
91 //};
92 
93 
94 }
95 
96 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
97 
98 // Change the vertices of the face whilst keeping everything else the same.
99 void Foam::removePoints::modifyFace
100 (
101  const label faceI,
102  const face& newFace,
103  polyTopoChange& meshMod
104 ) const
105 {
106  // Get other face data.
107  label patchI = -1;
108  label owner = mesh_.faceOwner()[faceI];
109  label neighbour = -1;
110 
111  if (mesh_.isInternalFace(faceI))
112  {
113  neighbour = mesh_.faceNeighbour()[faceI];
114  }
115  else
116  {
117  patchI = mesh_.boundaryMesh().whichPatch(faceI);
118  }
119 
120  label zoneID = mesh_.faceZones().whichZone(faceI);
121 
122  bool zoneFlip = false;
123 
124  if (zoneID >= 0)
125  {
126  const faceZone& fZone = mesh_.faceZones()[zoneID];
127 
128  zoneFlip = fZone.flipMap()[fZone.whichFace(faceI)];
129  }
130 
131  meshMod.setAction
132  (
134  (
135  newFace, // modified face
136  faceI, // label of face being modified
137  owner, // owner
138  neighbour, // neighbour
139  false, // face flip
140  patchI, // patch for face
141  false, // remove from zone
142  zoneID, // zone for face
143  zoneFlip // face flip in zone
144  )
145  );
146 }
147 
148 
149 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
150 
151 // Construct from mesh
152 Foam::removePoints::removePoints
153 (
154  const polyMesh& mesh,
155  const bool undoable
156 )
157 :
158  mesh_(mesh),
159  undoable_(undoable)
160 {}
161 
162 
163 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
164 
166 (
167  const scalar minCos,
168  boolList& pointCanBeDeleted
169 ) const
170 {
171  // Containers to store two edges per point:
172  // -1 : not filled
173  // >= 0 : edge label
174  // -2 : more than two edges for point
175  labelList edge0(mesh_.nPoints(), -1);
176  labelList edge1(mesh_.nPoints(), -1);
177 
178  const edgeList& edges = mesh_.edges();
179 
180  forAll(edges, edgeI)
181  {
182  const edge& e = edges[edgeI];
183 
184  forAll(e, eI)
185  {
186  label pointI = e[eI];
187 
188  if (edge0[pointI] == -2)
189  {
190  // Already too many edges
191  }
192  else if (edge0[pointI] == -1)
193  {
194  // Store first edge using point
195  edge0[pointI] = edgeI;
196  }
197  else
198  {
199  // Already one edge using point. Check second container.
200 
201  if (edge1[pointI] == -1)
202  {
203  // Store second edge using point
204  edge1[pointI] = edgeI;
205  }
206  else
207  {
208  // Third edge using point. Mark.
209  edge0[pointI] = -2;
210  edge1[pointI] = -2;
211  }
212  }
213  }
214  }
215 
216 
217  // Check the ones used by only 2 edges that these are sufficiently in line.
218  const pointField& points = mesh_.points();
219 
220  pointCanBeDeleted.setSize(mesh_.nPoints());
221  pointCanBeDeleted = false;
222  //label nDeleted = 0;
223 
224  forAll(edge0, pointI)
225  {
226  if (edge0[pointI] >= 0 && edge1[pointI] >= 0)
227  {
228  // Point used by two edges exactly
229 
230  const edge& e0 = edges[edge0[pointI]];
231  const edge& e1 = edges[edge1[pointI]];
232 
233  label common = e0.commonVertex(e1);
234  label vLeft = e0.otherVertex(common);
235  label vRight = e1.otherVertex(common);
236 
237  vector e0Vec = points[common] - points[vLeft];
238  e0Vec /= mag(e0Vec) + VSMALL;
239 
240  vector e1Vec = points[vRight] - points[common];
241  e1Vec /= mag(e1Vec) + VSMALL;
242 
243  if ((e0Vec & e1Vec) > minCos)
244  {
245  pointCanBeDeleted[pointI] = true;
246  //nDeleted++;
247  }
248  }
249  else if (edge0[pointI] == -1)
250  {
251  // point not used at all
252  pointCanBeDeleted[pointI] = true;
253  //nDeleted++;
254  }
255  }
256  edge0.clear();
257  edge1.clear();
258 
259 
260  // Protect any points on faces that would collapse down to nothing
261  // No particular intelligence so might protect too many points
262  forAll(mesh_.faces(), faceI)
263  {
264  const face& f = mesh_.faces()[faceI];
265 
266  label nCollapse = 0;
267  forAll(f, fp)
268  {
269  if (pointCanBeDeleted[f[fp]])
270  {
271  nCollapse++;
272  }
273  }
274 
275  if ((f.size() - nCollapse) < 3)
276  {
277  // Just unmark enough points
278  forAll(f, fp)
279  {
280  if (pointCanBeDeleted[f[fp]])
281  {
282  pointCanBeDeleted[f[fp]] = false;
283  --nCollapse;
284  if (nCollapse == 0)
285  {
286  break;
287  }
288  }
289  }
290  }
291  }
292 
293 
294  // Point can be deleted only if all processors want to delete it
296  (
297  mesh_,
298  pointCanBeDeleted,
299  andEqOp<bool>(),
300  true // null value
301  );
302 
303  label nDeleted = 0;
304  forAll(pointCanBeDeleted, pointI)
305  {
306  if (pointCanBeDeleted[pointI])
307  {
308  nDeleted++;
309  }
310  }
311 
312  return returnReduce(nDeleted, sumOp<label>());
313 }
314 
315 
317 (
318  const boolList& pointCanBeDeleted,
319  polyTopoChange& meshMod
320 )
321 {
322  // Count deleted points
323  label nDeleted = 0;
324  forAll(pointCanBeDeleted, pointI)
325  {
326  if (pointCanBeDeleted[pointI])
327  {
328  nDeleted++;
329  }
330  }
331 
332  // Faces (in mesh face labels) affected by points removed. Will hopefully
333  // be only a few.
334  labelHashSet facesAffected(4*nDeleted);
335 
336 
337  // Undo: from global mesh point to index in savedPoint_
338  Map<label> pointToSaved;
339 
340  // Size undo storage
341  if (undoable_)
342  {
343  savedPoints_.setSize(nDeleted);
344  pointToSaved.resize(2*nDeleted);
345  }
346 
347 
348  // Remove points
349  // ~~~~~~~~~~~~~
350 
351  nDeleted = 0;
352 
353  forAll(pointCanBeDeleted, pointI)
354  {
355  if (pointCanBeDeleted[pointI])
356  {
357  if (undoable_)
358  {
359  pointToSaved.insert(pointI, nDeleted);
360  savedPoints_[nDeleted++] = mesh_.points()[pointI];
361  }
362  meshMod.setAction(polyRemovePoint(pointI));
363 
364  // Store faces affected
365  const labelList& pFaces = mesh_.pointFaces()[pointI];
366 
367  forAll(pFaces, i)
368  {
369  facesAffected.insert(pFaces[i]);
370  }
371  }
372  }
373 
374 
375 
376  // Update faces
377  // ~~~~~~~~~~~~
378 
379 
380  if (undoable_)
381  {
382  savedFaceLabels_.setSize(facesAffected.size());
383  savedFaces_.setSize(facesAffected.size());
384  }
385  label nSaved = 0;
386 
387  forAllConstIter(labelHashSet, facesAffected, iter)
388  {
389  label faceI = iter.key();
390 
391  const face& f = mesh_.faces()[faceI];
392 
393  face newFace(f.size());
394 
395  label newI = 0;
396 
397  forAll(f, fp)
398  {
399  label pointI = f[fp];
400 
401  if (!pointCanBeDeleted[pointI])
402  {
403  newFace[newI++] = pointI;
404  }
405  }
406  newFace.setSize(newI);
407 
408  // Actually change the face to the new vertices
409  modifyFace(faceI, newFace, meshMod);
410 
411  // Save the face. Negative indices are into savedPoints_
412  if (undoable_)
413  {
414  savedFaceLabels_[nSaved] = faceI;
415 
416  face& savedFace = savedFaces_[nSaved++];
417  savedFace.setSize(f.size());
418 
419  forAll(f, fp)
420  {
421  label pointI = f[fp];
422 
423  if (pointCanBeDeleted[pointI])
424  {
425  savedFace[fp] = -pointToSaved[pointI]-1;
426  }
427  else
428  {
429  savedFace[fp] = pointI;
430  }
431  }
432  }
433  }
434 
435  if (undoable_)
436  {
437  // DEBUG: Compare the stored faces with the current ones.
438  if (debug)
439  {
440  forAll(savedFaceLabels_, saveI)
441  {
442  // Points from the mesh
443  List<point> meshPoints
444  (
446  (
447  mesh_.points(),
448  mesh_.faces()[savedFaceLabels_[saveI]] // mesh face
449  )
450  );
451 
452  // Points from the saved face
453  List<point> keptPoints
454  (
456  (
457  mesh_.points(),
458  savedPoints_,
459  savedFaces_[saveI] // saved face
460  )
461  );
462 
463  if (meshPoints != keptPoints)
464  {
465  FatalErrorIn("setRefinement")
466  << "faceI:" << savedFaceLabels_[saveI] << nl
467  << "meshPoints:" << meshPoints << nl
468  << "keptPoints:" << keptPoints << nl
469  << abort(FatalError);
470  }
471  }
472  }
473  }
474 }
475 
476 
478 {
479  if (undoable_)
480  {
481  forAll(savedFaceLabels_, localI)
482  {
483  if (savedFaceLabels_[localI] >= 0)
484  {
485  label newFaceI = map.reverseFaceMap()[savedFaceLabels_[localI]];
486 
487  if (newFaceI == -1)
488  {
490  (
491  "removePoints::updateMesh(const mapPolyMesh&)"
492  ) << "Old face " << savedFaceLabels_[localI]
493  << " seems to have dissapeared."
494  << abort(FatalError);
495  }
496  savedFaceLabels_[localI] = newFaceI;
497  }
498  }
499 
500 
501  // Renumber mesh vertices (indices >=0). Leave saved vertices
502  // (<0) intact.
503  forAll(savedFaces_, i)
504  {
505  face& f = savedFaces_[i];
506 
507  forAll(f, fp)
508  {
509  label pointI = f[fp];
510 
511  if (pointI >= 0)
512  {
513  f[fp] = map.reversePointMap()[pointI];
514 
515  if (f[fp] == -1)
516  {
518  (
519  "removePoints::updateMesh(const mapPolyMesh&)"
520  ) << "Old point " << pointI
521  << " seems to have dissapeared."
522  << abort(FatalError);
523  }
524  }
525  }
526  }
527 
528 
529  // DEBUG: Compare the stored faces with the current ones.
530  if (debug)
531  {
532  forAll(savedFaceLabels_, saveI)
533  {
534  if (savedFaceLabels_[saveI] >= 0)
535  {
536  const face& f = mesh_.faces()[savedFaceLabels_[saveI]];
537 
538  // Get kept points of saved faces.
539  const face& savedFace = savedFaces_[saveI];
540 
541  face keptFace(savedFace.size());
542  label keptFp = 0;
543 
544  forAll(savedFace, fp)
545  {
546  label pointI = savedFace[fp];
547 
548  if (pointI >= 0)
549  {
550  keptFace[keptFp++] = pointI;
551  }
552  }
553  keptFace.setSize(keptFp);
554 
555  // Compare as faces (since f might have rotated and
556  // face::operator== takes care of this)
557  if (keptFace != f)
558  {
559  FatalErrorIn("setRefinement")
560  << "faceI:" << savedFaceLabels_[saveI] << nl
561  << "face:" << f << nl
562  << "keptFace:" << keptFace << nl
563  << "saved points:"
565  (
566  mesh_.points(),
567  savedPoints_,
568  savedFace
569  )() << nl
570  << abort(FatalError);
571  }
572  }
573  }
574  }
575  }
576 }
577 
578 
579 // Given list of faces to undo picks up the local indices of the faces
580 // to restore. Additionally it also picks up all the faces that use
581 // any of the deleted points.
583 (
584  const labelList& undoFaces,
585  labelList& localFaces,
586  labelList& localPoints
587 ) const
588 {
589  if (!undoable_)
590  {
592  (
593  "removePoints::getUnrefimentSet(const labelList&"
594  ", labelList&, labelList&) const"
595  ) << "removePoints not constructed with"
596  << " unrefinement capability."
597  << abort(FatalError);
598  }
599 
600  if (debug)
601  {
602  // Check if synced.
603  faceSet undoFacesSet(mesh_, "undoFacesSet", undoFaces);
604  label sz = undoFacesSet.size();
605 
606  undoFacesSet.sync(mesh_);
607  if (sz != undoFacesSet.size())
608  {
610  (
611  "removePoints::getUnrefimentSet(const labelList&"
612  ", labelList&, labelList&) const"
613  ) << "undoFaces not synchronised across coupled faces." << endl
614  << "Size before sync:" << sz
615  << " after sync:" << undoFacesSet.size()
616  << abort(FatalError);
617  }
618  }
619 
620 
621  // Problem: input undoFaces are synced. But e.g.
622  // two faces, A (uncoupled) and B(coupled) share a deleted point. A gets
623  // passed in to be restored. Now picking up the deleted point and the
624  // original faces using it picks up face B. But not its coupled neighbour!
625  // Problem is that we cannot easily synchronise the deleted points
626  // (e.g. syncPointList) since they're not in the mesh anymore - only the
627  // faces are. So we have to collect the points-to-restore as indices
628  // in the faces (which is information we can synchronise)
629 
630 
631 
632  // Mark points-to-restore
633  labelHashSet localPointsSet(undoFaces.size());
634 
635  {
636  // Create set of faces to be restored
637  labelHashSet undoFacesSet(undoFaces);
638 
639  forAll(savedFaceLabels_, saveI)
640  {
641  if (savedFaceLabels_[saveI] < 0)
642  {
644  (
645  "removePoints::getUnrefimentSet(const labelList&"
646  ", labelList&, labelList&) const"
647  ) << "Illegal face label " << savedFaceLabels_[saveI]
648  << " at index " << saveI
649  << abort(FatalError);
650  }
651 
652  if (undoFacesSet.found(savedFaceLabels_[saveI]))
653  {
654  const face& savedFace = savedFaces_[saveI];
655 
656  forAll(savedFace, fp)
657  {
658  if (savedFace[fp] < 0)
659  {
660  label savedPointI = -savedFace[fp]-1;
661 
662  if (savedPoints_[savedPointI] == vector::max)
663  {
665  (
666  "removePoints::getUnrefimentSet"
667  "(const labelList&, labelList&, labelList&)"
668  " const"
669  ) << "Trying to restore point " << savedPointI
670  << " from mesh face " << savedFaceLabels_[saveI]
671  << " saved face:" << savedFace
672  << " which has already been undone."
673  << abort(FatalError);
674  }
675 
676  localPointsSet.insert(savedPointI);
677  }
678  }
679  }
680  }
681 
682 
683  // Per boundary face, per index in face whether the point needs
684  // restoring. Note that this is over all saved faces, not just over
685  // the ones in undoFaces.
686 
687  boolListList faceVertexRestore(mesh_.nFaces()-mesh_.nInternalFaces());
688 
689  // Populate with my local points-to-restore.
690  forAll(savedFaces_, saveI)
691  {
692  label bFaceI = savedFaceLabels_[saveI] - mesh_.nInternalFaces();
693 
694  if (bFaceI >= 0)
695  {
696  const face& savedFace = savedFaces_[saveI];
697 
698  boolList& fRestore = faceVertexRestore[bFaceI];
699 
700  fRestore.setSize(savedFace.size());
701  fRestore = false;
702 
703  forAll(savedFace, fp)
704  {
705  if (savedFace[fp] < 0)
706  {
707  label savedPointI = -savedFace[fp]-1;
708 
709  if (localPointsSet.found(savedPointI))
710  {
711  fRestore[fp] = true;
712  }
713  }
714  }
715  }
716  }
717 
719  (
720  mesh_,
721  faceVertexRestore,
722  faceEqOp<bool, orEqOp>(), // special operator to handle faces
723  Foam::dummyTransform() // no transformation
724  );
725 
726  // So now if any of the points-to-restore is used by any coupled face
727  // anywhere the corresponding index in faceVertexRestore will be set.
728 
729  // Now combine the localPointSet and the (sychronised)
730  // boundary-points-to-restore.
731 
732  forAll(savedFaces_, saveI)
733  {
734  label bFaceI = savedFaceLabels_[saveI] - mesh_.nInternalFaces();
735 
736  if (bFaceI >= 0)
737  {
738  const boolList& fRestore = faceVertexRestore[bFaceI];
739 
740  const face& savedFace = savedFaces_[saveI];
741 
742  forAll(fRestore, fp)
743  {
744  // Does neighbour require point restored?
745  if (fRestore[fp])
746  {
747  if (savedFace[fp] >= 0)
748  {
750  (
751  "removePoints::getUnrefimentSet"
752  "(const labelList&, labelList&, labelList&)"
753  " const"
754  ) << "Problem: on coupled face:"
755  << savedFaceLabels_[saveI]
756  << " fc:"
757  << mesh_.faceCentres()[savedFaceLabels_[saveI]]
758  << endl
759  << " my neighbour tries to restore the vertex"
760  << " at index " << fp
761  << " whereas my saved face:" << savedFace
762  << " does not indicate a deleted vertex"
763  << " at that position."
764  << abort(FatalError);
765  }
766 
767  label savedPointI = -savedFace[fp]-1;
768 
769  localPointsSet.insert(savedPointI);
770  }
771  }
772  }
773  }
774  }
775 
776  localPoints = localPointsSet.toc();
777 
778 
779  // Collect all saved faces using any localPointsSet
780  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
781 
782  labelHashSet localFacesSet(2*undoFaces.size());
783 
784  forAll(savedFaces_, saveI)
785  {
786  const face& savedFace = savedFaces_[saveI];
787 
788  forAll(savedFace, fp)
789  {
790  if (savedFace[fp] < 0)
791  {
792  label savedPointI = -savedFace[fp]-1;
793 
794  if (localPointsSet.found(savedPointI))
795  {
796  localFacesSet.insert(saveI);
797  }
798  }
799  }
800  }
801  localFaces = localFacesSet.toc();
802 
803 
804  // Note that at this point the localFaces to restore can contain points
805  // that are not going to be restored! The localFaces though will
806  // be guaranteed to be all the ones affected by the restoration of the
807  // localPoints.
808 }
809 
810 
812 (
813  const labelList& localFaces,
814  const labelList& localPoints,
815  polyTopoChange& meshMod
816 )
817 {
818  if (!undoable_)
819  {
821  (
822  "removePoints::setUnrefinement(const labelList&"
823  ", labelList&, polyTopoChange&)"
824  ) << "removePoints not constructed with"
825  << " unrefinement capability."
826  << abort(FatalError);
827  }
828 
829 
830  // Per savedPoint -1 or the restored point label
831  labelList addedPoints(savedPoints_.size(), -1);
832 
833  forAll(localPoints, i)
834  {
835  label localI = localPoints[i];
836 
837  if (savedPoints_[localI] == vector::max)
838  {
840  (
841  "removePoints::setUnrefinement(const labelList&"
842  ", labelList&, polyTopoChange&)"
843  ) << "Saved point " << localI << " already restored!"
844  << abort(FatalError);
845  }
846 
847  addedPoints[localI] = meshMod.setAction
848  (
850  (
851  savedPoints_[localI], // point
852  -1, // master point
853  -1, // zone for point
854  true // supports a cell
855  )
856  );
857 
858  // Mark the restored points so they are not restored again.
859  savedPoints_[localI] = vector::max;
860  }
861 
862  forAll(localFaces, i)
863  {
864  label saveI = localFaces[i];
865 
866  // Modify indices into saved points (so < 0) to point to the
867  // added points.
868  face& savedFace = savedFaces_[saveI];
869 
870  face newFace(savedFace.size());
871  label newFp = 0;
872 
873  bool hasSavedPoints = false;
874 
875  forAll(savedFace, fp)
876  {
877  if (savedFace[fp] < 0)
878  {
879  label addedPointI = addedPoints[-savedFace[fp]-1];
880 
881  if (addedPointI != -1)
882  {
883  savedFace[fp] = addedPointI;
884  newFace[newFp++] = addedPointI;
885  }
886  else
887  {
888  hasSavedPoints = true;
889  }
890  }
891  else
892  {
893  newFace[newFp++] = savedFace[fp];
894  }
895  }
896  newFace.setSize(newFp);
897 
898  modifyFace(savedFaceLabels_[saveI], newFace, meshMod);
899 
900  if (!hasSavedPoints)
901  {
902  // Face fully restored. Mark for compaction later on
903  savedFaceLabels_[saveI] = -1;
904  savedFaces_[saveI].clear();
905  }
906  }
907 
908 
909  // Compact face labels
910  label newSaveI = 0;
911 
912  forAll(savedFaceLabels_, saveI)
913  {
914  if (savedFaceLabels_[saveI] != -1)
915  {
916  if (newSaveI != saveI)
917  {
918  savedFaceLabels_[newSaveI] = savedFaceLabels_[saveI];
919  savedFaces_[newSaveI].transfer(savedFaces_[saveI]);
920  }
921  newSaveI++;
922  }
923  }
924 
925  savedFaceLabels_.setSize(newSaveI);
926  savedFaces_.setSize(newSaveI);
927 
928 
929  // Check that all faces have been restored that use any restored points
930  if (debug)
931  {
932  forAll(savedFaceLabels_, saveI)
933  {
934  const face& savedFace = savedFaces_[saveI];
935 
936  forAll(savedFace, fp)
937  {
938  if (savedFace[fp] < 0)
939  {
940  label addedPointI = addedPoints[-savedFace[fp]-1];
941 
942  if (addedPointI != -1)
943  {
944  FatalErrorIn("setUnrefinement")
945  << "Face:" << savedFaceLabels_[saveI]
946  << " savedVerts:" << savedFace
947  << " uses restored point:" << -savedFace[fp]-1
948  << " with new pointlabel:" << addedPointI
949  << abort(FatalError);
950  }
951  }
952  }
953  }
954  }
955 }
956 
957 
958 // ************************************************************************* //
const pointField & points
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:249
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.
const labelList & reverseFaceMap() const
Reverse face map.
Definition: mapPolyMesh.H:497
A list of face labels.
Definition: faceSet.H:48
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
label commonVertex(const edge &a) const
Return common vertex.
Definition: edgeI.H:121
dimensioned< scalar > mag(const dimensioned< Type > &)
labelList f(nPoints)
Class containing data for point removal.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
const labelList & reversePointMap() const
Reverse point map.
Definition: mapPolyMesh.H:465
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
void setUnrefinement(const labelList &localFaces, const labelList &localPoints, polyTopoChange &)
Restore selected faces and vertices.
Definition: removePoints.C:812
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
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
void resize(const label newSize)
Resize the hash table for efficiency.
Definition: HashTable.C:436
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
static const Vector max
Definition: Vector.H:82
static void syncBoundaryFaceList(const polyMesh &, UList< T > &, const CombineOp &cop, const TransformOp &top)
Synchronize values on boundary faces only.
void updateMesh(const mapPolyMesh &)
Force recalculation of locally stored data on topological change.
Definition: removePoints.C:477
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Namespace for OpenFOAM.
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:317
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
void getUnrefimentSet(const labelList &undoFaces, labelList &localFaces, labelList &localPoints) const
Given set of faces to restore calculates a consistent set of.
Definition: removePoints.C:583
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
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
void clear()
Clear the list, i.e. set size to zero.
Definition: List.C:379
static const char nl
Definition: Ostream.H:260
const double e
Elementary charge.
Definition: doubleFloat.H:78
void setSize(const label)
Reset size of List.
Definition: List.C:318
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual void sync(const polyMesh &mesh)
Sync faceSet across coupled patches.
Definition: faceSet.C:113
A List with indirect addressing.
Definition: fvMatrix.H:106
#define forAll(list, i)
Definition: UList.H:421
label otherVertex(const label a) const
Given one vertex, return the other.
Definition: edgeI.H:103
void operator()(List< T > &x, const List< T > &y) const
Definition: removePoints.C:53
Class describing modification of a face.
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
errorManip< error > abort(error &err)
Definition: errorManip.H:131
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
error FatalError
scalar y
Direct mesh changes based on v1.3 polyTopoChange syntax.
Dummy transform to be used with syncTools.
label rcIndex(const label i) const
Return the reverse circular index, i.e. the previous index.
Definition: UListI.H:65
Class containing data for point addition.
Definition: polyAddPoint.H:47
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints,-1);const cellModel &hex=*(cellModeller::lookup("hex"));labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){pointMap[i]=i;}for(label i=0;i< nPoints;i++){if(f[i] > 0.0){hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells]=cellShape(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei){if(edges[ei].mag(points)< SMALL){label start=pointMap[edges[ei].start()];while(start!=pointMap[start]){start=pointMap[start];}label end=pointMap[edges[ei].end()];while(end!=pointMap[end]){end=pointMap[end];}label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;}}cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){cellShape &cs=cellShapes[celli];forAll(cs, i){cs[i]=pointMap[cs[i]];}cs.collapse();}label bcIDs[11]={-1, 0, 2, 4,-1, 5,-1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={&wallPolyPatch::typeName,&wallPolyPatch::typeName,&wallPolyPatch::typeName,&wallPolyPatch::typeName,&symmetryPolyPatch::typeName,&wedgePolyPatch::typeName,&polyPatch::typeName,&polyPatch::typeName,&polyPatch::typeName,&polyPatch::typeName,&symmetryPolyPatch::typeName,&oldCyclicPolyPatch::typeName};enum patchTypeNames{PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={"piston","valve","liner","cylinderHead","axis","wedge","inflow","outflow","presin","presout","symmetryPlane","cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:235
Indexes into negList (negative index) or posList (zero or positive index).
void setRefinement(const boolList &, polyTopoChange &)
Play commands into polyTopoChange to remove points. Gets.
Definition: removePoints.C:317
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
defineTypeNameAndDebug(combustionModel, 0)
Combine-reduce operator to combine data on faces. Takes care.
Definition: removePoints.C:48
label countPointUsage(const scalar minCos, boolList &pointCanBeDeleted) const
Mark in pointCanBeDeleted the points that can be deleted.
Definition: removePoints.C:166
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116