isoSurface.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 | Copyright (C) 2011-2017 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "isoSurface.H"
27 #include "fvMesh.H"
28 #include "mergePoints.H"
30 #include "slicedVolFields.H"
31 #include "volFields.H"
32 #include "surfaceFields.H"
33 #include "OFstream.H"
34 #include "meshTools.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(isoSurface, 0);
41 }
42 
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 bool Foam::isoSurface::noTransform(const tensor& tt) const
47 {
48  return
49  (mag(tt.xx()-1) < mergeDistance_)
50  && (mag(tt.yy()-1) < mergeDistance_)
51  && (mag(tt.zz()-1) < mergeDistance_)
52  && (mag(tt.xy()) < mergeDistance_)
53  && (mag(tt.xz()) < mergeDistance_)
54  && (mag(tt.yx()) < mergeDistance_)
55  && (mag(tt.yz()) < mergeDistance_)
56  && (mag(tt.zx()) < mergeDistance_)
57  && (mag(tt.zy()) < mergeDistance_);
58 }
59 
60 
61 bool Foam::isoSurface::collocatedPatch(const polyPatch& pp)
62 {
63  const coupledPolyPatch& cpp = refCast<const coupledPolyPatch>(pp);
64  return cpp.parallel() && !cpp.separated();
65 }
66 
67 
68 Foam::PackedBoolList Foam::isoSurface::collocatedFaces
69 (
70  const coupledPolyPatch& pp
71 ) const
72 {
73  // Initialise to false
74  PackedBoolList collocated(pp.size());
75 
76  if (isA<processorPolyPatch>(pp))
77  {
78  if (collocatedPatch(pp))
79  {
80  forAll(pp, i)
81  {
82  collocated[i] = 1u;
83  }
84  }
85  }
86  else if (isA<cyclicPolyPatch>(pp))
87  {
88  if (collocatedPatch(pp))
89  {
90  forAll(pp, i)
91  {
92  collocated[i] = 1u;
93  }
94  }
95  }
96  else
97  {
99  << "Unhandled coupledPolyPatch type " << pp.type()
100  << abort(FatalError);
101  }
102  return collocated;
103 }
104 
105 
106 void Foam::isoSurface::syncUnseparatedPoints
107 (
108  pointField& pointValues,
109  const point& nullValue
110 ) const
111 {
112  // Until syncPointList handles separated coupled patches with multiple
113  // transforms do our own synchronisation of non-separated patches only
114  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
115 
116  if (Pstream::parRun())
117  {
118  // Send
119  forAll(patches, patchi)
120  {
121  if
122  (
123  isA<processorPolyPatch>(patches[patchi])
124  && patches[patchi].nPoints() > 0
125  && collocatedPatch(patches[patchi])
126  )
127  {
128  const processorPolyPatch& pp =
129  refCast<const processorPolyPatch>(patches[patchi]);
130 
131  const labelList& meshPts = pp.meshPoints();
132  const labelList& nbrPts = pp.neighbPoints();
133 
134  pointField patchInfo(meshPts.size());
135 
136  forAll(nbrPts, pointi)
137  {
138  label nbrPointi = nbrPts[pointi];
139  patchInfo[nbrPointi] = pointValues[meshPts[pointi]];
140  }
141 
142  OPstream toNbr
143  (
145  pp.neighbProcNo()
146  );
147  toNbr << patchInfo;
148  }
149  }
150 
151  // Receive and combine.
152 
153  forAll(patches, patchi)
154  {
155  if
156  (
157  isA<processorPolyPatch>(patches[patchi])
158  && patches[patchi].nPoints() > 0
159  && collocatedPatch(patches[patchi])
160  )
161  {
162  const processorPolyPatch& pp =
163  refCast<const processorPolyPatch>(patches[patchi]);
164 
165  pointField nbrPatchInfo(pp.nPoints());
166  {
167  // We do not know the number of points on the other side
168  // so cannot use Pstream::read.
169  IPstream fromNbr
170  (
172  pp.neighbProcNo()
173  );
174  fromNbr >> nbrPatchInfo;
175  }
176 
177  const labelList& meshPts = pp.meshPoints();
178 
179  forAll(meshPts, pointi)
180  {
181  label meshPointi = meshPts[pointi];
182  minEqOp<point>()
183  (
184  pointValues[meshPointi],
185  nbrPatchInfo[pointi]
186  );
187  }
188  }
189  }
190  }
191 
192  // Do the cyclics.
193  forAll(patches, patchi)
194  {
195  if (isA<cyclicPolyPatch>(patches[patchi]))
196  {
197  const cyclicPolyPatch& cycPatch =
198  refCast<const cyclicPolyPatch>(patches[patchi]);
199 
200  if (cycPatch.owner() && collocatedPatch(cycPatch))
201  {
202  // Owner does all.
203 
204  const edgeList& coupledPoints = cycPatch.coupledPoints();
205  const labelList& meshPts = cycPatch.meshPoints();
206  const cyclicPolyPatch& nbrPatch = cycPatch.neighbPatch();
207  const labelList& nbrMeshPoints = nbrPatch.meshPoints();
208 
209  pointField half0Values(coupledPoints.size());
210  pointField half1Values(coupledPoints.size());
211 
212  forAll(coupledPoints, i)
213  {
214  const edge& e = coupledPoints[i];
215  half0Values[i] = pointValues[meshPts[e[0]]];
216  half1Values[i] = pointValues[nbrMeshPoints[e[1]]];
217  }
218 
219  forAll(coupledPoints, i)
220  {
221  const edge& e = coupledPoints[i];
222  label p0 = meshPts[e[0]];
223  label p1 = nbrMeshPoints[e[1]];
224 
225  minEqOp<point>()(pointValues[p0], half1Values[i]);
226  minEqOp<point>()(pointValues[p1], half0Values[i]);
227  }
228  }
229  }
230  }
231 
232  // Synchronize multiple shared points.
233  const globalMeshData& pd = mesh_.globalData();
234 
235  if (pd.nGlobalPoints() > 0)
236  {
237  // Values on shared points.
238  pointField sharedPts(pd.nGlobalPoints(), nullValue);
239 
240  forAll(pd.sharedPointLabels(), i)
241  {
242  label meshPointi = pd.sharedPointLabels()[i];
243  // Fill my entries in the shared points
244  sharedPts[pd.sharedPointAddr()[i]] = pointValues[meshPointi];
245  }
246 
247  // Combine on master.
248  Pstream::listCombineGather(sharedPts, minEqOp<point>());
249  Pstream::listCombineScatter(sharedPts);
250 
251  // Now we will all have the same information. Merge it back with
252  // my local information.
253  forAll(pd.sharedPointLabels(), i)
254  {
255  label meshPointi = pd.sharedPointLabels()[i];
256  pointValues[meshPointi] = sharedPts[pd.sharedPointAddr()[i]];
257  }
258  }
259 }
260 
261 
262 Foam::scalar Foam::isoSurface::isoFraction
263 (
264  const scalar s0,
265  const scalar s1
266 ) const
267 {
268  scalar d = s1-s0;
269 
270  if (mag(d) > VSMALL)
271  {
272  return (iso_-s0)/d;
273  }
274  else
275  {
276  return -1.0;
277  }
278 }
279 
280 
281 bool Foam::isoSurface::isEdgeOfFaceCut
282 (
283  const scalarField& pVals,
284  const face& f,
285  const bool ownLower,
286  const bool neiLower
287 ) const
288 {
289  forAll(f, fp)
290  {
291  bool fpLower = (pVals[f[fp]] < iso_);
292  if
293  (
294  (fpLower != ownLower)
295  || (fpLower != neiLower)
296  || (fpLower != (pVals[f[f.fcIndex(fp)]] < iso_))
297  )
298  {
299  return true;
300  }
301  }
302  return false;
303 }
304 
305 
306 void Foam::isoSurface::getNeighbour
307 (
308  const labelList& boundaryRegion,
309  const volVectorField& meshC,
310  const volScalarField& cVals,
311  const label celli,
312  const label facei,
313  scalar& nbrValue,
314  point& nbrPoint
315 ) const
316 {
317  const labelList& own = mesh_.faceOwner();
318  const labelList& nei = mesh_.faceNeighbour();
319 
320  if (mesh_.isInternalFace(facei))
321  {
322  label nbr = (own[facei] == celli ? nei[facei] : own[facei]);
323  nbrValue = cVals[nbr];
324  nbrPoint = meshC[nbr];
325  }
326  else
327  {
328  label bFacei = facei-mesh_.nInternalFaces();
329  label patchi = boundaryRegion[bFacei];
330  const polyPatch& pp = mesh_.boundaryMesh()[patchi];
331  label patchFacei = facei-pp.start();
332 
333  nbrValue = cVals.boundaryField()[patchi][patchFacei];
334  nbrPoint = meshC.boundaryField()[patchi][patchFacei];
335  }
336 }
337 
338 
339 void Foam::isoSurface::calcCutTypes
340 (
341  const labelList& boundaryRegion,
342  const volVectorField& meshC,
343  const volScalarField& cVals,
344  const scalarField& pVals
345 )
346 {
347  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
348  const labelList& own = mesh_.faceOwner();
349  const labelList& nei = mesh_.faceNeighbour();
350 
351  faceCutType_.setSize(mesh_.nFaces());
352  faceCutType_ = NOTCUT;
353 
354  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
355  {
356  // CC edge.
357  bool ownLower = (cVals[own[facei]] < iso_);
358 
359  scalar nbrValue;
360  point nbrPoint;
361  getNeighbour
362  (
363  boundaryRegion,
364  meshC,
365  cVals,
366  own[facei],
367  facei,
368  nbrValue,
369  nbrPoint
370  );
371 
372  bool neiLower = (nbrValue < iso_);
373 
374  if (ownLower != neiLower)
375  {
376  faceCutType_[facei] = CUT;
377  }
378  else
379  {
380  // See if any mesh edge is cut by looping over all the edges of the
381  // face.
382  const face f = mesh_.faces()[facei];
383 
384  if (isEdgeOfFaceCut(pVals, f, ownLower, neiLower))
385  {
386  faceCutType_[facei] = CUT;
387  }
388  }
389  }
390 
391  forAll(patches, patchi)
392  {
393  const polyPatch& pp = patches[patchi];
394 
395  label facei = pp.start();
396 
397  forAll(pp, i)
398  {
399  bool ownLower = (cVals[own[facei]] < iso_);
400 
401  scalar nbrValue;
402  point nbrPoint;
403  getNeighbour
404  (
405  boundaryRegion,
406  meshC,
407  cVals,
408  own[facei],
409  facei,
410  nbrValue,
411  nbrPoint
412  );
413 
414  bool neiLower = (nbrValue < iso_);
415 
416  if (ownLower != neiLower)
417  {
418  faceCutType_[facei] = CUT;
419  }
420  else
421  {
422  // Mesh edge.
423  const face f = mesh_.faces()[facei];
424 
425  if (isEdgeOfFaceCut(pVals, f, ownLower, neiLower))
426  {
427  faceCutType_[facei] = CUT;
428  }
429  }
430 
431  facei++;
432  }
433  }
434 
435 
436 
437  nCutCells_ = 0;
438  cellCutType_.setSize(mesh_.nCells());
439  cellCutType_ = NOTCUT;
440 
441  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
442  {
443  if (faceCutType_[facei] != NOTCUT)
444  {
445  if (cellCutType_[own[facei]] == NOTCUT)
446  {
447  cellCutType_[own[facei]] = CUT;
448  nCutCells_++;
449  }
450  if (cellCutType_[nei[facei]] == NOTCUT)
451  {
452  cellCutType_[nei[facei]] = CUT;
453  nCutCells_++;
454  }
455  }
456  }
457  for (label facei = mesh_.nInternalFaces(); facei < mesh_.nFaces(); facei++)
458  {
459  if (faceCutType_[facei] != NOTCUT)
460  {
461  if (cellCutType_[own[facei]] == NOTCUT)
462  {
463  cellCutType_[own[facei]] = CUT;
464  nCutCells_++;
465  }
466  }
467  }
468 
469  if (debug)
470  {
471  Pout<< "isoSurface : detected " << nCutCells_
472  << " candidate cut cells (out of " << mesh_.nCells()
473  << ")." << endl;
474  }
475 }
476 
477 
478 // Caculate centre of surface.
479 Foam::point Foam::isoSurface::calcCentre(const triSurface& s)
480 {
481  vector sum = Zero;
482 
483  forAll(s, i)
484  {
485  sum += s[i].centre(s.points());
486  }
487  return sum/s.size();
488 }
489 
490 
491 void Foam::isoSurface::calcSnappedCc
492 (
493  const labelList& boundaryRegion,
494  const volVectorField& meshC,
495  const volScalarField& cVals,
496  const scalarField& pVals,
497 
498  DynamicList<point>& snappedPoints,
499  labelList& snappedCc
500 ) const
501 {
502  const pointField& pts = mesh_.points();
503  const pointField& cc = mesh_.cellCentres();
504 
505  snappedCc.setSize(mesh_.nCells());
506  snappedCc = -1;
507 
508  // Work arrays
509  DynamicList<point, 64> localTriPoints(64);
510 
511  forAll(mesh_.cells(), celli)
512  {
513  if (cellCutType_[celli] == CUT)
514  {
515  scalar cVal = cVals[celli];
516 
517  const cell& cFaces = mesh_.cells()[celli];
518 
519  localTriPoints.clear();
520  label nOther = 0;
521  point otherPointSum = Zero;
522 
523  // Create points for all intersections close to cell centre
524  // (i.e. from pyramid edges)
525 
526  forAll(cFaces, cFacei)
527  {
528  label facei = cFaces[cFacei];
529 
530  scalar nbrValue;
531  point nbrPoint;
532  getNeighbour
533  (
534  boundaryRegion,
535  meshC,
536  cVals,
537  celli,
538  facei,
539  nbrValue,
540  nbrPoint
541  );
542 
543  // Calculate intersection points of edges to cell centre
544  FixedList<scalar, 3> s;
545  FixedList<point, 3> pt;
546 
547  // From cc to neighbour cc.
548  s[2] = isoFraction(cVal, nbrValue);
549  pt[2] = (1.0-s[2])*cc[celli] + s[2]*nbrPoint;
550 
551  const face& f = mesh_.faces()[cFaces[cFacei]];
552 
553  forAll(f, fp)
554  {
555  // From cc to fp
556  label p0 = f[fp];
557  s[0] = isoFraction(cVal, pVals[p0]);
558  pt[0] = (1.0-s[0])*cc[celli] + s[0]*pts[p0];
559 
560  // From cc to fp+1
561  label p1 = f[f.fcIndex(fp)];
562  s[1] = isoFraction(cVal, pVals[p1]);
563  pt[1] = (1.0-s[1])*cc[celli] + s[1]*pts[p1];
564 
565  if
566  (
567  (s[0] >= 0.0 && s[0] <= 0.5)
568  && (s[1] >= 0.0 && s[1] <= 0.5)
569  && (s[2] >= 0.0 && s[2] <= 0.5)
570  )
571  {
572  localTriPoints.append(pt[0]);
573  localTriPoints.append(pt[1]);
574  localTriPoints.append(pt[2]);
575  }
576  else
577  {
578  // Get average of all other points
579  forAll(s, i)
580  {
581  if (s[i] >= 0.0 && s[i] <= 0.5)
582  {
583  otherPointSum += pt[i];
584  nOther++;
585  }
586  }
587  }
588  }
589  }
590 
591  if (localTriPoints.size() == 0)
592  {
593  // No complete triangles. Get average of all intersection
594  // points.
595  if (nOther > 0)
596  {
597  snappedCc[celli] = snappedPoints.size();
598  snappedPoints.append(otherPointSum/nOther);
599 
600  //Pout<< " point:" << pointi
601  // << " replacing coord:" << mesh_.points()[pointi]
602  // << " by average:" << collapsedPoint[pointi] << endl;
603  }
604  }
605  else if (localTriPoints.size() == 3)
606  {
607  // Single triangle. No need for any analysis. Average points.
609  points.transfer(localTriPoints);
610  snappedCc[celli] = snappedPoints.size();
611  snappedPoints.append(sum(points)/points.size());
612 
613  //Pout<< " point:" << pointi
614  // << " replacing coord:" << mesh_.points()[pointi]
615  // << " by average:" << collapsedPoint[pointi] << endl;
616  }
617  else
618  {
619  // Convert points into triSurface.
620 
621  // Merge points and compact out non-valid triangles
622  labelList triMap; // merged to unmerged triangle
623  labelList triPointReverseMap; // unmerged to merged point
624  triSurface surf
625  (
626  stitchTriPoints
627  (
628  false, // do not check for duplicate tris
629  localTriPoints,
630  triPointReverseMap,
631  triMap
632  )
633  );
634 
635  labelList faceZone;
636  label nZones = surf.markZones
637  (
638  boolList(surf.nEdges(), false),
639  faceZone
640  );
641 
642  if (nZones == 1)
643  {
644  snappedCc[celli] = snappedPoints.size();
645  snappedPoints.append(calcCentre(surf));
646  //Pout<< " point:" << pointi << " nZones:" << nZones
647  // << " replacing coord:" << mesh_.points()[pointi]
648  // << " by average:" << collapsedPoint[pointi] << endl;
649  }
650  }
651  }
652  }
653 }
654 
655 
656 void Foam::isoSurface::calcSnappedPoint
657 (
658  const PackedBoolList& isBoundaryPoint,
659  const labelList& boundaryRegion,
660  const volVectorField& meshC,
661  const volScalarField& cVals,
662  const scalarField& pVals,
663 
664  DynamicList<point>& snappedPoints,
665  labelList& snappedPoint
666 ) const
667 {
668  const pointField& pts = mesh_.points();
669  const pointField& cc = mesh_.cellCentres();
670 
671  pointField collapsedPoint(mesh_.nPoints(), point::max);
672 
673 
674  // Work arrays
675  DynamicList<point, 64> localTriPoints(100);
676 
677  forAll(mesh_.pointFaces(), pointi)
678  {
679  if (isBoundaryPoint.get(pointi) == 1)
680  {
681  continue;
682  }
683 
684  const labelList& pFaces = mesh_.pointFaces()[pointi];
685 
686  bool anyCut = false;
687 
688  forAll(pFaces, i)
689  {
690  label facei = pFaces[i];
691 
692  if (faceCutType_[facei] == CUT)
693  {
694  anyCut = true;
695  break;
696  }
697  }
698 
699  if (!anyCut)
700  {
701  continue;
702  }
703 
704 
705  localTriPoints.clear();
706  label nOther = 0;
707  point otherPointSum = Zero;
708 
709  forAll(pFaces, pFacei)
710  {
711  // Create points for all intersections close to point
712  // (i.e. from pyramid edges)
713 
714  label facei = pFaces[pFacei];
715  const face& f = mesh_.faces()[facei];
716  label own = mesh_.faceOwner()[facei];
717 
718  // Get neighbour value
719  scalar nbrValue;
720  point nbrPoint;
721  getNeighbour
722  (
723  boundaryRegion,
724  meshC,
725  cVals,
726  own,
727  facei,
728  nbrValue,
729  nbrPoint
730  );
731 
732  // Calculate intersection points of edges emanating from point
733  FixedList<scalar, 4> s;
734  FixedList<point, 4> pt;
735 
736  label fp = findIndex(f, pointi);
737  s[0] = isoFraction(pVals[pointi], cVals[own]);
738  pt[0] = (1.0-s[0])*pts[pointi] + s[0]*cc[own];
739 
740  s[1] = isoFraction(pVals[pointi], nbrValue);
741  pt[1] = (1.0-s[1])*pts[pointi] + s[1]*nbrPoint;
742 
743  label nextPointi = f[f.fcIndex(fp)];
744  s[2] = isoFraction(pVals[pointi], pVals[nextPointi]);
745  pt[2] = (1.0-s[2])*pts[pointi] + s[2]*pts[nextPointi];
746 
747  label prevPointi = f[f.rcIndex(fp)];
748  s[3] = isoFraction(pVals[pointi], pVals[prevPointi]);
749  pt[3] = (1.0-s[3])*pts[pointi] + s[3]*pts[prevPointi];
750 
751  if
752  (
753  (s[0] >= 0.0 && s[0] <= 0.5)
754  && (s[1] >= 0.0 && s[1] <= 0.5)
755  && (s[2] >= 0.0 && s[2] <= 0.5)
756  )
757  {
758  localTriPoints.append(pt[0]);
759  localTriPoints.append(pt[1]);
760  localTriPoints.append(pt[2]);
761  }
762  if
763  (
764  (s[0] >= 0.0 && s[0] <= 0.5)
765  && (s[1] >= 0.0 && s[1] <= 0.5)
766  && (s[3] >= 0.0 && s[3] <= 0.5)
767  )
768  {
769  localTriPoints.append(pt[3]);
770  localTriPoints.append(pt[0]);
771  localTriPoints.append(pt[1]);
772  }
773 
774  // Get average of points as fallback
775  forAll(s, i)
776  {
777  if (s[i] >= 0.0 && s[i] <= 0.5)
778  {
779  otherPointSum += pt[i];
780  nOther++;
781  }
782  }
783  }
784 
785  if (localTriPoints.size() == 0)
786  {
787  // No complete triangles. Get average of all intersection
788  // points.
789  if (nOther > 0)
790  {
791  collapsedPoint[pointi] = otherPointSum/nOther;
792  }
793  }
794  else if (localTriPoints.size() == 3)
795  {
796  // Single triangle. No need for any analysis. Average points.
798  points.transfer(localTriPoints);
799  collapsedPoint[pointi] = sum(points)/points.size();
800  }
801  else
802  {
803  // Convert points into triSurface.
804 
805  // Merge points and compact out non-valid triangles
806  labelList triMap; // merged to unmerged triangle
807  labelList triPointReverseMap; // unmerged to merged point
808  triSurface surf
809  (
810  stitchTriPoints
811  (
812  false, // do not check for duplicate tris
813  localTriPoints,
814  triPointReverseMap,
815  triMap
816  )
817  );
818 
819  labelList faceZone;
820  label nZones = surf.markZones
821  (
822  boolList(surf.nEdges(), false),
823  faceZone
824  );
825 
826  if (nZones == 1)
827  {
828  collapsedPoint[pointi] = calcCentre(surf);
829  }
830  }
831  }
832 
833 
834  // Synchronise snap location
835  syncUnseparatedPoints(collapsedPoint, point::max);
836 
837 
838  snappedPoint.setSize(mesh_.nPoints());
839  snappedPoint = -1;
840 
841  forAll(collapsedPoint, pointi)
842  {
843  if (collapsedPoint[pointi] != point::max)
844  {
845  snappedPoint[pointi] = snappedPoints.size();
846  snappedPoints.append(collapsedPoint[pointi]);
847  }
848  }
849 }
850 
851 
852 Foam::triSurface Foam::isoSurface::stitchTriPoints
853 (
854  const bool checkDuplicates,
855  const List<point>& triPoints,
856  labelList& triPointReverseMap, // unmerged to merged point
857  labelList& triMap // merged to unmerged triangle
858 ) const
859 {
860  label nTris = triPoints.size()/3;
861 
862  if ((triPoints.size() % 3) != 0)
863  {
865  << "Problem: number of points " << triPoints.size()
866  << " not a multiple of 3." << abort(FatalError);
867  }
868 
869  pointField newPoints;
871  (
872  triPoints,
873  mergeDistance_,
874  false,
875  triPointReverseMap,
876  newPoints
877  );
878 
879  // Check that enough merged.
880  if (debug)
881  {
882  pointField newNewPoints;
883  labelList oldToNew;
884  bool hasMerged = mergePoints
885  (
886  newPoints,
887  mergeDistance_,
888  true,
889  oldToNew,
890  newNewPoints
891  );
892 
893  if (hasMerged)
894  {
896  << "Merged points contain duplicates"
897  << " when merging with distance " << mergeDistance_ << endl
898  << "merged:" << newPoints.size() << " re-merged:"
899  << newNewPoints.size()
900  << abort(FatalError);
901  }
902  }
903 
904 
905  List<labelledTri> tris;
906  {
907  DynamicList<labelledTri> dynTris(nTris);
908  label rawPointi = 0;
909  DynamicList<label> newToOldTri(nTris);
910 
911  for (label oldTriI = 0; oldTriI < nTris; oldTriI++)
912  {
913  labelledTri tri
914  (
915  triPointReverseMap[rawPointi],
916  triPointReverseMap[rawPointi+1],
917  triPointReverseMap[rawPointi+2],
918  0
919  );
920  rawPointi += 3;
921 
922  if ((tri[0] != tri[1]) && (tri[0] != tri[2]) && (tri[1] != tri[2]))
923  {
924  newToOldTri.append(oldTriI);
925  dynTris.append(tri);
926  }
927  }
928 
929  triMap.transfer(newToOldTri);
930  tris.transfer(dynTris);
931  }
932 
933 
934 
935  // Determine 'flat hole' situation (see RMT paper).
936  // Two unconnected triangles get connected because (some of) the edges
937  // separating them get collapsed. Below only checks for duplicate triangles,
938  // not non-manifold edge connectivity.
939  if (checkDuplicates)
940  {
941  labelListList pointFaces;
942  invertManyToMany(newPoints.size(), tris, pointFaces);
943 
944  // Filter out duplicates.
945  DynamicList<label> newToOldTri(tris.size());
946 
947  forAll(tris, triI)
948  {
949  const labelledTri& tri = tris[triI];
950  const labelList& pFaces = pointFaces[tri[0]];
951 
952  // Find the maximum of any duplicates. Maximum since the tris
953  // below triI
954  // get overwritten so we cannot use them in a comparison.
955  label dupTriI = -1;
956  forAll(pFaces, i)
957  {
958  label nbrTriI = pFaces[i];
959 
960  if (nbrTriI > triI && (tris[nbrTriI] == tri))
961  {
962  //Pout<< "Duplicate : " << triI << " verts:" << tri
963  // << " to " << nbrTriI << " verts:" << tris[nbrTriI]
964  // << endl;
965  dupTriI = nbrTriI;
966  break;
967  }
968  }
969 
970  if (dupTriI == -1)
971  {
972  // There is no (higher numbered) duplicate triangle
973  label newTriI = newToOldTri.size();
974  newToOldTri.append(triMap[triI]);
975  tris[newTriI] = tris[triI];
976  }
977  }
978 
979  triMap.transfer(newToOldTri);
980  tris.setSize(triMap.size());
981 
982  if (debug)
983  {
984  Pout<< "isoSurface : merged from " << nTris
985  << " down to " << tris.size() << " unique triangles." << endl;
986  }
987 
988 
989  if (debug)
990  {
991  triSurface surf(tris, geometricSurfacePatchList(0), newPoints);
992 
993  forAll(surf, facei)
994  {
995  const labelledTri& f = surf[facei];
996  const labelList& fFaces = surf.faceFaces()[facei];
997 
998  forAll(fFaces, i)
999  {
1000  label nbrFacei = fFaces[i];
1001 
1002  if (nbrFacei <= facei)
1003  {
1004  // lower numbered faces already checked
1005  continue;
1006  }
1007 
1008  const labelledTri& nbrF = surf[nbrFacei];
1009 
1010  if (f == nbrF)
1011  {
1013  << "Check : "
1014  << " triangle " << facei << " vertices " << f
1015  << " fc:" << f.centre(surf.points())
1016  << " has the same vertices as triangle " << nbrFacei
1017  << " vertices " << nbrF
1018  << " fc:" << nbrF.centre(surf.points())
1019  << abort(FatalError);
1020  }
1021  }
1022  }
1023  }
1024  }
1025 
1026  return triSurface(tris, geometricSurfacePatchList(0), newPoints, true);
1027 }
1028 
1029 
1030 bool Foam::isoSurface::validTri(const triSurface& surf, const label facei)
1031 {
1032  // Simple check on indices ok.
1033 
1034  const labelledTri& f = surf[facei];
1035 
1036  if
1037  (
1038  (f[0] < 0) || (f[0] >= surf.points().size())
1039  || (f[1] < 0) || (f[1] >= surf.points().size())
1040  || (f[2] < 0) || (f[2] >= surf.points().size())
1041  )
1042  {
1044  << "triangle " << facei << " vertices " << f
1045  << " uses point indices outside point range 0.."
1046  << surf.points().size()-1 << endl;
1047 
1048  return false;
1049  }
1050 
1051  if ((f[0] == f[1]) || (f[0] == f[2]) || (f[1] == f[2]))
1052  {
1054  << "triangle " << facei
1055  << " uses non-unique vertices " << f
1056  << endl;
1057  return false;
1058  }
1059 
1060  // duplicate triangle check
1061 
1062  const labelList& fFaces = surf.faceFaces()[facei];
1063 
1064  // Check if faceNeighbours use same points as this face.
1065  // Note: discards normal information - sides of baffle are merged.
1066  forAll(fFaces, i)
1067  {
1068  label nbrFacei = fFaces[i];
1069 
1070  if (nbrFacei <= facei)
1071  {
1072  // lower numbered faces already checked
1073  continue;
1074  }
1075 
1076  const labelledTri& nbrF = surf[nbrFacei];
1077 
1078  if
1079  (
1080  ((f[0] == nbrF[0]) || (f[0] == nbrF[1]) || (f[0] == nbrF[2]))
1081  && ((f[1] == nbrF[0]) || (f[1] == nbrF[1]) || (f[1] == nbrF[2]))
1082  && ((f[2] == nbrF[0]) || (f[2] == nbrF[1]) || (f[2] == nbrF[2]))
1083  )
1084  {
1086  << "triangle " << facei << " vertices " << f
1087  << " fc:" << f.centre(surf.points())
1088  << " has the same vertices as triangle " << nbrFacei
1089  << " vertices " << nbrF
1090  << " fc:" << nbrF.centre(surf.points())
1091  << endl;
1092 
1093  return false;
1094  }
1095  }
1096  return true;
1097 }
1098 
1099 
1100 Foam::triSurface Foam::isoSurface::subsetMesh
1101 (
1102  const triSurface& s,
1103  const labelList& newToOldFaces,
1104  labelList& oldToNewPoints,
1105  labelList& newToOldPoints
1106 )
1107 {
1108  const boolList include
1109  (
1110  createWithValues<boolList>
1111  (
1112  s.size(),
1113  false,
1114  newToOldFaces,
1115  true
1116  )
1117  );
1118 
1119  newToOldPoints.setSize(s.points().size());
1120  oldToNewPoints.setSize(s.points().size());
1121  oldToNewPoints = -1;
1122  {
1123  label pointi = 0;
1124 
1125  forAll(include, oldFacei)
1126  {
1127  if (include[oldFacei])
1128  {
1129  // Renumber labels for triangle
1130  const labelledTri& tri = s[oldFacei];
1131 
1132  forAll(tri, fp)
1133  {
1134  label oldPointi = tri[fp];
1135 
1136  if (oldToNewPoints[oldPointi] == -1)
1137  {
1138  oldToNewPoints[oldPointi] = pointi;
1139  newToOldPoints[pointi++] = oldPointi;
1140  }
1141  }
1142  }
1143  }
1144  newToOldPoints.setSize(pointi);
1145  }
1146 
1147  // Extract points
1148  pointField newPoints(newToOldPoints.size());
1149  forAll(newToOldPoints, i)
1150  {
1151  newPoints[i] = s.points()[newToOldPoints[i]];
1152  }
1153  // Extract faces
1154  List<labelledTri> newTriangles(newToOldFaces.size());
1155 
1156  forAll(newToOldFaces, i)
1157  {
1158  // Get old vertex labels
1159  const labelledTri& tri = s[newToOldFaces[i]];
1160 
1161  newTriangles[i][0] = oldToNewPoints[tri[0]];
1162  newTriangles[i][1] = oldToNewPoints[tri[1]];
1163  newTriangles[i][2] = oldToNewPoints[tri[2]];
1164  newTriangles[i].region() = tri.region();
1165  }
1166 
1167  // Reuse storage.
1168  return triSurface(newTriangles, s.patches(), newPoints, true);
1169 }
1170 
1171 
1172 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1173 
1176  const volScalarField& cVals,
1177  const scalarField& pVals,
1178  const scalar iso,
1179  const bool regularise,
1180  const scalar mergeTol
1181 )
1182 :
1183  mesh_(cVals.mesh()),
1184  pVals_(pVals),
1185  iso_(iso),
1186  regularise_(regularise),
1187  mergeDistance_(mergeTol*mesh_.bounds().mag())
1188 {
1189  if (debug)
1190  {
1191  Pout<< "isoSurface:" << nl
1192  << " isoField : " << cVals.name() << nl
1193  << " cell min/max : "
1194  << min(cVals.primitiveField()) << " / "
1195  << max(cVals.primitiveField()) << nl
1196  << " point min/max : "
1197  << min(pVals_) << " / "
1198  << max(pVals_) << nl
1199  << " isoValue : " << iso << nl
1200  << " regularise : " << regularise_ << nl
1201  << " mergeTol : " << mergeTol << nl
1202  << endl;
1203  }
1204 
1205  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
1206 
1207 
1208  // Rewrite input field
1209  // ~~~~~~~~~~~~~~~~~~~
1210  // Rewrite input volScalarField to have interpolated values
1211  // on separated patches.
1212 
1213  cValsPtr_.reset(adaptPatchFields(cVals).ptr());
1214 
1215 
1216  // Construct cell centres field consistent with cVals
1217  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1218  // Generate field to interpolate. This is identical to the mesh.C()
1219  // except on separated coupled patches and on empty patches.
1220 
1221  slicedVolVectorField meshC
1222  (
1223  IOobject
1224  (
1225  "C",
1226  mesh_.pointsInstance(),
1227  mesh_.meshSubDir,
1228  mesh_,
1231  false
1232  ),
1233  mesh_,
1234  dimLength,
1235  mesh_.cellCentres(),
1236  mesh_.faceCentres()
1237  );
1238  forAll(patches, patchi)
1239  {
1240  const polyPatch& pp = patches[patchi];
1241 
1242  // Adapt separated coupled (proc and cyclic) patches
1243  if (pp.coupled())
1244  {
1245  fvPatchVectorField& pfld = const_cast<fvPatchVectorField&>
1246  (
1247  meshC.boundaryField()[patchi]
1248  );
1249 
1250  PackedBoolList isCollocated
1251  (
1252  collocatedFaces(refCast<const coupledPolyPatch>(pp))
1253  );
1254 
1255  forAll(isCollocated, i)
1256  {
1257  if (!isCollocated[i])
1258  {
1259  pfld[i] = mesh_.faceCentres()[pp.start()+i];
1260  }
1261  }
1262  }
1263  else if (isA<emptyPolyPatch>(pp))
1264  {
1265  typedef slicedVolVectorField::Boundary bType;
1266 
1267  bType& bfld = const_cast<bType&>(meshC.boundaryField());
1268 
1269  // Clear old value. Cannot resize it since is a slice.
1270  bfld.set(patchi, nullptr);
1271 
1272  // Set new value we can change
1273  bfld.set
1274  (
1275  patchi,
1277  (
1278  mesh_.boundary()[patchi],
1279  meshC
1280  )
1281  );
1282 
1283  // Change to face centres
1284  bfld[patchi] = pp.patchSlice(mesh_.faceCentres());
1285  }
1286  }
1287 
1288 
1289 
1290  // Pre-calculate patch-per-face to avoid whichPatch call.
1291  labelList boundaryRegion(mesh_.nFaces()-mesh_.nInternalFaces());
1292 
1293  forAll(patches, patchi)
1294  {
1295  const polyPatch& pp = patches[patchi];
1296 
1297  label facei = pp.start();
1298 
1299  forAll(pp, i)
1300  {
1301  boundaryRegion[facei-mesh_.nInternalFaces()] = patchi;
1302  facei++;
1303  }
1304  }
1305 
1306 
1307 
1308  // Determine if any cut through face/cell
1309  calcCutTypes(boundaryRegion, meshC, cValsPtr_(), pVals_);
1310 
1311 
1312  DynamicList<point> snappedPoints(nCutCells_);
1313 
1314  // Per cc -1 or a point inside snappedPoints.
1315  labelList snappedCc;
1316  if (regularise_)
1317  {
1318  calcSnappedCc
1319  (
1320  boundaryRegion,
1321  meshC,
1322  cValsPtr_(),
1323  pVals_,
1324 
1325  snappedPoints,
1326  snappedCc
1327  );
1328  }
1329  else
1330  {
1331  snappedCc.setSize(mesh_.nCells());
1332  snappedCc = -1;
1333  }
1334 
1335 
1336 
1337  if (debug)
1338  {
1339  Pout<< "isoSurface : shifted " << snappedPoints.size()
1340  << " cell centres to intersection." << endl;
1341  }
1342 
1343  label nCellSnaps = snappedPoints.size();
1344 
1345 
1346  // Per point -1 or a point inside snappedPoints.
1347  labelList snappedPoint;
1348  if (regularise_)
1349  {
1350  // Determine if point is on boundary.
1351  PackedBoolList isBoundaryPoint(mesh_.nPoints());
1352 
1353  forAll(patches, patchi)
1354  {
1355  // Mark all boundary points that are not physically coupled
1356  // (so anything but collocated coupled patches)
1357 
1358  if (patches[patchi].coupled())
1359  {
1360  const coupledPolyPatch& cpp =
1361  refCast<const coupledPolyPatch>
1362  (
1363  patches[patchi]
1364  );
1365 
1366  PackedBoolList isCollocated(collocatedFaces(cpp));
1367 
1368  forAll(isCollocated, i)
1369  {
1370  if (!isCollocated[i])
1371  {
1372  const face& f = mesh_.faces()[cpp.start()+i];
1373 
1374  forAll(f, fp)
1375  {
1376  isBoundaryPoint.set(f[fp], 1);
1377  }
1378  }
1379  }
1380  }
1381  else
1382  {
1383  const polyPatch& pp = patches[patchi];
1384 
1385  forAll(pp, i)
1386  {
1387  const face& f = mesh_.faces()[pp.start()+i];
1388 
1389  forAll(f, fp)
1390  {
1391  isBoundaryPoint.set(f[fp], 1);
1392  }
1393  }
1394  }
1395  }
1396 
1397  calcSnappedPoint
1398  (
1399  isBoundaryPoint,
1400  boundaryRegion,
1401  meshC,
1402  cValsPtr_(),
1403  pVals_,
1404 
1405  snappedPoints,
1406  snappedPoint
1407  );
1408  }
1409  else
1410  {
1411  snappedPoint.setSize(mesh_.nPoints());
1412  snappedPoint = -1;
1413  }
1414 
1415  if (debug)
1416  {
1417  Pout<< "isoSurface : shifted " << snappedPoints.size()-nCellSnaps
1418  << " vertices to intersection." << endl;
1419  }
1420 
1421 
1422  {
1423  DynamicList<point> triPoints(3*nCutCells_);
1424  DynamicList<label> triMeshCells(nCutCells_);
1425 
1426  generateTriPoints
1427  (
1428  cValsPtr_(),
1429  pVals_,
1430 
1431  meshC,
1432  mesh_.points(),
1433 
1434  snappedPoints,
1435  snappedCc,
1436  snappedPoint,
1437 
1438  triPoints, // 3 points of the triangle
1439  triMeshCells // per triangle the originating cell
1440  );
1441 
1442  if (debug)
1443  {
1444  Pout<< "isoSurface : generated " << triMeshCells.size()
1445  << " unmerged triangles from " << triPoints.size()
1446  << " unmerged points." << endl;
1447  }
1448 
1449 
1450  // Merge points and compact out non-valid triangles
1451  labelList triMap; // merged to unmerged triangle
1453  (
1454  stitchTriPoints
1455  (
1456  true, // check for duplicate tris
1457  triPoints,
1458  triPointMergeMap_, // unmerged to merged point
1459  triMap
1460  )
1461  );
1462 
1463  if (debug)
1464  {
1465  Pout<< "isoSurface : generated " << triMap.size()
1466  << " merged triangles." << endl;
1467  }
1468 
1469  meshCells_.setSize(triMap.size());
1470  forAll(triMap, i)
1471  {
1472  meshCells_[i] = triMeshCells[triMap[i]];
1473  }
1474  }
1475 
1476  if (debug)
1477  {
1478  Pout<< "isoSurface : checking " << size()
1479  << " triangles for validity." << endl;
1480 
1481  forAll(*this, triI)
1482  {
1483  // Copied from surfaceCheck
1484  validTri(*this, triI);
1485  }
1486  }
1487 
1488 
1489  if (debug)
1490  {
1491  fileName stlFile = mesh_.time().path() + ".stl";
1492  Pout<< "Dumping surface to " << stlFile << endl;
1493  triSurface::write(stlFile);
1494  }
1495 }
1496 
1497 
1498 // ************************************************************************* //
Foam::surfaceFields.
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const word & name() const
Return name.
Definition: IOobject.H:291
A class for handling file names.
Definition: fileName.H:69
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
friend Ostream & operator(Ostream &, const UList< T > &)
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
const double e
Elementary charge.
Definition: doubleFloat.H:78
List< geometricSurfacePatch > geometricSurfacePatchList
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
isoSurface(const volScalarField &cellIsoVals, const scalarField &pointIsoVals, const scalar iso, const bool regularise, const scalar mergeTol=1e-6)
Construct from cell values and point values. Uses boundaryField.
Definition: isoSurface.C:1175
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Macros for easy insertion into run-time selection tables.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().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))
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:310
const pointField & points
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
List< edge > edgeList
Definition: edgeList.H:38
label nPoints
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:91
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::polyBoundaryMesh.
static const char nl
Definition: Ostream.H:262
const Mesh & mesh() const
Return mesh.
defineTypeNameAndDebug(combustionModel, 0)
Merge points. See below.
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
void setSize(const label)
Reset size of List.
Definition: List.C:281
A bit-packed bool list.
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:394
label patchi
vector point
Point is a vector.
Definition: point.H:41
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:300
label mergePoints(const UList< Type > &points, const scalar mergeTol, const bool verbose, labelList &pointMap, const Type &origin=Type::zero)
Sorts and merges points. All points closer than/equal mergeTol get merged.
void invertManyToMany(const label len, const UList< InList > &, List< OutList > &)
Invert many-to-many.
dimensioned< scalar > mag(const dimensioned< Type > &)
fileName path() const
Return directory path name (part before last /)
Definition: fileName.C:250
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Triangulated surface description with patch information.
Definition: triSurface.H:65
Specialization of GeometricField which holds slices of given complete fields in a form that they act ...
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
const List< T >::subList patchSlice(const UList< T > &l) const
Slice list to patch.
Definition: polyPatch.H:333
Namespace for OpenFOAM.