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-2015 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.
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
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
295  syncTools::syncPointList
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  {
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  << "Old face " << savedFaceLabels_[localI]
491  << " seems to have dissapeared."
492  << abort(FatalError);
493  }
494  savedFaceLabels_[localI] = newFaceI;
495  }
496  }
497 
498 
499  // Renumber mesh vertices (indices >=0). Leave saved vertices
500  // (<0) intact.
501  forAll(savedFaces_, i)
502  {
503  face& f = savedFaces_[i];
504 
505  forAll(f, fp)
506  {
507  label pointI = f[fp];
508 
509  if (pointI >= 0)
510  {
511  f[fp] = map.reversePointMap()[pointI];
512 
513  if (f[fp] == -1)
514  {
516  << "Old point " << pointI
517  << " seems to have dissapeared."
518  << abort(FatalError);
519  }
520  }
521  }
522  }
523 
524 
525  // DEBUG: Compare the stored faces with the current ones.
526  if (debug)
527  {
528  forAll(savedFaceLabels_, saveI)
529  {
530  if (savedFaceLabels_[saveI] >= 0)
531  {
532  const face& f = mesh_.faces()[savedFaceLabels_[saveI]];
533 
534  // Get kept points of saved faces.
535  const face& savedFace = savedFaces_[saveI];
536 
537  face keptFace(savedFace.size());
538  label keptFp = 0;
539 
540  forAll(savedFace, fp)
541  {
542  label pointI = savedFace[fp];
543 
544  if (pointI >= 0)
545  {
546  keptFace[keptFp++] = pointI;
547  }
548  }
549  keptFace.setSize(keptFp);
550 
551  // Compare as faces (since f might have rotated and
552  // face::operator== takes care of this)
553  if (keptFace != f)
554  {
556  << "faceI:" << savedFaceLabels_[saveI] << nl
557  << "face:" << f << nl
558  << "keptFace:" << keptFace << nl
559  << "saved points:"
561  (
562  mesh_.points(),
563  savedPoints_,
564  savedFace
565  )() << nl
566  << abort(FatalError);
567  }
568  }
569  }
570  }
571  }
572 }
573 
574 
575 // Given list of faces to undo picks up the local indices of the faces
576 // to restore. Additionally it also picks up all the faces that use
577 // any of the deleted points.
579 (
580  const labelList& undoFaces,
581  labelList& localFaces,
582  labelList& localPoints
583 ) const
584 {
585  if (!undoable_)
586  {
588  << "removePoints not constructed with"
589  << " unrefinement capability."
590  << abort(FatalError);
591  }
592 
593  if (debug)
594  {
595  // Check if synced.
596  faceSet undoFacesSet(mesh_, "undoFacesSet", undoFaces);
597  label sz = undoFacesSet.size();
598 
599  undoFacesSet.sync(mesh_);
600  if (sz != undoFacesSet.size())
601  {
603  << "undoFaces not synchronised across coupled faces." << endl
604  << "Size before sync:" << sz
605  << " after sync:" << undoFacesSet.size()
606  << abort(FatalError);
607  }
608  }
609 
610 
611  // Problem: input undoFaces are synced. But e.g.
612  // two faces, A (uncoupled) and B(coupled) share a deleted point. A gets
613  // passed in to be restored. Now picking up the deleted point and the
614  // original faces using it picks up face B. But not its coupled neighbour!
615  // Problem is that we cannot easily synchronise the deleted points
616  // (e.g. syncPointList) since they're not in the mesh anymore - only the
617  // faces are. So we have to collect the points-to-restore as indices
618  // in the faces (which is information we can synchronise)
619 
620 
621 
622  // Mark points-to-restore
623  labelHashSet localPointsSet(undoFaces.size());
624 
625  {
626  // Create set of faces to be restored
627  labelHashSet undoFacesSet(undoFaces);
628 
629  forAll(savedFaceLabels_, saveI)
630  {
631  if (savedFaceLabels_[saveI] < 0)
632  {
634  << "Illegal face label " << savedFaceLabels_[saveI]
635  << " at index " << saveI
636  << abort(FatalError);
637  }
638 
639  if (undoFacesSet.found(savedFaceLabels_[saveI]))
640  {
641  const face& savedFace = savedFaces_[saveI];
642 
643  forAll(savedFace, fp)
644  {
645  if (savedFace[fp] < 0)
646  {
647  label savedPointI = -savedFace[fp]-1;
648 
649  if (savedPoints_[savedPointI] == vector::max)
650  {
652  << "Trying to restore point " << savedPointI
653  << " from mesh face " << savedFaceLabels_[saveI]
654  << " saved face:" << savedFace
655  << " which has already been undone."
656  << abort(FatalError);
657  }
658 
659  localPointsSet.insert(savedPointI);
660  }
661  }
662  }
663  }
664 
665 
666  // Per boundary face, per index in face whether the point needs
667  // restoring. Note that this is over all saved faces, not just over
668  // the ones in undoFaces.
669 
670  boolListList faceVertexRestore(mesh_.nFaces()-mesh_.nInternalFaces());
671 
672  // Populate with my local points-to-restore.
673  forAll(savedFaces_, saveI)
674  {
675  label bFaceI = savedFaceLabels_[saveI] - mesh_.nInternalFaces();
676 
677  if (bFaceI >= 0)
678  {
679  const face& savedFace = savedFaces_[saveI];
680 
681  boolList& fRestore = faceVertexRestore[bFaceI];
682 
683  fRestore.setSize(savedFace.size());
684  fRestore = false;
685 
686  forAll(savedFace, fp)
687  {
688  if (savedFace[fp] < 0)
689  {
690  label savedPointI = -savedFace[fp]-1;
691 
692  if (localPointsSet.found(savedPointI))
693  {
694  fRestore[fp] = true;
695  }
696  }
697  }
698  }
699  }
700 
702  (
703  mesh_,
704  faceVertexRestore,
705  faceEqOp<bool, orEqOp>(), // special operator to handle faces
706  Foam::dummyTransform() // no transformation
707  );
708 
709  // So now if any of the points-to-restore is used by any coupled face
710  // anywhere the corresponding index in faceVertexRestore will be set.
711 
712  // Now combine the localPointSet and the (sychronised)
713  // boundary-points-to-restore.
714 
715  forAll(savedFaces_, saveI)
716  {
717  label bFaceI = savedFaceLabels_[saveI] - mesh_.nInternalFaces();
718 
719  if (bFaceI >= 0)
720  {
721  const boolList& fRestore = faceVertexRestore[bFaceI];
722 
723  const face& savedFace = savedFaces_[saveI];
724 
725  forAll(fRestore, fp)
726  {
727  // Does neighbour require point restored?
728  if (fRestore[fp])
729  {
730  if (savedFace[fp] >= 0)
731  {
733  << "Problem: on coupled face:"
734  << savedFaceLabels_[saveI]
735  << " fc:"
736  << mesh_.faceCentres()[savedFaceLabels_[saveI]]
737  << endl
738  << " my neighbour tries to restore the vertex"
739  << " at index " << fp
740  << " whereas my saved face:" << savedFace
741  << " does not indicate a deleted vertex"
742  << " at that position."
743  << abort(FatalError);
744  }
745 
746  label savedPointI = -savedFace[fp]-1;
747 
748  localPointsSet.insert(savedPointI);
749  }
750  }
751  }
752  }
753  }
754 
755  localPoints = localPointsSet.toc();
756 
757 
758  // Collect all saved faces using any localPointsSet
759  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
760 
761  labelHashSet localFacesSet(2*undoFaces.size());
762 
763  forAll(savedFaces_, saveI)
764  {
765  const face& savedFace = savedFaces_[saveI];
766 
767  forAll(savedFace, fp)
768  {
769  if (savedFace[fp] < 0)
770  {
771  label savedPointI = -savedFace[fp]-1;
772 
773  if (localPointsSet.found(savedPointI))
774  {
775  localFacesSet.insert(saveI);
776  }
777  }
778  }
779  }
780  localFaces = localFacesSet.toc();
781 
782 
783  // Note that at this point the localFaces to restore can contain points
784  // that are not going to be restored! The localFaces though will
785  // be guaranteed to be all the ones affected by the restoration of the
786  // localPoints.
787 }
788 
789 
791 (
792  const labelList& localFaces,
793  const labelList& localPoints,
794  polyTopoChange& meshMod
795 )
796 {
797  if (!undoable_)
798  {
800  << "removePoints not constructed with"
801  << " unrefinement capability."
802  << abort(FatalError);
803  }
804 
805 
806  // Per savedPoint -1 or the restored point label
807  labelList addedPoints(savedPoints_.size(), -1);
808 
809  forAll(localPoints, i)
810  {
811  label localI = localPoints[i];
812 
813  if (savedPoints_[localI] == vector::max)
814  {
816  << "Saved point " << localI << " already restored!"
817  << abort(FatalError);
818  }
819 
820  addedPoints[localI] = meshMod.setAction
821  (
823  (
824  savedPoints_[localI], // point
825  -1, // master point
826  -1, // zone for point
827  true // supports a cell
828  )
829  );
830 
831  // Mark the restored points so they are not restored again.
832  savedPoints_[localI] = vector::max;
833  }
834 
835  forAll(localFaces, i)
836  {
837  label saveI = localFaces[i];
838 
839  // Modify indices into saved points (so < 0) to point to the
840  // added points.
841  face& savedFace = savedFaces_[saveI];
842 
843  face newFace(savedFace.size());
844  label newFp = 0;
845 
846  bool hasSavedPoints = false;
847 
848  forAll(savedFace, fp)
849  {
850  if (savedFace[fp] < 0)
851  {
852  label addedPointI = addedPoints[-savedFace[fp]-1];
853 
854  if (addedPointI != -1)
855  {
856  savedFace[fp] = addedPointI;
857  newFace[newFp++] = addedPointI;
858  }
859  else
860  {
861  hasSavedPoints = true;
862  }
863  }
864  else
865  {
866  newFace[newFp++] = savedFace[fp];
867  }
868  }
869  newFace.setSize(newFp);
870 
871  modifyFace(savedFaceLabels_[saveI], newFace, meshMod);
872 
873  if (!hasSavedPoints)
874  {
875  // Face fully restored. Mark for compaction later on
876  savedFaceLabels_[saveI] = -1;
877  savedFaces_[saveI].clear();
878  }
879  }
880 
881 
882  // Compact face labels
883  label newSaveI = 0;
884 
885  forAll(savedFaceLabels_, saveI)
886  {
887  if (savedFaceLabels_[saveI] != -1)
888  {
889  if (newSaveI != saveI)
890  {
891  savedFaceLabels_[newSaveI] = savedFaceLabels_[saveI];
892  savedFaces_[newSaveI].transfer(savedFaces_[saveI]);
893  }
894  newSaveI++;
895  }
896  }
897 
898  savedFaceLabels_.setSize(newSaveI);
899  savedFaces_.setSize(newSaveI);
900 
901 
902  // Check that all faces have been restored that use any restored points
903  if (debug)
904  {
905  forAll(savedFaceLabels_, saveI)
906  {
907  const face& savedFace = savedFaces_[saveI];
908 
909  forAll(savedFace, fp)
910  {
911  if (savedFace[fp] < 0)
912  {
913  label addedPointI = addedPoints[-savedFace[fp]-1];
914 
915  if (addedPointI != -1)
916  {
918  << "Face:" << savedFaceLabels_[saveI]
919  << " savedVerts:" << savedFace
920  << " uses restored point:" << -savedFace[fp]-1
921  << " with new pointlabel:" << addedPointI
922  << abort(FatalError);
923  }
924  }
925  }
926  }
927  }
928 }
929 
930 
931 // ************************************************************************* //
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
Foam::removePoints::removePoints
removePoints(const removePoints &)
Disallow default bitwise copy construct.
polyRemovePoint.H
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::polyRemovePoint
Class containing data for point removal.
Definition: polyRemovePoint.H:46
Foam::HashTable::toc
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:201
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::removePoints::savedFaces_
faceList savedFaces_
If undoable: per stored face the vertices. Negative indices.
Definition: removePoints.H:77
Foam::syncTools::syncBoundaryFaceList
static void syncBoundaryFaceList(const polyMesh &, UList< T > &, const CombineOp &cop, const TransformOp &top, const bool parRun=Pstream::parRun())
Synchronize values on boundary faces only.
Definition: syncToolsTemplates.C:1285
Foam::removePoints::mesh_
const polyMesh & mesh_
Reference to mesh.
Definition: removePoints.H:64
Foam::Vector< scalar >::max
static const Vector max
Definition: Vector.H:82
polyTopoChange.H
Foam::polyTopoChange
Direct mesh changes based on v1.3 polyTopoChange syntax.
Definition: polyTopoChange.H:97
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Foam::removePoints::getUnrefimentSet
void getUnrefimentSet(const labelList &undoFaces, labelList &localFaces, labelList &localPoints) const
Given set of faces to restore calculates a consistent set of.
Definition: removePoints.C:579
dummyTransform.H
Dummy transform to be used with syncTools.
Foam::Map< label >
Foam::faceSet
A list of face labels.
Definition: faceSet.H:48
Foam::polyModifyFace
Class describing modification of a face.
Definition: polyModifyFace.H:48
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::faceEqOp::operator()
void operator()(List< T > &x, const List< T > &y) const
Definition: removePoints.C:53
Foam::dummyTransform
Definition: dummyTransform.H:44
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
polyMesh.H
Foam::HashSet< label, Hash< label > >
syncTools.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
BiIndirectList.H
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
pFaces
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
Foam::faceSet::sync
virtual void sync(const polyMesh &mesh)
Sync faceSet across coupled patches.
Definition: faceSet.C:113
Foam::removePoints::savedFaceLabels_
labelList savedFaceLabels_
If undoable: per stored face the original mesh face label.
Definition: removePoints.H:73
removePoints.H
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::removePoints::modifyFace
void modifyFace(const label faceI, const face &, polyTopoChange &) const
Change the vertices of the face whilst keeping everything else.
Definition: removePoints.C:100
faceSet.H
Foam::BiIndirectList
Indexes into negList (negative index) or posList (zero or positive index).
Definition: BiIndirectList.H:49
polyAddPoint.H
Foam::removePoints::setUnrefinement
void setUnrefinement(const labelList &localFaces, const labelList &localPoints, polyTopoChange &)
Restore selected faces and vertices.
Definition: removePoints.C:791
Foam::FatalError
error FatalError
Foam::HashTable::size
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::faceEqOp
Combine-reduce operator to combine data on faces. Takes care.
Definition: removePoints.C:48
Foam::HashTable::found
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::faceZone::whichFace
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:314
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
PstreamReduceOps.H
polyModifyFace.H
Foam::edge::commonVertex
label commonVertex(const edge &a) const
Return common vertex.
Definition: edgeI.H:121
Foam::mapPolyMesh::reverseFaceMap
const labelList & reverseFaceMap() const
Reverse face map.
Definition: mapPolyMesh.H:497
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1004
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::polyTopoChange::setAction
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.
Definition: polyTopoChange.C:2530
Foam::sumOp
Definition: ops.H:162
Foam::removePoints::undoable_
const bool undoable_
Whether undoable.
Definition: removePoints.H:67
f
labelList f(nPoints)
Foam::mapPolyMesh::reversePointMap
const labelList & reversePointMap() const
Reverse point map.
Definition: mapPolyMesh.H:465
Foam::Vector< scalar >
Foam::List< T >
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:94
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: List.C:379
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::polyAddPoint
Class containing data for point addition.
Definition: polyAddPoint.H:47
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:106
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::removePoints::savedPoints_
pointField savedPoints_
If undoable: deleted points.
Definition: removePoints.H:70
Foam::removePoints::updateMesh
void updateMesh(const mapPolyMesh &)
Force recalculation of locally stored data on topological change.
Definition: removePoints.C:477
Foam::removePoints::setRefinement
void setRefinement(const boolList &, polyTopoChange &)
Play commands into polyTopoChange to remove points. Gets.
Definition: removePoints.C:317
Foam::edge::otherVertex
label otherVertex(const label a) const
Given one vertex, return the other.
Definition: edgeI.H:103
Foam::faceZone::flipMap
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:249
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::removePoints::countPointUsage
label countPointUsage(const scalar minCos, boolList &pointCanBeDeleted) const
Mark in pointCanBeDeleted the points that can be deleted.
Definition: removePoints.C:166
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::andEqOp
Definition: ops.H:81