isoSurfaceCell.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2016-2020 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 "isoSurfaceCell.H"
30#include "isoSurfacePoint.H"
31#include "dictionary.H"
32#include "polyMesh.H"
33#include "mergePoints.H"
34#include "tetMatcher.H"
35#include "syncTools.H"
36#include "triSurface.H"
37#include "triSurfaceTools.H"
38#include "Time.H"
39#include "triPoints.H"
40
41// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42
45
46
47// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
48
49namespace Foam
50{
52}
53
54
55// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
56
57Foam::scalar Foam::isoSurfaceCell::isoFraction
58(
59 const scalar s0,
60 const scalar s1
61) const
62{
63 const scalar d = s1-s0;
64
65 if (mag(d) > VSMALL)
66 {
67 return (iso_-s0)/d;
68 }
69
70 return -1.0;
71}
72
73
74Foam::labelPair Foam::isoSurfaceCell::findCommonPoints
75(
76 const labelledTri& tri0,
77 const labelledTri& tri1
78)
79{
80 labelPair common(-1, -1);
81
82 label fp0 = 0;
83 label fp1 = tri1.find(tri0[fp0]);
84
85 if (fp1 == -1)
86 {
87 fp0 = 1;
88 fp1 = tri1.find(tri0[fp0]);
89 }
90
91 if (fp1 != -1)
92 {
93 // So tri0[fp0] is tri1[fp1]
94
95 // Find next common point
96 label fp0p1 = tri0.fcIndex(fp0);
97 label fp1p1 = tri1.fcIndex(fp1);
98 label fp1m1 = tri1.rcIndex(fp1);
99
100 if (tri0[fp0p1] == tri1[fp1p1] || tri0[fp0p1] == tri1[fp1m1])
101 {
102 common[0] = tri0[fp0];
103 common[1] = tri0[fp0p1];
104 }
105 }
106 return common;
107}
108
109
110Foam::point Foam::isoSurfaceCell::calcCentre(const triSurface& s)
111{
112 vector sum = Zero;
113
114 forAll(s, i)
115 {
116 sum += s[i].centre(s.points());
117 }
118 return sum/s.size();
119}
120
121
122Foam::pointIndexHit Foam::isoSurfaceCell::collapseSurface
123(
124 const label celli,
125 pointField& localPoints,
126 DynamicList<labelledTri, 64>& localTris
127) const
128{
129 pointIndexHit info(false, Zero, localTris.size());
130
131 if (localTris.size() == 1)
132 {
133 const labelledTri& tri = localTris[0];
134 info.setPoint(tri.centre(localPoints));
135 info.setHit();
136 }
137 else if (localTris.size() == 2)
138 {
139 // Check if the two triangles share an edge.
140 const labelledTri& tri0 = localTris[0];
141 const labelledTri& tri1 = localTris[1];
142
143 labelPair shared = findCommonPoints(tri0, tri1);
144
145 if (shared[0] != -1)
146 {
147 const vector n0 = tri0.areaNormal(localPoints);
148 const vector n1 = tri1.areaNormal(localPoints);
149
150 // Merge any zero-sized triangles,
151 // or if they point in the same direction.
152
153 if
154 (
155 mag(n0) <= ROOTVSMALL
156 || mag(n1) <= ROOTVSMALL
157 || (n0 & n1) >= 0
158 )
159 {
160 info.setPoint
161 (
162 0.5
163 * (
164 tri0.centre(localPoints)
165 + tri1.centre(localPoints)
166 )
167 );
168 info.setHit();
169 }
170 }
171 }
172 else if (localTris.size())
173 {
174 // Check if single region. Rare situation.
175 triSurface surf
176 (
177 localTris,
179 localPoints,
180 true
181 );
182 localTris.clearStorage();
183
184 labelList faceZone;
185 label nZones = surf.markZones
186 (
187 boolList(surf.nEdges(), false),
188 faceZone
189 );
190
191 if (nZones == 1)
192 {
193 // Check that all normals make a decent angle
194 scalar minCos = GREAT;
195 const vector& n0 = surf.faceNormals()[0];
196 for (label i = 1; i < surf.size(); ++i)
197 {
198 scalar cosAngle = (n0 & surf.faceNormals()[i]);
199 if (cosAngle < minCos)
200 {
201 minCos = cosAngle;
202 }
203 }
204
205 if (minCos > 0)
206 {
207 info.setPoint(calcCentre(surf));
208 info.setHit();
209 }
210 }
211 }
212
213 return info;
214}
215
216
217void Foam::isoSurfaceCell::calcSnappedCc
218(
219 const bitSet& isTet,
220 const scalarField& cVals,
221 const scalarField& pVals,
222
223 DynamicList<point>& snappedPoints,
224 labelList& snappedCc
225) const
226{
227 const pointField& cc = mesh_.cellCentres();
228 const pointField& pts = mesh_.points();
229
230 snappedCc.setSize(mesh_.nCells());
231 snappedCc = -1;
232
233 // Work arrays
234 DynamicList<point, 64> localPoints(64);
235 DynamicList<labelledTri, 64> localTris(64);
236 Map<label> pointToLocal(64);
237
238 forAll(mesh_.cells(), celli)
239 {
240 if (cellCutType_[celli] == cutType::CUT) // No tet cuts
241 {
242 const scalar cVal = cVals[celli];
243
244 const cell& cFaces = mesh_.cells()[celli];
245
246 localPoints.clear();
247 localTris.clear();
248 pointToLocal.clear();
249
250 // Create points for all intersections close to cell centre
251 // (i.e. from pyramid edges)
252
253 for (const label facei : cFaces)
254 {
255 const face& f = mesh_.faces()[facei];
256
257 for (const label pointi : f)
258 {
259 scalar s = isoFraction(cVal, pVals[pointi]);
260
261 if (s >= 0.0 && s <= 0.5)
262 {
263 if (pointToLocal.insert(pointi, localPoints.size()))
264 {
265 localPoints.append((1.0-s)*cc[celli]+s*pts[pointi]);
266 }
267 }
268 }
269 }
270
271 if (localPoints.size() == 1)
272 {
273 // No need for any analysis.
274 snappedCc[celli] = snappedPoints.size();
275 snappedPoints.append(localPoints[0]);
276
277 //Pout<< "cell:" << celli
278 // << " at " << mesh_.cellCentres()[celli]
279 // << " collapsing " << localPoints
280 // << " intersections down to "
281 // << snappedPoints[snappedCc[celli]] << endl;
282 }
283 else if (localPoints.size() == 2)
284 {
285 //? No need for any analysis.???
286 snappedCc[celli] = snappedPoints.size();
287 snappedPoints.append(0.5*(localPoints[0]+localPoints[1]));
288
289 //Pout<< "cell:" << celli
290 // << " at " << mesh_.cellCentres()[celli]
291 // << " collapsing " << localPoints
292 // << " intersections down to "
293 // << snappedPoints[snappedCc[celli]] << endl;
294 }
295 else if (localPoints.size())
296 {
297 // Need to analyse
298 for (const label facei : cFaces)
299 {
300 const face& f = mesh_.faces()[facei];
301
302 // Do a tetrahedralisation. Each face to cc becomes pyr.
303 // Each pyr gets split into tets by diagonalisation
304 // of face.
305
306 const label fp0 = mesh_.tetBasePtIs()[facei];
307 label fp = f.fcIndex(fp0);
308 for (label i = 2; i < f.size(); ++i)
309 {
310 label nextFp = f.fcIndex(fp);
311 triFace tri(f[fp0], f[fp], f[nextFp]);
312
313 // Get fractions for the three edges to cell centre
314 FixedList<scalar, 3> s(3);
315 s[0] = isoFraction(cVal, pVals[tri[0]]);
316 s[1] = isoFraction(cVal, pVals[tri[1]]);
317 s[2] = isoFraction(cVal, pVals[tri[2]]);
318
319 if
320 (
321 (s[0] >= 0.0 && s[0] <= 0.5)
322 && (s[1] >= 0.0 && s[1] <= 0.5)
323 && (s[2] >= 0.0 && s[2] <= 0.5)
324 )
325 {
326 if
327 (
328 (mesh_.faceOwner()[facei] == celli)
329 == (cVal >= pVals[tri[0]])
330 )
331 {
332 localTris.append
333 (
334 labelledTri
335 (
336 pointToLocal[tri[1]],
337 pointToLocal[tri[0]],
338 pointToLocal[tri[2]],
339 0
340 )
341 );
342 }
343 else
344 {
345 localTris.append
346 (
347 labelledTri
348 (
349 pointToLocal[tri[0]],
350 pointToLocal[tri[1]],
351 pointToLocal[tri[2]],
352 0
353 )
354 );
355 }
356 }
357
358 fp = nextFp;
359 }
360 }
361
362 pointField surfPoints;
363 surfPoints.transfer(localPoints);
364 pointIndexHit info = collapseSurface
365 (
366 celli,
367 surfPoints,
368 localTris
369 );
370
371 if (info.hit())
372 {
373 snappedCc[celli] = snappedPoints.size();
374 snappedPoints.append(info.hitPoint());
375
376 //Pout<< "cell:" << celli
377 // << " at " << mesh_.cellCentres()[celli]
378 // << " collapsing " << surfPoints
379 // << " intersections down to "
380 // << snappedPoints[snappedCc[celli]] << endl;
381 }
382 }
383 }
384 }
385}
386
387
388void Foam::isoSurfaceCell::genPointTris
389(
390 const scalarField& cellValues,
391 const scalarField& pointValues,
392 const label pointi,
393 const label facei,
394 const label celli,
395 DynamicList<point, 64>& localTriPoints
396) const
397{
398 const pointField& cc = mesh_.cellCentres();
399 const pointField& pts = mesh_.points();
400 const face& f = mesh_.faces()[facei];
401
402 const label fp0 = mesh_.tetBasePtIs()[facei];
403 label fp = f.fcIndex(fp0);
404 for (label i = 2; i < f.size(); ++i)
405 {
406 label nextFp = f.fcIndex(fp);
407 triFace tri(f[fp0], f[fp], f[nextFp]);
408
409 label index = tri.find(pointi);
410
411 if (index == -1)
412 {
413 continue;
414 }
415
416 // Tet between index..index-1, index..index+1, index..cc
417 label b = tri[tri.fcIndex(index)];
418 label c = tri[tri.rcIndex(index)];
419
420 // Get fractions for the three edges emanating from point
421 FixedList<scalar, 3> s(3);
422 s[0] = isoFraction(pointValues[pointi], pointValues[b]);
423 s[1] = isoFraction(pointValues[pointi], pointValues[c]);
424 s[2] = isoFraction(pointValues[pointi], cellValues[celli]);
425
426 if
427 (
428 (s[0] >= 0.0 && s[0] <= 0.5)
429 && (s[1] >= 0.0 && s[1] <= 0.5)
430 && (s[2] >= 0.0 && s[2] <= 0.5)
431 )
432 {
433 point p0 = (1.0-s[0])*pts[pointi] + s[0]*pts[b];
434 point p1 = (1.0-s[1])*pts[pointi] + s[1]*pts[c];
435 point p2 = (1.0-s[2])*pts[pointi] + s[2]*cc[celli];
436
437 if
438 (
439 (mesh_.faceOwner()[facei] == celli)
440 == (pointValues[pointi] > cellValues[celli])
441 )
442 {
443 localTriPoints.append(p0);
444 localTriPoints.append(p1);
445 localTriPoints.append(p2);
446 }
447 else
448 {
449 localTriPoints.append(p1);
450 localTriPoints.append(p0);
451 localTriPoints.append(p2);
452 }
453 }
454
455 fp = nextFp;
456 }
457}
458
459
460void Foam::isoSurfaceCell::genPointTris
461(
462 const scalarField& pointValues,
463 const label pointi,
464 const label facei,
465 const label celli,
466 DynamicList<point, 64>& localTriPoints
467) const
468{
469 const pointField& pts = mesh_.points();
470 const cell& cFaces = mesh_.cells()[celli];
471
472 // Make tet from this face to the 4th point (same as cellcentre in
473 // non-tet cells)
474 const face& f = mesh_.faces()[facei];
475
476 // Find 4th point
477 label ccPointi = -1;
478 for (const label cfacei : cFaces)
479 {
480 const face& f1 = mesh_.faces()[cfacei];
481 for (const label p1 : f1)
482 {
483 if (!f.found(p1))
484 {
485 ccPointi = p1;
486 break;
487 }
488 }
489 if (ccPointi != -1)
490 {
491 break;
492 }
493 }
494
495
496 // Tet between index..index-1, index..index+1, index..cc
497 label index = f.find(pointi);
498 label b = f[f.fcIndex(index)];
499 label c = f[f.rcIndex(index)];
500
501 //Pout<< " p0:" << pointi << " b:" << b << " c:" << c
502 //<< " d:" << ccPointi << endl;
503
504 // Get fractions for the three edges emanating from point
505 FixedList<scalar, 3> s(3);
506 s[0] = isoFraction(pointValues[pointi], pointValues[b]);
507 s[1] = isoFraction(pointValues[pointi], pointValues[c]);
508 s[2] = isoFraction(pointValues[pointi], pointValues[ccPointi]);
509
510 if
511 (
512 (s[0] >= 0.0 && s[0] <= 0.5)
513 && (s[1] >= 0.0 && s[1] <= 0.5)
514 && (s[2] >= 0.0 && s[2] <= 0.5)
515 )
516 {
517 point p0 = (1.0-s[0])*pts[pointi] + s[0]*pts[b];
518 point p1 = (1.0-s[1])*pts[pointi] + s[1]*pts[c];
519 point p2 = (1.0-s[2])*pts[pointi] + s[2]*pts[ccPointi];
520
521 if (mesh_.faceOwner()[facei] != celli)
522 {
523 localTriPoints.append(p0);
524 localTriPoints.append(p1);
525 localTriPoints.append(p2);
526 }
527 else
528 {
529 localTriPoints.append(p1);
530 localTriPoints.append(p0);
531 localTriPoints.append(p2);
532 }
533 }
534}
535
536
537void Foam::isoSurfaceCell::calcSnappedPoint
538(
539 const bitSet& isTet,
540 const scalarField& cVals,
541 const scalarField& pVals,
542
543 DynamicList<point>& snappedPoints,
544 labelList& snappedPoint
545) const
546{
547 const labelList& faceOwn = mesh_.faceOwner();
548 const labelList& faceNei = mesh_.faceNeighbour();
549
550 // Determine if point is on boundary. Points on boundaries are never
551 // snapped. Coupled boundaries are handled explicitly so not marked here.
552 bitSet isBoundaryPoint(mesh_.nPoints());
553 const polyBoundaryMesh& patches = mesh_.boundaryMesh();
554
555 for (const polyPatch& pp : patches)
556 {
557 if (!pp.coupled())
558 {
559 for (const label facei : pp.range())
560 {
561 const face& f = mesh_.faces()[facei];
562
563 isBoundaryPoint.set(f);
564 }
565 }
566 }
567
568
569 const point greatPoint(GREAT, GREAT, GREAT);
570
571 pointField collapsedPoint(mesh_.nPoints(), greatPoint);
572
573
574 // Work arrays
575 DynamicList<point, 64> localTriPoints(100);
576 labelHashSet localPointCells(100);
577
578 forAll(mesh_.pointFaces(), pointi)
579 {
580 constexpr uint8_t realCut(cutType::CUT | cutType::TETCUT);
581
582 if (isBoundaryPoint.test(pointi))
583 {
584 continue;
585 }
586
587 const labelList& pFaces = mesh_.pointFaces()[pointi];
588
589 bool anyCut = false;
590
591 for (const label facei : pFaces)
592 {
593 if
594 (
595 (cellCutType_[faceOwn[facei]] & realCut) != 0
596 || (
597 mesh_.isInternalFace(facei)
598 && (cellCutType_[faceNei[facei]] & realCut) != 0
599 )
600 )
601 {
602 anyCut = true;
603 break;
604 }
605 }
606
607 if (!anyCut)
608 {
609 continue;
610 }
611
612
613 // Do a pointCells walk (by using pointFaces)
614
615 localPointCells.clear();
616 localTriPoints.clear();
617
618 for (const label facei : pFaces)
619 {
620 const label own = faceOwn[facei];
621
622 if (isTet.test(own))
623 {
624 // Since tets have no cell centre to include make sure
625 // we only generate a triangle once per point.
626 if (localPointCells.insert(own))
627 {
628 genPointTris(pVals, pointi, facei, own, localTriPoints);
629 }
630 }
631 else
632 {
633 genPointTris
634 (
635 cVals,
636 pVals,
637 pointi,
638 facei,
639 own,
640 localTriPoints
641 );
642 }
643
644 if (mesh_.isInternalFace(facei))
645 {
646 const label nei = faceNei[facei];
647
648 if (isTet.test(nei))
649 {
650 if (localPointCells.insert(nei))
651 {
652 genPointTris(pVals, pointi, facei, nei, localTriPoints);
653 }
654 }
655 else
656 {
657 genPointTris
658 (
659 cVals,
660 pVals,
661 pointi,
662 facei,
663 nei,
664 localTriPoints
665 );
666 }
667 }
668 }
669
670 if (localTriPoints.size() == 3)
671 {
672 // Single triangle. No need for any analysis. Average points.
674 points.transfer(localTriPoints);
675 collapsedPoint[pointi] = sum(points)/points.size();
676
677 //Pout<< " point:" << pointi
678 // << " replacing coord:" << mesh_.points()[pointi]
679 // << " by average:" << collapsedPoint[pointi] << endl;
680 }
681 else if (localTriPoints.size())
682 {
683 // Convert points into triSurface.
684
685 // Merge points and compact out non-valid triangles
686 labelList triMap; // merged to unmerged triangle
687 labelList triPointReverseMap; // unmerged to merged point
688 triSurface surf
689 (
690 stitchTriPoints
691 (
692 false, // do not check for duplicate tris
693 localTriPoints,
694 triPointReverseMap,
695 triMap
696 )
697 );
698
699 labelList faceZone;
700 label nZones = surf.markZones
701 (
702 boolList(surf.nEdges(), false),
703 faceZone
704 );
705
706 if (nZones == 1)
707 {
708 // Check that all normals make a decent angle
709 scalar minCos = GREAT;
710 const vector& n0 = surf.faceNormals()[0];
711 for (label i = 1; i < surf.size(); ++i)
712 {
713 const vector& n = surf.faceNormals()[i];
714 scalar cosAngle = (n0 & n);
715 if (cosAngle < minCos)
716 {
717 minCos = cosAngle;
718 }
719 }
720 if (minCos > 0)
721 {
722 collapsedPoint[pointi] = calcCentre(surf);
723 }
724 }
725 }
726 }
727
729 (
730 mesh_,
731 collapsedPoint,
732 minMagSqrEqOp<point>(),
733 greatPoint
734 );
735
736 snappedPoint.setSize(mesh_.nPoints());
737 snappedPoint = -1;
738
739 forAll(collapsedPoint, pointi)
740 {
741 // Cannot do == comparison since might be transformed so have
742 // truncation errors.
743 if (magSqr(collapsedPoint[pointi]) < 0.5*magSqr(greatPoint))
744 {
745 snappedPoint[pointi] = snappedPoints.size();
746 snappedPoints.append(collapsedPoint[pointi]);
747 }
748 }
749}
750
751
752Foam::triSurface Foam::isoSurfaceCell::stitchTriPoints
753(
754 const bool checkDuplicates,
755 const List<point>& triPoints,
756 labelList& triPointReverseMap, // unmerged to merged point
757 labelList& triMap // merged to unmerged triangle
758) const
759{
760 label nTris = triPoints.size()/3;
761
762 if ((triPoints.size() % 3) != 0)
763 {
765 << "Problem: number of points " << triPoints.size()
766 << " not a multiple of 3." << abort(FatalError);
767 }
768
769 pointField newPoints;
771 (
772 triPoints,
773 mergeDistance_,
774 false,
775 triPointReverseMap,
776 newPoints
777 );
778
779 // Check that enough merged.
780 if (debug)
781 {
782 Pout<< "isoSurfaceCell : merged from " << triPointReverseMap.size()
783 << " points down to " << newPoints.size() << endl;
784
785 labelList oldToNew;
786 const label nUnique = Foam::mergePoints
787 (
788 newPoints,
789 mergeDistance_,
790 true,
791 oldToNew
792 );
793
794 if (nUnique != newPoints.size())
795 {
797 << "Merged points contain duplicates"
798 << " when merging with distance " << mergeDistance_ << endl
799 << "merged:" << newPoints.size() << " re-merged:"
800 << nUnique
801 << abort(FatalError);
802 }
803 }
804
805
806 List<labelledTri> tris;
807 {
808 DynamicList<labelledTri> dynTris(nTris);
809 label rawPointi = 0;
810 DynamicList<label> newToOldTri(nTris);
811
812 for (label oldTriI = 0; oldTriI < nTris; ++oldTriI)
813 {
814 labelledTri tri
815 (
816 triPointReverseMap[rawPointi],
817 triPointReverseMap[rawPointi+1],
818 triPointReverseMap[rawPointi+2],
819 0
820 );
821 if (tri.valid())
822 {
823 newToOldTri.append(oldTriI);
824 dynTris.append(tri);
825 }
826
827 rawPointi += 3;
828 }
829
830 triMap.transfer(newToOldTri);
831 tris.transfer(dynTris);
832 }
833
834
835 // Use face centres to determine 'flat hole' situation (see RMT paper).
836 // Two unconnected triangles get connected because (some of) the edges
837 // separating them get collapsed. Below only checks for duplicate triangles,
838 // not non-manifold edge connectivity.
839 if (checkDuplicates)
840 {
841 if (debug)
842 {
843 Pout<< "isoSurfaceCell : merged from " << nTris
844 << " down to " << tris.size() << " triangles." << endl;
845 }
846
847 pointField centres(tris.size());
848 forAll(tris, triI)
849 {
850 centres[triI] = tris[triI].centre(newPoints);
851 }
852
853 labelList oldToMerged;
854 label nUnique = Foam::mergePoints
855 (
856 centres,
857 mergeDistance_,
858 false,
859 oldToMerged
860 );
861
862 if (debug)
863 {
864 Pout<< "isoSurfaceCell : detected "
865 << (oldToMerged.size() - nUnique)
866 << " duplicate triangles." << endl;
867 }
868
869 if (oldToMerged.size() != nUnique)
870 {
871 // Filter out duplicates.
872 label newTriI = 0;
873 DynamicList<label> newToOldTri(tris.size());
874 labelList newToMaster(nUnique, -1);
875 forAll(tris, triI)
876 {
877 label mergedI = oldToMerged[triI];
878
879 if (newToMaster[mergedI] == -1)
880 {
881 newToMaster[mergedI] = triI;
882 newToOldTri.append(triMap[triI]);
883 tris[newTriI++] = tris[triI];
884 }
885 }
886
887 triMap.transfer(newToOldTri);
888 tris.setSize(newTriI);
889 }
890 }
891
892 return triSurface(tris, geometricSurfacePatchList(), newPoints, true);
893}
894
895
896void Foam::isoSurfaceCell::calcAddressing
897(
898 const triSurface& surf,
899 List<FixedList<label, 3>>& faceEdges,
900 labelList& edgeFace0,
901 labelList& edgeFace1,
902 Map<labelList>& edgeFacesRest
903) const
904{
905 const pointField& points = surf.points();
906
907 pointField edgeCentres(3*surf.size());
908 label edgeI = 0;
909 forAll(surf, triI)
910 {
911 const labelledTri& tri = surf[triI];
912 edgeCentres[edgeI++] = 0.5*(points[tri[0]]+points[tri[1]]);
913 edgeCentres[edgeI++] = 0.5*(points[tri[1]]+points[tri[2]]);
914 edgeCentres[edgeI++] = 0.5*(points[tri[2]]+points[tri[0]]);
915 }
916
917 labelList oldToMerged;
918 const label nUnique = Foam::mergePoints
919 (
920 edgeCentres,
921 mergeDistance_,
922 false,
923 oldToMerged
924 );
925
926 if (debug)
927 {
928 Pout<< "isoSurfaceCell : detected "
929 << nUnique
930 << " edges on " << surf.size() << " triangles." << endl;
931 }
932
933 if (nUnique == edgeCentres.size())
934 {
935 // Nothing to do
936 return;
937 }
938
939
940 // Determine faceEdges
941 faceEdges.setSize(surf.size());
942 edgeI = 0;
943 forAll(surf, triI)
944 {
945 faceEdges[triI][0] = oldToMerged[edgeI++];
946 faceEdges[triI][1] = oldToMerged[edgeI++];
947 faceEdges[triI][2] = oldToMerged[edgeI++];
948 }
949
950
951 // Determine edgeFaces
952 edgeFace0.resize(nUnique);
953 edgeFace0 = -1;
954 edgeFace1.resize(nUnique);
955 edgeFace1 = -1;
956 edgeFacesRest.clear();
957
958 forAll(oldToMerged, oldEdgeI)
959 {
960 label triI = oldEdgeI / 3;
961 label edgeI = oldToMerged[oldEdgeI];
962
963 if (edgeFace0[edgeI] == -1)
964 {
965 edgeFace0[edgeI] = triI;
966 }
967 else if (edgeFace1[edgeI] == -1)
968 {
969 edgeFace1[edgeI] = triI;
970 }
971 else
972 {
973 //WarningInFunction
974 // << "Edge " << edgeI << " with centre "
975 // << " used by more than two triangles: " << edgeFace0[edgeI]
976 // << ", "
977 // << edgeFace1[edgeI] << " and " << triI << endl;
978 Map<labelList>::iterator iter = edgeFacesRest.find(edgeI);
979
980 if (iter != edgeFacesRest.end())
981 {
982 labelList& eFaces = iter();
983 label sz = eFaces.size();
984 eFaces.setSize(sz+1);
985 eFaces[sz] = triI;
986 }
987 else
988 {
989 edgeFacesRest.insert(edgeI, labelList(1, triI));
990 }
991 }
992 }
993}
994
995
996bool Foam::isoSurfaceCell::danglingTriangle
997(
998 const FixedList<label, 3>& fEdges,
999 const labelList& edgeFace1
1000)
1001{
1002 label nOpen = 0;
1003 for (const label edgei : fEdges)
1004 {
1005 if (edgeFace1[edgei] == -1)
1006 {
1007 ++nOpen;
1008 }
1009 }
1010
1011 return (nOpen == 1 || nOpen == 2 || nOpen == 3);
1012}
1013
1014
1015Foam::label Foam::isoSurfaceCell::markDanglingTriangles
1016(
1017 const List<FixedList<label, 3>>& faceEdges,
1018 const labelList& edgeFace0,
1019 const labelList& edgeFace1,
1020 const Map<labelList>& edgeFacesRest,
1021 boolList& keepTriangles
1022)
1023{
1024 keepTriangles.setSize(faceEdges.size());
1025 keepTriangles = true;
1026
1027 label nDangling = 0;
1028
1029 // Remove any dangling triangles
1030 forAllConstIters(edgeFacesRest, iter)
1031 {
1032 // These are all the non-manifold edges. Filter out all triangles
1033 // with only one connected edge (= this edge)
1034
1035 const label edgeI = iter.key();
1036 const labelList& otherEdgeFaces = iter.val();
1037
1038 // Remove all dangling triangles
1039 if (danglingTriangle(faceEdges[edgeFace0[edgeI]], edgeFace1))
1040 {
1041 keepTriangles[edgeFace0[edgeI]] = false;
1042 ++nDangling;
1043 }
1044 if (danglingTriangle(faceEdges[edgeFace1[edgeI]], edgeFace1))
1045 {
1046 keepTriangles[edgeFace1[edgeI]] = false;
1047 ++nDangling;
1048 }
1049 for (const label triI : otherEdgeFaces)
1050 {
1051 if (danglingTriangle(faceEdges[triI], edgeFace1))
1052 {
1053 keepTriangles[triI] = false;
1054 ++nDangling;
1055 }
1056 }
1057 }
1058 return nDangling;
1059}
1060
1061
1062Foam::triSurface Foam::isoSurfaceCell::subsetMesh
1063(
1064 const triSurface& s,
1065 const labelList& newToOldFaces,
1066 labelList& oldToNewPoints,
1067 labelList& newToOldPoints
1068)
1069{
1070 const boolList include
1071 (
1072 ListOps::createWithValue<bool>(s.size(), newToOldFaces, true, false)
1073 );
1074
1075 newToOldPoints.setSize(s.points().size());
1076 oldToNewPoints.setSize(s.points().size());
1077 oldToNewPoints = -1;
1078 {
1079 label pointi = 0;
1080
1081 forAll(include, oldFacei)
1082 {
1083 if (include[oldFacei])
1084 {
1085 // Renumber labels for face
1086 for (const label oldPointi : s[oldFacei])
1087 {
1088 if (oldToNewPoints[oldPointi] == -1)
1089 {
1090 oldToNewPoints[oldPointi] = pointi;
1091 newToOldPoints[pointi++] = oldPointi;
1092 }
1093 }
1094 }
1095 }
1096 newToOldPoints.setSize(pointi);
1097 }
1098
1099 // Extract points
1100 pointField newPoints(newToOldPoints.size());
1101 forAll(newToOldPoints, i)
1102 {
1103 newPoints[i] = s.points()[newToOldPoints[i]];
1104 }
1105 // Extract faces
1106 List<labelledTri> newTriangles(newToOldFaces.size());
1107
1108 forAll(newToOldFaces, i)
1109 {
1110 // Get old vertex labels
1111 const labelledTri& tri = s[newToOldFaces[i]];
1112
1113 newTriangles[i][0] = oldToNewPoints[tri[0]];
1114 newTriangles[i][1] = oldToNewPoints[tri[1]];
1115 newTriangles[i][2] = oldToNewPoints[tri[2]];
1116 newTriangles[i].region() = tri.region();
1117 }
1118
1119 // Reuse storage.
1120 return triSurface(newTriangles, s.patches(), newPoints, true);
1121}
1122
1123
1124// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1125
1127(
1128 const polyMesh& mesh,
1129 const scalarField& cellValues,
1130 const scalarField& pointValues,
1131 const scalar iso,
1132 const isoSurfaceParams& params,
1133 const bitSet& ignoreCells
1134)
1135:
1136 isoSurfaceBase(mesh, cellValues, pointValues, iso, params),
1137 mergeDistance_(params.mergeTol()*mesh.bounds().mag()),
1138 cellCutType_(mesh.nCells(), cutType::UNVISITED)
1139{
1140 const bool regularise = (params.filter() != filterType::NONE);
1141
1142 if (debug)
1143 {
1144 Pout<< "isoSurfaceCell:" << nl
1145 << " cell min/max : " << minMax(cVals_) << nl
1146 << " point min/max : " << minMax(pVals_) << nl
1147 << " isoValue : " << iso << nl
1148 << " filter : " << Switch(regularise) << nl
1149 << " mergeTol : " << params.mergeTol() << nl
1150 << " mesh span : " << mesh.bounds().mag() << nl
1151 << " mergeDistance : " << mergeDistance_ << nl
1152 << " ignoreCells : " << ignoreCells.count()
1153 << " / " << cVals_.size() << nl
1154 << endl;
1155 }
1156
1157
1158 label nBlockedCells = 0;
1159
1160 // Mark ignoreCells as blocked
1161 nBlockedCells += blockCells(cellCutType_, ignoreCells);
1162
1163
1164 // Some processor domains may require tetBasePtIs and others do not.
1165 // Do now to ensure things stay synchronized.
1166 (void)mesh_.tetBasePtIs();
1167
1168
1169 // Calculate a tet/non-tet filter
1170 bitSet isTet(mesh_.nCells());
1171 {
1172 for (label celli = 0; celli < mesh_.nCells(); ++celli)
1173 {
1174 if (tetMatcher::test(mesh_, celli))
1175 {
1176 isTet.set(celli);
1177 }
1178 }
1179 }
1180
1181 // Determine cell cuts
1182 nCutCells_ = calcCellCuts(cellCutType_);
1183
1184 if (debug)
1185 {
1186 Pout<< "isoSurfaceCell : candidate cells cut "
1187 << nCutCells_
1188 << " blocked " << nBlockedCells
1189 << " total " << mesh_.nCells() << endl;
1190 }
1191
1192 if (debug && isA<fvMesh>(mesh))
1193 {
1194 const auto& fvmesh = dynamicCast<const fvMesh>(mesh);
1195
1196 volScalarField debugField
1197 (
1198 IOobject
1199 (
1200 "isoSurfaceCell.cutType",
1201 fvmesh.time().timeName(),
1202 fvmesh.time(),
1205 false
1206 ),
1207 fvmesh,
1209 );
1210
1211 auto& debugFld = debugField.primitiveFieldRef();
1212
1213 forAll(cellCutType_, celli)
1214 {
1215 debugFld[celli] = cellCutType_[celli];
1216 }
1217
1218 Pout<< "Writing cut types:"
1219 << debugField.objectPath() << endl;
1220
1221 debugField.write();
1222 }
1223
1224
1225 DynamicList<point> snappedPoints(nCutCells_);
1226
1227 // Per cc -1 or a point inside snappedPoints.
1228 labelList snappedCc;
1229 if (regularise)
1230 {
1231 calcSnappedCc
1232 (
1233 isTet,
1234 cellValues,
1236 snappedPoints,
1237 snappedCc
1238 );
1239 }
1240 else
1241 {
1242 snappedCc.setSize(mesh_.nCells());
1243 snappedCc = -1;
1244 }
1245
1246 if (debug)
1247 {
1248 Pout<< "isoSurfaceCell : shifted " << snappedPoints.size()
1249 << " cell centres to intersection." << endl;
1250 }
1251
1252 snappedPoints.shrink();
1253 label nCellSnaps = snappedPoints.size();
1254
1255 // Per point -1 or a point inside snappedPoints.
1256 labelList snappedPoint;
1257 if (regularise)
1258 {
1259 calcSnappedPoint
1260 (
1261 isTet,
1262 cellValues,
1264 snappedPoints,
1265 snappedPoint
1266 );
1267 }
1268 else
1269 {
1270 snappedPoint.setSize(mesh_.nPoints());
1271 snappedPoint = -1;
1272 }
1273
1274 if (debug)
1275 {
1276 Pout<< "isoSurfaceCell : shifted " << snappedPoints.size()-nCellSnaps
1277 << " vertices to intersection." << endl;
1278 }
1279
1280
1281 // Use a triSurface as a temporary for various operations
1282 triSurface tmpsurf;
1283
1284 {
1285 DynamicList<point> triPoints(nCutCells_);
1286 DynamicList<label> triMeshCells(nCutCells_);
1287
1288 generateTriPoints
1289 (
1290 cellValues,
1292
1294 mesh_.points(),
1295
1296 snappedPoints,
1297 snappedCc,
1298 snappedPoint,
1299
1300 triPoints,
1301 triMeshCells
1302 );
1303
1304 if (debug)
1305 {
1306 Pout<< "isoSurfaceCell : generated " << triMeshCells.size()
1307 << " unmerged triangles." << endl;
1308 }
1309
1310
1311 label nOldPoints = triPoints.size();
1312
1313 // Trimmed to original triangle
1314 DynamicList<label> trimTriMap;
1315 // Trimmed to original point
1316 labelList trimTriPointMap;
1317 if (getClipBounds().valid())
1318 {
1319 isoSurfacePoint::trimToBox
1320 (
1322 triPoints, // new points
1323 trimTriMap, // map from (new) triangle to original
1324 trimTriPointMap, // map from (new) point to original
1325 interpolatedPoints_, // labels of newly introduced points
1326 interpolatedOldPoints_, // and their interpolation
1327 interpolationWeights_
1328 );
1329 triMeshCells = labelField(triMeshCells, trimTriMap);
1330 }
1331
1332
1333 // Merge points and compact out non-valid triangles
1334 labelList triMap; // merged to unmerged triangle
1335 tmpsurf = stitchTriPoints
1336 (
1337 regularise, // check for duplicate tris
1338 triPoints,
1339 triPointMergeMap_, // unmerged to merged point
1340 triMap // merged to unmerged triangle
1341 );
1342
1343 if (debug)
1344 {
1345 Pout<< "isoSurfaceCell : generated " << triMap.size()
1346 << " merged triangles." << endl;
1347 }
1348
1349 if (getClipBounds().valid())
1350 {
1351 // Adjust interpolatedPoints_
1352 inplaceRenumber(triPointMergeMap_, interpolatedPoints_);
1353
1354 // Adjust triPointMergeMap_
1355 labelList newTriPointMergeMap(nOldPoints, -1);
1356 forAll(trimTriPointMap, trimPointI)
1357 {
1358 label oldPointI = trimTriPointMap[trimPointI];
1359 if (oldPointI >= 0)
1360 {
1361 label pointI = triPointMergeMap_[trimPointI];
1362 if (pointI >= 0)
1363 {
1364 newTriPointMergeMap[oldPointI] = pointI;
1365 }
1366 }
1367 }
1368 triPointMergeMap_.transfer(newTriPointMergeMap);
1369 }
1370
1371 meshCells_.setSize(triMap.size());
1372 forAll(triMap, i)
1373 {
1374 meshCells_[i] = triMeshCells[triMap[i]];
1375 }
1376 }
1377
1378
1379 if (debug)
1380 {
1381 Pout<< "isoSurfaceCell : checking " << tmpsurf.size()
1382 << " triangles for validity." << endl;
1383
1384 forAll(tmpsurf, triI)
1385 {
1386 triSurfaceTools::validTri(tmpsurf, triI);
1387 }
1388 }
1389
1390
1391 if (regularise)
1392 {
1394 labelList edgeFace0, edgeFace1;
1395 Map<labelList> edgeFacesRest;
1396
1397
1398 while (true)
1399 {
1400 // Calculate addressing
1401 calcAddressing
1402 (
1403 tmpsurf,
1404 faceEdges,
1405 edgeFace0,
1406 edgeFace1,
1407 edgeFacesRest
1408 );
1409
1410 // See if any dangling triangles
1411 boolList keepTriangles;
1412 label nDangling = markDanglingTriangles
1413 (
1414 faceEdges,
1415 edgeFace0,
1416 edgeFace1,
1417 edgeFacesRest,
1418 keepTriangles
1419 );
1420
1421 if (debug)
1422 {
1423 Pout<< "isoSurfaceCell : detected " << nDangling
1424 << " dangling triangles." << endl;
1425 }
1426
1427 if (nDangling == 0)
1428 {
1429 break;
1430 }
1431
1432 // Create face map (new to old)
1433 labelList subsetTriMap(findIndices(keepTriangles, true));
1434
1435 labelList subsetPointMap;
1436 labelList reversePointMap;
1437 tmpsurf = subsetMesh
1438 (
1439 tmpsurf,
1440 subsetTriMap,
1441 reversePointMap,
1442 subsetPointMap
1443 );
1444 meshCells_ = labelField(meshCells_, subsetTriMap);
1445 inplaceRenumber(reversePointMap, triPointMergeMap_);
1446 }
1447 }
1448
1449
1450 // Transfer to mesh storage. Note, an iso-surface has no zones
1451 {
1452 // Recover the pointField
1453 pointField pts;
1454 tmpsurf.swapPoints(pts);
1455
1456 // Transcribe from triFace to face
1457 faceList faces;
1458 tmpsurf.triFaceFaces(faces);
1459
1460 tmpsurf.clearOut();
1461
1462 Mesh updated(std::move(pts), std::move(faces), surfZoneList());
1463
1464 this->Mesh::transfer(updated);
1465 }
1466}
1467
1468
1469// ************************************************************************* //
label n
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
DynamicList< T, SizeMin > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:434
static constexpr label size() noexcept
Return the number of elements in the FixedList.
Definition: FixedList.H:416
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
fileName objectPath() const
The complete path + object name.
Definition: IOobjectI.H:214
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:330
void transfer(List< T > &list)
Definition: List.C:447
void setSize(const label n)
Alias for resize()
Definition: List.H:218
A HashTable to objects of type <T> with a label key.
Definition: Map.H:60
typename parent_type::iterator iterator
Definition: Map.H:69
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:66
const labelListList & faceEdges() const
Return face-edge addressing.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:78
bool found(const T &val, label pos=0) const
True if the value if found in the list.
Definition: UListI.H:265
label find(const T &val, label pos=0) const
Find index of the first occurrence of the value.
Definition: UList.C:212
label rcIndex(const label i) const noexcept
Definition: UListI.H:67
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
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
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:133
Low-level components common to various iso-surface algorithms.
cutType
The type of cell/face cuts.
const scalarField & cellValues() const noexcept
The mesh cell values used for creating the iso-surface.
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.
const polyMesh & mesh() const noexcept
The mesh for which the iso-surface is associated.
const scalarField & pointValues() const noexcept
The mesh point values used for creating the iso-surface.
const polyMesh & mesh_
Reference to mesh.
label blockCells(UList< cutType > &cuts, const bitSet &ignoreCells) const
Mark ignoreCells as BLOCKED.
const scalarField & cVals_
Cell values.
A surface formed by the iso value. After "Polygonising A Scalar Field Using Tetrahedrons",...
Preferences for controlling iso-surface algorithms.
filterType filter() const noexcept
Get current filter type.
const boundBox & getClipBounds() const noexcept
Get optional clipping bounding box.
scalar mergeTol() const noexcept
Get current merge tolerance.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:462
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:906
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
const vectorField & cellCentres() const
label nPoints() const noexcept
Number of mesh points.
label nCells() const noexcept
Number of mesh cells.
virtual bool write(const bool valid=true) const
Write using setting from DB.
static void syncPointPositions(const polyMesh &mesh, List< point > &positions, const CombineOp &cop, const point &nullValue)
Synchronize locations on all mesh points.
Definition: syncTools.H:201
static bool test(const UList< face > &faces)
Definition: tetMatcher.C:89
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:89
Triangle storage. Null constructable (unfortunately triangle<point, point> is not)
Definition: triPoints.H:55
static bool validTri(const triSurface &, const label facei, const bool verbose=true)
Check single triangle for (topological) validity.
Triangulated surface description with patch information.
Definition: triSurface.H:79
void triFaceFaces(List< face > &plainFaceList) const
Create a list of faces from the triFaces.
Definition: triSurface.C:723
virtual void swapPoints(pointField &pts)
Swap points. Similar to movePoints, but returns the old points.
Definition: triSurface.C:625
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const volScalarField & p0
Definition: EEqn.H:36
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
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))
label nZones
Convenience macros for instantiating iso-surface interpolate methods.
#define defineIsoSurfaceInterpolateMethods(ThisClass)
Geometric merging of points. See below.
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:57
const dimensionSet dimless
Dimensionless.
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Inplace renumber the values (not the indices) of a list.
List< label > labelList
A List of labels.
Definition: List.H:66
PointIndexHit< point > pointIndexHit
A PointIndexHit for 3D points.
Definition: pointIndexHit.H:46
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
List< geometricSurfacePatch > geometricSurfacePatchList
A List of geometricSurfacePatch.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
vector point
Point is a vector.
Definition: point.H:43
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
label mergePoints(const PointList &points, labelList &pointToUnique, labelList &uniquePoints, const scalar mergeTol=SMALL, const bool verbose=false)
MinMax< label > minMax(const labelHashSet &set)
Find the min/max values of labelHashSet.
Definition: hashSets.C:61
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:52
List< bool > boolList
A List of bools.
Definition: List.H:64
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
error FatalError
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
List< surfZone > surfZoneList
Definition: surfZoneList.H:47
labelList findIndices(const ListType &input, typename ListType::const_reference val, label start=0)
Linear search to find all occurrences of given element.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::HEX);labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells].reset(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:235
face triFace(3)
labelList f(nPoints)
volScalarField & b
Definition: createFields.H:27
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278