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