Loading...
Searching...
No Matches
isoSurfacePoint.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-2017 OpenFOAM Foundation
9 Copyright (C) 2015-2025 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 "isoSurfacePoint.H"
30#include "mergePoints.H"
31#include "slicedVolFields.H"
32#include "volFields.H"
33#include "triSurfaceTools.H"
34#include "triSurface.H"
35#include "triangle.H"
36
37// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38
41
42
43// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44
45namespace Foam
46{
48}
49
51// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
53namespace Foam
55
56// Helper class for slicing triangles
57struct storeOp
58{
63 list(tris)
64 {}
65
66 void operator()(const triPoints& tri)
67 {
68 list.append(tri);
69 }
70
71 void operator()(triPoints&& tri)
72 {
73 list.append(std::move(tri));
74 }
75};
76
77
78// Is patch co-located (i.e. non-separated) coupled patch?
79static inline bool collocatedPatch(const polyPatch& pp)
80{
81 const auto* cpp = isA<coupledPolyPatch>(pp);
82 return cpp && cpp->parallel() && !cpp->separated();
83}
84
85
86// Collocated patch, with extra checks
87#if 0
88static bool isCollocatedPatch(const coupledPolyPatch& pp)
89{
91 {
92 return collocatedPatch(pp);
93 }
94
96 << "Unhandled coupledPolyPatch type " << pp.type() << nl
97 << abort(FatalError);
98
99 return false;
100}
101#endif
102
103} // End namespace Foam
104
105
106// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
107
108Foam::bitSet Foam::isoSurfacePoint::collocatedFaces
109(
110 const coupledPolyPatch& pp
111)
112{
113 // Initialise to false
114 bitSet collocated(pp.size());
115
117 {
118 if (collocatedPatch(pp))
119 {
120 // All on
121 collocated = true;
122 }
123 }
124 else
125 {
127 << "Unhandled coupledPolyPatch type " << pp.type()
128 << abort(FatalError);
129 }
130 return collocated;
131}
132
133
134void Foam::isoSurfacePoint::syncUnseparatedPoints
135(
136 pointField& pointValues,
137 const point& nullValue
138) const
139{
140 // Until syncPointList handles separated coupled patches with multiple
141 // transforms do our own synchronisation of non-separated patches only
142 const polyBoundaryMesh& patches = mesh_.boundaryMesh();
143
144 if (Pstream::parRun())
145 {
146 const labelList& procPatches = mesh_.globalData().processorPatches();
147
148 // Send
149 for (const label patchi : procPatches)
150 {
151 const polyPatch& pp = patches[patchi];
152 const auto& procPatch = refCast<const processorPolyPatch>(pp);
153
154 if (pp.nPoints() && collocatedPatch(pp))
155 {
156 const labelList& meshPts = procPatch.meshPoints();
157 const labelList& nbrPts = procPatch.neighbPoints();
158
159 pointField patchInfo(meshPts.size());
160
161 forAll(nbrPts, pointi)
162 {
163 const label nbrPointi = nbrPts[pointi];
164 patchInfo[nbrPointi] = pointValues[meshPts[pointi]];
165 }
166
167 OPstream toNbr
168 (
170 procPatch.neighbProcNo()
171 );
172 toNbr << patchInfo;
173 }
174 }
175
176 // Receive and combine.
177 for (const label patchi : procPatches)
178 {
179 const polyPatch& pp = patches[patchi];
180 const auto& procPatch = refCast<const processorPolyPatch>(pp);
181
182 if (pp.nPoints() && collocatedPatch(pp))
183 {
184 // We do not know the number of points on the other side
185 // so cannot use UIPstream::read
186
187 pointField nbrPatchInfo;
188 IPstream::recv(nbrPatchInfo, procPatch.neighbProcNo());
189
190 const labelList& meshPts = procPatch.meshPoints();
191
192 forAll(meshPts, pointi)
193 {
194 const label meshPointi = meshPts[pointi];
196 (
197 pointValues[meshPointi],
198 nbrPatchInfo[pointi]
199 );
200 }
201 }
202 }
203 }
204
205 // Do the cyclics.
206 for (const polyPatch& pp : patches)
207 {
209
210 if (cpp && cpp->owner() && collocatedPatch(*cpp))
211 {
212 // Owner does all.
213
214 const auto& cycPatch = *cpp;
215 const auto& nbrPatch = cycPatch.neighbPatch();
216
217 const edgeList& coupledPoints = cycPatch.coupledPoints();
218 const labelList& meshPts = cycPatch.meshPoints();
219 const labelList& nbrMeshPoints = nbrPatch.meshPoints();
220
221 pointField half0Values(coupledPoints.size());
222 pointField half1Values(coupledPoints.size());
223
224 forAll(coupledPoints, i)
225 {
226 const edge& e = coupledPoints[i];
227 half0Values[i] = pointValues[meshPts[e[0]]];
228 half1Values[i] = pointValues[nbrMeshPoints[e[1]]];
229 }
230
231 forAll(coupledPoints, i)
232 {
233 const edge& e = coupledPoints[i];
234 const label p0 = meshPts[e[0]];
235 const label p1 = nbrMeshPoints[e[1]];
236
237 minEqOp<point>()(pointValues[p0], half1Values[i]);
238 minEqOp<point>()(pointValues[p1], half0Values[i]);
239 }
240 }
241 }
242
243 // Synchronize multiple shared points.
244 const globalMeshData& pd = mesh_.globalData();
245
246 if (pd.nGlobalPoints() > 0)
247 {
248 // Values on shared points.
249 pointField sharedPts(pd.nGlobalPoints(), nullValue);
250
251 forAll(pd.sharedPointLabels(), i)
252 {
253 const label meshPointi = pd.sharedPointLabels()[i];
254 // Fill my entries in the shared points
255 sharedPts[pd.sharedPointAddr()[i]] = pointValues[meshPointi];
256 }
257
258 // Globally consistent
259 Pstream::listReduce(sharedPts, minOp<point>());
260
261 // Now we will all have the same information. Merge it back with
262 // my local information.
263 forAll(pd.sharedPointLabels(), i)
264 {
265 const label meshPointi = pd.sharedPointLabels()[i];
266 pointValues[meshPointi] = sharedPts[pd.sharedPointAddr()[i]];
267 }
268 }
269}
270
271
272Foam::scalar Foam::isoSurfacePoint::isoFraction
273(
274 const scalar s0,
275 const scalar s1
276) const
277{
278 const scalar d = s1-s0;
279
280 if (mag(d) > VSMALL)
281 {
282 return (iso_-s0)/d;
283 }
284
285 return -1.0;
286}
287
288
289bool Foam::isoSurfacePoint::isEdgeOfFaceCut
290(
291 const scalarField& pVals,
292 const face& f,
293 const bool ownLower,
294 const bool neiLower
295) const
296{
297 forAll(f, fp)
298 {
299 const bool fpLower = (pVals[f[fp]] < iso_);
300
301 if
302 (
303 fpLower != ownLower
304 || fpLower != neiLower
305 || fpLower != (pVals[f[f.fcIndex(fp)]] < iso_)
306 )
307 {
308 return true;
309 }
310 }
311
312 return false;
313}
314
315
316void Foam::isoSurfacePoint::getNeighbour
317(
319 const volVectorField& meshC,
320 const volScalarField& cVals,
321 const label celli,
322 const label facei,
323 scalar& nbrValue,
324 point& nbrPoint
325) const
326{
327 const labelList& own = mesh_.faceOwner();
328 const labelList& nei = mesh_.faceNeighbour();
329
330 if (mesh_.isInternalFace(facei))
331 {
332 const label nbr = (own[facei] == celli ? nei[facei] : own[facei]);
333 nbrValue = cVals[nbr];
334 nbrPoint = meshC[nbr];
335 }
336 else
337 {
338 const label bFacei = facei-mesh_.nInternalFaces();
339 const label patchi = boundaryRegion[bFacei];
340 const polyPatch& pp = mesh_.boundaryMesh()[patchi];
341 const label patchFacei = facei-pp.start();
342
343 nbrValue = cVals.boundaryField()[patchi][patchFacei];
344 nbrPoint = meshC.boundaryField()[patchi][patchFacei];
345 }
346}
347
348
349void Foam::isoSurfacePoint::calcCutTypes
350(
352 const volVectorField& meshC,
353 const volScalarField& cVals,
354 const scalarField& pVals
355)
356{
357 const polyBoundaryMesh& patches = mesh_.boundaryMesh();
358 const labelList& own = mesh_.faceOwner();
359 const labelList& nei = mesh_.faceNeighbour();
360
361 faceCutType_.resize(mesh_.nFaces());
362 faceCutType_ = cutType::NOTCUT;
363
364 for (label facei = 0; facei < mesh_.nInternalFaces(); ++facei)
365 {
366 const scalar ownValue = cVals[own[facei]];
367 const bool ownLower = (ownValue < iso_);
368
369 scalar nbrValue;
370 point nbrPoint;
371 getNeighbour
372 (
374 meshC,
375 cVals,
376 own[facei],
377 facei,
378 nbrValue,
379 nbrPoint
380 );
381
382 const bool neiLower = (nbrValue < iso_);
383
384 if (ownLower != neiLower)
385 {
386 faceCutType_[facei] = cutType::CUT;
387 }
388 else
389 {
390 // Is mesh edge cut?
391 // - by looping over all the edges of the face.
392 const face& f = mesh_.faces()[facei];
393
394 if (isEdgeOfFaceCut(pVals, f, ownLower, neiLower))
395 {
396 faceCutType_[facei] = cutType::CUT;
397 }
398 }
399 }
400
401 for (const polyPatch& pp : patches)
402 {
403 for (const label facei : pp.range())
404 {
405 const scalar ownValue = cVals[own[facei]];
406 const bool ownLower = (ownValue < iso_);
407
408 scalar nbrValue;
409 point nbrPoint;
410 getNeighbour
411 (
413 meshC,
414 cVals,
415 own[facei],
416 facei,
417 nbrValue,
418 nbrPoint
419 );
420
421 const bool neiLower = (nbrValue < iso_);
422
423 if (ownLower != neiLower)
424 {
425 faceCutType_[facei] = cutType::CUT;
426 }
427 else
428 {
429 // Is mesh edge cut?
430 // - by looping over all the edges of the face.
431 const face& f = mesh_.faces()[facei];
432
433 if (isEdgeOfFaceCut(pVals, f, ownLower, neiLower))
434 {
435 faceCutType_[facei] = cutType::CUT;
436 }
437 }
438 }
439 }
440
441 nCutCells_ = 0;
442 cellCutType_.resize(mesh_.nCells());
443 cellCutType_ = cutType::NOTCUT;
444
445
446 // Propagate internal face cuts into the cells.
447
448 for (label facei = 0; facei < mesh_.nInternalFaces(); ++facei)
449 {
450 if (faceCutType_[facei] == cutType::NOTCUT)
451 {
452 continue;
453 }
454
455 if (cellCutType_[own[facei]] == cutType::NOTCUT)
456 {
457 cellCutType_[own[facei]] = cutType::CUT;
458 ++nCutCells_;
459 }
460 if (cellCutType_[nei[facei]] == cutType::NOTCUT)
461 {
462 cellCutType_[nei[facei]] = cutType::CUT;
463 ++nCutCells_;
464 }
465 }
466
467
468 // Propagate boundary face cuts into the cells.
469
470 for (label facei = mesh_.nInternalFaces(); facei < mesh_.nFaces(); ++facei)
471 {
472 if (faceCutType_[facei] == cutType::NOTCUT)
473 {
474 continue;
475 }
476
477 if (cellCutType_[own[facei]] == cutType::NOTCUT)
478 {
479 cellCutType_[own[facei]] = cutType::CUT;
480 ++nCutCells_;
481 }
482 }
483
484 if (debug)
485 {
486 Pout<< "isoSurfacePoint : candidate cut cells "
487 << nCutCells_ << " / " << mesh_.nCells() << endl;
488 }
489}
490
491
492Foam::point Foam::isoSurfacePoint::calcCentre(const triSurface& s)
493{
494 // Calculate centre of surface.
495
496 vector sum = Zero;
497
498 forAll(s, i)
499 {
500 sum += s[i].centre(s.points());
501 }
502 return sum/s.size();
503}
504
505
506void Foam::isoSurfacePoint::calcSnappedCc
507(
509 const volVectorField& meshC,
510 const volScalarField& cVals,
511 const scalarField& pVals,
512
513 DynamicList<point>& snappedPoints,
514 labelList& snappedCc
515) const
516{
517 const pointField& pts = mesh_.points();
518 const pointField& cc = mesh_.cellCentres();
519
520 snappedCc.setSize(mesh_.nCells());
521 snappedCc = -1;
522
523 // Work arrays
524 DynamicList<point, 64> localTriPoints(64);
525
526 forAll(mesh_.cells(), celli)
527 {
528 if (cellCutType_[celli] == cutType::CUT)
529 {
530 const scalar cVal = cVals[celli];
531
532 localTriPoints.clear();
533 label nOther = 0;
534 point otherPointSum = Zero;
535
536 // Create points for all intersections close to cell centre
537 // (i.e. from pyramid edges)
538
539 for (const label facei : mesh_.cells()[celli])
540 {
541 const face& f = mesh_.faces()[facei];
542
543 scalar nbrValue;
544 point nbrPoint;
545 getNeighbour
546 (
548 meshC,
549 cVals,
550 celli,
551 facei,
552 nbrValue,
553 nbrPoint
554 );
555
556 // Calculate intersection points of edges to cell centre
559
560 // From cc to neighbour cc.
561 s[2] = isoFraction(cVal, nbrValue);
562 pt[2] = (1.0-s[2])*cc[celli] + s[2]*nbrPoint;
563
564 forAll(f, fp)
565 {
566 // From cc to fp
567 label p0 = f[fp];
568 s[0] = isoFraction(cVal, pVals[p0]);
569 pt[0] = (1.0-s[0])*cc[celli] + s[0]*pts[p0];
570
571 // From cc to fp+1
572 label p1 = f[f.fcIndex(fp)];
573 s[1] = isoFraction(cVal, pVals[p1]);
574 pt[1] = (1.0-s[1])*cc[celli] + s[1]*pts[p1];
575
576 if
577 (
578 (s[0] >= 0.0 && s[0] <= 0.5)
579 && (s[1] >= 0.0 && s[1] <= 0.5)
580 && (s[2] >= 0.0 && s[2] <= 0.5)
581 )
582 {
583 localTriPoints.append(pt[0]);
584 localTriPoints.append(pt[1]);
585 localTriPoints.append(pt[2]);
586 }
587 else
588 {
589 // Get average of all other points
590 forAll(s, i)
591 {
592 if (s[i] >= 0.0 && s[i] <= 0.5)
593 {
594 otherPointSum += pt[i];
595 ++nOther;
596 }
597 }
598 }
599 }
600 }
601
602 if (localTriPoints.size() == 0)
603 {
604 // No complete triangles. Get average of all intersection
605 // points.
606 if (nOther > 0)
607 {
608 snappedCc[celli] = snappedPoints.size();
609 snappedPoints.append(otherPointSum/nOther);
610
611 //Pout<< " point:" << pointi
612 // << " replacing coord:" << mesh_.points()[pointi]
613 // << " by average:" << collapsedPoint[pointi] << endl;
614 }
615 }
616 else if (localTriPoints.size() == 3)
617 {
618 // Single triangle. No need for any analysis. Average points.
619 snappedCc[celli] = snappedPoints.size();
620 snappedPoints.append(sum(localTriPoints)/3);
621 localTriPoints.clear();
622
623 //Pout<< " point:" << pointi
624 // << " replacing coord:" << mesh_.points()[pointi]
625 // << " by average:" << collapsedPoint[pointi] << endl;
626 }
627 else
628 {
629 // Convert points into triSurface.
630
631 // Merge points and compact out non-valid triangles
632 labelList triMap; // merged to unmerged triangle
633 labelList triPointReverseMap; // unmerged to merged point
634 triSurface surf
635 (
636 stitchTriPoints
637 (
638 false, // do not check for duplicate tris
639 localTriPoints,
640 triPointReverseMap,
641 triMap
642 )
643 );
644
646 label nZones = surf.markZones
647 (
648 boolList(surf.nEdges(), false),
650 );
651
652 if (nZones == 1)
653 {
654 snappedCc[celli] = snappedPoints.size();
655 snappedPoints.append(calcCentre(surf));
656 //Pout<< " point:" << pointi << " nZones:" << nZones
657 // << " replacing coord:" << mesh_.points()[pointi]
658 // << " by average:" << collapsedPoint[pointi] << endl;
659 }
660 }
661 }
662 }
663}
664
665
666void Foam::isoSurfacePoint::calcSnappedPoint
667(
668 const bitSet& isBoundaryPoint,
670 const volVectorField& meshC,
671 const volScalarField& cVals,
672 const scalarField& pVals,
673
674 DynamicList<point>& snappedPoints,
675 labelList& snappedPoint
676) const
677{
678 const pointField& pts = mesh_.points();
679 const pointField& cc = mesh_.cellCentres();
680
681 pointField collapsedPoint(mesh_.nPoints(), point::max);
682
683 // Work arrays
684 DynamicList<point, 64> localTriPoints(100);
685
686 forAll(mesh_.pointFaces(), pointi)
687 {
688 if (isBoundaryPoint.test(pointi))
689 {
690 continue;
691 }
692
693 const labelList& pFaces = mesh_.pointFaces()[pointi];
694
695 bool anyCut = false;
696
697 for (const label facei : pFaces)
698 {
699 if (faceCutType_[facei] == cutType::CUT)
700 {
701 anyCut = true;
702 break;
703 }
704 }
705
706 if (!anyCut)
707 {
708 continue;
709 }
710
711
712 localTriPoints.clear();
713 label nOther = 0;
714 point otherPointSum = Zero;
715
716 for (const label facei : pFaces)
717 {
718 // Create points for all intersections close to point
719 // (i.e. from pyramid edges)
720 const face& f = mesh_.faces()[facei];
721 const label own = mesh_.faceOwner()[facei];
722
723 // Get neighbour value
724 scalar nbrValue;
725 point nbrPoint;
726 getNeighbour
727 (
729 meshC,
730 cVals,
731 own,
732 facei,
733 nbrValue,
734 nbrPoint
735 );
736
737 // Calculate intersection points of edges emanating from point
740
741 label fp = f.find(pointi);
742 s[0] = isoFraction(pVals[pointi], cVals[own]);
743 pt[0] = (1.0-s[0])*pts[pointi] + s[0]*cc[own];
744
745 s[1] = isoFraction(pVals[pointi], nbrValue);
746 pt[1] = (1.0-s[1])*pts[pointi] + s[1]*nbrPoint;
747
748 label nextPointi = f[f.fcIndex(fp)];
749 s[2] = isoFraction(pVals[pointi], pVals[nextPointi]);
750 pt[2] = (1.0-s[2])*pts[pointi] + s[2]*pts[nextPointi];
751
752 label prevPointi = f[f.rcIndex(fp)];
753 s[3] = isoFraction(pVals[pointi], pVals[prevPointi]);
754 pt[3] = (1.0-s[3])*pts[pointi] + s[3]*pts[prevPointi];
755
756 if
757 (
758 (s[0] >= 0.0 && s[0] <= 0.5)
759 && (s[1] >= 0.0 && s[1] <= 0.5)
760 && (s[2] >= 0.0 && s[2] <= 0.5)
761 )
762 {
763 localTriPoints.append(pt[0]);
764 localTriPoints.append(pt[1]);
765 localTriPoints.append(pt[2]);
766 }
767 if
768 (
769 (s[0] >= 0.0 && s[0] <= 0.5)
770 && (s[1] >= 0.0 && s[1] <= 0.5)
771 && (s[3] >= 0.0 && s[3] <= 0.5)
772 )
773 {
774 localTriPoints.append(pt[3]);
775 localTriPoints.append(pt[0]);
776 localTriPoints.append(pt[1]);
777 }
778
779 // Get average of points as fallback
780 forAll(s, i)
781 {
782 if (s[i] >= 0.0 && s[i] <= 0.5)
783 {
784 otherPointSum += pt[i];
785 ++nOther;
786 }
787 }
788 }
789
790 if (localTriPoints.size() == 0)
791 {
792 // No complete triangles. Get average of all intersection
793 // points.
794 if (nOther > 0)
795 {
796 collapsedPoint[pointi] = otherPointSum/nOther;
797 }
798 }
799 else if (localTriPoints.size() == 3)
800 {
801 // Single triangle. No need for any analysis. Average points.
803 points.transfer(localTriPoints);
804 collapsedPoint[pointi] = sum(points)/points.size();
805 }
806 else
807 {
808 // Convert points into triSurface.
809
810 // Merge points and compact out non-valid triangles
811 labelList triMap; // merged to unmerged triangle
812 labelList triPointReverseMap; // unmerged to merged point
813 triSurface surf
814 (
815 stitchTriPoints
816 (
817 false, // do not check for duplicate tris
818 localTriPoints,
819 triPointReverseMap,
820 triMap
821 )
822 );
823
825 label nZones = surf.markZones
826 (
827 boolList(surf.nEdges(), false),
829 );
830
831 if (nZones == 1)
832 {
833 collapsedPoint[pointi] = calcCentre(surf);
834 }
835 }
836 }
837
838
839 // Synchronise snap location
840 syncUnseparatedPoints(collapsedPoint, point::max);
841
842
843 snappedPoint.setSize(mesh_.nPoints());
844 snappedPoint = -1;
845
846 forAll(collapsedPoint, pointi)
847 {
848 if (collapsedPoint[pointi] != point::max)
849 {
850 snappedPoint[pointi] = snappedPoints.size();
851 snappedPoints.append(collapsedPoint[pointi]);
852 }
853 }
854}
855
856
857Foam::triSurface Foam::isoSurfacePoint::stitchTriPoints
858(
859 const bool checkDuplicates,
860 const List<point>& triPoints,
861 labelList& triPointReverseMap, // unmerged to merged point
862 labelList& triMap // merged to unmerged triangle
863) const
864{
865 label nTris = triPoints.size()/3;
866
867 if ((triPoints.size() % 3) != 0)
868 {
870 << "Problem: number of points " << triPoints.size()
871 << " not a multiple of 3." << abort(FatalError);
872 }
873
874 pointField newPoints;
876 (
877 triPoints,
878 mergeDistance_,
879 false,
880 triPointReverseMap,
881 newPoints
882 );
883
884 // Check that enough merged.
885 if (debug)
886 {
887 Pout<< "isoSurfacePoint : merged from " << triPointReverseMap.size()
888 << " down to " << newPoints.size() << " unique points." << endl;
889
890 pointField newNewPoints;
891 labelList oldToNew;
892 bool hasMerged = mergePoints
893 (
894 newPoints,
895 mergeDistance_,
896 true,
897 oldToNew,
898 newNewPoints
899 );
900
901 if (hasMerged)
902 {
904 << "Merged points contain duplicates"
905 << " when merging with distance " << mergeDistance_ << endl
906 << "merged:" << oldToNew.size() << " re-merged:"
907 << newNewPoints.size()
908 << abort(FatalError);
909 }
910 }
911
912
914 {
915 DynamicList<labelledTri> dynTris(nTris);
916 label rawPointi = 0;
917 DynamicList<label> newToOldTri(nTris);
918
919 for (label oldTriI = 0; oldTriI < nTris; oldTriI++)
920 {
921 labelledTri tri
922 (
923 triPointReverseMap[rawPointi],
924 triPointReverseMap[rawPointi+1],
925 triPointReverseMap[rawPointi+2],
926 0
927 );
928 rawPointi += 3;
929
930 if (tri.good())
931 {
932 newToOldTri.append(oldTriI);
933 dynTris.append(tri);
934 }
935 }
936
937 triMap.transfer(newToOldTri);
938 tris.transfer(dynTris);
939 }
940
941
942
943 // Determine 'flat hole' situation (see RMT paper).
944 // Two unconnected triangles get connected because (some of) the edges
945 // separating them get collapsed. Below only checks for duplicate triangles,
946 // not non-manifold edge connectivity.
947 if (checkDuplicates)
948 {
949 labelListList pointFaces;
950 invertManyToMany(newPoints.size(), tris, pointFaces);
951
952 // Filter out duplicates.
953 DynamicList<label> newToOldTri(tris.size());
954
955 forAll(tris, triI)
956 {
957 const labelledTri& tri = tris[triI];
958 const labelList& pFaces = pointFaces[tri[0]];
959
960 // Find the maximum of any duplicates. Maximum since the tris
961 // below triI
962 // get overwritten so we cannot use them in a comparison.
963 label dupTriI = -1;
964 forAll(pFaces, i)
965 {
966 label nbrTriI = pFaces[i];
967
968 if (nbrTriI > triI && (tris[nbrTriI] == tri))
969 {
970 //Pout<< "Duplicate : " << triI << " verts:" << tri
971 // << " to " << nbrTriI << " verts:" << tris[nbrTriI]
972 // << endl;
973 dupTriI = nbrTriI;
974 break;
975 }
976 }
977
978 if (dupTriI == -1)
979 {
980 // There is no (higher numbered) duplicate triangle
981 label newTriI = newToOldTri.size();
982 newToOldTri.append(triMap[triI]);
983 tris[newTriI] = tris[triI];
984 }
985 }
986
987 triMap.transfer(newToOldTri);
988 tris.setSize(triMap.size());
989
990 if (debug)
991 {
992 Pout<< "isoSurfacePoint : merged from " << nTris
993 << " down to " << tris.size() << " unique triangles." << endl;
994 }
995
996
997 if (debug)
998 {
999 triSurface surf(tris, geometricSurfacePatchList(), newPoints);
1000
1001 forAll(surf, facei)
1002 {
1003 const labelledTri& f = surf[facei];
1004 const labelList& fFaces = surf.faceFaces()[facei];
1005
1006 forAll(fFaces, i)
1007 {
1008 label nbrFacei = fFaces[i];
1009
1010 if (nbrFacei <= facei)
1011 {
1012 // lower numbered faces already checked
1013 continue;
1014 }
1015
1016 const labelledTri& nbrF = surf[nbrFacei];
1017
1018 if (f == nbrF)
1019 {
1021 << "Check : "
1022 << " triangle " << facei << " vertices " << f
1023 << " fc:" << f.centre(surf.points())
1024 << " has the same vertices as triangle " << nbrFacei
1025 << " vertices " << nbrF
1026 << " fc:" << nbrF.centre(surf.points())
1027 << abort(FatalError);
1028 }
1029 }
1030 }
1031 }
1032 }
1033
1034 return triSurface(tris, geometricSurfacePatchList(), newPoints, true);
1035}
1036
1037
1038void Foam::isoSurfacePoint::trimToPlanes
1039(
1040 const PtrList<plane>& planes,
1041 const triPointRef& tri,
1042 DynamicList<point>& newTriPoints
1043)
1044{
1045 // Buffer for generated triangles
1046 DynamicList<triPoints> insideTrisA;
1047 storeOp insideOpA(insideTrisA);
1048
1049 // Buffer for generated triangles
1050 DynamicList<triPoints> insideTrisB;
1051 storeOp insideOpB(insideTrisB);
1052
1053 triPointRef::dummyOp dop;
1054
1055 // Store starting triangle in insideTrisA
1056 insideOpA(triPoints(tri.a(), tri.b(), tri.c()));
1057
1058
1059 bool useA = true;
1060
1061 forAll(planes, faceI)
1062 {
1063 const plane& pln = planes[faceI];
1064
1065 if (useA)
1066 {
1067 insideTrisB.clear();
1068 for (const triPoints& tri : insideTrisA)
1069 {
1070 triPointRef(tri).sliceWithPlane(pln, insideOpB, dop);
1071 }
1072 }
1073 else
1074 {
1075 insideTrisA.clear();
1076 for (const triPoints& tri : insideTrisB)
1077 {
1078 triPointRef(tri).sliceWithPlane(pln, insideOpA, dop);
1079 }
1080 }
1081 useA = !useA;
1082 }
1083
1084
1085 DynamicList<triPoints>& newTris = (useA ? insideTrisA : insideTrisB);
1086
1087 newTriPoints.reserve(newTriPoints.size() + 3*newTris.size());
1088
1089 // Transfer
1090 for (const triPoints& tri : newTris)
1091 {
1092 newTriPoints.append(tri[0]);
1093 newTriPoints.append(tri[1]);
1094 newTriPoints.append(tri[2]);
1095 }
1096}
1097
1098
1099void Foam::isoSurfacePoint::trimToBox
1100(
1101 const treeBoundBox& bb,
1103 DynamicList<label>& triMap // trimmed to original tris
1104)
1105{
1106 if (debug)
1107 {
1108 Pout<< "isoSurfacePoint : trimming to " << bb << endl;
1109 }
1110
1111 // Generate inwards pointing planes
1114 {
1115 const vector& n = treeBoundBox::faceNormals[faceI];
1116 planes.set(faceI, new plane(bb.faceCentre(faceI), -n));
1117 }
1118
1119 label nTris = triPoints.size()/3;
1120
1121 DynamicList<point> newTriPoints(triPoints.size()/16);
1122 triMap.setCapacity(nTris/16);
1123
1124 label vertI = 0;
1125 for (label triI = 0; triI < nTris; triI++)
1126 {
1127 const point& p0 = triPoints[vertI++];
1128 const point& p1 = triPoints[vertI++];
1129 const point& p2 = triPoints[vertI++];
1130
1131 label oldNPoints = newTriPoints.size();
1132 trimToPlanes
1133 (
1134 planes,
1135 triPointRef(p0, p1, p2),
1136 newTriPoints
1137 );
1138
1139 label nCells = (newTriPoints.size()-oldNPoints)/3;
1140 for (label i = 0; i < nCells; i++)
1141 {
1142 triMap.append(triI);
1143 }
1144 }
1145
1146 if (debug)
1147 {
1148 Pout<< "isoSurfacePoint : trimmed from " << nTris
1149 << " down to " << triMap.size()
1150 << " triangles." << endl;
1151 }
1152
1153 triPoints.transfer(newTriPoints);
1154}
1155
1156
1157void Foam::isoSurfacePoint::trimToBox
1158(
1159 const treeBoundBox& bb,
1160 DynamicList<point>& triPoints, // new points
1161 DynamicList<label>& triMap, // map from (new) triangle to original
1162 labelList& triPointMap, // map from (new) point to original
1163 labelList& interpolatedPoints, // labels of newly introduced points
1164 List<FixedList<label, 3>>& interpolatedOldPoints,// and their interpolation
1166)
1167{
1168 const List<point> oldTriPoints(triPoints);
1169
1170 // Trim triPoints, return map
1171 trimToBox(bb, triPoints, triMap);
1172
1173
1174 // Find point correspondence:
1175 // - one-to-one for preserved original points (triPointMap)
1176 // - interpolation for newly introduced points
1177 // (interpolatedOldPoints)
1178 label sz = oldTriPoints.size()/100;
1179 DynamicList<label> dynInterpolatedPoints(sz);
1180 DynamicList<FixedList<label, 3>> dynInterpolatedOldPoints(sz);
1181 DynamicList<FixedList<scalar, 3>> dynInterpolationWeights(sz);
1182
1183
1184 triPointMap.setSize(triPoints.size());
1185 forAll(triMap, triI)
1186 {
1187 label oldTriI = triMap[triI];
1188
1189 // Find point correspondence. Assumes coordinate is bit-copy.
1190 for (label i = 0; i < 3; i++)
1191 {
1192 label pointI = 3*triI+i;
1193 const point& pt = triPoints[pointI];
1194
1195 // Compare to old-triangle's points
1196 label matchPointI = -1;
1197 for (label j = 0; j < 3; j++)
1198 {
1199 label oldPointI = 3*oldTriI+j;
1200 if (pt == oldTriPoints[oldPointI])
1201 {
1202 matchPointI = oldPointI;
1203 break;
1204 }
1205 }
1206
1207 triPointMap[pointI] = matchPointI;
1208
1209 // If new point: calculate and store interpolation
1210 if (matchPointI == -1)
1211 {
1212 dynInterpolatedPoints.append(pointI);
1213
1214 FixedList<label, 3> oldPoints
1215 (
1216 {3*oldTriI, 3*oldTriI+1, 3*oldTriI+2}
1217 );
1218 dynInterpolatedOldPoints.append(oldPoints);
1219
1220 triPointRef tri(oldTriPoints, oldPoints);
1221 barycentric2D bary = tri.pointToBarycentric(pt);
1222 FixedList<scalar, 3> weights({bary.a(), bary.b(), bary.c()});
1223
1224 dynInterpolationWeights.append(weights);
1225 }
1226 }
1227 }
1228
1229 interpolatedPoints.transfer(dynInterpolatedPoints);
1230 interpolatedOldPoints.transfer(dynInterpolatedOldPoints);
1231 interpolationWeights.transfer(dynInterpolationWeights);
1232}
1233
1234
1235Foam::triSurface Foam::isoSurfacePoint::subsetMesh
1236(
1237 const triSurface& s,
1238 const labelList& newToOldFaces,
1239 labelList& oldToNewPoints,
1240 labelList& newToOldPoints
1241)
1242{
1243 const boolList include
1244 (
1245 ListOps::createWithValue<bool>(s.size(), newToOldFaces, true, false)
1246 );
1247
1248 newToOldPoints.setSize(s.points().size());
1249 oldToNewPoints.setSize(s.points().size());
1250 oldToNewPoints = -1;
1251 {
1252 label pointi = 0;
1253
1254 forAll(include, oldFacei)
1255 {
1256 if (include[oldFacei])
1257 {
1258 // Renumber labels for triangle
1259 const labelledTri& tri = s[oldFacei];
1260
1261 forAll(tri, fp)
1262 {
1263 label oldPointi = tri[fp];
1264
1265 if (oldToNewPoints[oldPointi] == -1)
1266 {
1267 oldToNewPoints[oldPointi] = pointi;
1268 newToOldPoints[pointi++] = oldPointi;
1269 }
1270 }
1271 }
1272 }
1273 newToOldPoints.setSize(pointi);
1274 }
1275
1276 // Extract points
1277 pointField newPoints(newToOldPoints.size());
1278 forAll(newToOldPoints, i)
1279 {
1280 newPoints[i] = s.points()[newToOldPoints[i]];
1281 }
1282 // Extract faces
1283 List<labelledTri> newTriangles(newToOldFaces.size());
1284
1285 forAll(newToOldFaces, i)
1286 {
1287 // Get old vertex labels
1288 const labelledTri& tri = s[newToOldFaces[i]];
1289
1290 newTriangles[i][0] = oldToNewPoints[tri[0]];
1291 newTriangles[i][1] = oldToNewPoints[tri[1]];
1292 newTriangles[i][2] = oldToNewPoints[tri[2]];
1293 newTriangles[i].region() = tri.region();
1294 }
1295
1296 // Reuse storage.
1297 return triSurface(newTriangles, s.patches(), newPoints, true);
1298}
1299
1300
1301// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1302
1304(
1305 const volScalarField& cellValues,
1306 const scalarField& pointValues,
1307 const scalar iso,
1308 const isoSurfaceParams& params,
1309 const bitSet& /*unused*/
1310)
1311:
1313 (
1314 cellValues.mesh(),
1315 cellValues.primitiveField(),
1316 pointValues,
1317 iso,
1318 params
1319 ),
1320 cValsPtr_(nullptr),
1321 mergeDistance_(params.mergeTol()*mesh_.bounds().mag())
1322{
1323 const bool regularise = (params.filter() != filterType::NONE);
1324 const fvMesh& fvmesh = cellValues.mesh();
1325
1326 if (debug)
1327 {
1328 Pout<< "isoSurfacePoint:" << nl
1329 << " isoField : " << cellValues.name() << nl
1330 << " cell min/max : " << minMax(cVals_) << nl
1331 << " point min/max : " << minMax(pVals_) << nl
1332 << " isoValue : " << iso << nl
1333 << " filter : " << Switch(regularise) << nl
1334 << " mergeTol : " << params.mergeTol() << nl
1335 << endl;
1336 }
1337
1338 const polyBoundaryMesh& patches = mesh_.boundaryMesh();
1339
1340 // Rewrite input field
1341 // ~~~~~~~~~~~~~~~~~~~
1342 // Rewrite input volScalarField to have interpolated values
1343 // on separated patches.
1344
1345 cValsPtr_.reset(adaptPatchFields(cellValues).ptr());
1346
1347
1348 // Construct cell centres field consistent with cVals
1349 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1350 // Generate field to interpolate. This is identical to the mesh.C()
1351 // except on separated coupled patches and on empty patches.
1352
1354 (
1355 IOobject
1356 (
1357 "C",
1358 fvmesh.pointsInstance(),
1360 fvmesh,
1364 ),
1365 fvmesh,
1366 dimLength,
1367 fvmesh.cellCentres(),
1368 fvmesh.faceCentres()
1369 );
1370 forAll(patches, patchi)
1371 {
1372 const polyPatch& pp = patches[patchi];
1373
1374 // Adapt separated coupled (proc and cyclic) patches
1375 if (pp.coupled())
1376 {
1377 auto& pfld = const_cast<fvPatchVectorField&>
1378 (
1379 meshC.boundaryField()[patchi]
1380 );
1381
1382 bitSet isCollocated
1383 (
1384 collocatedFaces(refCast<const coupledPolyPatch>(pp))
1385 );
1386
1387 forAll(isCollocated, i)
1388 {
1389 if (!isCollocated[i])
1390 {
1391 pfld[i] = mesh_.faceCentres()[pp.start()+i];
1392 }
1393 }
1394 }
1395 else if (isA<emptyPolyPatch>(pp))
1396 {
1397 auto& bfld = const_cast<slicedVolVectorField::Boundary&>
1398 (
1399 meshC.boundaryField()
1400 );
1401
1402 // Clear old value. Cannot resize it since is a slice.
1403 // Set new value we can change
1404
1405 bfld.release(patchi);
1406 bfld.set
1407 (
1408 patchi,
1409 new calculatedFvPatchField<vector>
1410 (
1411 fvmesh.boundary()[patchi],
1412 meshC
1413 )
1414 );
1415
1416 // Change to face centres
1417 bfld[patchi] = pp.patchSlice(mesh_.faceCentres());
1418 }
1419 }
1420
1421
1422 // Pre-calculate patch-per-face to avoid whichPatch call.
1423 labelList boundaryRegion(mesh_.nBoundaryFaces());
1424
1425 for (const polyPatch& pp : patches)
1426 {
1427 SubList<label>(boundaryRegion, pp.size(), pp.offset()) = pp.index();
1428 }
1429
1430
1431 // Determine if any cut through face/cell
1432 calcCutTypes(boundaryRegion, meshC, cValsPtr_(), pVals_);
1433
1434 if (debug)
1435 {
1436 volScalarField debugField
1437 (
1438 IOobject
1439 (
1440 "isoSurfacePoint.cutType",
1441 fvmesh.time().timeName(),
1442 fvmesh.time(),
1446 ),
1447 fvmesh,
1449 );
1450
1451 auto& debugFld = debugField.primitiveFieldRef();
1452
1453 forAll(cellCutType_, celli)
1454 {
1455 debugFld[celli] = cellCutType_[celli];
1456 }
1457
1458 Pout<< "Writing cut types:"
1459 << debugField.objectPath() << endl;
1460
1461 debugField.write();
1462 }
1463
1464
1465 DynamicList<point> snappedPoints(nCutCells_);
1466
1467 // Per cc -1 or a point inside snappedPoints.
1468 labelList snappedCc;
1469 if (regularise)
1470 {
1471 calcSnappedCc
1472 (
1473 boundaryRegion,
1474 meshC,
1475 cValsPtr_(),
1476 pVals_,
1477
1478 snappedPoints,
1479 snappedCc
1480 );
1481 }
1482 else
1483 {
1484 snappedCc.setSize(mesh_.nCells());
1485 snappedCc = -1;
1486 }
1487
1488
1489
1490 if (debug)
1491 {
1492 Pout<< "isoSurfacePoint : shifted " << snappedPoints.size()
1493 << " cell centres to intersection." << endl;
1494 }
1495
1496 label nCellSnaps = snappedPoints.size();
1497
1498
1499 // Per point -1 or a point inside snappedPoints.
1500 labelList snappedPoint;
1501 if (regularise)
1502 {
1503 // Determine if point is on boundary.
1504 bitSet isBoundaryPoint(mesh_.nPoints());
1505
1506 for (const polyPatch& pp : patches)
1507 {
1508 // Mark all boundary points that are not physically coupled
1509 // (so anything but collocated coupled patches)
1510
1511 if (pp.coupled())
1512 {
1513 const coupledPolyPatch& cpp =
1515
1516 bitSet isCollocated(collocatedFaces(cpp));
1517
1518 forAll(isCollocated, i)
1519 {
1520 if (!isCollocated[i])
1521 {
1522 const face& f = mesh_.faces()[cpp.start()+i];
1523
1524 isBoundaryPoint.set(f);
1525 }
1526 }
1527 }
1528 else
1529 {
1530 forAll(pp, i)
1531 {
1532 const face& f = mesh_.faces()[pp.start()+i];
1533
1534 isBoundaryPoint.set(f);
1535 }
1536 }
1537 }
1538
1539 calcSnappedPoint
1540 (
1541 isBoundaryPoint,
1542 boundaryRegion,
1543 meshC,
1544 cValsPtr_(),
1545 pVals_,
1546
1547 snappedPoints,
1548 snappedPoint
1549 );
1550 }
1551 else
1552 {
1553 snappedPoint.setSize(mesh_.nPoints());
1554 snappedPoint = -1;
1555 }
1556
1557 if (debug)
1558 {
1559 Pout<< "isoSurfacePoint : shifted " << snappedPoints.size()-nCellSnaps
1560 << " vertices to intersection." << endl;
1561 }
1562
1563
1564 // Use a triSurface as a temporary for various operations
1565 triSurface tmpsurf;
1566
1567 {
1568 DynamicList<point> triPoints(3*nCutCells_);
1569 DynamicList<label> triMeshCells(nCutCells_);
1570
1571 generateTriPoints
1572 (
1573 cValsPtr_(),
1574 pVals_,
1575
1576 meshC,
1577 mesh_.points(),
1578
1579 snappedPoints,
1580 snappedCc,
1581 snappedPoint,
1582
1583 triPoints, // 3 points of the triangle
1584 triMeshCells // per triangle the originating cell
1585 );
1586
1587 if (debug)
1588 {
1589 Pout<< "isoSurfacePoint : generated " << triMeshCells.size()
1590 << " unmerged triangles from " << triPoints.size()
1591 << " unmerged points." << endl;
1592 }
1593
1594 label nOldPoints = triPoints.size();
1595
1596 // Trimmed to original triangle
1597 DynamicList<label> trimTriMap;
1598 // Trimmed to original point
1599 labelList trimTriPointMap;
1600 if (getClipBounds().good())
1601 {
1602 trimToBox
1603 (
1604 treeBoundBox(getClipBounds()),
1605 triPoints, // new points
1606 trimTriMap, // map from (new) triangle to original
1607 trimTriPointMap, // map from (new) point to original
1608 interpolatedPoints_, // labels of newly introduced points
1609 interpolatedOldPoints_, // and their interpolation
1610 interpolationWeights_
1611 );
1612 triMeshCells = labelField(triMeshCells, trimTriMap);
1613 }
1614
1615
1616 // Merge points and compact out non-valid triangles
1617 labelList triMap; // merged to unmerged triangle
1618 tmpsurf = stitchTriPoints
1619 (
1620 true, // check for duplicate tris
1621 triPoints,
1622 triPointMergeMap_, // unmerged to merged point
1623 triMap
1624 );
1625
1626 if (debug)
1627 {
1628 Pout<< "isoSurfacePoint : generated " << triMap.size()
1629 << " merged triangles." << endl;
1630 }
1631
1632
1633 if (getClipBounds().good())
1634 {
1635 // Adjust interpolatedPoints_
1636 inplaceRenumber(triPointMergeMap_, interpolatedPoints_);
1637
1638 // Adjust triPointMergeMap_
1639 labelList newTriPointMergeMap(nOldPoints, -1);
1640 forAll(trimTriPointMap, trimPointI)
1641 {
1642 label oldPointI = trimTriPointMap[trimPointI];
1643 if (oldPointI >= 0)
1644 {
1645 label pointI = triPointMergeMap_[trimPointI];
1646 if (pointI >= 0)
1647 {
1648 newTriPointMergeMap[oldPointI] = pointI;
1649 }
1650 }
1651 }
1652 triPointMergeMap_.transfer(newTriPointMergeMap);
1653 }
1654
1655 meshCells_.setSize(triMap.size());
1656 forAll(triMap, i)
1657 {
1658 meshCells_[i] = triMeshCells[triMap[i]];
1659 }
1660 }
1661
1662 if (debug)
1663 {
1664 Pout<< "isoSurfacePoint : checking " << tmpsurf.size()
1665 << " triangles for validity." << endl;
1666
1667 forAll(tmpsurf, facei)
1668 {
1669 triSurfaceTools::validTri(tmpsurf, facei);
1670 }
1671
1672 fileName stlFile = mesh_.time().path() + ".stl";
1673 Pout<< "Dumping surface to " << stlFile << endl;
1674 tmpsurf.write(stlFile);
1675 }
1676
1677
1678 // Transfer to mesh storage. Note, an iso-surface has no zones
1679 {
1680 // Recover the pointField
1682 tmpsurf.swapPoints(pts);
1683
1684 // Transcribe from triFace to face
1685 faceList faces;
1686 tmpsurf.triFaceFaces(faces);
1687
1688 tmpsurf.clearOut();
1689
1690 Mesh updated(std::move(pts), std::move(faces), surfZoneList());
1691
1692 this->Mesh::transfer(updated);
1693 }
1694}
1695
1696
1697// ************************************************************************* //
label n
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void reserve(const label len)
Reserve allocation space for at least this size, allocating new space if required and retaining old c...
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
static constexpr label size() noexcept
Return the number of elements in the FixedList.
Definition FixedList.H:619
void transfer(FixedList< T, N > &list)
Transfer by swapping using a move assignment for the content of the individual list elements.
Definition FixedListI.H:384
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
GeometricBoundaryField< vector, fvPatchField, volMesh > Boundary
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
fileName objectPath() const
The complete path + object name.
Definition IOobjectI.H:313
static void recv(Type &value, const int fromProcNo, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm, IOstreamOption::streamFormat fmt=IOstreamOption::BINARY)
Receive and deserialize a value. Uses operator>> for de-serialization.
Definition IPstream.H:80
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition List.C:347
void setSize(label n)
Alias for resize().
Definition List.H:536
void transfer(pointField &pointLst, List< face > &faceLst)
Output inter-processor communications stream.
Definition OPstream.H:53
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
label nPoints() const
Number of points supporting patch faces.
static void listReduce(UList< T > &values, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce list elements (list must be equal size on all ranks), applying bop to each list element.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
A non-owning sub-view of a List (allocated or unallocated storage).
Definition SubList.H:61
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition Switch.H:81
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
label rcIndex(const label i) const noexcept
The reverse circular index. The previous index in the list which returns to the last at the beginning...
Definition UListI.H:106
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label find(const T &val) const
Find index of the first occurrence of the value.
Definition UList.C:160
label fcIndex(const label i) const noexcept
The forward circular index. The next index in the list which returns to the first at the end of the l...
Definition UListI.H:99
@ buffered
"buffered" : (MPI_Bsend, MPI_Recv)
Definition UPstream.H:82
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static const Form max
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition bitSetI.H:502
static constexpr label nFaces() noexcept
Number of faces for boundBox and HEX.
Definition boundBox.H:168
static const FixedList< vector, 6 > faceNormals
The unit normal per face.
Definition boundBox.H:136
The boundaryRegion persistent data saved as a Map<dictionary>.
This boundary condition is not designed to be evaluated; it is assumed that the value is assigned via...
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
Cyclic plane patch.
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
A subset of mesh faces organised as a primitive patch.
Definition faceZone.H:63
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
A class for handling file names.
Definition fileName.H:75
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition fvMesh.H:395
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
Abstract base class for interpolating in 1D.
Low-level components common to various iso-surface algorithms.
isoSurfaceBase(const isoSurfaceBase &)=delete
No copy construct.
const scalarField & cellValues() const noexcept
The mesh cell values used for creating the iso-surface.
const scalarField & pVals_
Point values.
labelList meshCells_
For every face, the original cell in mesh.
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.
const scalarField & cVals_
Cell values.
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.
A surface formed by the iso value. After "Regularised Marching Tetrahedra: improved iso-surface extra...
isoSurfacePoint(const volScalarField &cellValues, const scalarField &pointValues, const scalar iso, const isoSurfaceParams &params=isoSurfaceParams(), const bitSet &ignoreCells=bitSet())
Construct from cell values and point values.
A triFace with additional (region) index.
Definition labelledTri.H:56
label region() const noexcept
Return the region index.
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition plane.H:91
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition polyMesh.C:838
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh").
Definition polyMesh.H:411
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
label start() const noexcept
Return start label of this patch in the polyMesh face list.
Definition polyPatch.H:446
const vectorField & faceCentres() const
const vectorField & cellCentres() const
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
Standard boundBox with extra functionality for use in octree.
Triangle point storage. Default constructable (triangle is not).
Definition triangle.H:77
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:74
void triFaceFaces(List< face > &plainFaceList) const
Create a list of faces from the triFaces.
Definition triSurface.C:717
void write(Ostream &os) const
Write to Ostream in simple OpenFOAM format.
virtual void swapPoints(pointField &pts)
Swap points. Similar to movePoints, but returns the old points.
Definition triSurface.C:619
void sliceWithPlane(const plane &pln, AboveOp &aboveOp, BelowOp &belowOp) const
Decompose triangle into triangles above and below plane.
Definition triangle.C:107
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const volScalarField & p0
Definition EEqn.H:36
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
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)
Geometric merging of points. See below.
List< T > createWithValue(const label len, const labelUList &locations, const T &val, const T &deflt=T())
Create a List filled with default values and various locations with another specified value.
Namespace for bounding specifications. At the moment, mostly for tables.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
void invertManyToMany(const label len, const UList< InputIntListType > &input, List< OutputIntListType > &output)
Invert many-to-many.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Inplace renumber the values within a list.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
List< surfZone > surfZoneList
List of surfZone.
Barycentric2D< scalar > barycentric2D
A scalar version of the templated Barycentric2D.
List< face > faceList
List of faces.
Definition faceListFwd.H:41
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)
Calculate merge mapping, preserving the original point order. All points closer/equal mergeTol are to...
SlicedGeometricField< vector, fvPatchField, slicedFvPatchField, volMesh > slicedVolVectorField
MinMax< label > minMax(const labelHashSet &set)
Find the min/max values of labelHashSet.
Definition hashSets.C:54
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< label > labelField
Specialisation of Field<T> for label.
Definition labelField.H:48
vector point
Point is a vector.
Definition point.H:37
static bool collocatedPatch(const polyPatch &pp)
List< geometricSurfacePatch > geometricSurfacePatchList
List of geometricSurfacePatch.
List< bool > boolList
A List of bools.
Definition List.H:60
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition triangleFwd.H:46
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
fvPatchField< vector > fvPatchVectorField
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
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];} pointMap[start]=pointMap[end]=Foam::min(start, end);} } 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};constexpr 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< DynamicList< face > > pFaces[nBCs]
labelList f(nPoints)
const pointField & pts
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
void operator()(triPoints &&tri)
void operator()(const triPoints &tri)
storeOp(DynamicList< triPoints > &tris)
DynamicList< triPoints > & list