isoSurface.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-2019 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 "isoSurface.H"
27 #include "polyMesh.H"
28 #include "tetMatcher.H"
29 #include "tetPointRef.H"
30 #include "DynamicField.H"
31 #include "syncTools.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(isoSurface, 0);
39 }
40 
41 namespace Foam
42 {
43  template<>
45  {"none", "partial", "full"};
46 }
47 
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
53 bool Foam::isoSurface::isTriCut
54 (
55  const triFace& tri,
56  const scalarField& pointValues
57 ) const
58 {
59  bool aLower = (pointValues[tri[0]] < iso_);
60  bool bLower = (pointValues[tri[1]] < iso_);
61  bool cLower = (pointValues[tri[2]] < iso_);
62 
63  return !(aLower == bLower && aLower == cLower);
64 }
65 
66 
67 Foam::isoSurface::cellCutType Foam::isoSurface::calcCutType
68 (
69  const bool isTet,
70  const label celli
71 ) const
72 {
73  const cell& cFaces = mesh_.cells()[celli];
74 
75  if (isTet)
76  {
77  forAll(cFaces, cFacei)
78  {
79  const face& f = mesh_.faces()[cFaces[cFacei]];
80 
81  for (label fp = 1; fp < f.size() - 1; fp++)
82  {
83  triFace tri(f[0], f[fp], f[f.fcIndex(fp)]);
84 
85  if (isTriCut(tri, pVals_))
86  {
87  return cellCutType::cut;
88  }
89  }
90  }
91  return cellCutType::notCut;
92  }
93  else
94  {
95  bool cellLower = (cVals_[celli] < iso_);
96 
97  // First check if there is any cut in cell
98  bool edgeCut = false;
99 
100  forAll(cFaces, cFacei)
101  {
102  label facei = cFaces[cFacei];
103  const face& f = mesh_.faces()[facei];
104 
105  // Check pyramids cut
106  forAll(f, fp)
107  {
108  if ((pVals_[f[fp]] < iso_) != cellLower)
109  {
110  edgeCut = true;
111  break;
112  }
113  }
114 
115  if (edgeCut)
116  {
117  break;
118  }
119 
120  const label fp0 = tetBasePtIs_[facei] < 0 ? 0 : tetBasePtIs_[facei];
121 
122  label fp = f.fcIndex(fp0);
123  for (label i = 2; i < f.size(); i++)
124  {
125  label nextFp = f.fcIndex(fp);
126 
127  if (isTriCut(triFace(f[fp0], f[fp], f[nextFp]), pVals_))
128  {
129  edgeCut = true;
130  break;
131  }
132 
133  fp = nextFp;
134  }
135 
136  if (edgeCut)
137  {
138  break;
139  }
140  }
141 
142  if (edgeCut)
143  {
144  // Count actual cuts (expensive since addressing needed)
145  // Note: not needed if you don't want to preserve maxima/minima
146  // centred around cellcentre. In that case just always return cut
147 
148  const labelList& cPoints = mesh_.cellPoints(celli);
149 
150  label nPyrCuts = 0;
151 
152  forAll(cPoints, i)
153  {
154  if ((pVals_[cPoints[i]] < iso_) != cellLower)
155  {
156  nPyrCuts++;
157  }
158  }
159 
160  if (nPyrCuts == cPoints.size())
161  {
162  return cellCutType::sphere;
163  }
164  else
165  {
166  return cellCutType::cut;
167  }
168  }
169  else
170  {
171  return cellCutType::notCut;
172  }
173  }
174 }
175 
176 
177 Foam::label Foam::isoSurface::calcCutTypes
178 (
179  tetMatcher& tet,
180  List<cellCutType>& cellCutTypes
181 )
182 {
183  const cellList& cells = mesh_.cells();
184 
185  cellCutTypes.setSize(cells.size());
186  label nCutCells = 0;
187  forAll(cells, celli)
188  {
189  cellCutTypes[celli] = calcCutType(tet.isA(mesh_, celli), celli);
190 
191  if (cellCutTypes[celli] == cellCutType::cut)
192  {
193  nCutCells++;
194  }
195  }
196 
197  if (debug)
198  {
199  Pout<< "isoSurface : detected " << nCutCells
200  << " candidate cut cells." << endl;
201  }
202  return nCutCells;
203 }
204 
205 
206 Foam::scalar Foam::isoSurface::minTetQ
207 (
208  const label facei,
209  const label faceBasePtI
210 ) const
211 {
213  (
214  mesh_,
215  mesh_.cellCentres()[mesh_.faceOwner()[facei]],
216  facei,
217  true,
218  faceBasePtI
219  );
220 
221  if (mesh_.isInternalFace(facei))
222  {
223  q = min
224  (
225  q,
227  (
228  mesh_,
229  mesh_.cellCentres()[mesh_.faceNeighbour()[facei]],
230  facei,
231  false,
232  faceBasePtI
233  )
234  );
235  }
236  return q;
237 }
238 
239 
240 void Foam::isoSurface::fixTetBasePtIs()
241 {
242  // Determine points used by two faces on the same cell
243  const cellList& cells = mesh_.cells();
244  const faceList& faces = mesh_.faces();
245  const labelList& faceOwner = mesh_.faceOwner();
246  const labelList& faceNeighbour = mesh_.faceNeighbour();
247 
248 
249  // Get face triangulation base point
250  tetBasePtIs_ = mesh_.tetBasePtIs();
251 
252 
253  // Mark all cells with illegal base points as potentially problematic
254  PackedBoolList problemCells(cells.size(), false);
255  forAll(tetBasePtIs_, facei)
256  {
257  if (tetBasePtIs_[facei] == -1)
258  {
259  problemCells[faceOwner[facei]] = true;
260  if (mesh_.isInternalFace(facei))
261  {
262  problemCells[faceNeighbour[facei]] = true;
263  }
264  }
265  }
266 
267 
268  // Mark all points which are shared by just two faces within an adjacent
269  // problem cell as problematic
270  PackedBoolList problemPoints(mesh_.points().size(), false);
271  forAll(cells, celli)
272  {
273  if (problemCells[celli])
274  {
275  const cell& cFaces = cells[celli];
276 
277  Map<label> pointCount;
278  forAll(cFaces, i)
279  {
280  const label facei = cFaces[i];
281  const face& f = faces[facei];
282  forAll(f, fp)
283  {
284  const label pointi = f[fp];
285 
286  Map<label>::iterator pointFnd = pointCount.find(pointi);
287  if (pointFnd == pointCount.end())
288  {
289  pointCount.insert(pointi, 1);
290  }
291  else
292  {
293  ++ pointFnd();
294  }
295  }
296  }
297 
298  forAllConstIter(Map<label>, pointCount, iter)
299  {
300  if (iter() == 1)
301  {
302  FatalErrorInFunction << "point:" << iter.key()
303  << " at:" << mesh_.points()[iter.key()]
304  << " only used by one face" << exit(FatalError);
305  }
306  if (iter() == 2)
307  {
308  problemPoints[iter.key()] = true;
309  }
310  }
311  }
312  }
313 
314 
315  // For all faces which form a part of a problem-cell, check if the base
316  // point is adjacent to any problem points. If it is, re-calculate the base
317  // point so that it is not.
318  label nAdapted = 0;
319  forAll(tetBasePtIs_, facei)
320  {
321  if
322  (
323  problemCells[faceOwner[facei]]
324  || (mesh_.isInternalFace(facei) && problemCells[faceNeighbour[facei]])
325  )
326  {
327  const face& f = faces[facei];
328 
329  // Check if either of the points adjacent to the base point is a
330  // problem point. If not, the existing base point can be retained.
331  const label fp0 = tetBasePtIs_[facei] < 0 ? 0 : tetBasePtIs_[facei];
332 
333  const bool prevPointIsProblem = problemPoints[f[f.rcIndex(fp0)]];
334  const bool nextPointIsProblem = problemPoints[f[f.fcIndex(fp0)]];
335 
336  if (!prevPointIsProblem && !nextPointIsProblem)
337  {
338  continue;
339  }
340 
341  // A new base point is required. Pick the point that results in the
342  // least-worst-worst tet, and which is not adjacent to any problem
343  // points.
344  scalar maxQ = -GREAT;
345  label maxFp = -1;
346  forAll(f, fp)
347  {
348  const bool prevPointIsProblem = problemPoints[f[f.rcIndex(fp)]];
349  const bool nextPointIsProblem = problemPoints[f[f.fcIndex(fp)]];
350 
351  if (!prevPointIsProblem && !nextPointIsProblem)
352  {
353  const scalar q = minTetQ(facei, fp);
354  if (q > maxQ)
355  {
356  maxQ = q;
357  maxFp = fp;
358  }
359  }
360  }
361 
362  if (maxFp != -1)
363  {
364  // Success! Set the new base point
365  tetBasePtIs_[facei] = maxFp;
366  }
367  else
368  {
369  // No point was found on face that would not result in some
370  // duplicate triangle. Do what? Continue and hope? Spit an
371  // error? Silently or noisily reduce the filtering level?
372  }
373 
374  ++ nAdapted;
375  }
376  }
377 
378  if (debug)
379  {
380  Pout<< "isoSurface : adapted starting point of triangulation on "
381  << nAdapted << " faces." << endl;
382  }
383 }
384 
385 
386 Foam::label Foam::isoSurface::generatePoint
387 (
388  const label facei,
389  const bool edgeIsDiag,
390  const edge& vertices,
391 
392  DynamicList<edge>& pointToVerts,
393  DynamicList<label>& pointToFace,
394  DynamicList<bool>& pointFromDiag,
395  EdgeMap<label>& vertsToPoint
396 ) const
397 {
398  EdgeMap<label>::const_iterator edgeFnd = vertsToPoint.find(vertices);
399  if (edgeFnd != vertsToPoint.end())
400  {
401  return edgeFnd();
402  }
403  else
404  {
405  // Generate new point
406  label pointi = pointToVerts.size();
407 
408  pointToVerts.append(vertices);
409  pointToFace.append(facei);
410  pointFromDiag.append(edgeIsDiag);
411  vertsToPoint.insert(vertices, pointi);
412  return pointi;
413  }
414 }
415 
416 
417 void Foam::isoSurface::generateTriPoints
418 (
419  const label facei,
420  const FixedList<scalar, 4>& s,
421  const FixedList<point, 4>& p,
422  const FixedList<label, 4>& pIndex,
423  const FixedList<bool, 6>& edgeIsDiag,// Per tet edge whether is face diag
424 
425  DynamicList<edge>& pointToVerts,
426  DynamicList<label>& pointToFace,
427  DynamicList<bool>& pointFromDiag,
428 
429  EdgeMap<label>& vertsToPoint,
430  DynamicList<label>& verts // Every three verts is new triangle
431 ) const
432 {
433  int triIndex = 0;
434  if (s[0] < iso_)
435  {
436  triIndex |= 1;
437  }
438  if (s[1] < iso_)
439  {
440  triIndex |= 2;
441  }
442  if (s[2] < iso_)
443  {
444  triIndex |= 4;
445  }
446  if (s[3] < iso_)
447  {
448  triIndex |= 8;
449  }
450 
451  // Form the vertices of the triangles for each case
452  switch (triIndex)
453  {
454  case 0x00:
455  case 0x0F:
456  break;
457 
458  case 0x01:
459  case 0x0E:
460  {
461  verts.append
462  (
463  generatePoint
464  (
465  facei,
466  edgeIsDiag[0],
467  edge(pIndex[0], pIndex[1]),
468  pointToVerts, pointToFace, pointFromDiag, vertsToPoint
469  )
470  );
471  verts.append
472  (
473  generatePoint
474  (
475  facei,
476  edgeIsDiag[1],
477  edge(pIndex[0], pIndex[2]),
478  pointToVerts, pointToFace, pointFromDiag, vertsToPoint
479  )
480  );
481  verts.append
482  (
483  generatePoint
484  (
485  facei,
486  edgeIsDiag[2],
487  edge(pIndex[0], pIndex[3]),
488  pointToVerts, pointToFace, pointFromDiag, vertsToPoint
489  )
490  );
491 
492  if (triIndex == 0x0E)
493  {
494  // Flip normals
495  label sz = verts.size();
496  Swap(verts[sz-2], verts[sz-1]);
497  }
498  }
499  break;
500 
501  case 0x02:
502  case 0x0D:
503  {
504  verts.append
505  (
506  generatePoint
507  (
508  facei,
509  edgeIsDiag[0],
510  edge(pIndex[1], pIndex[0]),
511  pointToVerts, pointToFace, pointFromDiag, vertsToPoint
512  )
513  );
514  verts.append
515  (
516  generatePoint
517  (
518  facei,
519  edgeIsDiag[3],
520  edge(pIndex[1], pIndex[3]),
521  pointToVerts, pointToFace, pointFromDiag, vertsToPoint
522  )
523  );
524  verts.append
525  (
526  generatePoint
527  (
528  facei,
529  edgeIsDiag[4],
530  edge(pIndex[1], pIndex[2]),
531  pointToVerts, pointToFace, pointFromDiag, vertsToPoint
532  )
533  );
534 
535  if (triIndex == 0x0D)
536  {
537  // Flip normals
538  label sz = verts.size();
539  Swap(verts[sz-2], verts[sz-1]);
540  }
541  }
542  break;
543 
544  case 0x03:
545  case 0x0C:
546  {
547  label p0p2
548  (
549  generatePoint
550  (
551  facei,
552  edgeIsDiag[1],
553  edge(pIndex[0], pIndex[2]),
554  pointToVerts, pointToFace, pointFromDiag, vertsToPoint
555  )
556  );
557  label p1p3
558  (
559  generatePoint
560  (
561  facei,
562  edgeIsDiag[3],
563  edge(pIndex[1], pIndex[3]),
564  pointToVerts, pointToFace, pointFromDiag, vertsToPoint
565  )
566  );
567 
568  verts.append
569  (
570  generatePoint
571  (
572  facei,
573  edgeIsDiag[2],
574  edge(pIndex[0], pIndex[3]),
575  pointToVerts, pointToFace, pointFromDiag, vertsToPoint
576  )
577  );
578  verts.append(p1p3);
579  verts.append(p0p2);
580  verts.append(p1p3);
581  verts.append
582  (
583  generatePoint
584  (
585  facei,
586  edgeIsDiag[4],
587  edge(pIndex[1], pIndex[2]),
588  pointToVerts, pointToFace, pointFromDiag, vertsToPoint
589  )
590  );
591  verts.append(p0p2);
592 
593  if (triIndex == 0x0C)
594  {
595  // Flip normals
596  label sz = verts.size();
597  Swap(verts[sz-5], verts[sz-4]);
598  Swap(verts[sz-2], verts[sz-1]);
599  }
600  }
601  break;
602 
603  case 0x04:
604  case 0x0B:
605  {
606  verts.append
607  (
608  generatePoint
609  (
610  facei,
611  edgeIsDiag[1],
612  edge(pIndex[2], pIndex[0]),
613  pointToVerts, pointToFace, pointFromDiag, vertsToPoint
614  )
615  );
616  verts.append
617  (
618  generatePoint
619  (
620  facei,
621  edgeIsDiag[4],
622  edge(pIndex[2], pIndex[1]),
623  pointToVerts, pointToFace, pointFromDiag, vertsToPoint
624  )
625  );
626  verts.append
627  (
628  generatePoint
629  (
630  facei,
631  edgeIsDiag[5],
632  edge(pIndex[2], pIndex[3]),
633  pointToVerts, pointToFace, pointFromDiag, vertsToPoint
634  )
635  );
636 
637  if (triIndex == 0x0B)
638  {
639  // Flip normals
640  label sz = verts.size();
641  Swap(verts[sz-2], verts[sz-1]);
642  }
643  }
644  break;
645 
646  case 0x05:
647  case 0x0A:
648  {
649  label p0p1
650  (
651  generatePoint
652  (
653  facei,
654  edgeIsDiag[0],
655  edge(pIndex[0], pIndex[1]),
656  pointToVerts, pointToFace, pointFromDiag, vertsToPoint
657  )
658  );
659  label p2p3
660  (
661  generatePoint
662  (
663  facei,
664  edgeIsDiag[5],
665  edge(pIndex[2], pIndex[3]),
666  pointToVerts, pointToFace, pointFromDiag, vertsToPoint
667  )
668  );
669 
670  verts.append(p0p1);
671  verts.append(p2p3);
672  verts.append
673  (
674  generatePoint
675  (
676  facei,
677  edgeIsDiag[2],
678  edge(pIndex[0], pIndex[3]),
679  pointToVerts, pointToFace, pointFromDiag, vertsToPoint
680  )
681  );
682  verts.append(p0p1);
683  verts.append
684  (
685  generatePoint
686  (
687  facei,
688  edgeIsDiag[4],
689  edge(pIndex[1], pIndex[2]),
690  pointToVerts, pointToFace, pointFromDiag, vertsToPoint
691  )
692  );
693  verts.append(p2p3);
694 
695  if (triIndex == 0x0A)
696  {
697  // Flip normals
698  label sz = verts.size();
699  Swap(verts[sz-5], verts[sz-4]);
700  Swap(verts[sz-2], verts[sz-1]);
701  }
702  }
703  break;
704 
705  case 0x06:
706  case 0x09:
707  {
708  label p0p1
709  (
710  generatePoint
711  (
712  facei,
713  edgeIsDiag[0],
714  edge(pIndex[0], pIndex[1]),
715  pointToVerts, pointToFace, pointFromDiag, vertsToPoint
716  )
717  );
718  label p2p3
719  (
720  generatePoint
721  (
722  facei,
723  edgeIsDiag[5],
724  edge(pIndex[2], pIndex[3]),
725  pointToVerts, pointToFace, pointFromDiag, vertsToPoint
726  )
727  );
728 
729  verts.append(p0p1);
730  verts.append
731  (
732  generatePoint
733  (
734  facei,
735  edgeIsDiag[3],
736  edge(pIndex[1], pIndex[3]),
737  pointToVerts, pointToFace, pointFromDiag, vertsToPoint
738  )
739  );
740  verts.append(p2p3);
741  verts.append(p0p1);
742  verts.append(p2p3);
743  verts.append
744  (
745  generatePoint
746  (
747  facei,
748  edgeIsDiag[1],
749  edge(pIndex[0], pIndex[2]),
750  pointToVerts, pointToFace, pointFromDiag, vertsToPoint
751  )
752  );
753 
754  if (triIndex == 0x09)
755  {
756  // Flip normals
757  label sz = verts.size();
758  Swap(verts[sz-5], verts[sz-4]);
759  Swap(verts[sz-2], verts[sz-1]);
760  }
761  }
762  break;
763 
764  case 0x08:
765  case 0x07:
766  {
767  verts.append
768  (
769  generatePoint
770  (
771  facei,
772  edgeIsDiag[2],
773  edge(pIndex[3], pIndex[0]),
774  pointToVerts, pointToFace, pointFromDiag, vertsToPoint
775  )
776  );
777  verts.append
778  (
779  generatePoint
780  (
781  facei,
782  edgeIsDiag[5],
783  edge(pIndex[3], pIndex[2]),
784  pointToVerts, pointToFace, pointFromDiag, vertsToPoint
785  )
786  );
787  verts.append
788  (
789  generatePoint
790  (
791  facei,
792  edgeIsDiag[3],
793  edge(pIndex[3], pIndex[1]),
794  pointToVerts, pointToFace, pointFromDiag, vertsToPoint
795  )
796  );
797  if (triIndex == 0x07)
798  {
799  // Flip normals
800  label sz = verts.size();
801  Swap(verts[sz-2], verts[sz-1]);
802  }
803  }
804  break;
805  }
806 }
807 
808 
809 void Foam::isoSurface::generateTriPoints
810 (
811  const label celli,
812  const bool isTet,
813 
814  DynamicList<edge>& pointToVerts,
815  DynamicList<label>& pointToFace,
816  DynamicList<bool>& pointFromDiag,
817 
818  EdgeMap<label>& vertsToPoint,
819  DynamicList<label>& verts,
820  DynamicList<label>& faceLabels
821 ) const
822 {
823  const faceList& faces = mesh_.faces();
824  const labelList& faceOwner = mesh_.faceOwner();
825  const pointField& points = mesh_.points();
826  const cell& cFaces = mesh_.cells()[celli];
827 
828  if (isTet)
829  {
830  // For tets don't do cell-centre decomposition, just use the
831  // tet points and values
832 
833  label facei = cFaces[0];
834  const face& f0 = faces[facei];
835 
836  // Get the other point
837  const face& f1 = faces[cFaces[1]];
838  label oppositeI = -1;
839  forAll(f1, fp)
840  {
841  oppositeI = f1[fp];
842  if (findIndex(f0, oppositeI) == -1)
843  {
844  break;
845  }
846  }
847 
848 
849  label p0 = f0[0];
850  label p1 = f0[1];
851  label p2 = f0[2];
852  FixedList<bool, 6> edgeIsDiag(false);
853 
854  if (faceOwner[facei] == celli)
855  {
856  Swap(p1, p2);
857  }
858 
859  tetPointRef tet
860  (
861  points[p0],
862  points[p1],
863  points[p2],
864  points[oppositeI]
865  );
866 
867  label startTrii = verts.size();
868  generateTriPoints
869  (
870  facei,
872  ({
873  pVals_[p0],
874  pVals_[p1],
875  pVals_[p2],
876  pVals_[oppositeI]
877  }),
879  ({
880  points[p0],
881  points[p1],
882  points[p2],
883  points[oppositeI]
884  }),
886  ({
887  p0,
888  p1,
889  p2,
890  oppositeI
891  }),
892  edgeIsDiag,
893 
894  pointToVerts,
895  pointToFace,
896  pointFromDiag,
897  vertsToPoint,
898  verts // Every three verts is new triangle
899  );
900 
901  label nTris = (verts.size()-startTrii)/3;
902  for (label i = 0; i < nTris; i++)
903  {
904  faceLabels.append(facei);
905  }
906  }
907  else
908  {
909  const pointField& cellCentres = mesh_.cellCentres();
910 
911  forAll(cFaces, cFacei)
912  {
913  label facei = cFaces[cFacei];
914  const face& f = faces[facei];
915 
916  const label fp0 = tetBasePtIs_[facei] < 0 ? 0 : tetBasePtIs_[facei];
917 
918  label startTrii = verts.size();
919 
920  label fp = f.fcIndex(fp0);
921  for (label i = 2; i < f.size(); i++)
922  {
923  label nextFp = f.fcIndex(fp);
924 
925  FixedList<bool, 6> edgeIsDiag(false);
926 
927  label p0 = f[fp0];
928  label p1 = f[fp];
929  label p2 = f[nextFp];
930  if (faceOwner[facei] == celli)
931  {
932  Swap(p1, p2);
933  if (i != 2) edgeIsDiag[1] = true;
934  if (i != f.size()-1) edgeIsDiag[0] = true;
935  }
936  else
937  {
938  if (i != 2) edgeIsDiag[0] = true;
939  if (i != f.size()-1) edgeIsDiag[1] = true;
940  }
941 
942  tetPointRef tet
943  (
944  points[p0],
945  points[p1],
946  points[p2],
947  cellCentres[celli]
948  );
949 
950  generateTriPoints
951  (
952  facei,
954  ({
955  pVals_[p0],
956  pVals_[p1],
957  pVals_[p2],
958  cVals_[celli]
959  }),
961  ({
962  points[p0],
963  points[p1],
964  points[p2],
965  cellCentres[celli]
966  }),
968  ({
969  p0,
970  p1,
971  p2,
972  mesh_.nPoints()+celli
973  }),
974  edgeIsDiag,
975 
976  pointToVerts,
977  pointToFace,
978  pointFromDiag,
979  vertsToPoint,
980  verts // Every three verts is new triangle
981  );
982 
983  fp = nextFp;
984  }
985 
986  label nTris = (verts.size()-startTrii)/3;
987  for (label i = 0; i < nTris; i++)
988  {
989  faceLabels.append(facei);
990  }
991  }
992  }
993 }
994 
995 
996 void Foam::isoSurface::triangulateOutside
997 (
998  const bool filterDiag,
999  const primitivePatch& pp,
1000  const boolList& pointFromDiag,
1001  const labelList& pointToFace,
1002  const label cellID,
1003 
1004  DynamicList<face>& compactFaces,
1005  DynamicList<label>& compactCellIDs
1006 ) const
1007 {
1008  // We can form pockets:
1009  // - 1. triangle on face
1010  // - 2. multiple triangles on interior (from diag edges)
1011  // - the edge loop will be pocket since it is only the diag
1012  // edges that give it volume?
1013 
1014  // Retriangulate the exterior loops
1015  const labelListList& edgeLoops = pp.edgeLoops();
1016  const labelList& mp = pp.meshPoints();
1017 
1018  forAll(edgeLoops, loopi)
1019  {
1020  const labelList& loop = edgeLoops[loopi];
1021 
1022  if (loop.size() > 2)
1023  {
1024  compactFaces.append(face(0));
1025  face& f = compactFaces.last();
1026 
1027  f.setSize(loop.size());
1028  label fpi = 0;
1029  forAll(f, i)
1030  {
1031  label pointi = mp[loop[i]];
1032  if (filterDiag && pointFromDiag[pointi])
1033  {
1034  label prevPointi = mp[loop[loop.fcIndex(i)]];
1035  if
1036  (
1037  pointFromDiag[prevPointi]
1038  && (pointToFace[pointi] != pointToFace[prevPointi])
1039  )
1040  {
1041  f[fpi++] = pointi;
1042  }
1043  else
1044  {
1045  // Filter out diagonal point
1046  }
1047  }
1048  else
1049  {
1050  f[fpi++] = pointi;
1051  }
1052  }
1053 
1054  if (fpi > 2)
1055  {
1056  f.setSize(fpi);
1057  }
1058  else
1059  {
1060  // Keep original face
1061  forAll(f, i)
1062  {
1063  label pointi = mp[loop[i]];
1064  f[i] = pointi;
1065  }
1066  }
1067  compactCellIDs.append(cellID);
1068  }
1069  }
1070 }
1071 
1072 
1073 Foam::MeshedSurface<Foam::face> Foam::isoSurface::removeInsidePoints
1074 (
1075  const bool filterDiag,
1076  const MeshedSurface<face>& s,
1077  const boolList& pointFromDiag,
1078  const labelList& pointToFace,
1079  const labelList& start, // Per cell the starting triangle
1080  DynamicList<label>& pointCompactMap, // Per returned point the original
1081  DynamicList<label>& compactCellIDs // Per returned tri the cellID
1082 ) const
1083 {
1084  const pointField& points = s.points();
1085 
1086  pointCompactMap.clear();
1087  compactCellIDs.clear();
1088 
1089  DynamicList<face> compactFaces(s.size()/8);
1090 
1091  for (label celli = 0; celli < start.size()-1; celli++)
1092  {
1093  // All triangles of the current cell
1094 
1095  label nTris = start[celli+1]-start[celli];
1096 
1097  if (nTris)
1098  {
1099  const SubList<face> cellFaces(s, nTris, start[celli]);
1100  const primitivePatch pp(cellFaces, points);
1101 
1102  triangulateOutside
1103  (
1104  filterDiag,
1105  pp,
1106  pointFromDiag,
1107  pointToFace,
1108  //protectedFace,
1109  celli,
1110 
1111  compactFaces,
1112  compactCellIDs
1113  );
1114  }
1115  }
1116 
1117 
1118  // Compact out unused points
1119  // Pick up the used vertices
1120  labelList oldToCompact(points.size(), -1);
1121  DynamicField<point> compactPoints(points.size());
1122  pointCompactMap.clear();
1123 
1124  forAll(compactFaces, i)
1125  {
1126  face& f = compactFaces[i];
1127  forAll(f, fp)
1128  {
1129  label pointi = f[fp];
1130  label compacti = oldToCompact[pointi];
1131  if (compacti == -1)
1132  {
1133  compacti = compactPoints.size();
1134  oldToCompact[pointi] = compacti;
1135  compactPoints.append(points[pointi]);
1136  pointCompactMap.append(pointi);
1137  }
1138  f[fp] = compacti;
1139  }
1140  }
1141 
1142 
1143  MeshedSurface<face> cpSurface;
1144  cpSurface.reset
1145  (
1146  move(compactPoints),
1147  move(compactFaces),
1148  surfZoneList(s.surfZones())
1149  );
1150 
1151  return cpSurface;
1152 }
1153 
1154 
1155 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1156 
1159  const polyMesh& mesh,
1160  const scalarField& cVals,
1161  const scalarField& pVals,
1162  const scalar iso,
1163  const filterType filter
1164 )
1165 :
1166  mesh_(mesh),
1167  cVals_(cVals),
1168  pVals_(pVals),
1169  iso_(iso)
1170 {
1171  if (debug)
1172  {
1173  Pout<< "isoSurface : iso:" << iso_
1174  << " filter:" << filterTypeNames_[filter] << endl;
1175  }
1176 
1177  fixTetBasePtIs();
1178 
1179  tetMatcher tet;
1180 
1181  // Determine if any cut through cell
1182  List<cellCutType> cellCutTypes;
1183  const label nCutCells = calcCutTypes(tet, cellCutTypes);
1184 
1185  // Per cell: 5 pyramids cut, each generating 2 triangles
1186  // - pointToVerts : from generated iso point to originating mesh verts
1187  DynamicList<edge> pointToVerts(10*nCutCells);
1188  // - pointToFace : from generated iso point to originating mesh face
1189  DynamicList<label> pointToFace(10*nCutCells);
1190  // - pointToFace : from generated iso point whether is on face diagonal
1191  DynamicList<bool> pointFromDiag(10*nCutCells);
1192 
1193  // Per cell: number of intersected edges:
1194  // - four faces cut so 4 mesh edges + 4 face-diagonal edges
1195  // - 4 of the pyramid edges
1196  EdgeMap<label> vertsToPoint(12*nCutCells);
1197  DynamicList<label> verts(12*nCutCells);
1198  // Per cell: 5 pyramids cut (since only one pyramid not cut)
1199  DynamicList<label> faceLabels(5*nCutCells);
1200  DynamicList<label> cellLabels(5*nCutCells);
1201 
1202 
1203  labelList startTri(mesh_.nCells()+1, 0);
1204 
1205  for (label celli = 0; celli < mesh_.nCells(); celli++)
1206  {
1207  startTri[celli] = faceLabels.size();
1208  if (cellCutTypes[celli] != cellCutType::notCut)
1209  {
1210  generateTriPoints
1211  (
1212  celli,
1213  tet.isA(mesh_, celli),
1214 
1215  pointToVerts,
1216  pointToFace,
1217  pointFromDiag,
1218 
1219  vertsToPoint,
1220  verts,
1221  faceLabels
1222  );
1223 
1224  for (label i = startTri[celli]; i < faceLabels.size(); i++)
1225  {
1226  cellLabels.append(celli);
1227  }
1228  }
1229  }
1230  startTri[mesh_.nCells()] = faceLabels.size();
1231 
1232 
1233  pointToVerts_.transfer(pointToVerts);
1234  meshCells_.transfer(cellLabels);
1235  pointToFace_.transfer(pointToFace);
1236 
1238  (
1239  interpolate
1240  (
1241  mesh_.cellCentres(),
1242  mesh_.points()
1243  )
1244  );
1245 
1246 
1247  // Assign to MeshedSurface
1248  faceList allTris(faceLabels.size());
1249  label verti = 0;
1250  forAll(allTris, i)
1251  {
1252  allTris[i].setSize(3);
1253  allTris[i][0] = verts[verti++];
1254  allTris[i][1] = verts[verti++];
1255  allTris[i][2] = verts[verti++];
1256  }
1257 
1258 
1259  surfZoneList allZones(1);
1260  allZones[0] = surfZone
1261  (
1262  "allFaces",
1263  allTris.size(), // Size
1264  0, // Start
1265  0 // Index
1266  );
1267 
1269  (
1270  move(allPoints.ref()),
1271  move(allTris),
1272  move(allZones)
1273  );
1274 
1275 
1276  // Now:
1277  // - generated faces and points are assigned to *this
1278  // - per point we know:
1279  // - pointOnDiag: whether it is on a face-diagonal edge
1280  // - pointToFace_: from what pyramid (cell+face) it was produced
1281  // (note that the pyramid faces are shared between multiple mesh faces)
1282  // - pointToVerts_ : originating mesh vertex or cell centre
1283 
1284 
1285  if (debug)
1286  {
1287  Pout<< "isoSurface : generated " << size() << " faces." << endl;
1288  }
1289 
1290 
1291  if (filter != filterType::none)
1292  {
1293  // Triangulate outside (filter edges to cell centres and optionally
1294  // face diagonals)
1295  DynamicList<label> pointCompactMap; // Back to original point
1296  DynamicList<label> compactCellIDs; // Per returned tri the cellID
1298  (
1299  removeInsidePoints
1300  (
1301  (filter == filterType::full ? true : false),
1302  *this,
1303  pointFromDiag,
1304  pointToFace_,
1305  startTri,
1306  pointCompactMap,
1307  compactCellIDs
1308  )
1309  );
1310 
1311  pointToVerts_ = UIndirectList<edge>(pointToVerts_, pointCompactMap)();
1312  pointToFace_ = UIndirectList<label>(pointToFace_, pointCompactMap)();
1313  pointFromDiag = UIndirectList<bool>(pointFromDiag, pointCompactMap)();
1314  meshCells_.transfer(compactCellIDs);
1315 
1316  if (debug)
1317  {
1318  Pout<< "isoSurface :"
1319  << " after removing cell centre and face-diag triangles : "
1320  << size() << endl;
1321  }
1322 
1323 
1324  if (filter == filterType::full)
1325  {
1326  // We remove verts on face diagonals. This is in fact just
1327  // straightening the edges of the face through the cell. This can
1328  // close off 'pockets' of triangles and create open or
1329  // multiply-connected triangles
1330 
1331  // Solved by eroding open-edges
1332  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1333 
1334 
1335  // Mark points on mesh outside. Note that we extend with nCells
1336  // so we can easily index with pointToVerts_.
1337  PackedBoolList isBoundaryPoint(mesh.nPoints() + mesh.nCells());
1338  for
1339  (
1340  label facei = mesh.nInternalFaces();
1341  facei < mesh.nFaces();
1342  facei++
1343  )
1344  {
1345  isBoundaryPoint.set(mesh.faces()[facei]);
1346  }
1347 
1348 
1349  while (true)
1350  {
1351  const labelList& mp = meshPoints();
1352 
1353  PackedBoolList removeFace(this->size());
1354  label nFaces = 0;
1355  {
1356  const labelListList& edgeFaces =
1358  forAll(edgeFaces, edgei)
1359  {
1360  const labelList& eFaces = edgeFaces[edgei];
1361  if (eFaces.size() == 1)
1362  {
1363  // Open edge. Check that vertices do not originate
1364  // from a boundary face
1365  const edge& e = edges()[edgei];
1366  const edge& verts0 = pointToVerts_[mp[e[0]]];
1367  const edge& verts1 = pointToVerts_[mp[e[1]]];
1368  if
1369  (
1370  isBoundaryPoint[verts0[0]]
1371  && isBoundaryPoint[verts0[1]]
1372  && isBoundaryPoint[verts1[0]]
1373  && isBoundaryPoint[verts1[1]]
1374  )
1375  {
1376  // Open edge on boundary face. Keep
1377  }
1378  else
1379  {
1380  // Open edge. Mark for erosion
1381  if (removeFace.set(eFaces[0]))
1382  {
1383  nFaces++;
1384  }
1385  }
1386  }
1387  }
1388  }
1389 
1390  if (debug)
1391  {
1392  Pout<< "isoSurface :"
1393  << " removing " << nFaces
1394  << " faces since on open edges" << endl;
1395  }
1396 
1397  if (returnReduce(nFaces, sumOp<label>()) == 0)
1398  {
1399  break;
1400  }
1401 
1402  // Remove the faces
1403  labelHashSet keepFaces(2*size());
1404  forAll(removeFace, facei)
1405  {
1406  if (!removeFace[facei])
1407  {
1408  keepFaces.insert(facei);
1409  }
1410  }
1411 
1412  labelList pointMap;
1414  MeshedSurface<face> filteredSurf
1415  (
1417  (
1418  keepFaces,
1419  pointMap,
1420  faceMap
1421  )
1422  );
1423  MeshedSurface<face>::transfer(filteredSurf);
1424 
1425  pointToVerts_ = UIndirectList<edge>(pointToVerts_, pointMap)();
1426  pointToFace_ = UIndirectList<label>(pointToFace_, pointMap)();
1427  pointFromDiag = UIndirectList<bool>(pointFromDiag, pointMap)();
1428  meshCells_ = UIndirectList<label>(meshCells_, faceMap)();
1429  }
1430  }
1431  }
1432 }
1433 
1434 
1435 // ************************************************************************* //
A tetrahedron primitive.
Definition: tetrahedron.H:61
static const scalar GREAT
Definition: scalar.H:111
A cellMatcher for tet cells.
Definition: tetMatcher.H:51
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:54
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
label nInternalFaces() const
label rcIndex(const label i) const
Return the reverse circular index, i.e. the previous index.
Definition: UListI.H:65
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
label nFaces() const
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
A surface zone on a MeshedSurface.
Definition: surfZone.H:62
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: UListI.H:58
label nCells() const
static scalar minQuality(const polyMesh &mesh, const point &cC, const label fI, const bool isOwner, const label faceBasePtI)
Given a face and cc and starting index for triangulation determine.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const dimensionedScalar & mp
Proton mass.
void transfer(MeshedSurface< Face > &)
Transfer the contents of the argument and annul the argument.
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:111
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:51
scalar f1
Definition: createFields.H:28
const labelListList & edgeLoops() const
Return list of closed loops of boundary vertices.
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
A list of faces which address into the list of points.
A List obtained as a section of another List.
Definition: SubList.H:53
void Swap(T &a, T &b)
Definition: Swap.H:43
const cellShapeList & cells
void set(const PackedList< 1 > &)
Set specified bits.
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
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:68
face triFace(3)
static const NamedEnum< filterType, 3 > filterTypeNames_
Definition: isoSurface.H:73
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
const Field< PointType > & points() const
Return reference to global points.
const labelListList & edgeFaces() const
Return edge-face addressing.
Dynamically sized Field.
Definition: DynamicField.H:49
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1156
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
virtual void reset(pointField &&points, List< Face > &&faces, surfZoneList &&zones)
Reset primitive data (points, faces and zones)
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
List< surfZone > surfZoneList
Definition: surfZoneList.H:45
virtual bool isA(const primitiveMesh &mesh, const label celli)
Exact match. Uses faceSizeMatch.
Definition: tetMatcher.C:218
void setSize(const label)
Reset size of List.
Definition: List.C:281
A bit-packed bool list.
label size() const
The surface size is the number of faces.
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
A List with indirect addressing.
Definition: fvMatrix.H:106
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
label nPoints() const
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A class for managing temporary objects.
Definition: PtrList.H:53
T & last()
Return the last element of the list.
Definition: UListI.H:128
void clear()
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:236
isoSurface(const polyMesh &mesh, const scalarField &cellValues, const scalarField &pointValues, const scalar iso, const filterType filter)
Construct from dictionary.
Definition: isoSurface.C:1158
Namespace for OpenFOAM.
const List< surfZone > & surfZones() const
Const access to the surface zones.