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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2018 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "BiIndirectList.H"
30 #include "removePoints.H"
31 #include "PstreamReduceOps.H"
32 #include "polyMesh.H"
33 #include "polyTopoChange.H"
34 #include "polyRemovePoint.H"
35 #include "polyAddPoint.H"
36 #include "polyModifyFace.H"
37 #include "syncTools.H"
38 #include "faceSet.H"
39 #include "dummyTransform.H"
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 defineTypeNameAndDebug(removePoints, 0);
47 
48 // Combine-reduce operator to combine data on faces. Takes care
49 // of reverse orientation on coupled face.
50 template<class T, template<class> class CombineOp>
51 class faceEqOp
52 {
53 
54 public:
55 
56  void operator()(List<T>& x, const List<T>& y) const
57  {
58  if (y.size() > 0)
59  {
60  if (x.size() == 0)
61  {
62  x = y;
63  }
64  else
65  {
66  label j = 0;
67  forAll(x, i)
68  {
69  CombineOp<T>()(x[i], y[j]);
70  j = y.rcIndex(j);
71  }
72  }
73  }
74  }
75 };
76 
77 }
78 
79 
80 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
81 
82 void Foam::removePoints::modifyFace
83 (
84  const label facei,
85  const face& newFace,
86  polyTopoChange& meshMod
87 ) const
88 {
89  // Get other face data.
90  label patchi = -1;
91  label owner = mesh_.faceOwner()[facei];
92  label neighbour = -1;
93 
94  if (mesh_.isInternalFace(facei))
95  {
96  neighbour = mesh_.faceNeighbour()[facei];
97  }
98  else
99  {
100  patchi = mesh_.boundaryMesh().whichPatch(facei);
101  }
102 
103  label zoneID = mesh_.faceZones().whichZone(facei);
104 
105  bool zoneFlip = false;
106 
107  if (zoneID >= 0)
108  {
109  const faceZone& fZone = mesh_.faceZones()[zoneID];
110 
111  zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
112  }
113 
114  meshMod.setAction
115  (
116  polyModifyFace
117  (
118  newFace, // modified face
119  facei, // label of face being modified
120  owner, // owner
121  neighbour, // neighbour
122  false, // face flip
123  patchi, // patch for face
124  false, // remove from zone
125  zoneID, // zone for face
126  zoneFlip // face flip in zone
127  )
128  );
129 }
130 
131 
132 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
133 
134 Foam::removePoints::removePoints
135 (
136  const polyMesh& mesh,
137  const bool undoable
138 )
139 :
140  mesh_(mesh),
141  undoable_(undoable)
142 {}
143 
144 
145 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
146 
148 (
149  const scalar minCos,
150  boolList& pointCanBeDeleted
151 ) const
152 {
153  // Containers to store two edges per point:
154  // -1 : not filled
155  // >= 0 : edge label
156  // -2 : more than two edges for point
157  labelList edge0(mesh_.nPoints(), -1);
158  labelList edge1(mesh_.nPoints(), -1);
159 
160  const edgeList& edges = mesh_.edges();
161 
162  forAll(edges, edgeI)
163  {
164  const edge& e = edges[edgeI];
165 
166  forAll(e, eI)
167  {
168  label pointi = e[eI];
169 
170  if (edge0[pointi] == -2)
171  {
172  // Already too many edges
173  }
174  else if (edge0[pointi] == -1)
175  {
176  // Store first edge using point
177  edge0[pointi] = edgeI;
178  }
179  else
180  {
181  // Already one edge using point. Check second container.
182 
183  if (edge1[pointi] == -1)
184  {
185  // Store second edge using point
186  edge1[pointi] = edgeI;
187  }
188  else
189  {
190  // Third edge using point. Mark.
191  edge0[pointi] = -2;
192  edge1[pointi] = -2;
193  }
194  }
195  }
196  }
197 
198 
199  // Check the ones used by only 2 edges that these are sufficiently in line.
200  const pointField& points = mesh_.points();
201 
202  pointCanBeDeleted.setSize(mesh_.nPoints());
203  pointCanBeDeleted = false;
204  //label nDeleted = 0;
205 
206  forAll(edge0, pointi)
207  {
208  if (edge0[pointi] >= 0 && edge1[pointi] >= 0)
209  {
210  // Point used by two edges exactly
211 
212  const edge& e0 = edges[edge0[pointi]];
213  const edge& e1 = edges[edge1[pointi]];
214 
215  label common = e0.commonVertex(e1);
216  label vLeft = e0.otherVertex(common);
217  label vRight = e1.otherVertex(common);
218 
219  const vector e0Vec = normalised(points[common] - points[vLeft]);
220  const vector e1Vec = normalised(points[vRight] - points[common]);
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  facesAffected.insert(pFaces);
347  }
348  }
349 
350 
351 
352  // Update faces
353  // ~~~~~~~~~~~~
354 
355 
356  if (undoable_)
357  {
358  savedFaceLabels_.setSize(facesAffected.size());
359  savedFaces_.setSize(facesAffected.size());
360  }
361  label nSaved = 0;
362 
363  for (const label facei : facesAffected)
364  {
365  const face& f = mesh_.faces()[facei];
366 
367  face newFace(f.size());
368 
369  label newI = 0;
370 
371  forAll(f, fp)
372  {
373  label pointi = f[fp];
374 
375  if (!pointCanBeDeleted[pointi])
376  {
377  newFace[newI++] = pointi;
378  }
379  }
380  newFace.setSize(newI);
381 
382  // Actually change the face to the new vertices
383  modifyFace(facei, newFace, meshMod);
384 
385  // Save the face. Negative indices are into savedPoints_
386  if (undoable_)
387  {
388  savedFaceLabels_[nSaved] = facei;
389 
390  face& savedFace = savedFaces_[nSaved++];
391  savedFace.setSize(f.size());
392 
393  forAll(f, fp)
394  {
395  label pointi = f[fp];
396 
397  if (pointCanBeDeleted[pointi])
398  {
399  savedFace[fp] = -pointToSaved[pointi]-1;
400  }
401  else
402  {
403  savedFace[fp] = pointi;
404  }
405  }
406  }
407  }
408 
409  if (undoable_)
410  {
411  // DEBUG: Compare the stored faces with the current ones.
412  if (debug)
413  {
414  forAll(savedFaceLabels_, saveI)
415  {
416  // Points from the mesh
417  List<point> meshPoints
418  (
420  (
421  mesh_.points(),
422  mesh_.faces()[savedFaceLabels_[saveI]] // mesh face
423  )
424  );
425 
426  // Points from the saved face
427  List<point> keptPoints
428  (
430  (
431  mesh_.points(),
432  savedPoints_,
433  savedFaces_[saveI] // saved face
434  )()
435  );
436 
437  if (meshPoints != keptPoints)
438  {
440  << "facei:" << savedFaceLabels_[saveI] << nl
441  << "meshPoints:" << meshPoints << nl
442  << "keptPoints:" << keptPoints << nl
443  << abort(FatalError);
444  }
445  }
446  }
447  }
448 }
449 
450 
452 {
453  if (undoable_)
454  {
455  forAll(savedFaceLabels_, localI)
456  {
457  if (savedFaceLabels_[localI] >= 0)
458  {
459  label newFacei = map.reverseFaceMap()[savedFaceLabels_[localI]];
460 
461  if (newFacei == -1)
462  {
464  << "Old face " << savedFaceLabels_[localI]
465  << " seems to have disappeared."
466  << abort(FatalError);
467  }
468  savedFaceLabels_[localI] = newFacei;
469  }
470  }
471 
472 
473  // Renumber mesh vertices (indices >=0). Leave saved vertices
474  // (<0) intact.
475  forAll(savedFaces_, i)
476  {
477  face& f = savedFaces_[i];
478 
479  forAll(f, fp)
480  {
481  label pointi = f[fp];
482 
483  if (pointi >= 0)
484  {
485  f[fp] = map.reversePointMap()[pointi];
486 
487  if (f[fp] == -1)
488  {
490  << "Old point " << pointi
491  << " seems to have disappeared."
492  << abort(FatalError);
493  }
494  }
495  }
496  }
497 
498 
499  // DEBUG: Compare the stored faces with the current ones.
500  if (debug)
501  {
502  forAll(savedFaceLabels_, saveI)
503  {
504  if (savedFaceLabels_[saveI] >= 0)
505  {
506  const face& f = mesh_.faces()[savedFaceLabels_[saveI]];
507 
508  // Get kept points of saved faces.
509  const face& savedFace = savedFaces_[saveI];
510 
511  face keptFace(savedFace.size());
512  label keptFp = 0;
513 
514  forAll(savedFace, fp)
515  {
516  label pointi = savedFace[fp];
517 
518  if (pointi >= 0)
519  {
520  keptFace[keptFp++] = pointi;
521  }
522  }
523  keptFace.setSize(keptFp);
524 
525  // Compare as faces (since f might have rotated and
526  // face::operator== takes care of this)
527  if (keptFace != f)
528  {
530  << "facei:" << savedFaceLabels_[saveI] << nl
531  << "face:" << f << nl
532  << "keptFace:" << keptFace << nl
533  << "saved points:"
535  (
536  mesh_.points(),
537  savedPoints_,
538  savedFace
539  )() << nl
540  << abort(FatalError);
541  }
542  }
543  }
544  }
545  }
546 }
547 
548 
549 // Given list of faces to undo picks up the local indices of the faces
550 // to restore. Additionally it also picks up all the faces that use
551 // any of the deleted points.
553 (
554  const labelList& undoFaces,
555  labelList& localFaces,
556  labelList& localPoints
557 ) const
558 {
559  if (!undoable_)
560  {
562  << "removePoints not constructed with"
563  << " unrefinement capability."
564  << abort(FatalError);
565  }
566 
567  if (debug)
568  {
569  // Check if synced.
570  faceSet undoFacesSet(mesh_, "undoFacesSet", undoFaces);
571  label sz = undoFacesSet.size();
572 
573  undoFacesSet.sync(mesh_);
574  if (sz != undoFacesSet.size())
575  {
577  << "undoFaces not synchronised across coupled faces." << endl
578  << "Size before sync:" << sz
579  << " after sync:" << undoFacesSet.size()
580  << abort(FatalError);
581  }
582  }
583 
584 
585  // Problem: input undoFaces are synced. But e.g.
586  // two faces, A (uncoupled) and B(coupled) share a deleted point. A gets
587  // passed in to be restored. Now picking up the deleted point and the
588  // original faces using it picks up face B. But not its coupled neighbour!
589  // Problem is that we cannot easily synchronise the deleted points
590  // (e.g. syncPointList) since they're not in the mesh anymore - only the
591  // faces are. So we have to collect the points-to-restore as indices
592  // in the faces (which is information we can synchronise)
593 
594 
595 
596  // Mark points-to-restore
597  labelHashSet localPointsSet(undoFaces.size());
598 
599  {
600  // Create set of faces to be restored
601  labelHashSet undoFacesSet(undoFaces);
602 
603  forAll(savedFaceLabels_, saveI)
604  {
605  if (savedFaceLabels_[saveI] < 0)
606  {
608  << "Illegal face label " << savedFaceLabels_[saveI]
609  << " at index " << saveI
610  << abort(FatalError);
611  }
612 
613  if (undoFacesSet.found(savedFaceLabels_[saveI]))
614  {
615  const face& savedFace = savedFaces_[saveI];
616 
617  forAll(savedFace, fp)
618  {
619  if (savedFace[fp] < 0)
620  {
621  label savedPointi = -savedFace[fp]-1;
622 
623  if (savedPoints_[savedPointi] == vector::max)
624  {
626  << "Trying to restore point " << savedPointi
627  << " from mesh face " << savedFaceLabels_[saveI]
628  << " saved face:" << savedFace
629  << " which has already been undone."
630  << abort(FatalError);
631  }
632 
633  localPointsSet.insert(savedPointi);
634  }
635  }
636  }
637  }
638 
639 
640  // Per boundary face, per index in face whether the point needs
641  // restoring. Note that this is over all saved faces, not just over
642  // the ones in undoFaces.
643 
644  boolListList faceVertexRestore(mesh_.nBoundaryFaces());
645 
646  // Populate with my local points-to-restore.
647  forAll(savedFaces_, saveI)
648  {
649  label bFacei = savedFaceLabels_[saveI] - mesh_.nInternalFaces();
650 
651  if (bFacei >= 0)
652  {
653  const face& savedFace = savedFaces_[saveI];
654 
655  boolList& fRestore = faceVertexRestore[bFacei];
656 
657  fRestore.setSize(savedFace.size());
658  fRestore = false;
659 
660  forAll(savedFace, fp)
661  {
662  if (savedFace[fp] < 0)
663  {
664  label savedPointi = -savedFace[fp]-1;
665 
666  if (localPointsSet.found(savedPointi))
667  {
668  fRestore[fp] = true;
669  }
670  }
671  }
672  }
673  }
674 
676  (
677  mesh_,
678  faceVertexRestore,
679  faceEqOp<bool, orEqOp>(), // special operator to handle faces
680  Foam::dummyTransform() // no transformation
681  );
682 
683  // So now if any of the points-to-restore is used by any coupled face
684  // anywhere the corresponding index in faceVertexRestore will be set.
685 
686  // Now combine the localPointSet and the (sychronised)
687  // boundary-points-to-restore.
688 
689  forAll(savedFaces_, saveI)
690  {
691  label bFacei = savedFaceLabels_[saveI] - mesh_.nInternalFaces();
692 
693  if (bFacei >= 0)
694  {
695  const boolList& fRestore = faceVertexRestore[bFacei];
696 
697  const face& savedFace = savedFaces_[saveI];
698 
699  forAll(fRestore, fp)
700  {
701  // Does neighbour require point restored?
702  if (fRestore[fp])
703  {
704  if (savedFace[fp] >= 0)
705  {
707  << "Problem: on coupled face:"
708  << savedFaceLabels_[saveI]
709  << " fc:"
710  << mesh_.faceCentres()[savedFaceLabels_[saveI]]
711  << endl
712  << " my neighbour tries to restore the vertex"
713  << " at index " << fp
714  << " whereas my saved face:" << savedFace
715  << " does not indicate a deleted vertex"
716  << " at that position."
717  << abort(FatalError);
718  }
719 
720  label savedPointi = -savedFace[fp]-1;
721 
722  localPointsSet.insert(savedPointi);
723  }
724  }
725  }
726  }
727  }
728 
729  localPoints = localPointsSet.toc();
730 
731 
732  // Collect all saved faces using any localPointsSet
733  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
734 
735  labelHashSet localFacesSet(2*undoFaces.size());
736 
737  forAll(savedFaces_, saveI)
738  {
739  const face& savedFace = savedFaces_[saveI];
740 
741  forAll(savedFace, fp)
742  {
743  if (savedFace[fp] < 0)
744  {
745  label savedPointi = -savedFace[fp]-1;
746 
747  if (localPointsSet.found(savedPointi))
748  {
749  localFacesSet.insert(saveI);
750  }
751  }
752  }
753  }
754  localFaces = localFacesSet.toc();
755 
756 
757  // Note that at this point the localFaces to restore can contain points
758  // that are not going to be restored! The localFaces though will
759  // be guaranteed to be all the ones affected by the restoration of the
760  // localPoints.
761 }
762 
763 
765 (
766  const labelList& localFaces,
767  const labelList& localPoints,
768  polyTopoChange& meshMod
769 )
770 {
771  if (!undoable_)
772  {
774  << "removePoints not constructed with"
775  << " unrefinement capability."
776  << abort(FatalError);
777  }
778 
779 
780  // Per savedPoint -1 or the restored point label
781  labelList addedPoints(savedPoints_.size(), -1);
782 
783  forAll(localPoints, i)
784  {
785  label localI = localPoints[i];
786 
787  if (savedPoints_[localI] == vector::max)
788  {
790  << "Saved point " << localI << " already restored!"
791  << abort(FatalError);
792  }
793 
794  addedPoints[localI] = meshMod.setAction
795  (
797  (
798  savedPoints_[localI], // point
799  -1, // master point
800  -1, // zone for point
801  true // supports a cell
802  )
803  );
804 
805  // Mark the restored points so they are not restored again.
806  savedPoints_[localI] = vector::max;
807  }
808 
809  forAll(localFaces, i)
810  {
811  label saveI = localFaces[i];
812 
813  // Modify indices into saved points (so < 0) to point to the
814  // added points.
815  face& savedFace = savedFaces_[saveI];
816 
817  face newFace(savedFace.size());
818  label newFp = 0;
819 
820  bool hasSavedPoints = false;
821 
822  forAll(savedFace, fp)
823  {
824  if (savedFace[fp] < 0)
825  {
826  label addedPointi = addedPoints[-savedFace[fp]-1];
827 
828  if (addedPointi != -1)
829  {
830  savedFace[fp] = addedPointi;
831  newFace[newFp++] = addedPointi;
832  }
833  else
834  {
835  hasSavedPoints = true;
836  }
837  }
838  else
839  {
840  newFace[newFp++] = savedFace[fp];
841  }
842  }
843  newFace.setSize(newFp);
844 
845  modifyFace(savedFaceLabels_[saveI], newFace, meshMod);
846 
847  if (!hasSavedPoints)
848  {
849  // Face fully restored. Mark for compaction later on
850  savedFaceLabels_[saveI] = -1;
851  savedFaces_[saveI].clear();
852  }
853  }
854 
855 
856  // Compact face labels
857  label newSaveI = 0;
858 
859  forAll(savedFaceLabels_, saveI)
860  {
861  if (savedFaceLabels_[saveI] != -1)
862  {
863  if (newSaveI != saveI)
864  {
865  savedFaceLabels_[newSaveI] = savedFaceLabels_[saveI];
866  savedFaces_[newSaveI].transfer(savedFaces_[saveI]);
867  }
868  newSaveI++;
869  }
870  }
871 
872  savedFaceLabels_.setSize(newSaveI);
873  savedFaces_.setSize(newSaveI);
874 
875 
876  // Check that all faces have been restored that use any restored points
877  if (debug)
878  {
879  forAll(savedFaceLabels_, saveI)
880  {
881  const face& savedFace = savedFaces_[saveI];
882 
883  forAll(savedFace, fp)
884  {
885  if (savedFace[fp] < 0)
886  {
887  label addedPointi = addedPoints[-savedFace[fp]-1];
888 
889  if (addedPointi != -1)
890  {
892  << "Face:" << savedFaceLabels_[saveI]
893  << " savedVerts:" << savedFace
894  << " uses restored point:" << -savedFace[fp]-1
895  << " with new pointlabel:" << addedPointi
896  << abort(FatalError);
897  }
898  }
899  }
900  }
901  }
902 }
903 
904 
905 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
polyRemovePoint.H
Foam::edge::commonVertex
label commonVertex(const edge &other) const
Return vertex common with other edge or -1 on failure.
Definition: edgeI.H:172
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::polyRemovePoint
Class containing data for point removal.
Definition: polyRemovePoint.H:48
Foam::syncTools::syncBoundaryFaceList
static void syncBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop, const TransformOp &top, const bool parRun=Pstream::parRun())
Synchronize values on boundary faces only.
Definition: syncToolsTemplates.C:998
polyTopoChange.H
Foam::andEqOp
Definition: ops.H:85
Foam::polyTopoChange
Direct mesh changes based on v1.3 polyTopoChange syntax.
Definition: polyTopoChange.H:99
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
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:553
dummyTransform.H
Dummy transform to be used with syncTools.
Foam::Map< label >
Foam::faceSet
A list of face labels.
Definition: faceSet.H:51
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::faceEqOp::operator()
void operator()(List< T > &x, const List< T > &y) const
Definition: removePoints.C:56
Foam::dummyTransform
Definition: dummyTransform.H:47
polyMesh.H
Foam::HashSet< label, Hash< label > >
syncTools.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::syncTools::syncPointList
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
Definition: syncToolsTemplates.C:721
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
BiIndirectList.H
Foam::faceSet::sync
virtual void sync(const polyMesh &mesh)
Sync faceSet across coupled patches.
Definition: faceSet.C:130
removePoints.H
Foam::Field< vector >
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1107
faceSet.H
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:486
Foam::polyBoundaryMesh::whichPatch
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Definition: polyBoundaryMesh.C:812
Foam::BiIndirectList
Indexes into negList (negative index) or posList (zero or positive index).
Definition: BiIndirectList.H:52
polyAddPoint.H
Foam::removePoints::setUnrefinement
void setUnrefinement(const labelList &localFaces, const labelList &localPoints, polyTopoChange &)
Restore selected faces and vertices.
Definition: removePoints.C:765
Foam::ZoneMesh::whichZone
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: ZoneMesh.C:289
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::faceEqOp
Definition: removePoints.C:51
zoneID
const labelIOList & zoneID
Definition: interpolatedFaces.H:22
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
PstreamReduceOps.H
Inter-processor communication reduction functions.
polyModifyFace.H
Foam::edge::otherVertex
label otherVertex(const label pointLabel) const
Given the point label for one vertex, return the other one.
Definition: edgeI.H:188
Foam::mapPolyMesh::reverseFaceMap
const labelList & reverseFaceMap() const
Reverse face map.
Definition: mapPolyMesh.H:501
Foam::normalised
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:487
Foam::primitiveMesh::isInternalFace
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
Definition: primitiveMeshI.H:103
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::polyTopoChange::setAction
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.
Definition: polyTopoChange.C:2481
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::max
static const Vector< scalar > max
Definition: VectorSpace.H:117
f
labelList f(nPoints)
Foam::mapPolyMesh::reversePointMap
const labelList & reversePointMap() const
Reverse point map.
Definition: mapPolyMesh.H:469
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:116
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::polyAddPoint
Class containing data for point addition.
Definition: polyAddPoint.H:49
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::UIndirectList
A List with indirect addressing.
Definition: faMatrix.H:60
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::removePoints::updateMesh
void updateMesh(const mapPolyMesh &)
Force recalculation of locally stored data on topological change.
Definition: removePoints.C:451
Foam::removePoints::setRefinement
void setRefinement(const boolList &, polyTopoChange &)
Play commands into polyTopoChange to remove points. Gets.
Definition: removePoints.C:296
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:148
pFaces
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::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].reset(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::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1113
y
scalar y
Definition: LISASMDCalcMethod1.H:14