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-2023 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_.reserve(2*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  face& f = compactFaces.emplace_back(loop.size());
805 
806  label fpi = 0;
807  forAll(f, i)
808  {
809  const label pointi = mp[loop[i]];
810  if (filterDiag && pointFromDiag[pointi])
811  {
812  const label prevPointi = mp[loop[loop.fcIndex(i)]];
813  if
814  (
815  pointFromDiag[prevPointi]
816  && (pointToFace[pointi] != pointToFace[prevPointi])
817  )
818  {
819  f[fpi++] = pointi;
820  }
821  else
822  {
823  // Filter out diagonal point
824  }
825  }
826  else
827  {
828  f[fpi++] = pointi;
829  }
830  }
831 
832  if (fpi > 2)
833  {
834  f.resize(fpi);
835  }
836  else
837  {
838  // Keep original face
839  forAll(f, i)
840  {
841  const label pointi = mp[loop[i]];
842  f[i] = pointi;
843  }
844  }
845  compactCellIDs.append(cellID);
846  }
847  }
848 }
849 
850 
851 void Foam::isoSurfaceTopo::removeInsidePoints
852 (
853  Mesh& s,
854  const bool filterDiag,
855 
856  // inputs
857  const boolUList& pointFromDiag,
858  const labelUList& pointToFace,
859  const labelUList& start, // Per cell the starting triangle
860 
861  // outputs
862  DynamicList<label>& compactCellIDs // Per returned tri the cellID
863 )
864 {
865  const pointField& points = s.points();
866 
867  compactCellIDs.clear();
868  compactCellIDs.reserve(s.size()/4);
869 
870  DynamicList<face> compactFaces(s.size()/4);
871 
872  for (label celli = 0; celli < start.size()-1; ++celli)
873  {
874  // Triangles for the current cell
875 
876  const label nTris = start[celli+1]-start[celli];
877 
878  if (nTris)
879  {
880  const primitivePatch pp
881  (
882  SubList<face>(s, nTris, start[celli]),
883  points
884  );
885 
886  triangulateOutside
887  (
888  filterDiag,
889  pp,
890  pointFromDiag,
891  pointToFace,
892  celli,
893 
894  compactFaces,
895  compactCellIDs
896  );
897  }
898  }
899 
900  s.swapFaces(compactFaces); // Use new faces
901 }
902 
903 
904 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
905 
907 (
908  const polyMesh& mesh,
909  const scalarField& cellValues,
910  const scalarField& pointValues,
911  const scalar iso,
912  const isoSurfaceParams& params,
913  const bitSet& ignoreCells
914 )
915 :
916  isoSurfaceBase(mesh, cellValues, pointValues, iso, params)
917 {
918  // The cell cut type
919  List<cutType> cellCutType_(mesh.nCells(), cutType::UNVISITED);
920 
921  // Time description (for debug output)
922  const word timeDesc(word::printf("%08d", mesh_.time().timeIndex()));
923 
924  if (debug)
925  {
926  Pout<< "isoSurfaceTopo:" << nl
927  << " cell min/max : " << minMax(cVals_) << nl
928  << " point min/max : " << minMax(pVals_) << nl
929  << " isoValue : " << iso << nl
930  << " filter : "
931  << isoSurfaceParams::filterNames[params.filter()] << nl
932  << " mesh span : " << mesh.bounds().mag() << nl
933  << " ignoreCells : " << ignoreCells.count()
934  << " / " << cVals_.size() << nl
935  << endl;
936  }
937 
938  this->ignoreCyclics();
939 
940  label nBlockedCells = 0;
941 
942  // Mark ignoreCells as blocked
943  nBlockedCells += blockCells(cellCutType_, ignoreCells);
944 
945  // Mark cells outside bounding box as blocked
946  nBlockedCells +=
947  blockCells(cellCutType_, params.getClipBounds(), volumeType::OUTSIDE);
948 
949  // Adjusted tet base points to improve tet quality
950  labelList tetBasePtIs
951  (
953  );
954 
955 
956  // Determine cell cuts
957  const label nCutCells = calcCellCuts(cellCutType_);
958 
959  if (debug)
960  {
961  Pout<< "isoSurfaceTopo : candidate cells cut "
962  << nCutCells
963  << " blocked " << nBlockedCells
964  << " total " << mesh_.nCells() << endl;
965  }
966 
967  if (debug && isA<fvMesh>(mesh))
968  {
969  const auto& fvmesh = refCast<const fvMesh>(mesh);
970 
971  volScalarField debugField
972  (
973  IOobject
974  (
975  "isoSurfaceTopo.cutType",
976  fvmesh.time().timeName(),
977  fvmesh.time(),
981  ),
982  fvmesh,
984  );
985 
986  auto& debugFld = debugField.primitiveFieldRef();
987 
988  forAll(cellCutType_, celli)
989  {
990  debugFld[celli] = cellCutType_[celli];
991  }
992 
993  Info<< "Writing cut types: " << debugField.objectRelPath() << nl;
994  debugField.write();
995  }
996 
997  // Additional debugging
998  if (debug & 8)
999  {
1000  // Write debug cuts cells in VTK format
1001  {
1002  constexpr uint8_t realCut(cutType::CUT | cutType::TETCUT);
1003  labelList debugCutCells(nCutCells, Zero);
1004 
1005  label nout = 0;
1006  forAll(cellCutType_, celli)
1007  {
1008  if ((cellCutType_[celli] & realCut) != 0)
1009  {
1010  debugCutCells[nout] = celli;
1011  ++nout;
1012  if (nout >= nCutCells) break;
1013  }
1014  }
1015 
1016  // The mesh subset cut
1017  vtk::vtuCells vtuCells;
1018  vtuCells.reset(mesh_, debugCutCells);
1019 
1020  vtk::internalMeshWriter writer
1021  (
1022  mesh_,
1023  vtuCells,
1024  fileName
1025  (
1026  mesh_.time().globalPath()
1027  / ("isoSurfaceTopo." + timeDesc + "-cutCells")
1028  )
1029  );
1030 
1031  writer.writeGeometry();
1032 
1033  // CellData
1034  writer.beginCellData();
1035  writer.writeProcIDs();
1036  writer.writeCellData("cutField", cVals_);
1037 
1038  // PointData
1039  writer.beginPointData();
1040  writer.writePointData("cutField", pVals_);
1041 
1042  Info<< "isoSurfaceTopo : (debug) wrote "
1043  << returnReduce(nCutCells, sumOp<label>())
1044  << " cut cells: "
1045  << writer.output().name() << nl;
1046  }
1047  }
1048 
1049 
1050  tetCutAddressing tetCutAddr
1051  (
1052  nCutCells,
1053  this->snap(),
1054  (debug & 8) // Enable debug tets
1055  );
1056 
1057  labelList startTri(mesh_.nCells()+1, Zero);
1058  for (label celli = 0; celli < mesh_.nCells(); ++celli)
1059  {
1060  startTri[celli] = tetCutAddr.nFaces();
1061  if ((cellCutType_[celli] & cutType::ANYCUT) != 0)
1062  {
1063  generateTriPoints
1064  (
1065  celli,
1066  // Same as tetMatcher::test(mesh_, celli),
1067  bool(cellCutType_[celli] & cutType::TETCUT),
1068 
1069  tetBasePtIs,
1070  tetCutAddr
1071  );
1072  }
1073  }
1074  startTri.back() = tetCutAddr.nFaces();
1075 
1076  // Information not needed anymore:
1077  tetBasePtIs.clear();
1078  tetCutAddr.clearHashes();
1079 
1080 
1081  // From list of vertices -> triangular faces
1082  faceList allTriFaces(startTri.back());
1083  {
1084  auto& verts = tetCutAddr.cutPoints();
1085 
1086  label verti = 0;
1087  for (face& f : allTriFaces)
1088  {
1089  f.resize(3);
1090  f[0] = verts[verti++];
1091  f[1] = verts[verti++];
1092  f[2] = verts[verti++];
1093  }
1094  verts.clearStorage(); // Not needed anymore
1095  }
1096 
1097 
1098  // The cells cut by the triangular faces
1099  meshCells_.resize(startTri.back());
1100  for (label celli = 0; celli < startTri.size()-1; ++celli)
1101  {
1102  // All triangles for the current cell
1104  (
1105  meshCells_,
1106  (startTri[celli+1] - startTri[celli]),
1107  startTri[celli]
1108  ) = celli;
1109  }
1110 
1111 
1112  pointToVerts_.transfer(tetCutAddr.pointToVerts());
1113 
1114  pointField allTriPoints
1115  (
1116  this->interpolateTemplate
1117  (
1118  mesh_.cellCentres(),
1119  mesh_.points()
1120  )
1121  );
1122 
1123 
1124  // Assign to MeshedSurface
1125  static_cast<Mesh&>(*this) = Mesh
1126  (
1127  std::move(allTriPoints),
1128  std::move(allTriFaces),
1129  surfZoneList() // zones not required (one zone)
1130  );
1131 
1132  if (debug)
1133  {
1134  Pout<< "isoSurfaceTopo : generated "
1135  << Mesh::size() << " triangles "
1136  << Mesh::points().size() << " points" << endl;
1137  }
1138 
1139  // Write debug triangulated surface
1140  if ((debug & 8) && (params.filter() != filterType::NONE))
1141  {
1142  const Mesh& s = *this;
1143 
1144  vtk::surfaceWriter writer
1145  (
1146  s.points(),
1147  s,
1148  fileName
1149  (
1150  mesh_.time().globalPath()
1151  / ("isoSurfaceTopo." + timeDesc + "-triangles")
1152  )
1153  );
1154 
1155  writer.writeGeometry();
1156 
1157  // CellData
1158  writer.beginCellData();
1159  writer.writeProcIDs();
1160  writer.write("cellID", meshCells_);
1161 
1162  // PointData
1163  writer.beginPointData();
1164  {
1165  // NB: may have non-compact surface points
1166  // --> use points().size() not nPoints()!
1167 
1168  labelList pointStatus(s.points().size(), Zero);
1169 
1170  forAll(pointToVerts_, i)
1171  {
1172  const edge& verts = pointToVerts_[i];
1173  if (verts.first() == verts.second())
1174  {
1175  // Duplicate index (ie, snapped)
1176  pointStatus[i] = 1;
1177  }
1178  if (tetCutAddr.pointFromDiag().test(i))
1179  {
1180  // Point on triangulation diagonal
1181  pointStatus[i] = -1;
1182  }
1183  }
1184 
1185  writer.write("point-status", pointStatus);
1186  }
1187 
1188  Info<< "isoSurfaceTopo : (debug) wrote "
1189  << returnReduce(s.size(), sumOp<label>())
1190  << " triangles : "
1191  << writer.output().name() << nl;
1192  }
1193 
1194 
1195  // Now:
1196  // - generated faces and points are assigned to *this
1197  // - per point we know:
1198  // - pointOnDiag: whether it is on a face-diagonal edge
1199  // - pointToFace: from what pyramid (cell+face) it was produced
1200  // (note that the pyramid faces are shared between multiple mesh faces)
1201  // - pointToVerts_ : originating mesh vertex or cell centre
1202 
1203  if (params.filter() == filterType::NONE)
1204  {
1205  // Compact out unused (snapped) points
1206  if (this->snap())
1207  {
1208  Mesh& s = *this;
1209 
1210  labelList pointMap; // Back to original point
1211  s.compactPoints(pointMap); // Compact out unused points
1212  pointToVerts_ = UIndirectList<edge>(pointToVerts_, pointMap)();
1213  }
1214  }
1215  else
1216  {
1217  // Initial filtering
1218 
1219  Mesh& s = *this;
1220 
1221  // Triangulate outside
1222  // (filter edges to cell centres and optionally face diagonals)
1223  DynamicList<label> compactCellIDs; // Per tri the cell
1224 
1225  removeInsidePoints
1226  (
1227  *this,
1228  // Filter face diagonal
1229  (
1230  params.filter() == filterType::DIAGCELL
1231  || params.filter() == filterType::NONMANIFOLD
1232  ),
1233  tetCutAddr.pointFromDiag(),
1234  tetCutAddr.pointToFace(),
1235  startTri,
1236  compactCellIDs
1237  );
1238 
1239  labelList pointMap; // Back to original point
1240  s.compactPoints(pointMap); // Compact out unused points
1241 
1242  pointToVerts_ = UIndirectList<edge>(pointToVerts_, pointMap)();
1243  meshCells_.transfer(compactCellIDs);
1244 
1245  if (debug)
1246  {
1247  Pout<< "isoSurfaceTopo :"
1248  " after removing cell centre and face-diag triangles: "
1249  << Mesh::size() << " faces "
1250  << Mesh::points().size() << " points"
1251  << endl;
1252  }
1253  }
1254 
1255  // Diagonal filter information not needed anymore
1256  tetCutAddr.clearDiagonal();
1257 
1258 
1259  // For more advanced filtering (eg, removal of open edges)
1260  // need the boundary and other 'protected' points
1261 
1262  bitSet isProtectedPoint;
1263  if
1264  (
1265  (params.filter() == filterType::NONMANIFOLD)
1266  || tetCutAddr.debugCutTetsOn()
1267  )
1268  {
1269  // Mark points on mesh outside as 'protected'
1270  // - never erode these edges
1271 
1272  isProtectedPoint.resize(mesh_.nPoints());
1273 
1274  for
1275  (
1276  label facei = mesh_.nInternalFaces();
1277  facei < mesh_.nFaces();
1278  ++facei
1279  )
1280  {
1281  isProtectedPoint.set(mesh_.faces()[facei]);
1282  }
1283 
1284  // Include faces that would be exposed from mesh subset
1285  if (nBlockedCells)
1286  {
1287  const labelList& faceOwn = mesh_.faceOwner();
1288  const labelList& faceNei = mesh_.faceNeighbour();
1289 
1290  for (label facei = 0; facei < mesh.nInternalFaces(); ++facei)
1291  {
1292  // If only one cell is blocked, the face corresponds
1293  // to an exposed subMesh face
1294  if
1295  (
1296  (cellCutType_[faceOwn[facei]] == cutType::BLOCKED)
1297  != (cellCutType_[faceNei[facei]] == cutType::BLOCKED)
1298  )
1299  {
1300  isProtectedPoint.set(mesh_.faces()[facei]);
1301  }
1302  }
1303  }
1304  }
1305 
1306  // Initial cell cut information not needed anymore
1307  cellCutType_.clear();
1308 
1309 
1310  // Additional debugging
1311  if (tetCutAddr.debugCutTetsOn())
1312  {
1313  // Write debug cut tets in VTK format
1314  {
1315  const auto& debugCuts = tetCutAddr.debugCutTets();
1316 
1317  // The TET shapes, using the mesh_ points information
1318  vtk::vtuCells vtuCells;
1319  vtuCells.resetShapes(debugCuts);
1320 
1321  // Use all points and all cell-centres
1322  vtuCells.setNumPoints(mesh_.nPoints());
1323  vtuCells.addPointCellLabels(identity(mesh_.nCells()));
1324 
1325  vtk::internalMeshWriter writer
1326  (
1327  mesh_,
1328  vtuCells,
1329  fileName
1330  (
1331  mesh_.time().globalPath()
1332  / ("isoSurfaceTopo." + timeDesc + "-cutTets")
1333  )
1334  );
1335 
1336  writer.writeGeometry();
1337 
1338  // CellData
1339  writer.beginCellData();
1340  writer.writeProcIDs();
1341 
1342  // Quality of the cut tets
1343  {
1344  Field<scalar> cutTetQuality(debugCuts.size());
1345  forAll(cutTetQuality, teti)
1346  {
1347  cutTetQuality[teti] = tetPointRef
1348  (
1349  getMeshPointRef(mesh_, debugCuts[teti][0]),
1350  getMeshPointRef(mesh_, debugCuts[teti][1]),
1351  getMeshPointRef(mesh_, debugCuts[teti][2]),
1352  getMeshPointRef(mesh_, debugCuts[teti][3])
1353  ).quality();
1354  }
1355  writer.writeCellData("tetQuality", cutTetQuality);
1356  }
1357 
1358  // PointData
1359  if (this->snap())
1360  {
1361  writer.beginPointData();
1362 
1363  labelList pointStatus(vtuCells.nFieldPoints(), Zero);
1364 
1365  for (const edge& verts : pointToVerts_)
1366  {
1367  if (verts.first() == verts.second())
1368  {
1369  // Duplicate index (ie, snapped)
1370  pointStatus[verts.first()] = 1;
1371  }
1372  }
1373 
1374  writer.writePointData("point-status", pointStatus);
1375  }
1376 
1377  Info<< "isoSurfaceTopo : (debug) wrote "
1378  << returnReduce(debugCuts.size(), sumOp<label>())
1379  << " cut tets: "
1380  << writer.output().name() << nl;
1381  }
1382 
1383  // Determining open edges. Same logic as used later...
1384 
1385  labelHashSet openEdgeIds(0);
1386 
1387  {
1388  const Mesh& s = *this;
1389 
1390  const labelList& mp = s.meshPoints();
1391  const edgeList& surfEdges = s.edges();
1392  const labelListList& edgeFaces = s.edgeFaces();
1393  openEdgeIds.reserve(s.size());
1394 
1395  forAll(edgeFaces, edgei)
1396  {
1397  const labelList& eFaces = edgeFaces[edgei];
1398  if (eFaces.size() == 1)
1399  {
1400  // Open edge (not originating from a boundary face)
1401 
1402  const edge& e = surfEdges[edgei];
1403  const edge& verts0 = pointToVerts_[mp[e.first()]];
1404  const edge& verts1 = pointToVerts_[mp[e.second()]];
1405 
1406  if
1407  (
1408  isProtectedPoint.test(verts0.first())
1409  && isProtectedPoint.test(verts0.second())
1410  && isProtectedPoint.test(verts1.first())
1411  && isProtectedPoint.test(verts1.second())
1412  )
1413  {
1414  // Open edge on boundary face. Keep
1415  }
1416  else
1417  {
1418  // Open edge
1419  openEdgeIds.insert(edgei);
1420  }
1421  }
1422  }
1423 
1424  const label nOpenEdges
1425  (
1426  returnReduce(openEdgeIds.size(), sumOp<label>())
1427  );
1428 
1429  if (nOpenEdges)
1430  {
1431  const edgeList debugEdges
1432  (
1433  surfEdges,
1434  openEdgeIds.sortedToc()
1435  );
1436 
1437  vtk::lineWriter writer
1438  (
1439  s.points(),
1440  debugEdges,
1441  fileName
1442  (
1443  mesh_.time().globalPath()
1444  / ("isoSurfaceTopo." + timeDesc + "-openEdges")
1445  )
1446  );
1447 
1448  writer.writeGeometry();
1449 
1450  // CellData
1451  writer.beginCellData();
1452  writer.writeProcIDs();
1453 
1454  Info<< "isoSurfaceTopo : (debug) wrote "
1455  << nOpenEdges << " open edges: "
1456  << writer.output().name() << nl;
1457  }
1458  else
1459  {
1460  Info<< "isoSurfaceTopo : no open edges" << nl;
1461  }
1462  }
1463 
1464  // Write debug surface with snaps
1465  if (this->snap())
1466  {
1467  const Mesh& s = *this;
1468 
1469  vtk::surfaceWriter writer
1470  (
1471  s.points(),
1472  s,
1473  fileName
1474  (
1475  mesh_.time().globalPath()
1476  / ("isoSurfaceTopo." + timeDesc + "-surface")
1477  )
1478  );
1479 
1480  writer.writeGeometry();
1481 
1482  // CellData
1483  writer.beginCellData();
1484  writer.writeProcIDs();
1485  writer.write("cellID", meshCells_);
1486 
1487  // PointData
1488  writer.beginPointData();
1489  {
1490  // NB: may have non-compact surface points
1491  // --> use points().size() not nPoints()!
1492 
1493  labelList pointStatus(s.points().size(), Zero);
1494 
1495  forAll(pointToVerts_, i)
1496  {
1497  const edge& verts = pointToVerts_[i];
1498  if (verts.first() == verts.second())
1499  {
1500  // Duplicate index (ie, snapped)
1501  pointStatus[i] = 1;
1502  }
1503  }
1504 
1505  writer.write("point-status", pointStatus);
1506  }
1507 
1508  Info<< "isoSurfaceTopo : (debug) wrote "
1509  << returnReduce(s.size(), sumOp<label>())
1510  << " faces : "
1511  << writer.output().name() << nl;
1512  }
1513  }
1514  tetCutAddr.clearDebug();
1515 
1516 
1517  if (params.filter() == filterType::NONMANIFOLD)
1518  {
1519  // We remove verts on face diagonals. This is in fact just
1520  // straightening the edges of the face through the cell. This can
1521  // close off 'pockets' of triangles and create open or
1522  // multiply-connected triangles
1523 
1524  // Solved by eroding open-edges
1525  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1526 
1527  // The list of surface faces that should be retained after erosion
1528  Mesh& surf = *this;
1529  labelList faceAddr(identity(surf.size()));
1530 
1531  bitSet faceSelection;
1532 
1533  while (true)
1534  {
1535  // Shadow the surface for the purposes of erosion
1536  uindirectPrimitivePatch erosion
1537  (
1538  UIndirectList<face>(surf, faceAddr),
1539  surf.points()
1540  );
1541 
1542  faceSelection.clear();
1543  faceSelection.resize(erosion.size());
1544 
1545  const labelList& mp = erosion.meshPoints();
1546  const edgeList& surfEdges = erosion.edges();
1547  const labelListList& edgeFaces = erosion.edgeFaces();
1548 
1549  label nEdgeRemove = 0;
1550 
1551  forAll(edgeFaces, edgei)
1552  {
1553  const labelList& eFaces = edgeFaces[edgei];
1554  if (eFaces.size() == 1)
1555  {
1556  // Open edge (not originating from a boundary face)
1557 
1558  const edge& e = surfEdges[edgei];
1559  const edge& verts0 = pointToVerts_[mp[e.first()]];
1560  const edge& verts1 = pointToVerts_[mp[e.second()]];
1561 
1562  if
1563  (
1564  isProtectedPoint.test(verts0.first())
1565  && isProtectedPoint.test(verts0.second())
1566  && isProtectedPoint.test(verts1.first())
1567  && isProtectedPoint.test(verts1.second())
1568  )
1569  {
1570  // Open edge on boundary face. Keep
1571  }
1572  else
1573  {
1574  // Open edge. Mark for erosion
1575  faceSelection.set(eFaces[0]);
1576  ++nEdgeRemove;
1577  }
1578  }
1579  }
1580 
1581  if (debug)
1582  {
1583  Pout<< "isoSurfaceTopo :"
1584  << " removing " << faceSelection.count()
1585  << " / " << faceSelection.size()
1586  << " faces on " << nEdgeRemove << " open edges" << endl;
1587  }
1588 
1589  if (returnReduceAnd(faceSelection.none()))
1590  {
1591  break;
1592  }
1593 
1594  // Remove the faces from the addressing
1595  inplaceSubset(faceSelection, faceAddr, true); // True = remove
1596  }
1597 
1598 
1599  // Finished erosion (if any)
1600  // - retain the faces listed in the updated addressing
1601 
1602  if (surf.size() != faceAddr.size())
1603  {
1604  faceSelection.clear();
1605  faceSelection.resize(surf.size());
1606  faceSelection.set(faceAddr);
1607 
1608  inplaceSubsetMesh(faceSelection);
1609  }
1610  }
1611 }
1612 
1613 
1614 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1615 
1616 void Foam::isoSurfaceTopo::inplaceSubsetMesh(const bitSet& include)
1617 {
1618  labelList pointMap;
1620  Mesh filtered
1621  (
1622  Mesh::subsetMesh(include, pointMap, faceMap)
1623  );
1624  Mesh::transfer(filtered);
1625 
1626  meshCells_ = UIndirectList<label>(meshCells_, faceMap)();
1627 
1628  pointToVerts_ = UIndirectList<edge>(pointToVerts_, pointMap)();
1629 }
1630 
1631 
1632 // ************************************************************************* //
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:426
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:160
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition: List.C:326
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:517
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1122
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:489
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
T & first()
Access first element of the list, position [0].
Definition: UList.H:853
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const cellList & cells() const
SubList< label > subList
Declare type of subList.
Definition: List.H:144
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:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
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), works like std::iota() but returning a...
Definition: labelLists.C:44
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
fileName globalPath() const
Return global path for the case = rootPath/globalCaseName. Same as TimePaths::globalPath() ...
Definition: Time.H:512
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:1116
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 the current time index.
Definition: TimeStateI.H:43
const Field< point_type > & points() const noexcept
Return reference to global points.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
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:326
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 & emplace_back(Args &&... args)
Construct an element at the end of the list, return reference to the new list element.
Definition: ListI.H:212
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:412
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:616
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
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.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127