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