isoSurfaceCell.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) 2016-2022 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 "isoSurfaceCell.H"
30 #include "isoSurfacePoint.H"
31 #include "dictionary.H"
32 #include "polyMesh.H"
33 #include "mergePoints.H"
34 #include "tetMatcher.H"
35 #include "syncTools.H"
36 #include "triSurface.H"
37 #include "triSurfaceTools.H"
38 #include "Time.H"
39 #include "triangle.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 #include "isoSurfaceBaseMethods.H"
45 
46 
47 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51  defineTypeNameAndDebug(isoSurfaceCell, 0);
52 }
53 
54 
55 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
56 
57 Foam::scalar Foam::isoSurfaceCell::isoFraction
58 (
59  const scalar s0,
60  const scalar s1
61 ) const
62 {
63  const scalar d = s1-s0;
64 
65  if (mag(d) > VSMALL)
66  {
67  return (iso_-s0)/d;
68  }
69 
70  return -1.0;
71 }
72 
73 
74 Foam::labelPair Foam::isoSurfaceCell::findCommonPoints
75 (
76  const labelledTri& tri0,
77  const labelledTri& tri1
78 )
79 {
80  labelPair common(-1, -1);
81 
82  label fp0 = 0;
83  label fp1 = tri1.find(tri0[fp0]);
84 
85  if (fp1 == -1)
86  {
87  fp0 = 1;
88  fp1 = tri1.find(tri0[fp0]);
89  }
90 
91  if (fp1 != -1)
92  {
93  // So tri0[fp0] is tri1[fp1]
94 
95  // Find next common point
96  label fp0p1 = tri0.fcIndex(fp0);
97  label fp1p1 = tri1.fcIndex(fp1);
98  label fp1m1 = tri1.rcIndex(fp1);
99 
100  if (tri0[fp0p1] == tri1[fp1p1] || tri0[fp0p1] == tri1[fp1m1])
101  {
102  common[0] = tri0[fp0];
103  common[1] = tri0[fp0p1];
104  }
105  }
106  return common;
107 }
108 
109 
110 Foam::point Foam::isoSurfaceCell::calcCentre(const triSurface& s)
111 {
112  vector sum = Zero;
113 
114  forAll(s, i)
115  {
116  sum += s[i].centre(s.points());
117  }
118  return sum/s.size();
119 }
120 
121 
122 Foam::pointIndexHit Foam::isoSurfaceCell::collapseSurface
123 (
124  const label celli,
125  pointField& localPoints,
126  DynamicList<labelledTri, 64>& localTris
127 ) const
128 {
129  pointIndexHit info(false, Zero, localTris.size());
130 
131  if (localTris.size() == 1)
132  {
133  const labelledTri& tri = localTris[0];
134  info.hitPoint(tri.centre(localPoints));
135  }
136  else if (localTris.size() == 2)
137  {
138  // Check if the two triangles share an edge.
139  const labelledTri& tri0 = localTris[0];
140  const labelledTri& tri1 = localTris[1];
141 
142  labelPair shared = findCommonPoints(tri0, tri1);
143 
144  if (shared[0] != -1)
145  {
146  const vector n0 = tri0.areaNormal(localPoints);
147  const vector n1 = tri1.areaNormal(localPoints);
148 
149  // Merge any zero-sized triangles,
150  // or if they point in the same direction.
151 
152  if
153  (
154  mag(n0) <= ROOTVSMALL
155  || mag(n1) <= ROOTVSMALL
156  || (n0 & n1) >= 0
157  )
158  {
159  info.hitPoint
160  (
161  0.5
162  * (
163  tri0.centre(localPoints)
164  + tri1.centre(localPoints)
165  )
166  );
167  }
168  }
169  }
170  else if (localTris.size())
171  {
172  // Check if single region. Rare situation.
173  triSurface surf
174  (
175  localTris,
177  localPoints,
178  true
179  );
180  localTris.clearStorage();
181 
182  labelList faceZone;
183  label nZones = surf.markZones
184  (
185  boolList(surf.nEdges(), false),
186  faceZone
187  );
188 
189  if (nZones == 1)
190  {
191  // Check that all normals make a decent angle
192  scalar minCos = GREAT;
193  const vector& n0 = surf.faceNormals()[0];
194  for (label i = 1; i < surf.size(); ++i)
195  {
196  scalar cosAngle = (n0 & surf.faceNormals()[i]);
197  if (cosAngle < minCos)
198  {
199  minCos = cosAngle;
200  }
201  }
202 
203  if (minCos > 0)
204  {
205  info.hitPoint(calcCentre(surf));
206  }
207  }
208  }
209 
210  return info;
211 }
212 
213 
214 void Foam::isoSurfaceCell::calcSnappedCc
215 (
216  const bitSet& isTet,
217  const scalarField& cVals,
218  const scalarField& pVals,
219 
220  DynamicList<point>& snappedPoints,
221  labelList& snappedCc
222 ) const
223 {
224  const pointField& cc = mesh_.cellCentres();
225  const pointField& pts = mesh_.points();
226 
227  snappedCc.setSize(mesh_.nCells());
228  snappedCc = -1;
229 
230  // Work arrays
231  DynamicList<point, 64> localPoints(64);
232  DynamicList<labelledTri, 64> localTris(64);
233  Map<label> pointToLocal(64);
234 
235  forAll(mesh_.cells(), celli)
236  {
237  if (cellCutType_[celli] == cutType::CUT) // No tet cuts
238  {
239  const scalar cVal = cVals[celli];
240 
241  const cell& cFaces = mesh_.cells()[celli];
242 
243  localPoints.clear();
244  localTris.clear();
245  pointToLocal.clear();
246 
247  // Create points for all intersections close to cell centre
248  // (i.e. from pyramid edges)
249 
250  for (const label facei : cFaces)
251  {
252  const face& f = mesh_.faces()[facei];
253 
254  for (const label pointi : f)
255  {
256  scalar s = isoFraction(cVal, pVals[pointi]);
257 
258  if (s >= 0.0 && s <= 0.5)
259  {
260  if (pointToLocal.insert(pointi, localPoints.size()))
261  {
262  localPoints.append((1.0-s)*cc[celli]+s*pts[pointi]);
263  }
264  }
265  }
266  }
267 
268  if (localPoints.size() == 1)
269  {
270  // No need for any analysis.
271  snappedCc[celli] = snappedPoints.size();
272  snappedPoints.append(localPoints[0]);
273 
274  //Pout<< "cell:" << celli
275  // << " at " << mesh_.cellCentres()[celli]
276  // << " collapsing " << localPoints
277  // << " intersections down to "
278  // << snappedPoints[snappedCc[celli]] << endl;
279  }
280  else if (localPoints.size() == 2)
281  {
282  //? No need for any analysis.???
283  snappedCc[celli] = snappedPoints.size();
284  snappedPoints.append(0.5*(localPoints[0]+localPoints[1]));
285 
286  //Pout<< "cell:" << celli
287  // << " at " << mesh_.cellCentres()[celli]
288  // << " collapsing " << localPoints
289  // << " intersections down to "
290  // << snappedPoints[snappedCc[celli]] << endl;
291  }
292  else if (localPoints.size())
293  {
294  // Need to analyse
295  for (const label facei : cFaces)
296  {
297  const face& f = mesh_.faces()[facei];
298 
299  // Do a tetrahedralisation. Each face to cc becomes pyr.
300  // Each pyr gets split into tets by diagonalisation
301  // of face.
302 
303  const label fp0 = mesh_.tetBasePtIs()[facei];
304  label fp = f.fcIndex(fp0);
305  for (label i = 2; i < f.size(); ++i)
306  {
307  label nextFp = f.fcIndex(fp);
308  triFace tri(f[fp0], f[fp], f[nextFp]);
309 
310  // Get fractions for the three edges to cell centre
311  FixedList<scalar, 3> s(3);
312  s[0] = isoFraction(cVal, pVals[tri[0]]);
313  s[1] = isoFraction(cVal, pVals[tri[1]]);
314  s[2] = isoFraction(cVal, pVals[tri[2]]);
315 
316  if
317  (
318  (s[0] >= 0.0 && s[0] <= 0.5)
319  && (s[1] >= 0.0 && s[1] <= 0.5)
320  && (s[2] >= 0.0 && s[2] <= 0.5)
321  )
322  {
323  if
324  (
325  (mesh_.faceOwner()[facei] == celli)
326  == (cVal >= pVals[tri[0]])
327  )
328  {
329  localTris.append
330  (
331  labelledTri
332  (
333  pointToLocal[tri[1]],
334  pointToLocal[tri[0]],
335  pointToLocal[tri[2]],
336  0
337  )
338  );
339  }
340  else
341  {
342  localTris.append
343  (
344  labelledTri
345  (
346  pointToLocal[tri[0]],
347  pointToLocal[tri[1]],
348  pointToLocal[tri[2]],
349  0
350  )
351  );
352  }
353  }
354 
355  fp = nextFp;
356  }
357  }
358 
359  pointField surfPoints;
360  surfPoints.transfer(localPoints);
361  pointIndexHit info = collapseSurface
362  (
363  celli,
364  surfPoints,
365  localTris
366  );
367 
368  if (info.hit())
369  {
370  snappedCc[celli] = snappedPoints.size();
371  snappedPoints.append(info.point());
372 
373  //Pout<< "cell:" << celli
374  // << " at " << mesh_.cellCentres()[celli]
375  // << " collapsing " << surfPoints
376  // << " intersections down to "
377  // << snappedPoints[snappedCc[celli]] << endl;
378  }
379  }
380  }
381  }
382 }
383 
384 
385 void Foam::isoSurfaceCell::genPointTris
386 (
387  const scalarField& cellValues,
388  const scalarField& pointValues,
389  const label pointi,
390  const label facei,
391  const label celli,
392  DynamicList<point, 64>& localTriPoints
393 ) const
394 {
395  const pointField& cc = mesh_.cellCentres();
396  const pointField& pts = mesh_.points();
397  const face& f = mesh_.faces()[facei];
398 
399  const label fp0 = mesh_.tetBasePtIs()[facei];
400  label fp = f.fcIndex(fp0);
401  for (label i = 2; i < f.size(); ++i)
402  {
403  label nextFp = f.fcIndex(fp);
404  triFace tri(f[fp0], f[fp], f[nextFp]);
405 
406  label index = tri.find(pointi);
407 
408  if (index == -1)
409  {
410  continue;
411  }
412 
413  // Tet between index..index-1, index..index+1, index..cc
414  label b = tri[tri.fcIndex(index)];
415  label c = tri[tri.rcIndex(index)];
416 
417  // Get fractions for the three edges emanating from point
418  FixedList<scalar, 3> s(3);
419  s[0] = isoFraction(pointValues[pointi], pointValues[b]);
420  s[1] = isoFraction(pointValues[pointi], pointValues[c]);
421  s[2] = isoFraction(pointValues[pointi], cellValues[celli]);
422 
423  if
424  (
425  (s[0] >= 0.0 && s[0] <= 0.5)
426  && (s[1] >= 0.0 && s[1] <= 0.5)
427  && (s[2] >= 0.0 && s[2] <= 0.5)
428  )
429  {
430  point p0 = (1.0-s[0])*pts[pointi] + s[0]*pts[b];
431  point p1 = (1.0-s[1])*pts[pointi] + s[1]*pts[c];
432  point p2 = (1.0-s[2])*pts[pointi] + s[2]*cc[celli];
433 
434  if
435  (
436  (mesh_.faceOwner()[facei] == celli)
437  == (pointValues[pointi] > cellValues[celli])
438  )
439  {
440  localTriPoints.append(p0);
441  localTriPoints.append(p1);
442  localTriPoints.append(p2);
443  }
444  else
445  {
446  localTriPoints.append(p1);
447  localTriPoints.append(p0);
448  localTriPoints.append(p2);
449  }
450  }
451 
452  fp = nextFp;
453  }
454 }
455 
456 
457 void Foam::isoSurfaceCell::genPointTris
458 (
459  const scalarField& pointValues,
460  const label pointi,
461  const label facei,
462  const label celli,
463  DynamicList<point, 64>& localTriPoints
464 ) const
465 {
466  const pointField& pts = mesh_.points();
467  const cell& cFaces = mesh_.cells()[celli];
468 
469  // Make tet from this face to the 4th point (same as cellcentre in
470  // non-tet cells)
471  const face& f = mesh_.faces()[facei];
472 
473  // Find 4th point
474  label ccPointi = -1;
475  for (const label cfacei : cFaces)
476  {
477  const face& f1 = mesh_.faces()[cfacei];
478  for (const label p1 : f1)
479  {
480  if (!f.found(p1))
481  {
482  ccPointi = p1;
483  break;
484  }
485  }
486  if (ccPointi != -1)
487  {
488  break;
489  }
490  }
491 
492 
493  // Tet between index..index-1, index..index+1, index..cc
494  label index = f.find(pointi);
495  label b = f[f.fcIndex(index)];
496  label c = f[f.rcIndex(index)];
497 
498  //Pout<< " p0:" << pointi << " b:" << b << " c:" << c
499  //<< " d:" << ccPointi << endl;
500 
501  // Get fractions for the three edges emanating from point
502  FixedList<scalar, 3> s(3);
503  s[0] = isoFraction(pointValues[pointi], pointValues[b]);
504  s[1] = isoFraction(pointValues[pointi], pointValues[c]);
505  s[2] = isoFraction(pointValues[pointi], pointValues[ccPointi]);
506 
507  if
508  (
509  (s[0] >= 0.0 && s[0] <= 0.5)
510  && (s[1] >= 0.0 && s[1] <= 0.5)
511  && (s[2] >= 0.0 && s[2] <= 0.5)
512  )
513  {
514  point p0 = (1.0-s[0])*pts[pointi] + s[0]*pts[b];
515  point p1 = (1.0-s[1])*pts[pointi] + s[1]*pts[c];
516  point p2 = (1.0-s[2])*pts[pointi] + s[2]*pts[ccPointi];
517 
518  if (mesh_.faceOwner()[facei] != celli)
519  {
520  localTriPoints.append(p0);
521  localTriPoints.append(p1);
522  localTriPoints.append(p2);
523  }
524  else
525  {
526  localTriPoints.append(p1);
527  localTriPoints.append(p0);
528  localTriPoints.append(p2);
529  }
530  }
531 }
532 
533 
534 void Foam::isoSurfaceCell::calcSnappedPoint
535 (
536  const bitSet& isTet,
537  const scalarField& cVals,
538  const scalarField& pVals,
539 
540  DynamicList<point>& snappedPoints,
541  labelList& snappedPoint
542 ) const
543 {
544  const labelList& faceOwn = mesh_.faceOwner();
545  const labelList& faceNei = mesh_.faceNeighbour();
546 
547  // Determine if point is on boundary. Points on boundaries are never
548  // snapped. Coupled boundaries are handled explicitly so not marked here.
549  bitSet isBoundaryPoint(mesh_.nPoints());
550  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
551 
552  for (const polyPatch& pp : patches)
553  {
554  if (!pp.coupled())
555  {
556  for (const label facei : pp.range())
557  {
558  const face& f = mesh_.faces()[facei];
559 
560  isBoundaryPoint.set(f);
561  }
562  }
563  }
564 
565 
566  const point greatPoint(GREAT, GREAT, GREAT);
567 
568  pointField collapsedPoint(mesh_.nPoints(), greatPoint);
569 
570 
571  // Work arrays
572  DynamicList<point, 64> localTriPoints(100);
573  labelHashSet localPointCells(100);
574 
575  forAll(mesh_.pointFaces(), pointi)
576  {
577  constexpr uint8_t realCut(cutType::CUT | cutType::TETCUT);
578 
579  if (isBoundaryPoint.test(pointi))
580  {
581  continue;
582  }
583 
584  const labelList& pFaces = mesh_.pointFaces()[pointi];
585 
586  bool anyCut = false;
587 
588  for (const label facei : pFaces)
589  {
590  if
591  (
592  (cellCutType_[faceOwn[facei]] & realCut) != 0
593  || (
594  mesh_.isInternalFace(facei)
595  && (cellCutType_[faceNei[facei]] & realCut) != 0
596  )
597  )
598  {
599  anyCut = true;
600  break;
601  }
602  }
603 
604  if (!anyCut)
605  {
606  continue;
607  }
608 
609 
610  // Do a pointCells walk (by using pointFaces)
611 
612  localPointCells.clear();
613  localTriPoints.clear();
614 
615  for (const label facei : pFaces)
616  {
617  const label own = faceOwn[facei];
618 
619  if (isTet.test(own))
620  {
621  // Since tets have no cell centre to include make sure
622  // we only generate a triangle once per point.
623  if (localPointCells.insert(own))
624  {
625  genPointTris(pVals, pointi, facei, own, localTriPoints);
626  }
627  }
628  else
629  {
630  genPointTris
631  (
632  cVals,
633  pVals,
634  pointi,
635  facei,
636  own,
637  localTriPoints
638  );
639  }
640 
641  if (mesh_.isInternalFace(facei))
642  {
643  const label nei = faceNei[facei];
644 
645  if (isTet.test(nei))
646  {
647  if (localPointCells.insert(nei))
648  {
649  genPointTris(pVals, pointi, facei, nei, localTriPoints);
650  }
651  }
652  else
653  {
654  genPointTris
655  (
656  cVals,
657  pVals,
658  pointi,
659  facei,
660  nei,
661  localTriPoints
662  );
663  }
664  }
665  }
666 
667  if (localTriPoints.size() == 3)
668  {
669  // Single triangle. No need for any analysis. Average points.
671  points.transfer(localTriPoints);
672  collapsedPoint[pointi] = sum(points)/points.size();
673 
674  //Pout<< " point:" << pointi
675  // << " replacing coord:" << mesh_.points()[pointi]
676  // << " by average:" << collapsedPoint[pointi] << endl;
677  }
678  else if (localTriPoints.size())
679  {
680  // Convert points into triSurface.
681 
682  // Merge points and compact out non-valid triangles
683  labelList triMap; // merged to unmerged triangle
684  labelList triPointReverseMap; // unmerged to merged point
685  triSurface surf
686  (
687  stitchTriPoints
688  (
689  false, // do not check for duplicate tris
690  localTriPoints,
691  triPointReverseMap,
692  triMap
693  )
694  );
695 
696  labelList faceZone;
697  label nZones = surf.markZones
698  (
699  boolList(surf.nEdges(), false),
700  faceZone
701  );
702 
703  if (nZones == 1)
704  {
705  // Check that all normals make a decent angle
706  scalar minCos = GREAT;
707  const vector& n0 = surf.faceNormals()[0];
708  for (label i = 1; i < surf.size(); ++i)
709  {
710  const vector& n = surf.faceNormals()[i];
711  scalar cosAngle = (n0 & n);
712  if (cosAngle < minCos)
713  {
714  minCos = cosAngle;
715  }
716  }
717  if (minCos > 0)
718  {
719  collapsedPoint[pointi] = calcCentre(surf);
720  }
721  }
722  }
723  }
724 
726  (
727  mesh_,
728  collapsedPoint,
729  minMagSqrEqOp<point>(),
730  greatPoint
731  );
732 
733  snappedPoint.setSize(mesh_.nPoints());
734  snappedPoint = -1;
735 
736  forAll(collapsedPoint, pointi)
737  {
738  // Cannot do == comparison since might be transformed so have
739  // truncation errors.
740  if (magSqr(collapsedPoint[pointi]) < 0.5*magSqr(greatPoint))
741  {
742  snappedPoint[pointi] = snappedPoints.size();
743  snappedPoints.append(collapsedPoint[pointi]);
744  }
745  }
746 }
747 
748 
749 Foam::triSurface Foam::isoSurfaceCell::stitchTriPoints
750 (
751  const bool checkDuplicates,
752  const List<point>& triPoints,
753  labelList& triPointReverseMap, // unmerged to merged point
754  labelList& triMap // merged to unmerged triangle
755 ) const
756 {
757  label nTris = triPoints.size()/3;
758 
759  if ((triPoints.size() % 3) != 0)
760  {
762  << "Problem: number of points " << triPoints.size()
763  << " not a multiple of 3." << abort(FatalError);
764  }
765 
766  pointField newPoints;
768  (
769  triPoints,
770  mergeDistance_,
771  false,
772  triPointReverseMap,
773  newPoints
774  );
775 
776  // Check that enough merged.
777  if (debug)
778  {
779  Pout<< "isoSurfaceCell : merged from " << triPointReverseMap.size()
780  << " points down to " << newPoints.size() << endl;
781 
782  labelList oldToNew;
783  const label nUnique = Foam::mergePoints
784  (
785  newPoints,
786  mergeDistance_,
787  true,
788  oldToNew
789  );
790 
791  if (nUnique != newPoints.size())
792  {
794  << "Merged points contain duplicates"
795  << " when merging with distance " << mergeDistance_ << endl
796  << "merged:" << newPoints.size() << " re-merged:"
797  << nUnique
798  << abort(FatalError);
799  }
800  }
801 
802 
803  List<labelledTri> tris;
804  {
805  DynamicList<labelledTri> dynTris(nTris);
806  label rawPointi = 0;
807  DynamicList<label> newToOldTri(nTris);
808 
809  for (label oldTriI = 0; oldTriI < nTris; ++oldTriI)
810  {
811  labelledTri tri
812  (
813  triPointReverseMap[rawPointi],
814  triPointReverseMap[rawPointi+1],
815  triPointReverseMap[rawPointi+2],
816  0
817  );
818  if (tri.good())
819  {
820  newToOldTri.append(oldTriI);
821  dynTris.append(tri);
822  }
823 
824  rawPointi += 3;
825  }
826 
827  triMap.transfer(newToOldTri);
828  tris.transfer(dynTris);
829  }
830 
831 
832  // Use face centres to determine 'flat hole' situation (see RMT paper).
833  // Two unconnected triangles get connected because (some of) the edges
834  // separating them get collapsed. Below only checks for duplicate triangles,
835  // not non-manifold edge connectivity.
836  if (checkDuplicates)
837  {
838  if (debug)
839  {
840  Pout<< "isoSurfaceCell : merged from " << nTris
841  << " down to " << tris.size() << " triangles." << endl;
842  }
843 
844  pointField centres(tris.size());
845  forAll(tris, triI)
846  {
847  centres[triI] = tris[triI].centre(newPoints);
848  }
849 
850  labelList oldToMerged;
851  label nUnique = Foam::mergePoints
852  (
853  centres,
854  mergeDistance_,
855  false,
856  oldToMerged
857  );
858 
859  if (debug)
860  {
861  Pout<< "isoSurfaceCell : detected "
862  << (oldToMerged.size() - nUnique)
863  << " duplicate triangles." << endl;
864  }
865 
866  if (oldToMerged.size() != nUnique)
867  {
868  // Filter out duplicates.
869  label newTriI = 0;
870  DynamicList<label> newToOldTri(tris.size());
871  labelList newToMaster(nUnique, -1);
872  forAll(tris, triI)
873  {
874  label mergedI = oldToMerged[triI];
875 
876  if (newToMaster[mergedI] == -1)
877  {
878  newToMaster[mergedI] = triI;
879  newToOldTri.append(triMap[triI]);
880  tris[newTriI++] = tris[triI];
881  }
882  }
883 
884  triMap.transfer(newToOldTri);
885  tris.setSize(newTriI);
886  }
887  }
888 
889  return triSurface(tris, geometricSurfacePatchList(), newPoints, true);
890 }
891 
892 
893 void Foam::isoSurfaceCell::calcAddressing
894 (
895  const triSurface& surf,
896  List<FixedList<label, 3>>& faceEdges,
897  labelList& edgeFace0,
898  labelList& edgeFace1,
899  Map<labelList>& edgeFacesRest
900 ) const
901 {
902  const pointField& points = surf.points();
903 
904  pointField edgeCentres(3*surf.size());
905  label edgeI = 0;
906  forAll(surf, triI)
907  {
908  const labelledTri& tri = surf[triI];
909  edgeCentres[edgeI++] = 0.5*(points[tri[0]]+points[tri[1]]);
910  edgeCentres[edgeI++] = 0.5*(points[tri[1]]+points[tri[2]]);
911  edgeCentres[edgeI++] = 0.5*(points[tri[2]]+points[tri[0]]);
912  }
913 
914  labelList oldToMerged;
915  const label nUnique = Foam::mergePoints
916  (
917  edgeCentres,
918  mergeDistance_,
919  false,
920  oldToMerged
921  );
922 
923  if (debug)
924  {
925  Pout<< "isoSurfaceCell : detected "
926  << nUnique
927  << " edges on " << surf.size() << " triangles." << endl;
928  }
929 
930  if (nUnique == edgeCentres.size())
931  {
932  // Nothing to do
933  return;
934  }
935 
936 
937  // Determine faceEdges
938  faceEdges.setSize(surf.size());
939  edgeI = 0;
940  forAll(surf, triI)
941  {
942  faceEdges[triI][0] = oldToMerged[edgeI++];
943  faceEdges[triI][1] = oldToMerged[edgeI++];
944  faceEdges[triI][2] = oldToMerged[edgeI++];
945  }
946 
947 
948  // Determine edgeFaces
949  edgeFace0.resize_nocopy(nUnique);
950  edgeFace0 = -1;
951  edgeFace1.resize_nocopy(nUnique);
952  edgeFace1 = -1;
953  edgeFacesRest.clear();
954 
955  forAll(oldToMerged, oldEdgeI)
956  {
957  label triI = oldEdgeI / 3;
958  label edgeI = oldToMerged[oldEdgeI];
959 
960  if (edgeFace0[edgeI] == -1)
961  {
962  edgeFace0[edgeI] = triI;
963  }
964  else if (edgeFace1[edgeI] == -1)
965  {
966  edgeFace1[edgeI] = triI;
967  }
968  else
969  {
970  //WarningInFunction
971  // << "Edge " << edgeI << " with centre "
972  // << " used by more than two triangles: " << edgeFace0[edgeI]
973  // << ", "
974  // << edgeFace1[edgeI] << " and " << triI << endl;
975 
976  edgeFacesRest(edgeI).push_back(triI);
977  }
978  }
979 }
980 
981 
982 bool Foam::isoSurfaceCell::danglingTriangle
983 (
984  const FixedList<label, 3>& fEdges,
985  const labelList& edgeFace1
986 )
987 {
988  label nOpen = 0;
989  for (const label edgei : fEdges)
990  {
991  if (edgeFace1[edgei] == -1)
992  {
993  ++nOpen;
994  }
995  }
996 
997  return (nOpen == 1 || nOpen == 2 || nOpen == 3);
998 }
999 
1000 
1001 Foam::label Foam::isoSurfaceCell::markDanglingTriangles
1002 (
1003  const List<FixedList<label, 3>>& faceEdges,
1004  const labelList& edgeFace0,
1005  const labelList& edgeFace1,
1006  const Map<labelList>& edgeFacesRest,
1007  boolList& keepTriangles
1008 )
1009 {
1010  keepTriangles.setSize(faceEdges.size());
1011  keepTriangles = true;
1012 
1013  label nDangling = 0;
1014 
1015  // Remove any dangling triangles
1016  forAllConstIters(edgeFacesRest, iter)
1017  {
1018  // These are all the non-manifold edges. Filter out all triangles
1019  // with only one connected edge (= this edge)
1020 
1021  const label edgeI = iter.key();
1022  const labelList& otherEdgeFaces = iter.val();
1023 
1024  // Remove all dangling triangles
1025  if (danglingTriangle(faceEdges[edgeFace0[edgeI]], edgeFace1))
1026  {
1027  keepTriangles[edgeFace0[edgeI]] = false;
1028  ++nDangling;
1029  }
1030  if (danglingTriangle(faceEdges[edgeFace1[edgeI]], edgeFace1))
1031  {
1032  keepTriangles[edgeFace1[edgeI]] = false;
1033  ++nDangling;
1034  }
1035  for (const label triI : otherEdgeFaces)
1036  {
1037  if (danglingTriangle(faceEdges[triI], edgeFace1))
1038  {
1039  keepTriangles[triI] = false;
1040  ++nDangling;
1041  }
1042  }
1043  }
1044  return nDangling;
1045 }
1046 
1047 
1048 Foam::triSurface Foam::isoSurfaceCell::subsetMesh
1049 (
1050  const triSurface& s,
1051  const labelList& newToOldFaces,
1052  labelList& oldToNewPoints,
1053  labelList& newToOldPoints
1054 )
1055 {
1056  const boolList include
1057  (
1058  ListOps::createWithValue<bool>(s.size(), newToOldFaces, true, false)
1059  );
1060 
1061  newToOldPoints.setSize(s.points().size());
1062  oldToNewPoints.setSize(s.points().size());
1063  oldToNewPoints = -1;
1064  {
1065  label pointi = 0;
1066 
1067  forAll(include, oldFacei)
1068  {
1069  if (include[oldFacei])
1070  {
1071  // Renumber labels for face
1072  for (const label oldPointi : s[oldFacei])
1073  {
1074  if (oldToNewPoints[oldPointi] == -1)
1075  {
1076  oldToNewPoints[oldPointi] = pointi;
1077  newToOldPoints[pointi++] = oldPointi;
1078  }
1079  }
1080  }
1081  }
1082  newToOldPoints.setSize(pointi);
1083  }
1084 
1085  // Extract points
1086  pointField newPoints(newToOldPoints.size());
1087  forAll(newToOldPoints, i)
1088  {
1089  newPoints[i] = s.points()[newToOldPoints[i]];
1090  }
1091  // Extract faces
1092  List<labelledTri> newTriangles(newToOldFaces.size());
1093 
1094  forAll(newToOldFaces, i)
1095  {
1096  // Get old vertex labels
1097  const labelledTri& tri = s[newToOldFaces[i]];
1098 
1099  newTriangles[i][0] = oldToNewPoints[tri[0]];
1100  newTriangles[i][1] = oldToNewPoints[tri[1]];
1101  newTriangles[i][2] = oldToNewPoints[tri[2]];
1102  newTriangles[i].region() = tri.region();
1103  }
1104 
1105  // Reuse storage.
1106  return triSurface(newTriangles, s.patches(), newPoints, true);
1107 }
1108 
1109 
1110 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1111 
1113 (
1114  const polyMesh& mesh,
1115  const scalarField& cellValues,
1116  const scalarField& pointValues,
1117  const scalar iso,
1118  const isoSurfaceParams& params,
1119  const bitSet& ignoreCells
1120 )
1121 :
1122  isoSurfaceBase(mesh, cellValues, pointValues, iso, params),
1123  mergeDistance_(params.mergeTol()*mesh.bounds().mag()),
1124  cellCutType_(mesh.nCells(), cutType::UNVISITED)
1125 {
1126  const bool regularise = (params.filter() != filterType::NONE);
1127 
1128  if (debug)
1129  {
1130  Pout<< "isoSurfaceCell:" << nl
1131  << " cell min/max : " << minMax(cVals_) << nl
1132  << " point min/max : " << minMax(pVals_) << nl
1133  << " isoValue : " << iso << nl
1134  << " filter : " << Switch(regularise) << nl
1135  << " mergeTol : " << params.mergeTol() << nl
1136  << " mesh span : " << mesh.bounds().mag() << nl
1137  << " mergeDistance : " << mergeDistance_ << nl
1138  << " ignoreCells : " << ignoreCells.count()
1139  << " / " << cVals_.size() << nl
1140  << endl;
1141  }
1142 
1143 
1144  label nBlockedCells = 0;
1145 
1146  // Mark ignoreCells as blocked
1147  nBlockedCells += blockCells(cellCutType_, ignoreCells);
1148 
1149 
1150  // Some processor domains may require tetBasePtIs and others do not.
1151  // Do now to ensure things stay synchronized.
1152  (void)mesh_.tetBasePtIs();
1153 
1154 
1155  // Calculate a tet/non-tet filter
1156  bitSet isTet(mesh_.nCells());
1157  {
1158  for (label celli = 0; celli < mesh_.nCells(); ++celli)
1159  {
1160  if (tetMatcher::test(mesh_, celli))
1161  {
1162  isTet.set(celli);
1163  }
1164  }
1165  }
1166 
1167  // Determine cell cuts
1168  nCutCells_ = calcCellCuts(cellCutType_);
1169 
1170  if (debug)
1171  {
1172  Pout<< "isoSurfaceCell : candidate cells cut "
1173  << nCutCells_
1174  << " blocked " << nBlockedCells
1175  << " total " << mesh_.nCells() << endl;
1176  }
1177 
1178  if (debug && isA<fvMesh>(mesh))
1179  {
1180  const auto& fvmesh = refCast<const fvMesh>(mesh);
1181 
1182  volScalarField debugField
1183  (
1184  IOobject
1185  (
1186  "isoSurfaceCell.cutType",
1187  fvmesh.time().timeName(),
1188  fvmesh.time(),
1192  ),
1193  fvmesh,
1195  );
1196 
1197  auto& debugFld = debugField.primitiveFieldRef();
1198 
1199  forAll(cellCutType_, celli)
1200  {
1201  debugFld[celli] = cellCutType_[celli];
1202  }
1203 
1204  Pout<< "Writing cut types:"
1205  << debugField.objectPath() << endl;
1206 
1207  debugField.write();
1208  }
1209 
1210 
1211  DynamicList<point> snappedPoints(nCutCells_);
1212 
1213  // Per cc -1 or a point inside snappedPoints.
1214  labelList snappedCc;
1215  if (regularise)
1216  {
1217  calcSnappedCc
1218  (
1219  isTet,
1220  cellValues,
1221  pointValues,
1222  snappedPoints,
1223  snappedCc
1224  );
1225  }
1226  else
1227  {
1228  snappedCc.setSize(mesh_.nCells());
1229  snappedCc = -1;
1230  }
1231 
1232  if (debug)
1233  {
1234  Pout<< "isoSurfaceCell : shifted " << snappedPoints.size()
1235  << " cell centres to intersection." << endl;
1236  }
1237 
1238  snappedPoints.shrink();
1239  label nCellSnaps = snappedPoints.size();
1240 
1241  // Per point -1 or a point inside snappedPoints.
1242  labelList snappedPoint;
1243  if (regularise)
1244  {
1245  calcSnappedPoint
1246  (
1247  isTet,
1248  cellValues,
1249  pointValues,
1250  snappedPoints,
1251  snappedPoint
1252  );
1253  }
1254  else
1255  {
1256  snappedPoint.setSize(mesh_.nPoints());
1257  snappedPoint = -1;
1258  }
1259 
1260  if (debug)
1261  {
1262  Pout<< "isoSurfaceCell : shifted " << snappedPoints.size()-nCellSnaps
1263  << " vertices to intersection." << endl;
1264  }
1265 
1266 
1267  // Use a triSurface as a temporary for various operations
1268  triSurface tmpsurf;
1269 
1270  {
1271  DynamicList<point> triPoints(nCutCells_);
1272  DynamicList<label> triMeshCells(nCutCells_);
1273 
1274  generateTriPoints
1275  (
1276  cellValues,
1277  pointValues,
1278 
1279  mesh_.cellCentres(),
1280  mesh_.points(),
1281 
1282  snappedPoints,
1283  snappedCc,
1284  snappedPoint,
1285 
1286  triPoints,
1287  triMeshCells
1288  );
1289 
1290  if (debug)
1291  {
1292  Pout<< "isoSurfaceCell : generated " << triMeshCells.size()
1293  << " unmerged triangles." << endl;
1294  }
1295 
1296 
1297  label nOldPoints = triPoints.size();
1298 
1299  // Trimmed to original triangle
1300  DynamicList<label> trimTriMap;
1301  // Trimmed to original point
1302  labelList trimTriPointMap;
1303  if (getClipBounds().good())
1304  {
1305  isoSurfacePoint::trimToBox
1306  (
1307  treeBoundBox(getClipBounds()),
1308  triPoints, // new points
1309  trimTriMap, // map from (new) triangle to original
1310  trimTriPointMap, // map from (new) point to original
1311  interpolatedPoints_, // labels of newly introduced points
1312  interpolatedOldPoints_, // and their interpolation
1313  interpolationWeights_
1314  );
1315  triMeshCells = labelField(triMeshCells, trimTriMap);
1316  }
1317 
1318 
1319  // Merge points and compact out non-valid triangles
1320  labelList triMap; // merged to unmerged triangle
1321  tmpsurf = stitchTriPoints
1322  (
1323  regularise, // check for duplicate tris
1324  triPoints,
1325  triPointMergeMap_, // unmerged to merged point
1326  triMap // merged to unmerged triangle
1327  );
1328 
1329  if (debug)
1330  {
1331  Pout<< "isoSurfaceCell : generated " << triMap.size()
1332  << " merged triangles." << endl;
1333  }
1334 
1335  if (getClipBounds().good())
1336  {
1337  // Adjust interpolatedPoints_
1338  inplaceRenumber(triPointMergeMap_, interpolatedPoints_);
1339 
1340  // Adjust triPointMergeMap_
1341  labelList newTriPointMergeMap(nOldPoints, -1);
1342  forAll(trimTriPointMap, trimPointI)
1343  {
1344  label oldPointI = trimTriPointMap[trimPointI];
1345  if (oldPointI >= 0)
1346  {
1347  label pointI = triPointMergeMap_[trimPointI];
1348  if (pointI >= 0)
1349  {
1350  newTriPointMergeMap[oldPointI] = pointI;
1351  }
1352  }
1353  }
1354  triPointMergeMap_.transfer(newTriPointMergeMap);
1355  }
1356 
1357  meshCells_.setSize(triMap.size());
1358  forAll(triMap, i)
1359  {
1360  meshCells_[i] = triMeshCells[triMap[i]];
1361  }
1362  }
1363 
1364 
1365  if (debug)
1366  {
1367  Pout<< "isoSurfaceCell : checking " << tmpsurf.size()
1368  << " triangles for validity." << endl;
1369 
1370  forAll(tmpsurf, triI)
1371  {
1372  triSurfaceTools::validTri(tmpsurf, triI);
1373  }
1374  }
1375 
1376 
1377  if (regularise)
1378  {
1379  List<FixedList<label, 3>> faceEdges;
1380  labelList edgeFace0, edgeFace1;
1381  Map<labelList> edgeFacesRest;
1382 
1383 
1384  while (true)
1385  {
1386  // Calculate addressing
1387  calcAddressing
1388  (
1389  tmpsurf,
1390  faceEdges,
1391  edgeFace0,
1392  edgeFace1,
1393  edgeFacesRest
1394  );
1395 
1396  // See if any dangling triangles
1397  boolList keepTriangles;
1398  label nDangling = markDanglingTriangles
1399  (
1400  faceEdges,
1401  edgeFace0,
1402  edgeFace1,
1403  edgeFacesRest,
1404  keepTriangles
1405  );
1406 
1407  if (debug)
1408  {
1409  Pout<< "isoSurfaceCell : detected " << nDangling
1410  << " dangling triangles." << endl;
1411  }
1412 
1413  if (nDangling == 0)
1414  {
1415  break;
1416  }
1417 
1418  // Create face map (new to old)
1419  labelList subsetTriMap(findIndices(keepTriangles, true));
1420 
1421  labelList subsetPointMap;
1422  labelList reversePointMap;
1423  tmpsurf = subsetMesh
1424  (
1425  tmpsurf,
1426  subsetTriMap,
1427  reversePointMap,
1428  subsetPointMap
1429  );
1430  meshCells_ = labelField(meshCells_, subsetTriMap);
1431  inplaceRenumber(reversePointMap, triPointMergeMap_);
1432  }
1433  }
1434 
1435 
1436  // Transfer to mesh storage. Note, an iso-surface has no zones
1437  {
1438  // Recover the pointField
1439  pointField pts;
1440  tmpsurf.swapPoints(pts);
1441 
1442  // Transcribe from triFace to face
1443  faceList faces;
1444  tmpsurf.triFaceFaces(faces);
1445 
1446  tmpsurf.clearOut();
1447 
1448  Mesh updated(std::move(pts), std::move(faces), surfZoneList());
1449 
1450  this->Mesh::transfer(updated);
1451  }
1452 }
1453 
1454 
1455 // ************************************************************************* //
const scalarField & pointValues() const noexcept
The mesh point values used for creating the iso-surface.
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:48
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
unsigned int count(const bool on=true) const
Count number of bits set.
Definition: bitSetI.H:426
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
label blockCells(UList< cutType > &cuts, const bitSet &ignoreCells) const
Mark ignoreCells as BLOCKED.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static void syncPointPositions(const polyMesh &mesh, List< point > &positions, const CombineOp &cop, const point &nullValue)
Synchronize locations on all mesh points.
Definition: syncTools.H:219
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition: List.C:326
label nPoints() const noexcept
Number of mesh points.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
const scalarField & cVals_
Cell values.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:517
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:899
Preferences for controlling iso-surface algorithms.
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:489
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
bool found(const T &val, label pos=0) const
Same as contains()
Definition: UList.H:879
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
Definition: pointIndexHit.H:58
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, any/none. Also accepts 0/1 as a string and shortcuts t/f, y/n.
Definition: Switch.H:77
labelList findIndices(const ListType &input, typename ListType::const_reference val, label start=0)
Linear search to find all occurrences of given element.
This class describes the interaction of an object (often a face) and a point. It carries the info of ...
Definition: pointIndexHit.H:44
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
Low-level components common to various iso-surface algorithms.
const scalar iso_
Iso value.
MinMax< label > minMax(const labelHashSet &set)
Find the min/max values of labelHashSet.
Definition: hashSets.C:54
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
label fcIndex(const label i) const noexcept
The forward circular index. The next index in the list which returns to the first at the end of the l...
Definition: UListI.H:97
label mergePoints(const PointList &points, labelList &pointToUnique, labelList &uniquePoints, const scalar mergeTol=SMALL, const bool verbose=false)
Calculate merge mapping, preserving the original point order. All points closer/equal mergeTol are to...
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
#define defineIsoSurfaceInterpolateMethods(ThisClass)
const boundBox & getClipBounds() const noexcept
Get optional clipping bounding box.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
const pointField & points
const polyMesh & mesh_
Reference to mesh.
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ const Field< vector > edgeCentres(faMeshTools::flattenEdgeField(aMesh.edgeCentres(), true))
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
label calcCellCuts(List< cutType > &cuts) const
Populate a list of candidate cell cuts using getCellCutType()
face triFace(3)
scalar mag() const
The magnitude/length of the bounding box diagonal.
Definition: boundBoxI.H:198
labelList meshCells_
For every face, the original cell in mesh.
const scalarField & cellValues() const noexcept
The mesh cell values used for creating the iso-surface.
const polyMesh & mesh() const noexcept
The mesh for which the iso-surface is associated.
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
Vector< scalar > vector
Definition: vector.H:57
A surface formed by the iso value. After "Polygonising A Scalar Field Using Tetrahedrons", Paul Bourke (http://paulbourke.net/geometry/polygonise) and "Regularised Marching Tetrahedra: improved iso-surface extraction", G.M. Treece, R.W. Prager and A.H. Gee.
const vectorField & cellCentres() const
errorManip< error > abort(error &err)
Definition: errorManip.H:139
label find(const T &val) const
Find index of the first occurrence of the value.
Definition: UList.C:173
meshedSurface Mesh
isoSurfaceCell(const polyMesh &mesh, const scalarField &cellValues, const scalarField &pointValues, const scalar iso, const isoSurfaceParams &params=isoSurfaceParams(), const bitSet &ignoreCells=bitSet())
Construct from cell and point values.
const scalarField & pVals_
Point values.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
int debug
Static debugging option.
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:51
List< geometricSurfacePatch > geometricSurfacePatchList
List of geometricSurfacePatch.
defineTypeNameAndDebug(combustionModel, 0)
Geometric merging of points. See below.
labelList f(nPoints)
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
static bool validTri(const triSurface &, const label facei, const bool verbose=true)
Check single triangle for (topological) validity.
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Inplace renumber the values (not the indices) of a list.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
vector point
Point is a vector.
Definition: point.H:37
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label nCells() const noexcept
Number of mesh cells.
const dimensionedScalar c
Speed of light in a vacuum.
const polyBoundaryMesh & patches
List< surfZone > surfZoneList
List of surfZone.
Definition: surfZoneList.H:32
Nothing to be read.
void transfer(pointField &pointLst, List< face > &faceLst)
Transfer the components.
const labelListList & faceEdges() const
Return face-edge addressing.
const boundBox & bounds() const noexcept
Return mesh bounding box.
Definition: polyMesh.H:616
label n
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
List< label > labelList
A List of labels.
Definition: List.H:62
Triangulated surface description with patch information.
Definition: triSurface.H:71
label rcIndex(const label i) const noexcept
The reverse circular index. The previous index in the list which returns to the last at the beginning...
Definition: UListI.H:104
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
static bool test(const UList< face > &faces)
Test if given list of faces satisfies criteria for TET. (4 tri)
Definition: tetMatcher.C:82
Convenience macros for instantiating iso-surface interpolate methods.
filterType filter() const noexcept
Get current filter type.
List< bool > boolList
A List of bools.
Definition: List.H:60
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Do not request registration (bool: false)
const volScalarField & p0
Definition: EEqn.H:36
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
scalar mergeTol() const noexcept
Get current merge tolerance.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
const pointField & pts
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127