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