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