isoSurfaceCell.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "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  const vector n0 = tri0.normal(localPoints);
289  const vector n1 = tri1.normal(localPoints);
290 
291  if ((n0 & n1) < 0)
292  {
293  // Too big an angle. Do not simplify.
294  }
295  else
296  {
297  info.setPoint
298  (
299  0.5
300  * (
301  tri0.centre(localPoints)
302  + tri1.centre(localPoints)
303  )
304  );
305  info.setHit();
306  }
307  }
308  }
309  else if (localTris.size())
310  {
311  // Check if single region. Rare situation.
312  triSurface surf
313  (
314  localTris,
316  localPoints,
317  true
318  );
319  localTris.clearStorage();
320 
321  labelList faceZone;
322  label nZones = surf.markZones
323  (
324  boolList(surf.nEdges(), false),
325  faceZone
326  );
327 
328  if (nZones == 1)
329  {
330  // Check that all normals make a decent angle
331  scalar minCos = great;
332  const vector& n0 = surf.faceNormals()[0];
333  for (label i = 1; i < surf.size(); i++)
334  {
335  scalar cosAngle = (n0 & surf.faceNormals()[i]);
336  if (cosAngle < minCos)
337  {
338  minCos = cosAngle;
339  }
340  }
341 
342  if (minCos > 0)
343  {
344  info.setPoint(calcCentre(surf));
345  info.setHit();
346  }
347  }
348  }
349 
350  return info;
351 }
352 
353 
354 void Foam::isoSurfaceCell::calcSnappedCc
355 (
356  const PackedBoolList& isTet,
357  const scalarField& cVals,
358  const scalarField& pVals,
359 
360  DynamicList<point>& snappedPoints,
361  labelList& snappedCc
362 ) const
363 {
364  const pointField& cc = mesh_.cellCentres();
365  const pointField& pts = mesh_.points();
366 
367  snappedCc.setSize(mesh_.nCells());
368  snappedCc = -1;
369 
370  // Work arrays
371  DynamicList<point, 64> localPoints(64);
372  DynamicList<labelledTri, 64> localTris(64);
373  Map<label> pointToLocal(64);
374 
375  forAll(mesh_.cells(), celli)
376  {
377  if (cellCutType_[celli] == CUT && isTet.get(celli) == 0)
378  {
379  scalar cVal = cVals[celli];
380 
381  const cell& cFaces = mesh_.cells()[celli];
382 
383  localPoints.clear();
384  localTris.clear();
385  pointToLocal.clear();
386 
387  // Create points for all intersections close to cell centre
388  // (i.e. from pyramid edges)
389 
390  forAll(cFaces, cFacei)
391  {
392  const face& f = mesh_.faces()[cFaces[cFacei]];
393 
394  forAll(f, fp)
395  {
396  label pointi = f[fp];
397 
398  scalar s = isoFraction(cVal, pVals[pointi]);
399 
400  if (s >= 0.0 && s <= 0.5)
401  {
402  if (pointToLocal.insert(pointi, localPoints.size()))
403  {
404  localPoints.append((1.0-s)*cc[celli]+s*pts[pointi]);
405  }
406  }
407  }
408  }
409 
410  if (localPoints.size() == 1)
411  {
412  // No need for any analysis.
413  snappedCc[celli] = snappedPoints.size();
414  snappedPoints.append(localPoints[0]);
415 
416  // Pout<< "cell:" << celli
417  // << " at " << mesh_.cellCentres()[celli]
418  // << " collapsing " << localPoints
419  // << " intersections down to "
420  // << snappedPoints[snappedCc[celli]] << endl;
421  }
422  else if (localPoints.size() == 2)
423  {
424  //? No need for any analysis.???
425  snappedCc[celli] = snappedPoints.size();
426  snappedPoints.append(0.5*(localPoints[0]+localPoints[1]));
427 
428  // Pout<< "cell:" << celli
429  // << " at " << mesh_.cellCentres()[celli]
430  // << " collapsing " << localPoints
431  // << " intersections down to "
432  // << snappedPoints[snappedCc[celli]] << endl;
433  }
434  else if (localPoints.size())
435  {
436  // Need to analyse
437  forAll(cFaces, cFacei)
438  {
439  label facei = cFaces[cFacei];
440  const face& f = mesh_.faces()[facei];
441 
442  // Do a tetrahedralisation. Each face to cc becomes pyr.
443  // Each pyr gets split into tets by diagonalisation
444  // of face.
445 
446  const label fp0 = mesh_.tetBasePtIs()[facei];
447  label fp = f.fcIndex(fp0);
448  for (label i = 2; i < f.size(); i++)
449  {
450  label nextFp = f.fcIndex(fp);
451  triFace tri(f[fp0], f[fp], f[nextFp]);
452 
453  // Get fractions for the three edges to cell centre
454  FixedList<scalar, 3> s(3);
455  s[0] = isoFraction(cVal, pVals[tri[0]]);
456  s[1] = isoFraction(cVal, pVals[tri[1]]);
457  s[2] = isoFraction(cVal, pVals[tri[2]]);
458 
459  if
460  (
461  (s[0] >= 0.0 && s[0] <= 0.5)
462  && (s[1] >= 0.0 && s[1] <= 0.5)
463  && (s[2] >= 0.0 && s[2] <= 0.5)
464  )
465  {
466  if
467  (
468  (mesh_.faceOwner()[facei] == celli)
469  == (cVal >= pVals[tri[0]])
470  )
471  {
472  localTris.append
473  (
474  labelledTri
475  (
476  pointToLocal[tri[1]],
477  pointToLocal[tri[0]],
478  pointToLocal[tri[2]],
479  0
480  )
481  );
482  }
483  else
484  {
485  localTris.append
486  (
487  labelledTri
488  (
489  pointToLocal[tri[0]],
490  pointToLocal[tri[1]],
491  pointToLocal[tri[2]],
492  0
493  )
494  );
495  }
496  }
497 
498  fp = nextFp;
499  }
500  }
501 
502  pointField surfPoints;
503  surfPoints.transfer(localPoints);
504  pointIndexHit info = collapseSurface
505  (
506  celli,
507  surfPoints,
508  localTris
509  );
510 
511  if (info.hit())
512  {
513  snappedCc[celli] = snappedPoints.size();
514  snappedPoints.append(info.hitPoint());
515 
516  // Pout<< "cell:" << celli
517  // << " at " << mesh_.cellCentres()[celli]
518  // << " collapsing " << surfPoints
519  // << " intersections down to "
520  // << snappedPoints[snappedCc[celli]] << endl;
521  }
522  }
523  }
524  }
525 }
526 
527 
528 void Foam::isoSurfaceCell::genPointTris
529 (
530  const scalarField& cellValues,
531  const scalarField& pointValues,
532  const label pointi,
533  const label facei,
534  const label celli,
535  DynamicList<point, 64>& localTriPoints
536 ) const
537 {
538  const pointField& cc = mesh_.cellCentres();
539  const pointField& pts = mesh_.points();
540  const face& f = mesh_.faces()[facei];
541 
542  const label fp0 = mesh_.tetBasePtIs()[facei];
543  label fp = f.fcIndex(fp0);
544  for (label i = 2; i < f.size(); i++)
545  {
546  label nextFp = f.fcIndex(fp);
547  triFace tri(f[fp0], f[fp], f[nextFp]);
548 
549  label index = findIndex(tri, pointi);
550 
551  if (index == -1)
552  {
553  continue;
554  }
555 
556  // Tet between index..index-1, index..index+1, index..cc
557  label b = tri[tri.fcIndex(index)];
558  label c = tri[tri.rcIndex(index)];
559 
560  // Get fractions for the three edges emanating from point
561  FixedList<scalar, 3> s(3);
562  s[0] = isoFraction(pointValues[pointi], pointValues[b]);
563  s[1] = isoFraction(pointValues[pointi], pointValues[c]);
564  s[2] = isoFraction(pointValues[pointi], cellValues[celli]);
565 
566  if
567  (
568  (s[0] >= 0.0 && s[0] <= 0.5)
569  && (s[1] >= 0.0 && s[1] <= 0.5)
570  && (s[2] >= 0.0 && s[2] <= 0.5)
571  )
572  {
573  point p0 = (1.0-s[0])*pts[pointi] + s[0]*pts[b];
574  point p1 = (1.0-s[1])*pts[pointi] + s[1]*pts[c];
575  point p2 = (1.0-s[2])*pts[pointi] + s[2]*cc[celli];
576 
577  if
578  (
579  (mesh_.faceOwner()[facei] == celli)
580  == (pointValues[pointi] > cellValues[celli])
581  )
582  {
583  localTriPoints.append(p0);
584  localTriPoints.append(p1);
585  localTriPoints.append(p2);
586  }
587  else
588  {
589  localTriPoints.append(p1);
590  localTriPoints.append(p0);
591  localTriPoints.append(p2);
592  }
593  }
594 
595  fp = nextFp;
596  }
597 }
598 
599 
600 void Foam::isoSurfaceCell::genPointTris
601 (
602  const scalarField& pointValues,
603  const label pointi,
604  const label facei,
605  const label celli,
606  DynamicList<point, 64>& localTriPoints
607 ) const
608 {
609  const pointField& pts = mesh_.points();
610  const cell& cFaces = mesh_.cells()[celli];
611 
612  FixedList<label, 4> tet;
613 
614  // Make tet from this face to the 4th point (same as cellcentre in
615  // non-tet cells)
616  const face& f = mesh_.faces()[facei];
617 
618  // Find 4th point
619  label ccPointi = -1;
620  forAll(cFaces, cFacei)
621  {
622  const face& f1 = mesh_.faces()[cFaces[cFacei]];
623  forAll(f1, fp)
624  {
625  label p1 = f1[fp];
626 
627  if (findIndex(f, p1) == -1)
628  {
629  ccPointi = p1;
630  break;
631  }
632  }
633  if (ccPointi != -1)
634  {
635  break;
636  }
637  }
638 
639 
640  // Tet between index..index-1, index..index+1, index..cc
641  label index = findIndex(f, pointi);
642  label b = f[f.fcIndex(index)];
643  label c = f[f.rcIndex(index)];
644 
645  // Pout<< " p0:" << pointi << " b:" << b << " c:" << c
646  //<< " d:" << ccPointi << endl;
647 
648  // Get fractions for the three edges emanating from point
649  FixedList<scalar, 3> s(3);
650  s[0] = isoFraction(pointValues[pointi], pointValues[b]);
651  s[1] = isoFraction(pointValues[pointi], pointValues[c]);
652  s[2] = isoFraction(pointValues[pointi], pointValues[ccPointi]);
653 
654  if
655  (
656  (s[0] >= 0.0 && s[0] <= 0.5)
657  && (s[1] >= 0.0 && s[1] <= 0.5)
658  && (s[2] >= 0.0 && s[2] <= 0.5)
659  )
660  {
661  point p0 = (1.0-s[0])*pts[pointi] + s[0]*pts[b];
662  point p1 = (1.0-s[1])*pts[pointi] + s[1]*pts[c];
663  point p2 = (1.0-s[2])*pts[pointi] + s[2]*pts[ccPointi];
664 
665  if (mesh_.faceOwner()[facei] != celli)
666  {
667  localTriPoints.append(p0);
668  localTriPoints.append(p1);
669  localTriPoints.append(p2);
670  }
671  else
672  {
673  localTriPoints.append(p1);
674  localTriPoints.append(p0);
675  localTriPoints.append(p2);
676  }
677  }
678 }
679 
680 
681 void Foam::isoSurfaceCell::calcSnappedPoint
682 (
683  const PackedBoolList& isTet,
684  const scalarField& cVals,
685  const scalarField& pVals,
686 
687  DynamicList<point>& snappedPoints,
688  labelList& snappedPoint
689 ) const
690 {
691  // Determine if point is on boundary. Points on boundaries are never
692  // snapped. Coupled boundaries are handled explicitly so not marked here.
693  PackedBoolList isBoundaryPoint(mesh_.nPoints());
694  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
695  forAll(patches, patchi)
696  {
697  const polyPatch& pp = patches[patchi];
698 
699  if (!pp.coupled())
700  {
701  label facei = pp.start();
702  forAll(pp, i)
703  {
704  const face& f = mesh_.faces()[facei++];
705 
706  forAll(f, fp)
707  {
708  isBoundaryPoint.set(f[fp], 1);
709  }
710  }
711  }
712  }
713 
714 
715  const point greatPoint(great, great, great);
716 
717  pointField collapsedPoint(mesh_.nPoints(), greatPoint);
718 
719 
720  // Work arrays
721  DynamicList<point, 64> localTriPoints(100);
722  labelHashSet localPointCells(100);
723 
724  forAll(mesh_.pointFaces(), pointi)
725  {
726  if (isBoundaryPoint.get(pointi) == 1)
727  {
728  continue;
729  }
730 
731  const labelList& pFaces = mesh_.pointFaces()[pointi];
732 
733  bool anyCut = false;
734 
735  forAll(pFaces, i)
736  {
737  label facei = pFaces[i];
738 
739  if
740  (
741  cellCutType_[mesh_.faceOwner()[facei]] == CUT
742  || (
743  mesh_.isInternalFace(facei)
744  && cellCutType_[mesh_.faceNeighbour()[facei]] == CUT
745  )
746  )
747  {
748  anyCut = true;
749  break;
750  }
751  }
752 
753  if (!anyCut)
754  {
755  continue;
756  }
757 
758 
759  // Do a pointCells walk (by using pointFaces)
760 
761  localPointCells.clear();
762  localTriPoints.clear();
763 
764  forAll(pFaces, pFacei)
765  {
766  label facei = pFaces[pFacei];
767  label own = mesh_.faceOwner()[facei];
768 
769  if (isTet.get(own) == 1)
770  {
771  // Since tets have no cell centre to include make sure
772  // we only generate a triangle once per point.
773  if (localPointCells.insert(own))
774  {
775  genPointTris(pVals, pointi, facei, own, localTriPoints);
776  }
777  }
778  else
779  {
780  genPointTris
781  (
782  cVals,
783  pVals,
784  pointi,
785  facei,
786  own,
787  localTriPoints
788  );
789  }
790 
791  if (mesh_.isInternalFace(facei))
792  {
793  label nei = mesh_.faceNeighbour()[facei];
794 
795  if (isTet.get(nei) == 1)
796  {
797  if (localPointCells.insert(nei))
798  {
799  genPointTris(pVals, pointi, facei, nei, localTriPoints);
800  }
801  }
802  else
803  {
804  genPointTris
805  (
806  cVals,
807  pVals,
808  pointi,
809  facei,
810  nei,
811  localTriPoints
812  );
813  }
814  }
815  }
816 
817  if (localTriPoints.size() == 3)
818  {
819  // Single triangle. No need for any analysis. Average points.
821  points.transfer(localTriPoints);
822  collapsedPoint[pointi] = sum(points)/points.size();
823 
824  // Pout<< " point:" << pointi
825  // << " replacing coord:" << mesh_.points()[pointi]
826  // << " by average:" << collapsedPoint[pointi] << endl;
827  }
828  else if (localTriPoints.size())
829  {
830  // Convert points into triSurface.
831 
832  // Merge points and compact out non-valid triangles
833  labelList triMap; // merged to unmerged triangle
834  labelList triPointReverseMap; // unmerged to merged point
835  triSurface surf
836  (
837  stitchTriPoints
838  (
839  false, // do not check for duplicate tris
840  localTriPoints,
841  triPointReverseMap,
842  triMap
843  )
844  );
845 
846  labelList faceZone;
847  label nZones = surf.markZones
848  (
849  boolList(surf.nEdges(), false),
850  faceZone
851  );
852 
853  if (nZones == 1)
854  {
855  // Check that all normals make a decent angle
856  scalar minCos = great;
857  const vector& n0 = surf.faceNormals()[0];
858  for (label i = 1; i < surf.size(); i++)
859  {
860  const vector& n = surf.faceNormals()[i];
861  scalar cosAngle = (n0 & n);
862  if (cosAngle < minCos)
863  {
864  minCos = cosAngle;
865  }
866  }
867  if (minCos > 0)
868  {
869  collapsedPoint[pointi] = calcCentre(surf);
870  }
871  }
872  }
873  }
874 
876  (
877  mesh_,
878  collapsedPoint,
879  minMagSqrEqOp<point>(),
880  greatPoint
881  );
882 
883  snappedPoint.setSize(mesh_.nPoints());
884  snappedPoint = -1;
885 
886  forAll(collapsedPoint, pointi)
887  {
888  // Cannot do == comparison since might be transformed so have
889  // truncation errors.
890  if (magSqr(collapsedPoint[pointi]) < 0.5*magSqr(greatPoint))
891  {
892  snappedPoint[pointi] = snappedPoints.size();
893  snappedPoints.append(collapsedPoint[pointi]);
894  }
895  }
896 }
897 
898 
899 Foam::triSurface Foam::isoSurfaceCell::stitchTriPoints
900 (
901  const bool checkDuplicates,
902  const List<point>& triPoints,
903  labelList& triPointReverseMap, // unmerged to merged point
904  labelList& triMap // merged to unmerged triangle
905 ) const
906 {
907  label nTris = triPoints.size()/3;
908 
909  if ((triPoints.size() % 3) != 0)
910  {
912  << "Problem: number of points " << triPoints.size()
913  << " not a multiple of 3." << abort(FatalError);
914  }
915 
916  pointField newPoints;
918  (
919  triPoints,
920  mergeDistance_,
921  false,
922  triPointReverseMap,
923  newPoints
924  );
925 
926  // Check that enough merged.
927  if (debug)
928  {
929  Pout<< "isoSurfaceCell : merged from " << triPoints.size()
930  << " points down to " << newPoints.size() << endl;
931 
932  pointField newNewPoints;
933  labelList oldToNew;
934  bool hasMerged = mergePoints
935  (
936  newPoints,
937  mergeDistance_,
938  true,
939  oldToNew,
940  newNewPoints
941  );
942 
943  if (hasMerged)
944  {
946  << "Merged points contain duplicates"
947  << " when merging with distance " << mergeDistance_ << endl
948  << "merged:" << newPoints.size() << " re-merged:"
949  << newNewPoints.size()
950  << abort(FatalError);
951  }
952  }
953 
954 
955  List<labelledTri> tris;
956  {
957  DynamicList<labelledTri> dynTris(nTris);
958  label rawPointi = 0;
959  DynamicList<label> newToOldTri(nTris);
960 
961  for (label oldTriI = 0; oldTriI < nTris; oldTriI++)
962  {
963  labelledTri tri
964  (
965  triPointReverseMap[rawPointi],
966  triPointReverseMap[rawPointi+1],
967  triPointReverseMap[rawPointi+2],
968  0
969  );
970  if ((tri[0] != tri[1]) && (tri[0] != tri[2]) && (tri[1] != tri[2]))
971  {
972  newToOldTri.append(oldTriI);
973  dynTris.append(tri);
974  }
975 
976  rawPointi += 3;
977  }
978 
979  triMap.transfer(newToOldTri);
980  tris.transfer(dynTris);
981  }
982 
983 
984  // Use face centres to determine 'flat hole' situation (see RMT paper).
985  // Two unconnected triangles get connected because (some of) the edges
986  // separating them get collapsed. Below only checks for duplicate triangles,
987  // not non-manifold edge connectivity.
988  if (checkDuplicates)
989  {
990  if (debug)
991  {
992  Pout<< "isoSurfaceCell : merged from " << nTris
993  << " down to " << tris.size() << " triangles." << endl;
994  }
995 
996  pointField centres(tris.size());
997  forAll(tris, triI)
998  {
999  centres[triI] = tris[triI].centre(newPoints);
1000  }
1001 
1002  pointField mergedCentres;
1003  labelList oldToMerged;
1004  bool hasMerged = mergePoints
1005  (
1006  centres,
1007  mergeDistance_,
1008  false,
1009  oldToMerged,
1010  mergedCentres
1011  );
1012 
1013  if (debug)
1014  {
1015  Pout<< "isoSurfaceCell : detected "
1016  << centres.size()-mergedCentres.size()
1017  << " duplicate triangles." << endl;
1018  }
1019 
1020  if (hasMerged)
1021  {
1022  // Filter out duplicates.
1023  label newTriI = 0;
1024  DynamicList<label> newToOldTri(tris.size());
1025  labelList newToMaster(mergedCentres.size(), -1);
1026  forAll(tris, triI)
1027  {
1028  label mergedI = oldToMerged[triI];
1029 
1030  if (newToMaster[mergedI] == -1)
1031  {
1032  newToMaster[mergedI] = triI;
1033  newToOldTri.append(triMap[triI]);
1034  tris[newTriI++] = tris[triI];
1035  }
1036  }
1037 
1038  triMap.transfer(newToOldTri);
1039  tris.setSize(newTriI);
1040  }
1041  }
1042 
1043  return triSurface(tris, geometricSurfacePatchList(0), newPoints, true);
1044 }
1045 
1046 
1047 bool Foam::isoSurfaceCell::validTri(const triSurface& surf, const label facei)
1048 {
1049  // Simple check on indices ok.
1050 
1051  const labelledTri& f = surf[facei];
1052 
1053  forAll(f, fp)
1054  {
1055  if (f[fp] < 0 || f[fp] >= surf.points().size())
1056  {
1058  << "triangle " << facei << " vertices " << f
1059  << " uses point indices outside point range 0.."
1060  << surf.points().size()-1 << endl;
1061 
1062  return false;
1063  }
1064  }
1065 
1066  if ((f[0] == f[1]) || (f[0] == f[2]) || (f[1] == f[2]))
1067  {
1069  << "triangle " << facei
1070  << " uses non-unique vertices " << f
1071  << endl;
1072  return false;
1073  }
1074 
1075  // duplicate triangle check
1076 
1077  const labelList& fFaces = surf.faceFaces()[facei];
1078 
1079  // Check if faceNeighbours use same points as this face.
1080  // Note: discards normal information - sides of baffle are merged.
1081  forAll(fFaces, i)
1082  {
1083  label nbrFacei = fFaces[i];
1084 
1085  if (nbrFacei <= facei)
1086  {
1087  // lower numbered faces already checked
1088  continue;
1089  }
1090 
1091  const labelledTri& nbrF = surf[nbrFacei];
1092 
1093  if
1094  (
1095  ((f[0] == nbrF[0]) || (f[0] == nbrF[1]) || (f[0] == nbrF[2]))
1096  && ((f[1] == nbrF[0]) || (f[1] == nbrF[1]) || (f[1] == nbrF[2]))
1097  && ((f[2] == nbrF[0]) || (f[2] == nbrF[1]) || (f[2] == nbrF[2]))
1098  )
1099  {
1101  << "triangle " << facei << " vertices " << f
1102  << " coords:" << f.points(surf.points())
1103  << " has the same vertices as triangle " << nbrFacei
1104  << " vertices " << nbrF
1105  << endl;
1106 
1107  return false;
1108  }
1109  }
1110  return true;
1111 }
1112 
1113 
1114 void Foam::isoSurfaceCell::calcAddressing
1115 (
1116  const triSurface& surf,
1117  List<FixedList<label, 3>>& faceEdges,
1118  labelList& edgeFace0,
1119  labelList& edgeFace1,
1120  Map<labelList>& edgeFacesRest
1121 ) const
1122 {
1123  const pointField& points = surf.points();
1124 
1125  pointField edgeCentres(3*surf.size());
1126  label edgeI = 0;
1127  forAll(surf, triI)
1128  {
1129  const labelledTri& tri = surf[triI];
1130  edgeCentres[edgeI++] = 0.5*(points[tri[0]]+points[tri[1]]);
1131  edgeCentres[edgeI++] = 0.5*(points[tri[1]]+points[tri[2]]);
1132  edgeCentres[edgeI++] = 0.5*(points[tri[2]]+points[tri[0]]);
1133  }
1134 
1135  pointField mergedCentres;
1136  labelList oldToMerged;
1137  bool hasMerged = mergePoints
1138  (
1139  edgeCentres,
1140  mergeDistance_,
1141  false,
1142  oldToMerged,
1143  mergedCentres
1144  );
1145 
1146  if (debug)
1147  {
1148  Pout<< "isoSurfaceCell : detected "
1149  << mergedCentres.size()
1150  << " edges on " << surf.size() << " triangles." << endl;
1151  }
1152 
1153  if (!hasMerged)
1154  {
1155  return;
1156  }
1157 
1158 
1159  // Determine faceEdges
1160  faceEdges.setSize(surf.size());
1161  edgeI = 0;
1162  forAll(surf, triI)
1163  {
1164  faceEdges[triI][0] = oldToMerged[edgeI++];
1165  faceEdges[triI][1] = oldToMerged[edgeI++];
1166  faceEdges[triI][2] = oldToMerged[edgeI++];
1167  }
1168 
1169 
1170  // Determine edgeFaces
1171  edgeFace0.setSize(mergedCentres.size());
1172  edgeFace0 = -1;
1173  edgeFace1.setSize(mergedCentres.size());
1174  edgeFace1 = -1;
1175  edgeFacesRest.clear();
1176 
1177  forAll(oldToMerged, oldEdgeI)
1178  {
1179  label triI = oldEdgeI / 3;
1180  label edgeI = oldToMerged[oldEdgeI];
1181 
1182  if (edgeFace0[edgeI] == -1)
1183  {
1184  edgeFace0[edgeI] = triI;
1185  }
1186  else if (edgeFace1[edgeI] == -1)
1187  {
1188  edgeFace1[edgeI] = triI;
1189  }
1190  else
1191  {
1192  // WarningInFunction
1193  // << "Edge " << edgeI << " with centre " << mergedCentres[edgeI]
1194  // << " used by more than two triangles: " << edgeFace0[edgeI]
1195  // << ", "
1196  // << edgeFace1[edgeI] << " and " << triI << endl;
1197  Map<labelList>::iterator iter = edgeFacesRest.find(edgeI);
1198 
1199  if (iter != edgeFacesRest.end())
1200  {
1201  labelList& eFaces = iter();
1202  label sz = eFaces.size();
1203  eFaces.setSize(sz+1);
1204  eFaces[sz] = triI;
1205  }
1206  else
1207  {
1208  edgeFacesRest.insert(edgeI, labelList(1, triI));
1209  }
1210  }
1211  }
1212 }
1213 
1214 
1215 bool Foam::isoSurfaceCell::danglingTriangle
1216 (
1217  const FixedList<label, 3>& fEdges,
1218  const labelList& edgeFace1
1219 )
1220 {
1221  label nOpen = 0;
1222  forAll(fEdges, i)
1223  {
1224  if (edgeFace1[fEdges[i]] == -1)
1225  {
1226  nOpen++;
1227  }
1228  }
1229 
1230  if (nOpen == 1 || nOpen == 2 || nOpen == 3)
1231  {
1232  return true;
1233  }
1234  else
1235  {
1236  return false;
1237  }
1238 }
1239 
1240 
1241 Foam::label Foam::isoSurfaceCell::markDanglingTriangles
1242 (
1243  const List<FixedList<label, 3>>& faceEdges,
1244  const labelList& edgeFace0,
1245  const labelList& edgeFace1,
1246  const Map<labelList>& edgeFacesRest,
1247  boolList& keepTriangles
1248 )
1249 {
1250  keepTriangles.setSize(faceEdges.size());
1251  keepTriangles = true;
1252 
1253  label nDangling = 0;
1254 
1255  // Remove any dangling triangles
1256  forAllConstIter(Map<labelList>, edgeFacesRest, iter)
1257  {
1258  // These are all the non-manifold edges. Filter out all triangles
1259  // with only one connected edge (= this edge)
1260 
1261  label edgeI = iter.key();
1262  const labelList& otherEdgeFaces = iter();
1263 
1264  // Remove all dangling triangles
1265  if (danglingTriangle(faceEdges[edgeFace0[edgeI]], edgeFace1))
1266  {
1267  keepTriangles[edgeFace0[edgeI]] = false;
1268  nDangling++;
1269  }
1270  if (danglingTriangle(faceEdges[edgeFace1[edgeI]], edgeFace1))
1271  {
1272  keepTriangles[edgeFace1[edgeI]] = false;
1273  nDangling++;
1274  }
1275  forAll(otherEdgeFaces, i)
1276  {
1277  label triI = otherEdgeFaces[i];
1278  if (danglingTriangle(faceEdges[triI], edgeFace1))
1279  {
1280  keepTriangles[triI] = false;
1281  nDangling++;
1282  }
1283  }
1284  }
1285  return nDangling;
1286 }
1287 
1288 
1289 Foam::triSurface Foam::isoSurfaceCell::subsetMesh
1290 (
1291  const triSurface& s,
1292  const labelList& newToOldFaces,
1293  labelList& oldToNewPoints,
1294  labelList& newToOldPoints
1295 )
1296 {
1297  const boolList include
1298  (
1299  createWithValues<boolList>
1300  (
1301  s.size(),
1302  false,
1303  newToOldFaces,
1304  true
1305  )
1306  );
1307 
1308  newToOldPoints.setSize(s.points().size());
1309  oldToNewPoints.setSize(s.points().size());
1310  oldToNewPoints = -1;
1311  {
1312  label pointi = 0;
1313 
1314  forAll(include, oldFacei)
1315  {
1316  if (include[oldFacei])
1317  {
1318  // Renumber labels for face
1319  const triSurface::FaceType& f = s[oldFacei];
1320 
1321  forAll(f, fp)
1322  {
1323  label oldPointi = f[fp];
1324 
1325  if (oldToNewPoints[oldPointi] == -1)
1326  {
1327  oldToNewPoints[oldPointi] = pointi;
1328  newToOldPoints[pointi++] = oldPointi;
1329  }
1330  }
1331  }
1332  }
1333  newToOldPoints.setSize(pointi);
1334  }
1335 
1336  // Extract points
1337  pointField newPoints(newToOldPoints.size());
1338  forAll(newToOldPoints, i)
1339  {
1340  newPoints[i] = s.points()[newToOldPoints[i]];
1341  }
1342  // Extract faces
1343  List<labelledTri> newTriangles(newToOldFaces.size());
1344 
1345  forAll(newToOldFaces, i)
1346  {
1347  // Get old vertex labels
1348  const labelledTri& tri = s[newToOldFaces[i]];
1349 
1350  newTriangles[i][0] = oldToNewPoints[tri[0]];
1351  newTriangles[i][1] = oldToNewPoints[tri[1]];
1352  newTriangles[i][2] = oldToNewPoints[tri[2]];
1353  newTriangles[i].region() = tri.region();
1354  }
1355 
1356  // Reuse storage.
1357  return triSurface(newTriangles, s.patches(), newPoints, true);
1358 }
1359 
1360 
1361 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1362 
1365  const polyMesh& mesh,
1366  const scalarField& cVals,
1367  const scalarField& pVals,
1368  const scalar iso,
1369  const bool regularise,
1370  const scalar mergeTol
1371 )
1372 :
1373  mesh_(mesh),
1374  cVals_(cVals),
1375  pVals_(pVals),
1376  iso_(iso),
1377  mergeDistance_(mergeTol*mesh.bounds().mag())
1378 {
1379  if (debug)
1380  {
1381  Pout<< "isoSurfaceCell : mergeTol:" << mergeTol
1382  << " mesh span:" << mesh.bounds().mag()
1383  << " mergeDistance:" << mergeDistance_ << endl;
1384  }
1385 
1386  // Determine if cell is tet
1387  PackedBoolList isTet(mesh_.nCells());
1388  {
1389  tetMatcher tet;
1390 
1391  forAll(isTet, celli)
1392  {
1393  if (tet.isA(mesh_, celli))
1394  {
1395  isTet.set(celli, 1);
1396  }
1397  }
1398  }
1399 
1400 
1401  // Determine if any cut through cell
1402  calcCutTypes(isTet, cVals, pVals);
1403 
1404  DynamicList<point> snappedPoints(nCutCells_);
1405 
1406  // Per cc -1 or a point inside snappedPoints.
1407  labelList snappedCc;
1408  if (regularise)
1409  {
1410  calcSnappedCc
1411  (
1412  isTet,
1413  cVals,
1414  pVals,
1415  snappedPoints,
1416  snappedCc
1417  );
1418  }
1419  else
1420  {
1421  snappedCc.setSize(mesh_.nCells());
1422  snappedCc = -1;
1423  }
1424 
1425  if (debug)
1426  {
1427  Pout<< "isoSurfaceCell : shifted " << snappedPoints.size()
1428  << " cell centres to intersection." << endl;
1429  }
1430 
1431  snappedPoints.shrink();
1432  label nCellSnaps = snappedPoints.size();
1433 
1434  // Per point -1 or a point inside snappedPoints.
1435  labelList snappedPoint;
1436  if (regularise)
1437  {
1438  calcSnappedPoint
1439  (
1440  isTet,
1441  cVals,
1442  pVals,
1443  snappedPoints,
1444  snappedPoint
1445  );
1446  }
1447  else
1448  {
1449  snappedPoint.setSize(mesh_.nPoints());
1450  snappedPoint = -1;
1451  }
1452 
1453  if (debug)
1454  {
1455  Pout<< "isoSurfaceCell : shifted " << snappedPoints.size()-nCellSnaps
1456  << " vertices to intersection." << endl;
1457  }
1458 
1459 
1460  {
1461  DynamicList<point> triPoints(nCutCells_);
1462  DynamicList<label> triMeshCells(nCutCells_);
1463 
1464  generateTriPoints
1465  (
1466  cVals,
1467  pVals,
1468 
1469  mesh_.cellCentres(),
1470  mesh_.points(),
1471 
1472  snappedPoints,
1473  snappedCc,
1474  snappedPoint,
1475 
1476  triPoints,
1477  triMeshCells
1478  );
1479 
1480  if (debug)
1481  {
1482  Pout<< "isoSurfaceCell : generated " << triMeshCells.size()
1483  << " unmerged triangles." << endl;
1484  }
1485 
1486  // Merge points and compact out non-valid triangles
1487  labelList triMap;
1489  (
1490  stitchTriPoints
1491  (
1492  regularise, // check for duplicate tris
1493  triPoints,
1494  triPointMergeMap_, // unmerged to merged point
1495  triMap // merged to unmerged triangle
1496  )
1497  );
1498 
1499  if (debug)
1500  {
1501  Pout<< "isoSurfaceCell : generated " << triMap.size()
1502  << " merged triangles." << endl;
1503  }
1504 
1505  meshCells_.setSize(triMap.size());
1506  forAll(triMap, i)
1507  {
1508  meshCells_[i] = triMeshCells[triMap[i]];
1509  }
1510  }
1511 
1512 
1513  if (debug)
1514  {
1515  Pout<< "isoSurfaceCell : checking " << size()
1516  << " triangles for validity." << endl;
1517 
1518  forAll(*this, triI)
1519  {
1520  // Copied from surfaceCheck
1521  validTri(*this, triI);
1522  }
1523  }
1524 
1525 
1526  if (regularise)
1527  {
1528  List<FixedList<label, 3>> faceEdges;
1529  labelList edgeFace0, edgeFace1;
1530  Map<labelList> edgeFacesRest;
1531 
1532 
1533  while (true)
1534  {
1535  // Calculate addressing
1536  calcAddressing
1537  (
1538  *this,
1539  faceEdges,
1540  edgeFace0,
1541  edgeFace1,
1542  edgeFacesRest
1543  );
1544 
1545  // See if any dangling triangles
1546  boolList keepTriangles;
1547  label nDangling = markDanglingTriangles
1548  (
1549  faceEdges,
1550  edgeFace0,
1551  edgeFace1,
1552  edgeFacesRest,
1553  keepTriangles
1554  );
1555 
1556  if (debug)
1557  {
1558  Pout<< "isoSurfaceCell : detected " << nDangling
1559  << " dangling triangles." << endl;
1560  }
1561 
1562  if (nDangling == 0)
1563  {
1564  break;
1565  }
1566 
1567  // Create face map (new to old)
1568  labelList subsetTriMap(findIndices(keepTriangles, true));
1569 
1570  labelList subsetPointMap;
1571  labelList reversePointMap;
1573  (
1574  subsetMesh
1575  (
1576  *this,
1577  subsetTriMap,
1578  reversePointMap,
1579  subsetPointMap
1580  )
1581  );
1582  meshCells_ = labelField(meshCells_, subsetTriMap);
1583  inplaceRenumber(reversePointMap, triPointMergeMap_);
1584  }
1585  }
1586 }
1587 
1588 
1589 // ************************************************************************* //
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:256
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 occurrences 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:97
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Templated 3D Vector derived from VectorSpace adding construction from 3 components, element access using x(), y() and z() member functions and the inner-product (dot-product) and cross product operators.
Definition: Vector.H:57
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 occurrence 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:66
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