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