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