isoSurfaceTopo.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2019 OpenFOAM Foundation
9  Copyright (C) 2019-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "isoSurfaceTopo.H"
30 #include "polyMesh.H"
31 #include "volFields.H"
32 #include "edgeHashes.H"
33 #include "tetCell.H"
34 #include "DynamicField.H"
35 #include "syncTools.H"
36 #include "indirectPrimitivePatch.H"
39 #include "foamVtkLineWriter.H"
40 #include "foamVtkSurfaceWriter.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 #include "isoSurfaceBaseMethods.H"
46 
47 
48 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
49 
50 namespace Foam
51 {
52  defineTypeNameAndDebug(isoSurfaceTopo, 0);
53 }
54 
55 
56 // Get/set snapIndex (0, 1 or 2) at given position
57 // 0 = no snap
58 // 1 = snap to first edge end
59 // 2 = snap to second edge end
60 // NB: 4 lower bits left free for regular tet-cut information
61 
62 #undef SNAP_END_VALUE
63 #undef SNAP_END_ENCODE
64 #define SNAP_END_ENCODE(pos, val) (((val) << (4 + 2 * pos)))
65 #define SNAP_END_VALUE(pos, val) (((val) >> (4 + 2 * pos)) & 0x3)
66 
67 
68 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
69 
70 namespace Foam
71 {
72 
73 // Check for tet values above/below given (iso) value
74 // Result encoded as an integer, with possible snapping information too
75 inline static int getTetCutIndex
76 (
77  scalar p0,
78  scalar p1,
79  scalar p2,
80  scalar p3,
81  const scalar val,
82  const bool doSnap
83 ) noexcept
84 {
85  int cutIndex
86  (
87  (p0 < val ? 1 : 0) // point 0
88  | (p1 < val ? 2 : 0) // point 1
89  | (p2 < val ? 4 : 0) // point 2
90  | (p3 < val ? 8 : 0) // point 3
91  );
92 
93  if (doSnap && cutIndex && cutIndex != 0xF)
94  {
95  // Not all below or all
96 
97  // Calculate distances (for snapping)
98  p0 -= val; if (cutIndex & 1) p0 *= -1;
99  p1 -= val; if (cutIndex & 2) p1 *= -1;
100  p2 -= val; if (cutIndex & 4) p2 *= -1;
101  p3 -= val; if (cutIndex & 8) p3 *= -1;
102 
103  // Add snap index into regular edge cut index
104  // Snap to end if less than approx 1% of the distance.
105  // - only valid if there is also a corresponding sign change
106  #undef ADD_SNAP_INDEX
107  #define ADD_SNAP_INDEX(pos, d1, d2, idx1, idx2) \
108  switch (cutIndex & (idx1 | idx2)) \
109  { \
110  case idx1 : /* first below, second above */ \
111  case idx2 : /* first above, second below */ \
112  cutIndex |= SNAP_END_ENCODE \
113  ( \
114  pos, \
115  ((d1 * 100 < d2) ? 1 : (d2 * 100 < d1) ? 2 : 0) \
116  ); \
117  break; \
118  }
119 
120  ADD_SNAP_INDEX(0, p0, p1, 1, 2); // Edge 0: 0 -> 1
121  ADD_SNAP_INDEX(1, p0, p2, 1, 4); // Edge 1: 0 -> 2
122  ADD_SNAP_INDEX(2, p0, p3, 1, 8); // Edge 2: 0 -> 3
123  ADD_SNAP_INDEX(3, p3, p1, 8, 2); // Edge 3: 3 -> 1
124  ADD_SNAP_INDEX(4, p1, p2, 2, 4); // Edge 4: 1 -> 2
125  ADD_SNAP_INDEX(5, p3, p2, 8, 4); // Edge 5: 3 -> 2
126  #undef ADD_SNAP_INDEX
127  }
128 
129  return cutIndex;
130 }
131 
132 
133 // Append three labels to list.
134 // Filter out degenerate (eg, snapped) tris. Flip face as requested
135 inline static void appendTriLabels
136 (
137  DynamicList<label>& verts,
138  const label a,
139  const label b,
140  const label c,
141  const bool flip // Flip normals
142 )
143 {
144  if (a != b && b != c && c != a)
145  {
146  verts.append(a);
147  if (flip)
148  {
149  verts.append(c);
150  verts.append(b);
151  }
152  else
153  {
154  verts.append(b);
155  verts.append(c);
156  }
157  }
158 }
159 
160 
161 // Return point reference to mesh points or cell-centres
162 inline static const point& getMeshPointRef
163 (
164  const polyMesh& mesh,
165  const label pointi
166 )
167 {
168  return
169  (
170  pointi < mesh.nPoints()
171  ? mesh.points()[pointi]
172  : mesh.cellCentres()[pointi - mesh.nPoints()]
173  );
174 }
175 
176 } // End namespace Foam
177 
178 
179 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
180 
181 Foam::isoSurfaceTopo::tetCutAddressing::tetCutAddressing
182 (
183  const label nCutCells,
184  const bool useSnap,
185  const bool useDebugCuts
186 )
187 :
188  vertsToPointLookup_(12*nCutCells),
189  snapVertsLookup_(0),
190 
191  pointToFace_(10*nCutCells),
192  pointFromDiag_(10*nCutCells),
193 
194  pointToVerts_(10*nCutCells),
195  cutPoints_(12*nCutCells),
196 
197  debugCutTets_(),
198  debugCutTetsOn_(useDebugCuts)
199 {
200  // Per cell: 5 pyramids cut, each generating 2 triangles
201 
202  // Per cell: number of intersected edges:
203  // - four faces cut so 4 mesh edges + 4 face-diagonal edges
204  // - 4 of the pyramid edges
205 
206  if (useSnap)
207  {
208  // Some, but not all, cells may have point snapping
209  snapVertsLookup_.resize(4*nCutCells);
210  }
211  if (debugCutTetsOn_)
212  {
213  debugCutTets_.reserve(6*nCutCells);
214  }
215 }
216 
217 
218 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
219 
220 void Foam::isoSurfaceTopo::tetCutAddressing::clearDebug()
221 {
222  debugCutTets_.clearStorage();
223 }
224 
225 
226 void Foam::isoSurfaceTopo::tetCutAddressing::clearDiagonal()
227 {
228  pointToFace_.clearStorage();
229  pointFromDiag_.clearStorage();
230 }
231 
232 
233 void Foam::isoSurfaceTopo::tetCutAddressing::clearHashes()
234 {
235  vertsToPointLookup_.clear();
236  snapVertsLookup_.clear();
237 }
238 
239 
240 Foam::label Foam::isoSurfaceTopo::tetCutAddressing::generatePoint
241 (
242  label facei,
243  bool edgeIsDiagonal,
244  const int snapEnd,
245  const edge& vertices
246 )
247 {
248  // Generate new point, unless it already exists for edge
249  // or corresponds to a snapped point (from another edge)
250 
251  label pointi = vertsToPointLookup_.lookup(vertices, -1);
252  if (pointi == -1)
253  {
254  bool addNewPoint(true);
255 
256  const label snapPointi =
257  (
258  (snapEnd == 1) ? vertices.first()
259  : (snapEnd == 2) ? vertices.second()
260  : -1
261  );
262 
263  if (snapPointi == -1)
264  {
265  // No snapped point
266  pointi = pointToVerts_.size();
267  pointToVerts_.append(vertices);
268  }
269  else
270  {
271  // Snapped point. No corresponding face or diagonal
272  facei = -1;
273  edgeIsDiagonal = false;
274 
275  pointi = snapVertsLookup_.lookup(snapPointi, -1);
276  addNewPoint = (pointi == -1);
277  if (addNewPoint)
278  {
279  pointi = pointToVerts_.size();
280  snapVertsLookup_.insert(snapPointi, pointi);
281  pointToVerts_.append(edge(snapPointi, snapPointi));
282  }
283  }
284 
285  if (addNewPoint)
286  {
287  pointToFace_.append(facei);
288  pointFromDiag_.append(edgeIsDiagonal);
289  }
290 
291  vertsToPointLookup_.insert(vertices, pointi);
292  }
293 
294  return pointi;
295 }
296 
297 
298 bool Foam::isoSurfaceTopo::tetCutAddressing::generatePoints
299 (
300  const label facei,
301  const int tetCutIndex,
302  const tetCell& tetLabels,
303 
304  // Per tet edge whether is face diag etc
305  const FixedList<bool, 6>& edgeIsDiagonal
306 )
307 {
308  bool flip(false);
309  const label nCutPointsOld(cutPoints_.size());
310 
311  // Form the vertices of the triangles for each case
312  switch (tetCutIndex & 0x0F)
313  {
314  case 0x00:
315  case 0x0F:
316  break;
317 
318  // Cut point 0
319  case 0x0E: flip = true; [[fallthrough]]; // Point 0 above cut
320  case 0x01: // Point 0 below cut
321  {
322  const label cutA
323  (
324  generatePoint
325  (
326  facei,
327  edgeIsDiagonal[0],
328  SNAP_END_VALUE(0, tetCutIndex),
329  tetLabels.edge(0) // 0 -> 1
330  )
331  );
332  const label cutB
333  (
334  generatePoint
335  (
336  facei,
337  edgeIsDiagonal[1],
338  SNAP_END_VALUE(1, tetCutIndex),
339  tetLabels.edge(1) // 0 -> 2
340  )
341  );
342  const label cutC
343  (
344  generatePoint
345  (
346  facei,
347  edgeIsDiagonal[2],
348  SNAP_END_VALUE(2, tetCutIndex),
349  tetLabels.edge(2) // 0 -> 3
350  )
351  );
352 
353  appendTriLabels(cutPoints_, cutA, cutB, cutC, flip);
354  }
355  break;
356 
357  // Cut point 1
358  case 0x0D: flip = true; [[fallthrough]]; // Point 1 above cut
359  case 0x02: // Point 1 below cut
360  {
361  const label cutA
362  (
363  generatePoint
364  (
365  facei,
366  edgeIsDiagonal[0],
367  SNAP_END_VALUE(0, tetCutIndex),
368  tetLabels.edge(0) // 0 -> 1
369  )
370  );
371  const label cutB
372  (
373  generatePoint
374  (
375  facei,
376  edgeIsDiagonal[3],
377  SNAP_END_VALUE(3, tetCutIndex),
378  tetLabels.edge(3) // 3 -> 1
379  )
380  );
381  const label cutC
382  (
383  generatePoint
384  (
385  facei,
386  edgeIsDiagonal[4],
387  SNAP_END_VALUE(4, tetCutIndex),
388  tetLabels.edge(4) // 1 -> 2
389  )
390  );
391 
392  appendTriLabels(cutPoints_, cutA, cutB, cutC, flip);
393  }
394  break;
395 
396  // Cut point 0/1 | 2/3
397  case 0x0C: flip = true; [[fallthrough]]; // Point 0/1 above cut
398  case 0x03: // Point 0/1 below cut
399  {
400  const label cutA
401  (
402  generatePoint
403  (
404  facei,
405  edgeIsDiagonal[1],
406  SNAP_END_VALUE(1, tetCutIndex),
407  tetLabels.edge(1) // 0 -> 2
408  )
409  );
410  const label cutB
411  (
412  generatePoint
413  (
414  facei,
415  edgeIsDiagonal[2],
416  SNAP_END_VALUE(2, tetCutIndex),
417  tetLabels.edge(2) // 0 -> 3
418  )
419  );
420  const label cutC
421  (
422  generatePoint
423  (
424  facei,
425  edgeIsDiagonal[3],
426  SNAP_END_VALUE(3, tetCutIndex),
427  tetLabels.edge(3) // 3 -> 1
428  )
429  );
430  const label cutD
431  (
432  generatePoint
433  (
434  facei,
435  edgeIsDiagonal[4],
436  SNAP_END_VALUE(4, tetCutIndex),
437  tetLabels.edge(4) // 1 -> 2
438  )
439  );
440 
441  appendTriLabels(cutPoints_, cutA, cutB, cutC, flip);
442  appendTriLabels(cutPoints_, cutA, cutC, cutD, flip);
443  }
444  break;
445 
446  // Cut point 2
447  case 0x0B: flip = true; [[fallthrough]]; // Point 2 above cut
448  case 0x04: // Point 2 below cut
449  {
450  const label cutA
451  (
452  generatePoint
453  (
454  facei,
455  edgeIsDiagonal[1],
456  SNAP_END_VALUE(1, tetCutIndex),
457  tetLabels.edge(1) // 0 -> 2
458  )
459  );
460  const label cutB
461  (
462  generatePoint
463  (
464  facei,
465  edgeIsDiagonal[4],
466  SNAP_END_VALUE(4, tetCutIndex),
467  tetLabels.edge(4) // 1 -> 2
468  )
469  );
470  const label cutC
471  (
472  generatePoint
473  (
474  facei,
475  edgeIsDiagonal[5],
476  SNAP_END_VALUE(5, tetCutIndex),
477  tetLabels.edge(5) // 3 -> 2
478  )
479  );
480 
481  appendTriLabels(cutPoints_, cutA, cutB, cutC, flip);
482  }
483  break;
484 
485  // Cut point 0/2 | 1/3
486  case 0x0A: flip = true; [[fallthrough]]; // Point 0/2 above cut
487  case 0x05: // Point 0/2 below cut
488  {
489  const label cutA
490  (
491  generatePoint
492  (
493  facei,
494  edgeIsDiagonal[0],
495  SNAP_END_VALUE(0, tetCutIndex),
496  tetLabels.edge(0) // 0 -> 1
497  )
498  );
499  const label cutB
500  (
501  generatePoint
502  (
503  facei,
504  edgeIsDiagonal[4],
505  SNAP_END_VALUE(4, tetCutIndex),
506  tetLabels.edge(4) // 1 -> 2
507  )
508  );
509  const label cutC
510  (
511  generatePoint
512  (
513  facei,
514  edgeIsDiagonal[5],
515  SNAP_END_VALUE(5, tetCutIndex),
516  tetLabels.edge(5) // 3 -> 2
517  )
518  );
519  const label cutD
520  (
521  generatePoint
522  (
523  facei,
524  edgeIsDiagonal[2],
525  SNAP_END_VALUE(2, tetCutIndex),
526  tetLabels.edge(2) // 0 -> 3
527  )
528  );
529 
530  appendTriLabels(cutPoints_, cutA, cutB, cutC, flip);
531  appendTriLabels(cutPoints_, cutA, cutC, cutD, flip);
532  }
533  break;
534 
535  // Cut point 1/2 | 0/3
536  case 0x09: flip = true; [[fallthrough]]; // Point 1/2 above cut
537  case 0x06: // Point 1/2 below cut
538  {
539  const label cutA
540  (
541  generatePoint
542  (
543  facei,
544  edgeIsDiagonal[0],
545  SNAP_END_VALUE(0, tetCutIndex),
546  tetLabels.edge(0) // 0 -> 1
547  )
548  );
549  const label cutB
550  (
551  generatePoint
552  (
553  facei,
554  edgeIsDiagonal[3],
555  SNAP_END_VALUE(3, tetCutIndex),
556  tetLabels.edge(3) // 3 -> 1
557  )
558  );
559  const label cutC
560  (
561  generatePoint
562  (
563  facei,
564  edgeIsDiagonal[5],
565  SNAP_END_VALUE(5, tetCutIndex),
566  tetLabels.edge(5) // 3 -> 2
567  )
568  );
569  const label cutD
570  (
571  generatePoint
572  (
573  facei,
574  edgeIsDiagonal[1],
575  SNAP_END_VALUE(1, tetCutIndex),
576  tetLabels.edge(1) // 0 -> 2
577  )
578  );
579 
580  appendTriLabels(cutPoints_, cutA, cutB, cutC, flip);
581  appendTriLabels(cutPoints_, cutA, cutC, cutD, flip);
582  }
583  break;
584 
585  // Cut point 3
586  case 0x07: flip = true; [[fallthrough]]; // Point 3 above cut
587  case 0x08: // Point 3 below cut
588  {
589  const label cutA
590  (
591  generatePoint
592  (
593  facei,
594  edgeIsDiagonal[2],
595  SNAP_END_VALUE(2, tetCutIndex),
596  tetLabels.edge(2) // 0 -> 3
597  )
598  );
599  const label cutB
600  (
601  generatePoint
602  (
603  facei,
604  edgeIsDiagonal[5],
605  SNAP_END_VALUE(5, tetCutIndex),
606  tetLabels.edge(5) // 3 -> 2
607  )
608  );
609  const label cutC
610  (
611  generatePoint
612  (
613  facei,
614  edgeIsDiagonal[3],
615  SNAP_END_VALUE(3, tetCutIndex),
616  tetLabels.edge(3) // 3 -> 1
617  )
618  );
619 
620  appendTriLabels(cutPoints_, cutA, cutB, cutC, flip);
621  }
622  break;
623  }
624 
625  const bool added(nCutPointsOld != cutPoints_.size());
626 
627  if (added && debugCutTetsOn_)
628  {
629  debugCutTets_.append(tetLabels.shape());
630  }
631 
632  return added;
633 }
634 
635 
636 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
637 
638 // Requires mesh_, tetBasePtIs
639 void Foam::isoSurfaceTopo::generateTriPoints
640 (
641  const label celli,
642  const bool isTet,
643  const labelList& tetBasePtIs,
644  tetCutAddressing& tetCutAddr
645 ) const
646 {
647  const faceList& faces = mesh_.faces();
648  const labelList& faceOwner = mesh_.faceOwner();
649  const cell& cFaces = mesh_.cells()[celli];
650  const bool doSnap = this->snap();
651 
652  if (isTet)
653  {
654  // For tets don't do cell-centre decomposition, just use the
655  // tet points and values
656 
657  const label facei = cFaces[0];
658  const face& f0 = faces[facei];
659 
660  // Get the other point from f1. Tbd: check if not duplicate face
661  // (ACMI / ignoreBoundaryFaces_).
662  const face& f1 = faces[cFaces[1]];
663  label apexi = -1;
664  forAll(f1, fp)
665  {
666  apexi = f1[fp];
667  if (!f0.found(apexi))
668  {
669  break;
670  }
671  }
672 
673  const label p0 = f0[0];
674  label p1 = f0[1];
675  label p2 = f0[2];
676 
677  if (faceOwner[facei] == celli)
678  {
679  std::swap(p1, p2);
680  }
681 
682  const tetCell tetLabels(p0, p1, p2, apexi);
683  const int tetCutIndex
684  (
686  (
687  pVals_[p0],
688  pVals_[p1],
689  pVals_[p2],
690  pVals_[apexi],
691  iso_,
692  doSnap
693  )
694  );
695 
696  tetCutAddr.generatePoints
697  (
698  facei,
699  tetCutIndex,
700  tetLabels,
701  FixedList<bool, 6>(false) // Not face diagonal
702  );
703  }
704  else
705  {
706  for (const label facei : cFaces)
707  {
708  if
709  (
710  !mesh_.isInternalFace(facei)
712  )
713  {
714  continue;
715  }
716 
717  const face& f = faces[facei];
718 
719  label fp0 = tetBasePtIs[facei];
720 
721  // Fallback
722  if (fp0 < 0)
723  {
724  fp0 = 0;
725  }
726 
727  const label p0 = f[fp0];
728  label fp = f.fcIndex(fp0);
729  for (label i = 2; i < f.size(); ++i)
730  {
731  label p1 = f[fp];
732  fp = f.fcIndex(fp);
733  label p2 = f[fp];
734 
735  FixedList<bool, 6> edgeIsDiagonal(false);
736  if (faceOwner[facei] == celli)
737  {
738  std::swap(p1, p2);
739  if (i != 2) edgeIsDiagonal[1] = true;
740  if (i != f.size()-1) edgeIsDiagonal[0] = true;
741  }
742  else
743  {
744  if (i != 2) edgeIsDiagonal[0] = true;
745  if (i != f.size()-1) edgeIsDiagonal[1] = true;
746  }
747 
748  const tetCell tetLabels(p0, p1, p2, mesh_.nPoints()+celli);
749  const int tetCutIndex
750  (
752  (
753  pVals_[p0],
754  pVals_[p1],
755  pVals_[p2],
756  cVals_[celli],
757  iso_,
758  doSnap
759  )
760  );
761 
762  tetCutAddr.generatePoints
763  (
764  facei,
765  tetCutIndex,
766  tetLabels,
767  edgeIsDiagonal
768  );
769  }
770  }
771  }
772 }
773 
774 
775 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
776 
777 void Foam::isoSurfaceTopo::triangulateOutside
778 (
779  const bool filterDiag,
780  const primitivePatch& pp,
781  const boolUList& pointFromDiag,
782  const labelUList& pointToFace,
783  const label cellID,
784 
785  // outputs
786  DynamicList<face>& compactFaces,
787  DynamicList<label>& compactCellIDs
788 )
789 {
790  // We can form pockets:
791  // - 1. triangle on face
792  // - 2. multiple triangles on interior (from diag edges)
793  // - the edge loop will be pocket since it is only the diag
794  // edges that give it volume?
795 
796  // Retriangulate the exterior loops
797  const labelListList& edgeLoops = pp.edgeLoops();
798  const labelList& mp = pp.meshPoints();
799 
800  for (const labelList& loop : edgeLoops)
801  {
802  if (loop.size() > 2)
803  {
804  compactFaces.append(face(loop.size()));
805  face& f = compactFaces.last();
806 
807  label fpi = 0;
808  forAll(f, i)
809  {
810  const label pointi = mp[loop[i]];
811  if (filterDiag && pointFromDiag[pointi])
812  {
813  const label prevPointi = mp[loop[loop.fcIndex(i)]];
814  if
815  (
816  pointFromDiag[prevPointi]
817  && (pointToFace[pointi] != pointToFace[prevPointi])
818  )
819  {
820  f[fpi++] = pointi;
821  }
822  else
823  {
824  // Filter out diagonal point
825  }
826  }
827  else
828  {
829  f[fpi++] = pointi;
830  }
831  }
832 
833  if (fpi > 2)
834  {
835  f.resize(fpi);
836  }
837  else
838  {
839  // Keep original face
840  forAll(f, i)
841  {
842  const label pointi = mp[loop[i]];
843  f[i] = pointi;
844  }
845  }
846  compactCellIDs.append(cellID);
847  }
848  }
849 }
850 
851 
852 void Foam::isoSurfaceTopo::removeInsidePoints
853 (
854  Mesh& s,
855  const bool filterDiag,
856 
857  // inputs
858  const boolUList& pointFromDiag,
859  const labelUList& pointToFace,
860  const labelUList& start, // Per cell the starting triangle
861 
862  // outputs
863  DynamicList<label>& compactCellIDs // Per returned tri the cellID
864 )
865 {
866  const pointField& points = s.points();
867 
868  compactCellIDs.clear();
869  compactCellIDs.reserve(s.size()/4);
870 
871  DynamicList<face> compactFaces(s.size()/4);
872 
873  for (label celli = 0; celli < start.size()-1; ++celli)
874  {
875  // Triangles for the current cell
876 
877  const label nTris = start[celli+1]-start[celli];
878 
879  if (nTris)
880  {
881  const primitivePatch pp
882  (
883  SubList<face>(s, nTris, start[celli]),
884  points
885  );
886 
887  triangulateOutside
888  (
889  filterDiag,
890  pp,
891  pointFromDiag,
892  pointToFace,
893  celli,
894 
895  compactFaces,
896  compactCellIDs
897  );
898  }
899  }
900 
901  s.swapFaces(compactFaces); // Use new faces
902 }
903 
904 
905 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
906 
908 (
909  const polyMesh& mesh,
910  const scalarField& cellValues,
911  const scalarField& pointValues,
912  const scalar iso,
913  const isoSurfaceParams& params,
914  const bitSet& ignoreCells
915 )
916 :
917  isoSurfaceBase(mesh, cellValues, pointValues, iso, params)
918 {
919  // The cell cut type
920  List<cutType> cellCutType_(mesh.nCells(), cutType::UNVISITED);
921 
922  // Time description (for debug output)
923  const word timeDesc(word::printf("%08d", mesh_.time().timeIndex()));
924 
925  if (debug)
926  {
927  Pout<< "isoSurfaceTopo:" << nl
928  << " cell min/max : " << minMax(cVals_) << nl
929  << " point min/max : " << minMax(pVals_) << nl
930  << " isoValue : " << iso << nl
931  << " filter : "
932  << isoSurfaceParams::filterNames[params.filter()] << nl
933  << " mesh span : " << mesh.bounds().mag() << nl
934  << " ignoreCells : " << ignoreCells.count()
935  << " / " << cVals_.size() << nl
936  << endl;
937  }
938 
939  this->ignoreCyclics();
940 
941  label nBlockedCells = 0;
942 
943  // Mark ignoreCells as blocked
944  nBlockedCells += blockCells(cellCutType_, ignoreCells);
945 
946  // Mark cells outside bounding box as blocked
947  nBlockedCells +=
948  blockCells(cellCutType_, params.getClipBounds(), volumeType::OUTSIDE);
949 
950  // Adjusted tet base points to improve tet quality
951  labelList tetBasePtIs
952  (
954  );
955 
956 
957  // Determine cell cuts
958  const label nCutCells = calcCellCuts(cellCutType_);
959 
960  if (debug)
961  {
962  Pout<< "isoSurfaceTopo : candidate cells cut "
963  << nCutCells
964  << " blocked " << nBlockedCells
965  << " total " << mesh_.nCells() << endl;
966  }
967 
968  if (debug && isA<fvMesh>(mesh))
969  {
970  const auto& fvmesh = refCast<const fvMesh>(mesh);
971 
972  volScalarField debugField
973  (
974  IOobject
975  (
976  "isoSurfaceTopo.cutType",
977  fvmesh.time().timeName(),
978  fvmesh.time(),
982  ),
983  fvmesh,
985  );
986 
987  auto& debugFld = debugField.primitiveFieldRef();
988 
989  forAll(cellCutType_, celli)
990  {
991  debugFld[celli] = cellCutType_[celli];
992  }
993 
994  Info<< "Writing cut types: " << debugField.objectRelPath() << nl;
995  debugField.write();
996  }
997 
998  // Additional debugging
999  if (debug & 8)
1000  {
1001  // Write debug cuts cells in VTK format
1002  {
1003  constexpr uint8_t realCut(cutType::CUT | cutType::TETCUT);
1004  labelList debugCutCells(nCutCells, Zero);
1005 
1006  label nout = 0;
1007  forAll(cellCutType_, celli)
1008  {
1009  if ((cellCutType_[celli] & realCut) != 0)
1010  {
1011  debugCutCells[nout] = celli;
1012  ++nout;
1013  if (nout >= nCutCells) break;
1014  }
1015  }
1016 
1017  // The mesh subset cut
1018  vtk::vtuCells vtuCells;
1019  vtuCells.reset(mesh_, debugCutCells);
1020 
1021  vtk::internalMeshWriter writer
1022  (
1023  mesh_,
1024  vtuCells,
1025  fileName
1026  (
1027  mesh_.time().globalPath()
1028  / ("isoSurfaceTopo." + timeDesc + "-cutCells")
1029  )
1030  );
1031 
1032  writer.writeGeometry();
1033 
1034  // CellData
1035  writer.beginCellData();
1036  writer.writeProcIDs();
1037  writer.writeCellData("cutField", cVals_);
1038 
1039  // PointData
1040  writer.beginPointData();
1041  writer.writePointData("cutField", pVals_);
1042 
1043  Info<< "isoSurfaceTopo : (debug) wrote "
1044  << returnReduce(nCutCells, sumOp<label>())
1045  << " cut cells: "
1046  << writer.output().name() << nl;
1047  }
1048  }
1049 
1050 
1051  tetCutAddressing tetCutAddr
1052  (
1053  nCutCells,
1054  this->snap(),
1055  (debug & 8) // Enable debug tets
1056  );
1057 
1058  labelList startTri(mesh_.nCells()+1, Zero);
1059  for (label celli = 0; celli < mesh_.nCells(); ++celli)
1060  {
1061  startTri[celli] = tetCutAddr.nFaces();
1062  if ((cellCutType_[celli] & cutType::ANYCUT) != 0)
1063  {
1064  generateTriPoints
1065  (
1066  celli,
1067  // Same as tetMatcher::test(mesh_, celli),
1068  bool(cellCutType_[celli] & cutType::TETCUT),
1069 
1070  tetBasePtIs,
1071  tetCutAddr
1072  );
1073  }
1074  }
1075  startTri.last() = tetCutAddr.nFaces();
1076 
1077  // Information not needed anymore:
1078  tetBasePtIs.clear();
1079  tetCutAddr.clearHashes();
1080 
1081 
1082  // From list of vertices -> triangular faces
1083  faceList allTriFaces(startTri.last());
1084  {
1085  auto& verts = tetCutAddr.cutPoints();
1086 
1087  label verti = 0;
1088  for (face& f : allTriFaces)
1089  {
1090  f.resize(3);
1091  f[0] = verts[verti++];
1092  f[1] = verts[verti++];
1093  f[2] = verts[verti++];
1094  }
1095  verts.clearStorage(); // Not needed anymore
1096  }
1097 
1098 
1099  // The cells cut by the triangular faces
1100  meshCells_.resize(startTri.last());
1101  for (label celli = 0; celli < startTri.size()-1; ++celli)
1102  {
1103  // All triangles for the current cell
1105  (
1106  meshCells_,
1107  (startTri[celli+1] - startTri[celli]),
1108  startTri[celli]
1109  ) = celli;
1110  }
1111 
1112 
1113  pointToVerts_.transfer(tetCutAddr.pointToVerts());
1114 
1115  pointField allTriPoints
1116  (
1117  this->interpolateTemplate
1118  (
1119  mesh_.cellCentres(),
1120  mesh_.points()
1121  )
1122  );
1123 
1124 
1125  // Assign to MeshedSurface
1126  static_cast<Mesh&>(*this) = Mesh
1127  (
1128  std::move(allTriPoints),
1129  std::move(allTriFaces),
1130  surfZoneList() // zones not required (one zone)
1131  );
1132 
1133  if (debug)
1134  {
1135  Pout<< "isoSurfaceTopo : generated "
1136  << Mesh::size() << " triangles "
1137  << Mesh::points().size() << " points" << endl;
1138  }
1139 
1140  // Write debug triangulated surface
1141  if ((debug & 8) && (params.filter() != filterType::NONE))
1142  {
1143  const Mesh& s = *this;
1144 
1145  vtk::surfaceWriter writer
1146  (
1147  s.points(),
1148  s,
1149  fileName
1150  (
1151  mesh_.time().globalPath()
1152  / ("isoSurfaceTopo." + timeDesc + "-triangles")
1153  )
1154  );
1155 
1156  writer.writeGeometry();
1157 
1158  // CellData
1159  writer.beginCellData();
1160  writer.writeProcIDs();
1161  writer.write("cellID", meshCells_);
1162 
1163  // PointData
1164  writer.beginPointData();
1165  {
1166  // NB: may have non-compact surface points
1167  // --> use points().size() not nPoints()!
1168 
1169  labelList pointStatus(s.points().size(), Zero);
1170 
1171  forAll(pointToVerts_, i)
1172  {
1173  const edge& verts = pointToVerts_[i];
1174  if (verts.first() == verts.second())
1175  {
1176  // Duplicate index (ie, snapped)
1177  pointStatus[i] = 1;
1178  }
1179  if (tetCutAddr.pointFromDiag().test(i))
1180  {
1181  // Point on triangulation diagonal
1182  pointStatus[i] = -1;
1183  }
1184  }
1185 
1186  writer.write("point-status", pointStatus);
1187  }
1188 
1189  Info<< "isoSurfaceTopo : (debug) wrote "
1190  << returnReduce(s.size(), sumOp<label>())
1191  << " triangles : "
1192  << writer.output().name() << nl;
1193  }
1194 
1195 
1196  // Now:
1197  // - generated faces and points are assigned to *this
1198  // - per point we know:
1199  // - pointOnDiag: whether it is on a face-diagonal edge
1200  // - pointToFace: from what pyramid (cell+face) it was produced
1201  // (note that the pyramid faces are shared between multiple mesh faces)
1202  // - pointToVerts_ : originating mesh vertex or cell centre
1203 
1204  if (params.filter() == filterType::NONE)
1205  {
1206  // Compact out unused (snapped) points
1207  if (this->snap())
1208  {
1209  Mesh& s = *this;
1210 
1211  labelList pointMap; // Back to original point
1212  s.compactPoints(pointMap); // Compact out unused points
1213  pointToVerts_ = UIndirectList<edge>(pointToVerts_, pointMap)();
1214  }
1215  }
1216  else
1217  {
1218  // Initial filtering
1219 
1220  Mesh& s = *this;
1221 
1222  // Triangulate outside
1223  // (filter edges to cell centres and optionally face diagonals)
1224  DynamicList<label> compactCellIDs; // Per tri the cell
1225 
1226  removeInsidePoints
1227  (
1228  *this,
1229  // Filter face diagonal
1230  (
1231  params.filter() == filterType::DIAGCELL
1232  || params.filter() == filterType::NONMANIFOLD
1233  ),
1234  tetCutAddr.pointFromDiag(),
1235  tetCutAddr.pointToFace(),
1236  startTri,
1237  compactCellIDs
1238  );
1239 
1240  labelList pointMap; // Back to original point
1241  s.compactPoints(pointMap); // Compact out unused points
1242 
1243  pointToVerts_ = UIndirectList<edge>(pointToVerts_, pointMap)();
1244  meshCells_.transfer(compactCellIDs);
1245 
1246  if (debug)
1247  {
1248  Pout<< "isoSurfaceTopo :"
1249  " after removing cell centre and face-diag triangles: "
1250  << Mesh::size() << " faces "
1251  << Mesh::points().size() << " points"
1252  << endl;
1253  }
1254  }
1255 
1256  // Diagonal filter information not needed anymore
1257  tetCutAddr.clearDiagonal();
1258 
1259 
1260  // For more advanced filtering (eg, removal of open edges)
1261  // need the boundary and other 'protected' points
1262 
1263  bitSet isProtectedPoint;
1264  if
1265  (
1266  (params.filter() == filterType::NONMANIFOLD)
1267  || tetCutAddr.debugCutTetsOn()
1268  )
1269  {
1270  // Mark points on mesh outside as 'protected'
1271  // - never erode these edges
1272 
1273  isProtectedPoint.resize(mesh_.nPoints());
1274 
1275  for
1276  (
1277  label facei = mesh_.nInternalFaces();
1278  facei < mesh_.nFaces();
1279  ++facei
1280  )
1281  {
1282  isProtectedPoint.set(mesh_.faces()[facei]);
1283  }
1284 
1285  // Include faces that would be exposed from mesh subset
1286  if (nBlockedCells)
1287  {
1288  const labelList& faceOwn = mesh_.faceOwner();
1289  const labelList& faceNei = mesh_.faceNeighbour();
1290 
1291  for (label facei = 0; facei < mesh.nInternalFaces(); ++facei)
1292  {
1293  // If only one cell is blocked, the face corresponds
1294  // to an exposed subMesh face
1295  if
1296  (
1297  (cellCutType_[faceOwn[facei]] == cutType::BLOCKED)
1298  != (cellCutType_[faceNei[facei]] == cutType::BLOCKED)
1299  )
1300  {
1301  isProtectedPoint.set(mesh_.faces()[facei]);
1302  }
1303  }
1304  }
1305  }
1306 
1307  // Initial cell cut information not needed anymore
1308  cellCutType_.clear();
1309 
1310 
1311  // Additional debugging
1312  if (tetCutAddr.debugCutTetsOn())
1313  {
1314  // Write debug cut tets in VTK format
1315  {
1316  const auto& debugCuts = tetCutAddr.debugCutTets();
1317 
1318  // The TET shapes, using the mesh_ points information
1319  vtk::vtuCells vtuCells;
1320  vtuCells.resetShapes(debugCuts);
1321 
1322  // Use all points and all cell-centres
1323  vtuCells.setNumPoints(mesh_.nPoints());
1324  vtuCells.addPointCellLabels(identity(mesh_.nCells()));
1325 
1326  vtk::internalMeshWriter writer
1327  (
1328  mesh_,
1329  vtuCells,
1330  fileName
1331  (
1332  mesh_.time().globalPath()
1333  / ("isoSurfaceTopo." + timeDesc + "-cutTets")
1334  )
1335  );
1336 
1337  writer.writeGeometry();
1338 
1339  // CellData
1340  writer.beginCellData();
1341  writer.writeProcIDs();
1342 
1343  // Quality of the cut tets
1344  {
1345  Field<scalar> cutTetQuality(debugCuts.size());
1346  forAll(cutTetQuality, teti)
1347  {
1348  cutTetQuality[teti] = tetPointRef
1349  (
1350  getMeshPointRef(mesh_, debugCuts[teti][0]),
1351  getMeshPointRef(mesh_, debugCuts[teti][1]),
1352  getMeshPointRef(mesh_, debugCuts[teti][2]),
1353  getMeshPointRef(mesh_, debugCuts[teti][3])
1354  ).quality();
1355  }
1356  writer.writeCellData("tetQuality", cutTetQuality);
1357  }
1358 
1359  // PointData
1360  if (this->snap())
1361  {
1362  writer.beginPointData();
1363 
1364  labelList pointStatus(vtuCells.nFieldPoints(), Zero);
1365 
1366  for (const edge& verts : pointToVerts_)
1367  {
1368  if (verts.first() == verts.second())
1369  {
1370  // Duplicate index (ie, snapped)
1371  pointStatus[verts.first()] = 1;
1372  }
1373  }
1374 
1375  writer.writePointData("point-status", pointStatus);
1376  }
1377 
1378  Info<< "isoSurfaceTopo : (debug) wrote "
1379  << returnReduce(debugCuts.size(), sumOp<label>())
1380  << " cut tets: "
1381  << writer.output().name() << nl;
1382  }
1383 
1384  // Determining open edges. Same logic as used later...
1385 
1386  labelHashSet openEdgeIds(0);
1387 
1388  {
1389  const Mesh& s = *this;
1390 
1391  const labelList& mp = s.meshPoints();
1392  const edgeList& surfEdges = s.edges();
1393  const labelListList& edgeFaces = s.edgeFaces();
1394  openEdgeIds.resize(2*s.size());
1395 
1396  forAll(edgeFaces, edgei)
1397  {
1398  const labelList& eFaces = edgeFaces[edgei];
1399  if (eFaces.size() == 1)
1400  {
1401  // Open edge (not originating from a boundary face)
1402 
1403  const edge& e = surfEdges[edgei];
1404  const edge& verts0 = pointToVerts_[mp[e.first()]];
1405  const edge& verts1 = pointToVerts_[mp[e.second()]];
1406 
1407  if
1408  (
1409  isProtectedPoint.test(verts0.first())
1410  && isProtectedPoint.test(verts0.second())
1411  && isProtectedPoint.test(verts1.first())
1412  && isProtectedPoint.test(verts1.second())
1413  )
1414  {
1415  // Open edge on boundary face. Keep
1416  }
1417  else
1418  {
1419  // Open edge
1420  openEdgeIds.insert(edgei);
1421  }
1422  }
1423  }
1424 
1425  const label nOpenEdges
1426  (
1427  returnReduce(openEdgeIds.size(), sumOp<label>())
1428  );
1429 
1430  if (nOpenEdges)
1431  {
1432  const edgeList debugEdges
1433  (
1434  surfEdges,
1435  openEdgeIds.sortedToc()
1436  );
1437 
1438  vtk::lineWriter writer
1439  (
1440  s.points(),
1441  debugEdges,
1442  fileName
1443  (
1444  mesh_.time().globalPath()
1445  / ("isoSurfaceTopo." + timeDesc + "-openEdges")
1446  )
1447  );
1448 
1449  writer.writeGeometry();
1450 
1451  // CellData
1452  writer.beginCellData();
1453  writer.writeProcIDs();
1454 
1455  Info<< "isoSurfaceTopo : (debug) wrote "
1456  << nOpenEdges << " open edges: "
1457  << writer.output().name() << nl;
1458  }
1459  else
1460  {
1461  Info<< "isoSurfaceTopo : no open edges" << nl;
1462  }
1463  }
1464 
1465  // Write debug surface with snaps
1466  if (this->snap())
1467  {
1468  const Mesh& s = *this;
1469 
1470  vtk::surfaceWriter writer
1471  (
1472  s.points(),
1473  s,
1474  fileName
1475  (
1476  mesh_.time().globalPath()
1477  / ("isoSurfaceTopo." + timeDesc + "-surface")
1478  )
1479  );
1480 
1481  writer.writeGeometry();
1482 
1483  // CellData
1484  writer.beginCellData();
1485  writer.writeProcIDs();
1486  writer.write("cellID", meshCells_);
1487 
1488  // PointData
1489  writer.beginPointData();
1490  {
1491  // NB: may have non-compact surface points
1492  // --> use points().size() not nPoints()!
1493 
1494  labelList pointStatus(s.points().size(), Zero);
1495 
1496  forAll(pointToVerts_, i)
1497  {
1498  const edge& verts = pointToVerts_[i];
1499  if (verts.first() == verts.second())
1500  {
1501  // Duplicate index (ie, snapped)
1502  pointStatus[i] = 1;
1503  }
1504  }
1505 
1506  writer.write("point-status", pointStatus);
1507  }
1508 
1509  Info<< "isoSurfaceTopo : (debug) wrote "
1510  << returnReduce(s.size(), sumOp<label>())
1511  << " faces : "
1512  << writer.output().name() << nl;
1513  }
1514  }
1515  tetCutAddr.clearDebug();
1516 
1517 
1518  if (params.filter() == filterType::NONMANIFOLD)
1519  {
1520  // We remove verts on face diagonals. This is in fact just
1521  // straightening the edges of the face through the cell. This can
1522  // close off 'pockets' of triangles and create open or
1523  // multiply-connected triangles
1524 
1525  // Solved by eroding open-edges
1526  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1527 
1528  // The list of surface faces that should be retained after erosion
1529  Mesh& surf = *this;
1530  labelList faceAddr(identity(surf.size()));
1531 
1532  bitSet faceSelection;
1533 
1534  while (true)
1535  {
1536  // Shadow the surface for the purposes of erosion
1537  uindirectPrimitivePatch erosion
1538  (
1539  UIndirectList<face>(surf, faceAddr),
1540  surf.points()
1541  );
1542 
1543  faceSelection.clear();
1544  faceSelection.resize(erosion.size());
1545 
1546  const labelList& mp = erosion.meshPoints();
1547  const edgeList& surfEdges = erosion.edges();
1548  const labelListList& edgeFaces = erosion.edgeFaces();
1549 
1550  label nEdgeRemove = 0;
1551 
1552  forAll(edgeFaces, edgei)
1553  {
1554  const labelList& eFaces = edgeFaces[edgei];
1555  if (eFaces.size() == 1)
1556  {
1557  // Open edge (not originating from a boundary face)
1558 
1559  const edge& e = surfEdges[edgei];
1560  const edge& verts0 = pointToVerts_[mp[e.first()]];
1561  const edge& verts1 = pointToVerts_[mp[e.second()]];
1562 
1563  if
1564  (
1565  isProtectedPoint.test(verts0.first())
1566  && isProtectedPoint.test(verts0.second())
1567  && isProtectedPoint.test(verts1.first())
1568  && isProtectedPoint.test(verts1.second())
1569  )
1570  {
1571  // Open edge on boundary face. Keep
1572  }
1573  else
1574  {
1575  // Open edge. Mark for erosion
1576  faceSelection.set(eFaces[0]);
1577  ++nEdgeRemove;
1578  }
1579  }
1580  }
1581 
1582  if (debug)
1583  {
1584  Pout<< "isoSurfaceTopo :"
1585  << " removing " << faceSelection.count()
1586  << " / " << faceSelection.size()
1587  << " faces on " << nEdgeRemove << " open edges" << endl;
1588  }
1589 
1590  if (returnReduceAnd(faceSelection.none()))
1591  {
1592  break;
1593  }
1594 
1595  // Remove the faces from the addressing
1596  inplaceSubset(faceSelection, faceAddr, true); // True = remove
1597  }
1598 
1599 
1600  // Finished erosion (if any)
1601  // - retain the faces listed in the updated addressing
1602 
1603  if (surf.size() != faceAddr.size())
1604  {
1605  faceSelection.clear();
1606  faceSelection.resize(surf.size());
1607  faceSelection.set(faceAddr);
1608 
1609  inplaceSubsetMesh(faceSelection);
1610  }
1611  }
1612 }
1613 
1614 
1615 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1616 
1617 void Foam::isoSurfaceTopo::inplaceSubsetMesh(const bitSet& include)
1618 {
1619  labelList pointMap;
1621  Mesh filtered
1622  (
1623  Mesh::subsetMesh(include, pointMap, faceMap)
1624  );
1625  Mesh::transfer(filtered);
1626 
1627  meshCells_ = UIndirectList<label>(meshCells_, faceMap)();
1628 
1629  pointToVerts_ = UIndirectList<edge>(pointToVerts_, pointMap)();
1630 }
1631 
1632 
1633 // ************************************************************************* //
void ignoreCyclics()
Set ignoreBoundaryFaces to ignore cyclics (cyclicACMI)
Remove pyramid edge points, face-diagonals.
bitSet ignoreBoundaryFaces_
Optional boundary faces to ignore.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
unsigned int count(const bool on=true) const
Count number of bits set.
Definition: bitSetI.H:493
isoSurfaceTopo(const polyMesh &mesh, const scalarField &cellValues, const scalarField &pointValues, const scalar iso, const isoSurfaceParams &params=isoSurfaceParams(), const bitSet &ignoreCells=bitSet())
Construct from cell and point values.
#define ADD_SNAP_INDEX(pos, d1, d2, idx1, idx2)
label blockCells(UList< cutType > &cuts, const bitSet &ignoreCells) const
Mark ignoreCells as BLOCKED.
Marching tet iso surface algorithm with optional filtering to keep only points originating from mesh ...
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/"finiteArea-edgesCentres"))
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:132
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition: List.C:452
label nPoints() const noexcept
Number of mesh points.
const scalarField & cVals_
Cell values.
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:500
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1117
List< edge > edgeList
List of edge.
Definition: edgeList.H:32
tetrahedron< point, const point & > tetPointRef
A tetrahedron using referred points.
Definition: tetrahedron.H:72
Preferences for controlling iso-surface algorithms.
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:472
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
T & first()
Access first element of the list, position [0].
Definition: UList.H:814
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
const cellList & cells() const
SubList< label > subList
Declare type of subList.
Definition: List.H:122
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
label nFaces() const noexcept
Number of mesh faces.
void inplaceSubset(const BoolListType &select, ListType &input, const bool invert=false)
Inplace extract elements of the input list when select is true.
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
UList< bool > boolUList
A UList of bools.
Definition: UList.H:76
const labelListList & edgeLoops() const
Return list of closed loops of boundary vertices.
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
tmp< Field< Type > > interpolateTemplate(const Field< Type > &cellData, const Field< Type > &pointData) const
Interpolates cellData and pointData fields.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Low-level components common to various iso-surface algorithms.
const scalar iso_
Iso value.
MinMax< label > minMax(const labelHashSet &set)
Find the min/max values of labelHashSet.
Definition: hashSets.C:54
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1073
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:414
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
void inplaceSubsetMesh(const bitSet &include)
Subset the surface using the selected faces.
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
pointField vertices(const blockVertexList &bvl)
#define defineIsoSurfaceInterpolateMethods(ThisClass)
const boundBox & getClipBounds() const noexcept
Get optional clipping bounding box.
bool returnReduceAnd(const bool value, const label comm=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
const pointField & points
const polyMesh & mesh_
Reference to mesh.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:31
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
A class for handling words, derived from Foam::string.
Definition: word.H:63
const Time & time() const noexcept
Return time registry.
static void appendTriLabels(DynamicList< label > &verts, const label a, const label b, const label c, const bool flip)
label calcCellCuts(List< cutType > &cuts) const
Populate a list of candidate cell cuts using getCellCutType()
static labelList adjustTetBasePtIs(const polyMesh &mesh, const bool report=false)
Return an adjusted list of tet base points.
scalar mag() const
The magnitude/length of the bounding box diagonal.
Definition: boundBoxI.H:198
#define SNAP_END_VALUE(pos, val)
labelList meshCells_
For every face, the original cell in mesh.
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1111
const polyMesh & mesh() const noexcept
The mesh for which the iso-surface is associated.
const labelListList & edgeFaces() const
Return edge-face addressing.
label nInternalFaces() const noexcept
Number of internal faces.
label timeIndex() const noexcept
Return current time index.
Definition: TimeStateI.H:30
const Field< point_type > & points() const noexcept
Return reference to global points.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1098
A location outside the volume.
Definition: volumeType.H:66
const vectorField & cellCentres() const
bool test(const label pos) const
Test for True value at specified position, never auto-vivify entries.
Definition: bitSet.H:341
meshedSurface Mesh
static word printf(const char *fmt, const PrimitiveType &val)
Use a printf-style formatter for a primitive.
const direction noexcept
Definition: Scalar.H:258
const scalarField & pVals_
Point values.
int debug
Static debugging option.
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field...
T & last()
Access last element of the list, position [size()-1].
Definition: UList.H:828
static int getTetCutIndex(scalar p0, scalar p1, scalar p2, scalar p3, const scalar val, const bool doSnap) noexcept
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
static const Enum< filterType > filterNames
Names for the filtering types.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
label size() const
The surface size is the number of faces.
vector point
Point is a vector.
Definition: point.H:37
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label nCells() const noexcept
Number of mesh cells.
static constexpr Foam::label BLOCKED
Definition: regionSplit.C:37
const dimensionedScalar c
Speed of light in a vacuum.
List< surfZone > surfZoneList
List of surfZone.
Definition: surfZoneList.H:32
Nothing to be read.
void clearStorage()
Clear the list and delete storage.
Definition: DynamicListI.H:396
static const point & getMeshPointRef(const polyMesh &mesh, const label pointi)
messageStream Info
Information stream (stdout output on master, null elsewhere)
const boundBox & bounds() const noexcept
Return mesh bounding box.
Definition: polyMesh.H:592
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:73
List< label > labelList
A List of labels.
Definition: List.H:62
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
bool snap() const noexcept
Get point snapping flag.
Convenience macros for instantiating iso-surface interpolate methods.
filterType filter() const noexcept
Get current filter type.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Do not request registration (bool: false)
const volScalarField & p0
Definition: EEqn.H:36
const dimensionedScalar mp
Proton mass.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
fileName globalPath() const
Return global path for the case.
Definition: TimePathsI.H:73
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133