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 | Copyright (C) 2011-2016 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 "isoSurfaceCell.H"
27 #include "dictionary.H"
28 #include "polyMesh.H"
29 #include "mergePoints.H"
30 #include "tetMatcher.H"
31 #include "syncTools.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(isoSurfaceCell, 0);
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 Foam::scalar Foam::isoSurfaceCell::isoFraction
45 (
46  const scalar s0,
47  const scalar s1
48 ) const
49 {
50  scalar d = s1-s0;
51 
52  if (mag(d) > VSMALL)
53  {
54  return (iso_-s0)/d;
55  }
56  else
57  {
58  return -1.0;
59  }
60 }
61 
62 
63 bool Foam::isoSurfaceCell::isTriCut
64 (
65  const triFace& tri,
66  const scalarField& pointValues
67 ) const
68 {
69  bool aLower = (pointValues[tri[0]] < iso_);
70  bool bLower = (pointValues[tri[1]] < iso_);
71  bool cLower = (pointValues[tri[2]] < iso_);
72 
73  return !(aLower == bLower && aLower == cLower);
74 }
75 
76 
77 Foam::isoSurfaceCell::cellCutType Foam::isoSurfaceCell::calcCutType
78 (
79  const PackedBoolList& isTet,
80  const scalarField& cellValues,
81  const scalarField& pointValues,
82  const label celli
83 ) const
84 {
85  const cell& cFaces = mesh_.cells()[celli];
86 
87  if (isTet.get(celli) == 1)
88  {
89  forAll(cFaces, cFacei)
90  {
91  const face& f = mesh_.faces()[cFaces[cFacei]];
92 
93  for (label fp = 1; fp < f.size() - 1; fp++)
94  {
95  triFace tri(f[0], f[fp], f[f.fcIndex(fp)]);
96 
97  if (isTriCut(tri, pointValues))
98  {
99  return CUT;
100  }
101  }
102  }
103  return NOTCUT;
104  }
105  else
106  {
107  bool cellLower = (cellValues[celli] < iso_);
108 
109  // First check if there is any cut in cell
110  bool edgeCut = false;
111 
112  forAll(cFaces, cFacei)
113  {
114  label facei = cFaces[cFacei];
115  const face& f = mesh_.faces()[facei];
116 
117  // Check pyramids cut
118  forAll(f, fp)
119  {
120  if ((pointValues[f[fp]] < iso_) != cellLower)
121  {
122  edgeCut = true;
123  break;
124  }
125  }
126 
127  if (edgeCut)
128  {
129  break;
130  }
131 
132  const label fp0 = mesh_.tetBasePtIs()[facei];
133  label fp = f.fcIndex(fp0);
134  for (label i = 2; i < f.size(); i++)
135  {
136  label nextFp = f.fcIndex(fp);
137 
138  if (isTriCut(triFace(f[fp0], f[fp], f[nextFp]), pointValues))
139  {
140  edgeCut = true;
141  break;
142  }
143 
144  fp = nextFp;
145  }
146 
147  if (edgeCut)
148  {
149  break;
150  }
151  }
152 
153  if (edgeCut)
154  {
155  // Count actual cuts (expensive since addressing needed)
156  // Note: not needed if you don't want to preserve maxima/minima
157  // centred around cellcentre. In that case just always return CUT
158 
159  const labelList& cPoints = mesh_.cellPoints(celli);
160 
161  label nPyrCuts = 0;
162 
163  forAll(cPoints, i)
164  {
165  if ((pointValues[cPoints[i]] < iso_) != cellLower)
166  {
167  nPyrCuts++;
168  }
169  }
170 
171  if (nPyrCuts == cPoints.size())
172  {
173  return SPHERE;
174  }
175  else
176  {
177  return CUT;
178  }
179  }
180  else
181  {
182  return NOTCUT;
183  }
184  }
185 }
186 
187 
188 void Foam::isoSurfaceCell::calcCutTypes
189 (
190  const PackedBoolList& isTet,
191  const scalarField& cVals,
192  const scalarField& pVals
193 )
194 {
195  cellCutType_.setSize(mesh_.nCells());
196  nCutCells_ = 0;
197  forAll(mesh_.cells(), celli)
198  {
199  cellCutType_[celli] = calcCutType(isTet, cVals, pVals, celli);
200 
201  if (cellCutType_[celli] == CUT)
202  {
203  nCutCells_++;
204  }
205  }
206 
207  if (debug)
208  {
209  Pout<< "isoSurfaceCell : detected " << nCutCells_
210  << " candidate cut cells." << endl;
211  }
212 }
213 
214 
215 Foam::labelPair Foam::isoSurfaceCell::findCommonPoints
216 (
217  const labelledTri& tri0,
218  const labelledTri& tri1
219 )
220 {
221  labelPair common(-1, -1);
222 
223  label fp0 = 0;
224  label fp1 = findIndex(tri1, tri0[fp0]);
225 
226  if (fp1 == -1)
227  {
228  fp0 = 1;
229  fp1 = findIndex(tri1, tri0[fp0]);
230  }
231 
232  if (fp1 != -1)
233  {
234  // So tri0[fp0] is tri1[fp1]
235 
236  // Find next common point
237  label fp0p1 = tri0.fcIndex(fp0);
238  label fp1p1 = tri1.fcIndex(fp1);
239  label fp1m1 = tri1.rcIndex(fp1);
240 
241  if (tri0[fp0p1] == tri1[fp1p1] || tri0[fp0p1] == tri1[fp1m1])
242  {
243  common[0] = tri0[fp0];
244  common[1] = tri0[fp0p1];
245  }
246  }
247  return common;
248 }
249 
250 
251 Foam::point Foam::isoSurfaceCell::calcCentre(const triSurface& s)
252 {
253  vector sum = Zero;
254 
255  forAll(s, i)
256  {
257  sum += s[i].centre(s.points());
258  }
259  return sum/s.size();
260 }
261 
262 
263 Foam::pointIndexHit Foam::isoSurfaceCell::collapseSurface
264 (
265  const label celli,
266  pointField& localPoints,
267  DynamicList<labelledTri, 64>& localTris
268 ) const
269 {
270  pointIndexHit info(false, Zero, localTris.size());
271 
272  if (localTris.size() == 1)
273  {
274  const labelledTri& tri = localTris[0];
275  info.setPoint(tri.centre(localPoints));
276  info.setHit();
277  }
278  else if (localTris.size() == 2)
279  {
280  // Check if the two triangles share an edge.
281  const labelledTri& tri0 = localTris[0];
282  const labelledTri& tri1 = localTris[1];
283 
284  labelPair shared = findCommonPoints(tri0, tri1);
285 
286  if (shared[0] != -1)
287  {
288  vector n0 = tri0.normal(localPoints);
289  n0 /= mag(n0);
290  vector n1 = tri1.normal(localPoints);
291  n1 /= mag(n1);
292 
293  if ((n0 & n1) < 0)
294  {
295  // Too big an angle. Do not simplify.
296  }
297  else
298  {
299  info.setPoint
300  (
301  0.5
302  * (
303  tri0.centre(localPoints)
304  + tri1.centre(localPoints)
305  )
306  );
307  info.setHit();
308  }
309  }
310  }
311  else if (localTris.size())
312  {
313  // Check if single region. Rare situation.
314  triSurface surf
315  (
316  localTris,
318  localPoints,
319  true
320  );
321  localTris.clearStorage();
322 
323  labelList faceZone;
324  label nZones = surf.markZones
325  (
326  boolList(surf.nEdges(), false),
327  faceZone
328  );
329 
330  if (nZones == 1)
331  {
332  // Check that all normals make a decent angle
333  scalar minCos = GREAT;
334  const vector& n0 = surf.faceNormals()[0];
335  for (label i = 1; i < surf.size(); i++)
336  {
337  scalar cosAngle = (n0 & surf.faceNormals()[i]);
338  if (cosAngle < minCos)
339  {
340  minCos = cosAngle;
341  }
342  }
343 
344  if (minCos > 0)
345  {
346  info.setPoint(calcCentre(surf));
347  info.setHit();
348  }
349  }
350  }
351 
352  return info;
353 }
354 
355 
356 void Foam::isoSurfaceCell::calcSnappedCc
357 (
358  const PackedBoolList& isTet,
359  const scalarField& cVals,
360  const scalarField& pVals,
361 
362  DynamicList<point>& snappedPoints,
363  labelList& snappedCc
364 ) const
365 {
366  const pointField& cc = mesh_.cellCentres();
367  const pointField& pts = mesh_.points();
368 
369  snappedCc.setSize(mesh_.nCells());
370  snappedCc = -1;
371 
372  // Work arrays
373  DynamicList<point, 64> localPoints(64);
374  DynamicList<labelledTri, 64> localTris(64);
375  Map<label> pointToLocal(64);
376 
377  forAll(mesh_.cells(), celli)
378  {
379  if (cellCutType_[celli] == CUT && isTet.get(celli) == 0)
380  {
381  scalar cVal = cVals[celli];
382 
383  const cell& cFaces = mesh_.cells()[celli];
384 
385  localPoints.clear();
386  localTris.clear();
387  pointToLocal.clear();
388 
389  // Create points for all intersections close to cell centre
390  // (i.e. from pyramid edges)
391 
392  forAll(cFaces, cFacei)
393  {
394  const face& f = mesh_.faces()[cFaces[cFacei]];
395 
396  forAll(f, fp)
397  {
398  label pointi = f[fp];
399 
400  scalar s = isoFraction(cVal, pVals[pointi]);
401 
402  if (s >= 0.0 && s <= 0.5)
403  {
404  if (pointToLocal.insert(pointi, localPoints.size()))
405  {
406  localPoints.append((1.0-s)*cc[celli]+s*pts[pointi]);
407  }
408  }
409  }
410  }
411 
412  if (localPoints.size() == 1)
413  {
414  // No need for any analysis.
415  snappedCc[celli] = snappedPoints.size();
416  snappedPoints.append(localPoints[0]);
417 
418  //Pout<< "cell:" << celli
419  // << " at " << mesh_.cellCentres()[celli]
420  // << " collapsing " << localPoints
421  // << " intersections down to "
422  // << snappedPoints[snappedCc[celli]] << endl;
423  }
424  else if (localPoints.size() == 2)
425  {
426  //? No need for any analysis.???
427  snappedCc[celli] = snappedPoints.size();
428  snappedPoints.append(0.5*(localPoints[0]+localPoints[1]));
429 
430  //Pout<< "cell:" << celli
431  // << " at " << mesh_.cellCentres()[celli]
432  // << " collapsing " << localPoints
433  // << " intersections down to "
434  // << snappedPoints[snappedCc[celli]] << endl;
435  }
436  else if (localPoints.size())
437  {
438  // Need to analyse
439  forAll(cFaces, cFacei)
440  {
441  label facei = cFaces[cFacei];
442  const face& f = mesh_.faces()[facei];
443 
444  // Do a tetrahedralisation. Each face to cc becomes pyr.
445  // Each pyr gets split into tets by diagonalisation
446  // of face.
447 
448  const label fp0 = mesh_.tetBasePtIs()[facei];
449  label fp = f.fcIndex(fp0);
450  for (label i = 2; i < f.size(); i++)
451  {
452  label nextFp = f.fcIndex(fp);
453  triFace tri(f[fp0], f[fp], f[nextFp]);
454 
455  // Get fractions for the three edges to cell centre
456  FixedList<scalar, 3> s(3);
457  s[0] = isoFraction(cVal, pVals[tri[0]]);
458  s[1] = isoFraction(cVal, pVals[tri[1]]);
459  s[2] = isoFraction(cVal, pVals[tri[2]]);
460 
461  if
462  (
463  (s[0] >= 0.0 && s[0] <= 0.5)
464  && (s[1] >= 0.0 && s[1] <= 0.5)
465  && (s[2] >= 0.0 && s[2] <= 0.5)
466  )
467  {
468  if
469  (
470  (mesh_.faceOwner()[facei] == celli)
471  == (cVal >= pVals[tri[0]])
472  )
473  {
474  localTris.append
475  (
476  labelledTri
477  (
478  pointToLocal[tri[1]],
479  pointToLocal[tri[0]],
480  pointToLocal[tri[2]],
481  0
482  )
483  );
484  }
485  else
486  {
487  localTris.append
488  (
489  labelledTri
490  (
491  pointToLocal[tri[0]],
492  pointToLocal[tri[1]],
493  pointToLocal[tri[2]],
494  0
495  )
496  );
497  }
498  }
499 
500  fp = nextFp;
501  }
502  }
503 
504  pointField surfPoints;
505  surfPoints.transfer(localPoints);
506  pointIndexHit info = collapseSurface
507  (
508  celli,
509  surfPoints,
510  localTris
511  );
512 
513  if (info.hit())
514  {
515  snappedCc[celli] = snappedPoints.size();
516  snappedPoints.append(info.hitPoint());
517 
518  //Pout<< "cell:" << celli
519  // << " at " << mesh_.cellCentres()[celli]
520  // << " collapsing " << surfPoints
521  // << " intersections down to "
522  // << snappedPoints[snappedCc[celli]] << endl;
523  }
524  }
525  }
526  }
527 }
528 
529 
530 void Foam::isoSurfaceCell::genPointTris
531 (
532  const scalarField& cellValues,
533  const scalarField& pointValues,
534  const label pointi,
535  const label facei,
536  const label celli,
537  DynamicList<point, 64>& localTriPoints
538 ) const
539 {
540  const pointField& cc = mesh_.cellCentres();
541  const pointField& pts = mesh_.points();
542  const face& f = mesh_.faces()[facei];
543 
544  const label fp0 = mesh_.tetBasePtIs()[facei];
545  label fp = f.fcIndex(fp0);
546  for (label i = 2; i < f.size(); i++)
547  {
548  label nextFp = f.fcIndex(fp);
549  triFace tri(f[fp0], f[fp], f[nextFp]);
550 
551  label index = findIndex(tri, pointi);
552 
553  if (index == -1)
554  {
555  continue;
556  }
557 
558  // Tet between index..index-1, index..index+1, index..cc
559  label b = tri[tri.fcIndex(index)];
560  label c = tri[tri.rcIndex(index)];
561 
562  // Get fractions for the three edges emanating from point
563  FixedList<scalar, 3> s(3);
564  s[0] = isoFraction(pointValues[pointi], pointValues[b]);
565  s[1] = isoFraction(pointValues[pointi], pointValues[c]);
566  s[2] = isoFraction(pointValues[pointi], cellValues[celli]);
567 
568  if
569  (
570  (s[0] >= 0.0 && s[0] <= 0.5)
571  && (s[1] >= 0.0 && s[1] <= 0.5)
572  && (s[2] >= 0.0 && s[2] <= 0.5)
573  )
574  {
575  point p0 = (1.0-s[0])*pts[pointi] + s[0]*pts[b];
576  point p1 = (1.0-s[1])*pts[pointi] + s[1]*pts[c];
577  point p2 = (1.0-s[2])*pts[pointi] + s[2]*cc[celli];
578 
579  if
580  (
581  (mesh_.faceOwner()[facei] == celli)
582  == (pointValues[pointi] > cellValues[celli])
583  )
584  {
585  localTriPoints.append(p0);
586  localTriPoints.append(p1);
587  localTriPoints.append(p2);
588  }
589  else
590  {
591  localTriPoints.append(p1);
592  localTriPoints.append(p0);
593  localTriPoints.append(p2);
594  }
595  }
596 
597  fp = nextFp;
598  }
599 }
600 
601 
602 void Foam::isoSurfaceCell::genPointTris
603 (
604  const scalarField& pointValues,
605  const label pointi,
606  const label facei,
607  const label celli,
608  DynamicList<point, 64>& localTriPoints
609 ) const
610 {
611  const pointField& pts = mesh_.points();
612  const cell& cFaces = mesh_.cells()[celli];
613 
614  FixedList<label, 4> tet;
615 
616  // Make tet from this face to the 4th point (same as cellcentre in
617  // non-tet cells)
618  const face& f = mesh_.faces()[facei];
619 
620  // Find 4th point
621  label ccPointi = -1;
622  forAll(cFaces, cFacei)
623  {
624  const face& f1 = mesh_.faces()[cFaces[cFacei]];
625  forAll(f1, fp)
626  {
627  label p1 = f1[fp];
628 
629  if (findIndex(f, p1) == -1)
630  {
631  ccPointi = p1;
632  break;
633  }
634  }
635  if (ccPointi != -1)
636  {
637  break;
638  }
639  }
640 
641 
642  // Tet between index..index-1, index..index+1, index..cc
643  label index = findIndex(f, pointi);
644  label b = f[f.fcIndex(index)];
645  label c = f[f.rcIndex(index)];
646 
647  //Pout<< " p0:" << pointi << " b:" << b << " c:" << c
648  //<< " d:" << ccPointi << endl;
649 
650  // Get fractions for the three edges emanating from point
651  FixedList<scalar, 3> s(3);
652  s[0] = isoFraction(pointValues[pointi], pointValues[b]);
653  s[1] = isoFraction(pointValues[pointi], pointValues[c]);
654  s[2] = isoFraction(pointValues[pointi], pointValues[ccPointi]);
655 
656  if
657  (
658  (s[0] >= 0.0 && s[0] <= 0.5)
659  && (s[1] >= 0.0 && s[1] <= 0.5)
660  && (s[2] >= 0.0 && s[2] <= 0.5)
661  )
662  {
663  point p0 = (1.0-s[0])*pts[pointi] + s[0]*pts[b];
664  point p1 = (1.0-s[1])*pts[pointi] + s[1]*pts[c];
665  point p2 = (1.0-s[2])*pts[pointi] + s[2]*pts[ccPointi];
666 
667  if (mesh_.faceOwner()[facei] != celli)
668  {
669  localTriPoints.append(p0);
670  localTriPoints.append(p1);
671  localTriPoints.append(p2);
672  }
673  else
674  {
675  localTriPoints.append(p1);
676  localTriPoints.append(p0);
677  localTriPoints.append(p2);
678  }
679  }
680 }
681 
682 
683 void Foam::isoSurfaceCell::calcSnappedPoint
684 (
685  const PackedBoolList& isTet,
686  const scalarField& cVals,
687  const scalarField& pVals,
688 
689  DynamicList<point>& snappedPoints,
690  labelList& snappedPoint
691 ) const
692 {
693  // Determine if point is on boundary. Points on boundaries are never
694  // snapped. Coupled boundaries are handled explicitly so not marked here.
695  PackedBoolList isBoundaryPoint(mesh_.nPoints());
696  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
697  forAll(patches, patchi)
698  {
699  const polyPatch& pp = patches[patchi];
700 
701  if (!pp.coupled())
702  {
703  label facei = pp.start();
704  forAll(pp, i)
705  {
706  const face& f = mesh_.faces()[facei++];
707 
708  forAll(f, fp)
709  {
710  isBoundaryPoint.set(f[fp], 1);
711  }
712  }
713  }
714  }
715 
716 
717  const point greatPoint(GREAT, GREAT, GREAT);
718 
719  pointField collapsedPoint(mesh_.nPoints(), greatPoint);
720 
721 
722  // Work arrays
723  DynamicList<point, 64> localTriPoints(100);
724  labelHashSet localPointCells(100);
725 
726  forAll(mesh_.pointFaces(), pointi)
727  {
728  if (isBoundaryPoint.get(pointi) == 1)
729  {
730  continue;
731  }
732 
733  const labelList& pFaces = mesh_.pointFaces()[pointi];
734 
735  bool anyCut = false;
736 
737  forAll(pFaces, i)
738  {
739  label facei = pFaces[i];
740 
741  if
742  (
743  cellCutType_[mesh_.faceOwner()[facei]] == CUT
744  || (
745  mesh_.isInternalFace(facei)
746  && cellCutType_[mesh_.faceNeighbour()[facei]] == CUT
747  )
748  )
749  {
750  anyCut = true;
751  break;
752  }
753  }
754 
755  if (!anyCut)
756  {
757  continue;
758  }
759 
760 
761  // Do a pointCells walk (by using pointFaces)
762 
763  localPointCells.clear();
764  localTriPoints.clear();
765 
766  forAll(pFaces, pFacei)
767  {
768  label facei = pFaces[pFacei];
769  label own = mesh_.faceOwner()[facei];
770 
771  if (isTet.get(own) == 1)
772  {
773  // Since tets have no cell centre to include make sure
774  // we only generate a triangle once per point.
775  if (localPointCells.insert(own))
776  {
777  genPointTris(pVals, pointi, facei, own, localTriPoints);
778  }
779  }
780  else
781  {
782  genPointTris
783  (
784  cVals,
785  pVals,
786  pointi,
787  facei,
788  own,
789  localTriPoints
790  );
791  }
792 
793  if (mesh_.isInternalFace(facei))
794  {
795  label nei = mesh_.faceNeighbour()[facei];
796 
797  if (isTet.get(nei) == 1)
798  {
799  if (localPointCells.insert(nei))
800  {
801  genPointTris(pVals, pointi, facei, nei, localTriPoints);
802  }
803  }
804  else
805  {
806  genPointTris
807  (
808  cVals,
809  pVals,
810  pointi,
811  facei,
812  nei,
813  localTriPoints
814  );
815  }
816  }
817  }
818 
819  if (localTriPoints.size() == 3)
820  {
821  // Single triangle. No need for any analysis. Average points.
823  points.transfer(localTriPoints);
824  collapsedPoint[pointi] = sum(points)/points.size();
825 
826  //Pout<< " point:" << pointi
827  // << " replacing coord:" << mesh_.points()[pointi]
828  // << " by average:" << collapsedPoint[pointi] << endl;
829  }
830  else if (localTriPoints.size())
831  {
832  // Convert points into triSurface.
833 
834  // Merge points and compact out non-valid triangles
835  labelList triMap; // merged to unmerged triangle
836  labelList triPointReverseMap; // unmerged to merged point
837  triSurface surf
838  (
839  stitchTriPoints
840  (
841  false, // do not check for duplicate tris
842  localTriPoints,
843  triPointReverseMap,
844  triMap
845  )
846  );
847 
848  labelList faceZone;
849  label nZones = surf.markZones
850  (
851  boolList(surf.nEdges(), false),
852  faceZone
853  );
854 
855  if (nZones == 1)
856  {
857  // Check that all normals make a decent angle
858  scalar minCos = GREAT;
859  const vector& n0 = surf.faceNormals()[0];
860  for (label i = 1; i < surf.size(); i++)
861  {
862  const vector& n = surf.faceNormals()[i];
863  scalar cosAngle = (n0 & n);
864  if (cosAngle < minCos)
865  {
866  minCos = cosAngle;
867  }
868  }
869  if (minCos > 0)
870  {
871  collapsedPoint[pointi] = calcCentre(surf);
872  }
873  }
874  }
875  }
876 
878  (
879  mesh_,
880  collapsedPoint,
881  minMagSqrEqOp<point>(),
882  greatPoint
883  );
884 
885  snappedPoint.setSize(mesh_.nPoints());
886  snappedPoint = -1;
887 
888  forAll(collapsedPoint, pointi)
889  {
890  // Cannot do == comparison since might be transformed so have
891  // truncation errors.
892  if (magSqr(collapsedPoint[pointi]) < 0.5*magSqr(greatPoint))
893  {
894  snappedPoint[pointi] = snappedPoints.size();
895  snappedPoints.append(collapsedPoint[pointi]);
896  }
897  }
898 }
899 
900 
901 Foam::triSurface Foam::isoSurfaceCell::stitchTriPoints
902 (
903  const bool checkDuplicates,
904  const List<point>& triPoints,
905  labelList& triPointReverseMap, // unmerged to merged point
906  labelList& triMap // merged to unmerged triangle
907 ) const
908 {
909  label nTris = triPoints.size()/3;
910 
911  if ((triPoints.size() % 3) != 0)
912  {
914  << "Problem: number of points " << triPoints.size()
915  << " not a multiple of 3." << abort(FatalError);
916  }
917 
918  pointField newPoints;
920  (
921  triPoints,
922  mergeDistance_,
923  false,
924  triPointReverseMap,
925  newPoints
926  );
927 
928  // Check that enough merged.
929  if (debug)
930  {
931  Pout<< "isoSurfaceCell : merged from " << triPoints.size()
932  << " points down to " << newPoints.size() << endl;
933 
934  pointField newNewPoints;
935  labelList oldToNew;
936  bool hasMerged = mergePoints
937  (
938  newPoints,
939  mergeDistance_,
940  true,
941  oldToNew,
942  newNewPoints
943  );
944 
945  if (hasMerged)
946  {
948  << "Merged points contain duplicates"
949  << " when merging with distance " << mergeDistance_ << endl
950  << "merged:" << newPoints.size() << " re-merged:"
951  << newNewPoints.size()
952  << abort(FatalError);
953  }
954  }
955 
956 
957  List<labelledTri> tris;
958  {
959  DynamicList<labelledTri> dynTris(nTris);
960  label rawPointi = 0;
961  DynamicList<label> newToOldTri(nTris);
962 
963  for (label oldTriI = 0; oldTriI < nTris; oldTriI++)
964  {
965  labelledTri tri
966  (
967  triPointReverseMap[rawPointi],
968  triPointReverseMap[rawPointi+1],
969  triPointReverseMap[rawPointi+2],
970  0
971  );
972  if ((tri[0] != tri[1]) && (tri[0] != tri[2]) && (tri[1] != tri[2]))
973  {
974  newToOldTri.append(oldTriI);
975  dynTris.append(tri);
976  }
977 
978  rawPointi += 3;
979  }
980 
981  triMap.transfer(newToOldTri);
982  tris.transfer(dynTris);
983  }
984 
985 
986  // Use face centres to determine 'flat hole' situation (see RMT paper).
987  // Two unconnected triangles get connected because (some of) the edges
988  // separating them get collapsed. Below only checks for duplicate triangles,
989  // not non-manifold edge connectivity.
990  if (checkDuplicates)
991  {
992  if (debug)
993  {
994  Pout<< "isoSurfaceCell : merged from " << nTris
995  << " down to " << tris.size() << " triangles." << endl;
996  }
997 
998  pointField centres(tris.size());
999  forAll(tris, triI)
1000  {
1001  centres[triI] = tris[triI].centre(newPoints);
1002  }
1003 
1004  pointField mergedCentres;
1005  labelList oldToMerged;
1006  bool hasMerged = mergePoints
1007  (
1008  centres,
1009  mergeDistance_,
1010  false,
1011  oldToMerged,
1012  mergedCentres
1013  );
1014 
1015  if (debug)
1016  {
1017  Pout<< "isoSurfaceCell : detected "
1018  << centres.size()-mergedCentres.size()
1019  << " duplicate triangles." << endl;
1020  }
1021 
1022  if (hasMerged)
1023  {
1024  // Filter out duplicates.
1025  label newTriI = 0;
1026  DynamicList<label> newToOldTri(tris.size());
1027  labelList newToMaster(mergedCentres.size(), -1);
1028  forAll(tris, triI)
1029  {
1030  label mergedI = oldToMerged[triI];
1031 
1032  if (newToMaster[mergedI] == -1)
1033  {
1034  newToMaster[mergedI] = triI;
1035  newToOldTri.append(triMap[triI]);
1036  tris[newTriI++] = tris[triI];
1037  }
1038  }
1039 
1040  triMap.transfer(newToOldTri);
1041  tris.setSize(newTriI);
1042  }
1043  }
1044 
1045  return triSurface(tris, geometricSurfacePatchList(0), newPoints, true);
1046 }
1047 
1048 
1049 bool Foam::isoSurfaceCell::validTri(const triSurface& surf, const label facei)
1050 {
1051  // Simple check on indices ok.
1052 
1053  const labelledTri& f = surf[facei];
1054 
1055  forAll(f, fp)
1056  {
1057  if (f[fp] < 0 || f[fp] >= surf.points().size())
1058  {
1060  << "triangle " << facei << " vertices " << f
1061  << " uses point indices outside point range 0.."
1062  << surf.points().size()-1 << endl;
1063 
1064  return false;
1065  }
1066  }
1067 
1068  if ((f[0] == f[1]) || (f[0] == f[2]) || (f[1] == f[2]))
1069  {
1071  << "triangle " << facei
1072  << " uses non-unique vertices " << f
1073  << endl;
1074  return false;
1075  }
1076 
1077  // duplicate triangle check
1078 
1079  const labelList& fFaces = surf.faceFaces()[facei];
1080 
1081  // Check if faceNeighbours use same points as this face.
1082  // Note: discards normal information - sides of baffle are merged.
1083  forAll(fFaces, i)
1084  {
1085  label nbrFacei = fFaces[i];
1086 
1087  if (nbrFacei <= facei)
1088  {
1089  // lower numbered faces already checked
1090  continue;
1091  }
1092 
1093  const labelledTri& nbrF = surf[nbrFacei];
1094 
1095  if
1096  (
1097  ((f[0] == nbrF[0]) || (f[0] == nbrF[1]) || (f[0] == nbrF[2]))
1098  && ((f[1] == nbrF[0]) || (f[1] == nbrF[1]) || (f[1] == nbrF[2]))
1099  && ((f[2] == nbrF[0]) || (f[2] == nbrF[1]) || (f[2] == nbrF[2]))
1100  )
1101  {
1103  << "triangle " << facei << " vertices " << f
1104  << " coords:" << f.points(surf.points())
1105  << " has the same vertices as triangle " << nbrFacei
1106  << " vertices " << nbrF
1107  << endl;
1108 
1109  return false;
1110  }
1111  }
1112  return true;
1113 }
1114 
1115 
1116 void Foam::isoSurfaceCell::calcAddressing
1117 (
1118  const triSurface& surf,
1119  List<FixedList<label, 3>>& faceEdges,
1120  labelList& edgeFace0,
1121  labelList& edgeFace1,
1122  Map<labelList>& edgeFacesRest
1123 ) const
1124 {
1125  const pointField& points = surf.points();
1126 
1127  pointField edgeCentres(3*surf.size());
1128  label edgeI = 0;
1129  forAll(surf, triI)
1130  {
1131  const labelledTri& tri = surf[triI];
1132  edgeCentres[edgeI++] = 0.5*(points[tri[0]]+points[tri[1]]);
1133  edgeCentres[edgeI++] = 0.5*(points[tri[1]]+points[tri[2]]);
1134  edgeCentres[edgeI++] = 0.5*(points[tri[2]]+points[tri[0]]);
1135  }
1136 
1137  pointField mergedCentres;
1138  labelList oldToMerged;
1139  bool hasMerged = mergePoints
1140  (
1141  edgeCentres,
1142  mergeDistance_,
1143  false,
1144  oldToMerged,
1145  mergedCentres
1146  );
1147 
1148  if (debug)
1149  {
1150  Pout<< "isoSurfaceCell : detected "
1151  << mergedCentres.size()
1152  << " edges on " << surf.size() << " triangles." << endl;
1153  }
1154 
1155  if (!hasMerged)
1156  {
1157  return;
1158  }
1159 
1160 
1161  // Determine faceEdges
1162  faceEdges.setSize(surf.size());
1163  edgeI = 0;
1164  forAll(surf, triI)
1165  {
1166  faceEdges[triI][0] = oldToMerged[edgeI++];
1167  faceEdges[triI][1] = oldToMerged[edgeI++];
1168  faceEdges[triI][2] = oldToMerged[edgeI++];
1169  }
1170 
1171 
1172  // Determine edgeFaces
1173  edgeFace0.setSize(mergedCentres.size());
1174  edgeFace0 = -1;
1175  edgeFace1.setSize(mergedCentres.size());
1176  edgeFace1 = -1;
1177  edgeFacesRest.clear();
1178 
1179  forAll(oldToMerged, oldEdgeI)
1180  {
1181  label triI = oldEdgeI / 3;
1182  label edgeI = oldToMerged[oldEdgeI];
1183 
1184  if (edgeFace0[edgeI] == -1)
1185  {
1186  edgeFace0[edgeI] = triI;
1187  }
1188  else if (edgeFace1[edgeI] == -1)
1189  {
1190  edgeFace1[edgeI] = triI;
1191  }
1192  else
1193  {
1194  //WarningInFunction
1195  // << "Edge " << edgeI << " with centre " << mergedCentres[edgeI]
1196  // << " used by more than two triangles: " << edgeFace0[edgeI]
1197  // << ", "
1198  // << edgeFace1[edgeI] << " and " << triI << endl;
1199  Map<labelList>::iterator iter = edgeFacesRest.find(edgeI);
1200 
1201  if (iter != edgeFacesRest.end())
1202  {
1203  labelList& eFaces = iter();
1204  label sz = eFaces.size();
1205  eFaces.setSize(sz+1);
1206  eFaces[sz] = triI;
1207  }
1208  else
1209  {
1210  edgeFacesRest.insert(edgeI, labelList(1, triI));
1211  }
1212  }
1213  }
1214 }
1215 
1216 
1217 bool Foam::isoSurfaceCell::danglingTriangle
1218 (
1219  const FixedList<label, 3>& fEdges,
1220  const labelList& edgeFace1
1221 )
1222 {
1223  label nOpen = 0;
1224  forAll(fEdges, i)
1225  {
1226  if (edgeFace1[fEdges[i]] == -1)
1227  {
1228  nOpen++;
1229  }
1230  }
1231 
1232  if (nOpen == 1 || nOpen == 2 || nOpen == 3)
1233  {
1234  return true;
1235  }
1236  else
1237  {
1238  return false;
1239  }
1240 }
1241 
1242 
1243 Foam::label Foam::isoSurfaceCell::markDanglingTriangles
1244 (
1245  const List<FixedList<label, 3>>& faceEdges,
1246  const labelList& edgeFace0,
1247  const labelList& edgeFace1,
1248  const Map<labelList>& edgeFacesRest,
1249  boolList& keepTriangles
1250 )
1251 {
1252  keepTriangles.setSize(faceEdges.size());
1253  keepTriangles = true;
1254 
1255  label nDangling = 0;
1256 
1257  // Remove any dangling triangles
1258  forAllConstIter(Map<labelList>, edgeFacesRest, iter)
1259  {
1260  // These are all the non-manifold edges. Filter out all triangles
1261  // with only one connected edge (= this edge)
1262 
1263  label edgeI = iter.key();
1264  const labelList& otherEdgeFaces = iter();
1265 
1266  // Remove all dangling triangles
1267  if (danglingTriangle(faceEdges[edgeFace0[edgeI]], edgeFace1))
1268  {
1269  keepTriangles[edgeFace0[edgeI]] = false;
1270  nDangling++;
1271  }
1272  if (danglingTriangle(faceEdges[edgeFace1[edgeI]], edgeFace1))
1273  {
1274  keepTriangles[edgeFace1[edgeI]] = false;
1275  nDangling++;
1276  }
1277  forAll(otherEdgeFaces, i)
1278  {
1279  label triI = otherEdgeFaces[i];
1280  if (danglingTriangle(faceEdges[triI], edgeFace1))
1281  {
1282  keepTriangles[triI] = false;
1283  nDangling++;
1284  }
1285  }
1286  }
1287  return nDangling;
1288 }
1289 
1290 
1291 Foam::triSurface Foam::isoSurfaceCell::subsetMesh
1292 (
1293  const triSurface& s,
1294  const labelList& newToOldFaces,
1295  labelList& oldToNewPoints,
1296  labelList& newToOldPoints
1297 )
1298 {
1299  const boolList include
1300  (
1301  createWithValues<boolList>
1302  (
1303  s.size(),
1304  false,
1305  newToOldFaces,
1306  true
1307  )
1308  );
1309 
1310  newToOldPoints.setSize(s.points().size());
1311  oldToNewPoints.setSize(s.points().size());
1312  oldToNewPoints = -1;
1313  {
1314  label pointi = 0;
1315 
1316  forAll(include, oldFacei)
1317  {
1318  if (include[oldFacei])
1319  {
1320  // Renumber labels for face
1321  const triSurface::FaceType& f = s[oldFacei];
1322 
1323  forAll(f, fp)
1324  {
1325  label oldPointi = f[fp];
1326 
1327  if (oldToNewPoints[oldPointi] == -1)
1328  {
1329  oldToNewPoints[oldPointi] = pointi;
1330  newToOldPoints[pointi++] = oldPointi;
1331  }
1332  }
1333  }
1334  }
1335  newToOldPoints.setSize(pointi);
1336  }
1337 
1338  // Extract points
1339  pointField newPoints(newToOldPoints.size());
1340  forAll(newToOldPoints, i)
1341  {
1342  newPoints[i] = s.points()[newToOldPoints[i]];
1343  }
1344  // Extract faces
1345  List<labelledTri> newTriangles(newToOldFaces.size());
1346 
1347  forAll(newToOldFaces, i)
1348  {
1349  // Get old vertex labels
1350  const labelledTri& tri = s[newToOldFaces[i]];
1351 
1352  newTriangles[i][0] = oldToNewPoints[tri[0]];
1353  newTriangles[i][1] = oldToNewPoints[tri[1]];
1354  newTriangles[i][2] = oldToNewPoints[tri[2]];
1355  newTriangles[i].region() = tri.region();
1356  }
1357 
1358  // Reuse storage.
1359  return triSurface(newTriangles, s.patches(), newPoints, true);
1360 }
1361 
1362 
1363 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1364 
1367  const polyMesh& mesh,
1368  const scalarField& cVals,
1369  const scalarField& pVals,
1370  const scalar iso,
1371  const bool regularise,
1372  const scalar mergeTol
1373 )
1374 :
1375  mesh_(mesh),
1376  cVals_(cVals),
1377  pVals_(pVals),
1378  iso_(iso),
1379  mergeDistance_(mergeTol*mesh.bounds().mag())
1380 {
1381  if (debug)
1382  {
1383  Pout<< "isoSurfaceCell : mergeTol:" << mergeTol
1384  << " mesh span:" << mesh.bounds().mag()
1385  << " mergeDistance:" << mergeDistance_ << endl;
1386  }
1387 
1388  // Determine if cell is tet
1389  PackedBoolList isTet(mesh_.nCells());
1390  {
1391  tetMatcher tet;
1392 
1393  forAll(isTet, celli)
1394  {
1395  if (tet.isA(mesh_, celli))
1396  {
1397  isTet.set(celli, 1);
1398  }
1399  }
1400  }
1401 
1402 
1403  // Determine if any cut through cell
1404  calcCutTypes(isTet, cVals, pVals);
1405 
1406  DynamicList<point> snappedPoints(nCutCells_);
1407 
1408  // Per cc -1 or a point inside snappedPoints.
1409  labelList snappedCc;
1410  if (regularise)
1411  {
1412  calcSnappedCc
1413  (
1414  isTet,
1415  cVals,
1416  pVals,
1417  snappedPoints,
1418  snappedCc
1419  );
1420  }
1421  else
1422  {
1423  snappedCc.setSize(mesh_.nCells());
1424  snappedCc = -1;
1425  }
1426 
1427  if (debug)
1428  {
1429  Pout<< "isoSurfaceCell : shifted " << snappedPoints.size()
1430  << " cell centres to intersection." << endl;
1431  }
1432 
1433  snappedPoints.shrink();
1434  label nCellSnaps = snappedPoints.size();
1435 
1436  // Per point -1 or a point inside snappedPoints.
1437  labelList snappedPoint;
1438  if (regularise)
1439  {
1440  calcSnappedPoint
1441  (
1442  isTet,
1443  cVals,
1444  pVals,
1445  snappedPoints,
1446  snappedPoint
1447  );
1448  }
1449  else
1450  {
1451  snappedPoint.setSize(mesh_.nPoints());
1452  snappedPoint = -1;
1453  }
1454 
1455  if (debug)
1456  {
1457  Pout<< "isoSurfaceCell : shifted " << snappedPoints.size()-nCellSnaps
1458  << " vertices to intersection." << endl;
1459  }
1460 
1461 
1462  {
1463  DynamicList<point> triPoints(nCutCells_);
1464  DynamicList<label> triMeshCells(nCutCells_);
1465 
1466  generateTriPoints
1467  (
1468  cVals,
1469  pVals,
1470 
1471  mesh_.cellCentres(),
1472  mesh_.points(),
1473 
1474  snappedPoints,
1475  snappedCc,
1476  snappedPoint,
1477 
1478  triPoints,
1479  triMeshCells
1480  );
1481 
1482  if (debug)
1483  {
1484  Pout<< "isoSurfaceCell : generated " << triMeshCells.size()
1485  << " unmerged triangles." << endl;
1486  }
1487 
1488  // Merge points and compact out non-valid triangles
1489  labelList triMap;
1491  (
1492  stitchTriPoints
1493  (
1494  regularise, // check for duplicate tris
1495  triPoints,
1496  triPointMergeMap_, // unmerged to merged point
1497  triMap // merged to unmerged triangle
1498  )
1499  );
1500 
1501  if (debug)
1502  {
1503  Pout<< "isoSurfaceCell : generated " << triMap.size()
1504  << " merged triangles." << endl;
1505  }
1506 
1507  meshCells_.setSize(triMap.size());
1508  forAll(triMap, i)
1509  {
1510  meshCells_[i] = triMeshCells[triMap[i]];
1511  }
1512  }
1513 
1514 
1515  if (debug)
1516  {
1517  Pout<< "isoSurfaceCell : checking " << size()
1518  << " triangles for validity." << endl;
1519 
1520  forAll(*this, triI)
1521  {
1522  // Copied from surfaceCheck
1523  validTri(*this, triI);
1524  }
1525  }
1526 
1527 
1528  if (regularise)
1529  {
1530  List<FixedList<label, 3>> faceEdges;
1531  labelList edgeFace0, edgeFace1;
1532  Map<labelList> edgeFacesRest;
1533 
1534 
1535  while (true)
1536  {
1537  // Calculate addressing
1538  calcAddressing
1539  (
1540  *this,
1541  faceEdges,
1542  edgeFace0,
1543  edgeFace1,
1544  edgeFacesRest
1545  );
1546 
1547  // See if any dangling triangles
1548  boolList keepTriangles;
1549  label nDangling = markDanglingTriangles
1550  (
1551  faceEdges,
1552  edgeFace0,
1553  edgeFace1,
1554  edgeFacesRest,
1555  keepTriangles
1556  );
1557 
1558  if (debug)
1559  {
1560  Pout<< "isoSurfaceCell : detected " << nDangling
1561  << " dangling triangles." << endl;
1562  }
1563 
1564  if (nDangling == 0)
1565  {
1566  break;
1567  }
1568 
1569  // Create face map (new to old)
1570  labelList subsetTriMap(findIndices(keepTriangles, true));
1571 
1572  labelList subsetPointMap;
1573  labelList reversePointMap;
1575  (
1576  subsetMesh
1577  (
1578  *this,
1579  subsetTriMap,
1580  reversePointMap,
1581  subsetPointMap
1582  )
1583  );
1584  meshCells_ = labelField(meshCells_, subsetTriMap);
1585  inplaceRenumber(reversePointMap, triPointMergeMap_);
1586  }
1587  }
1588 }
1589 
1590 
1591 // ************************************************************************* //
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:49
A cellMatcher for tet cells.
Definition: tetMatcher.H:51
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
static void syncPointPositions(const polyMesh &mesh, List< point > &l, const CombineOp &cop, const point &nullValue)
Synchronize locations on all mesh points.
Definition: syncTools.H:196
friend Ostream & operator(Ostream &, const UList< T > &)
List< geometricSurfacePatch > geometricSurfacePatchList
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
HashTable< T, label, Hash< label > >::iterator iterator
Definition: Map.H:56
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
labelList findIndices(const ListType &, typename ListType::const_reference, const label start=0)
Find all occurences of given element. Linear search.
Macros for easy insertion into run-time selection tables.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().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))
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
const pointField & points
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
face triFace(3)
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:90
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:91
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
dimensioned< scalar > magSqr(const dimensioned< Type > &)
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:240
defineTypeNameAndDebug(combustionModel, 0)
Merge points. See below.
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
virtual bool isA(const primitiveMesh &mesh, const label celli)
Exact match. Uses faceSizeMatch.
Definition: tetMatcher.C:218
isoSurfaceCell(const polyMesh &mesh, const scalarField &cellValues, const scalarField &pointValues, const scalar iso, const bool regularise, const scalar mergeTol=1e-6)
Construct from dictionary.
void setSize(const label)
Reset size of List.
Definition: List.C:281
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:430
A bit-packed bool list.
label patchi
vector point
Point is a vector.
Definition: point.H:41
#define WarningInFunction
Report a warning using Foam::Warning.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
label mergePoints(const UList< Type > &points, const scalar mergeTol, const bool verbose, labelList &pointMap, const Type &origin=Type::zero)
Sorts and merges points. All points closer than/equal mergeTol get merged.
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Triangulated surface description with patch information.
Definition: triSurface.H:65
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
Namespace for OpenFOAM.
A HashTable to objects of type <T> with a label key.
Definition: Map.H:49