surfaceFeatures.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2017-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 "surfaceFeatures.H"
30 #include "triSurface.H"
31 #include "indexedOctree.H"
32 #include "treeDataEdge.H"
33 #include "treeDataPoint.H"
34 #include "meshTools.H"
35 #include "Fstream.H"
36 #include "unitConversion.H"
37 #include "edgeHashes.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43  defineTypeNameAndDebug(surfaceFeatures, 0);
44 
45  const scalar surfaceFeatures::parallelTolerance = sin(degToRad(1.0));
46 
47 
49 // Check if the point is on the line
50 static bool onLine(const Foam::point& p, const linePointRef& line)
51 {
52  const point& a = line.start();
53  const point& b = line.end();
54 
55  if
56  (
57  (p.x() < min(a.x(), b.x()) || p.x() > max(a.x(), b.x()))
58  || (p.y() < min(a.y(), b.y()) || p.y() > max(a.y(), b.y()))
59  || (p.z() < min(a.z(), b.z()) || p.z() > max(a.z(), b.z()))
60  )
61  {
62  return false;
63  }
64 
65  return true;
66 }
68 
69 } // End namespace Foam
70 
71 
72 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
73 
74 Foam::pointIndexHit Foam::surfaceFeatures::edgeNearest
75 (
76  const linePointRef& line,
77  const point& sample
78 )
79 {
80  pointHit eHit = line.nearestDist(sample);
81 
82  // Classification of position on edge.
83  label endPoint;
84 
85  if (eHit.hit())
86  {
87  // Nearest point is on edge itself.
88  // Note: might be at or very close to endpoint. We should use tolerance
89  // here.
90  endPoint = -1;
91  }
92  else
93  {
94  // Nearest point has to be one of the end points. Find out
95  // which one.
96  if
97  (
98  eHit.point().distSqr(line.start())
99  < eHit.point().distSqr(line.end())
100  )
101  {
102  endPoint = 0;
103  }
104  else
105  {
106  endPoint = 1;
107  }
108  }
109 
110  return pointIndexHit(eHit, endPoint);
111 }
112 
113 
114 // Go from selected edges only to a value for every edge
116  const
117 {
118  List<edgeStatus> edgeStat(surf_.nEdges(), NONE);
119 
120  // Region edges first
121  for (label i = 0; i < externalStart_; i++)
122  {
123  edgeStat[featureEdges_[i]] = REGION;
124  }
125 
126  // External edges
127  for (label i = externalStart_; i < internalStart_; i++)
128  {
129  edgeStat[featureEdges_[i]] = EXTERNAL;
130  }
131 
132  // Internal edges
133  for (label i = internalStart_; i < featureEdges_.size(); i++)
134  {
135  edgeStat[featureEdges_[i]] = INTERNAL;
136  }
138  return edgeStat;
139 }
140 
141 
142 // Set from value for every edge
144 (
145  const List<edgeStatus>& edgeStat,
146  const scalar includedAngle
147 )
148 {
149  // Count
150 
151  label nRegion = 0;
152  label nExternal = 0;
153  label nInternal = 0;
154 
155  forAll(edgeStat, edgeI)
156  {
157  if (edgeStat[edgeI] == REGION)
158  {
159  nRegion++;
160  }
161  else if (edgeStat[edgeI] == EXTERNAL)
162  {
163  nExternal++;
164  }
165  else if (edgeStat[edgeI] == INTERNAL)
166  {
167  nInternal++;
168  }
169  }
170 
171  externalStart_ = nRegion;
172  internalStart_ = externalStart_ + nExternal;
173 
174 
175  // Copy
176 
177  featureEdges_.setSize(internalStart_ + nInternal);
178 
179  label regionI = 0;
180  label externalI = externalStart_;
181  label internalI = internalStart_;
182 
183  forAll(edgeStat, edgeI)
184  {
185  if (edgeStat[edgeI] == REGION)
186  {
187  featureEdges_[regionI++] = edgeI;
188  }
189  else if (edgeStat[edgeI] == EXTERNAL)
190  {
191  featureEdges_[externalI++] = edgeI;
192  }
193  else if (edgeStat[edgeI] == INTERNAL)
194  {
195  featureEdges_[internalI++] = edgeI;
196  }
197  }
198 
199  const scalar minCos = Foam::cos(degToRad(180.0 - includedAngle));
200 
201  calcFeatPoints(edgeStat, minCos);
202 }
203 
204 
205 //construct feature points where more than 2 feature edges meet
206 void Foam::surfaceFeatures::calcFeatPoints
207 (
208  const List<edgeStatus>& edgeStat,
209  const scalar minCos
210 )
211 {
212  DynamicList<label> featurePoints(surf_.nPoints()/1000);
213 
214  const labelListList& pointEdges = surf_.pointEdges();
215  const edgeList& edges = surf_.edges();
216  const pointField& localPoints = surf_.localPoints();
217 
218  forAll(pointEdges, pointi)
219  {
220  const labelList& pEdges = pointEdges[pointi];
221 
222  label nFeatEdges = 0;
223 
224  forAll(pEdges, i)
225  {
226  if (edgeStat[pEdges[i]] != NONE)
227  {
228  nFeatEdges++;
229  }
230  }
231 
232  if (nFeatEdges > 2)
233  {
234  featurePoints.append(pointi);
235  }
236  else if (nFeatEdges == 2)
237  {
238  // Check the angle between the two edges
239  DynamicList<vector> edgeVecs(2);
240 
241  forAll(pEdges, i)
242  {
243  const label edgeI = pEdges[i];
244 
245  if (edgeStat[edgeI] != NONE)
246  {
247  vector vec = edges[edgeI].vec(localPoints);
248  scalar magVec = mag(vec);
249  if (magVec > SMALL)
250  {
251  edgeVecs.append(vec/magVec);
252  }
253  }
254  }
255 
256  if (edgeVecs.size() == 2 && mag(edgeVecs[0] & edgeVecs[1]) < minCos)
257  {
258  featurePoints.append(pointi);
259  }
260  }
261  }
262 
263  featurePoints_.transfer(featurePoints);
264 }
265 
266 
267 void Foam::surfaceFeatures::classifyFeatureAngles
268 (
269  const labelListList& edgeFaces,
270  List<edgeStatus>& edgeStat,
271  const scalar minCos,
272  const bool geometricTestOnly
273 ) const
274 {
275  const vectorField& faceNormals = surf_.faceNormals();
276  const pointField& points = surf_.points();
277 
278  // Special case: minCos=1
279  bool selectAll = (mag(minCos-1.0) < SMALL);
280 
281  forAll(edgeFaces, edgeI)
282  {
283  const labelList& eFaces = edgeFaces[edgeI];
284 
285  if (eFaces.size() != 2)
286  {
287  // Non-manifold. What to do here? Is region edge? external edge?
288  edgeStat[edgeI] = REGION;
289  }
290  else
291  {
292  label face0 = eFaces[0];
293  label face1 = eFaces[1];
294 
295  if
296  (
297  !geometricTestOnly
298  && surf_[face0].region() != surf_[face1].region()
299  )
300  {
301  edgeStat[edgeI] = REGION;
302  }
303  else if
304  (
305  selectAll
306  || ((faceNormals[face0] & faceNormals[face1]) < minCos)
307  )
308  {
309  // Check if convex or concave by looking at angle
310  // between face centres and normal
311  vector f0Tof1 =
312  surf_[face1].centre(points)
313  - surf_[face0].centre(points);
314 
315  if ((f0Tof1 & faceNormals[face0]) >= 0.0)
316  {
317  edgeStat[edgeI] = INTERNAL;
318  }
319  else
320  {
321  edgeStat[edgeI] = EXTERNAL;
322  }
323  }
324  }
325  }
326 }
327 
328 
329 // Returns next feature edge connected to pointi with correct value.
330 Foam::label Foam::surfaceFeatures::nextFeatEdge
331 (
332  const List<edgeStatus>& edgeStat,
333  const labelList& featVisited,
334  const label unsetVal,
335  const label prevEdgeI,
336  const label vertI
337 ) const
338 {
339  const labelList& pEdges = surf_.pointEdges()[vertI];
340 
341  label nextEdgeI = -1;
342 
343  forAll(pEdges, i)
344  {
345  label edgeI = pEdges[i];
346 
347  if
348  (
349  edgeI != prevEdgeI
350  && edgeStat[edgeI] != NONE
351  && featVisited[edgeI] == unsetVal
352  )
353  {
354  if (nextEdgeI == -1)
355  {
356  nextEdgeI = edgeI;
357  }
358  else
359  {
360  // More than one feature edge to choose from. End of segment.
361  return -1;
362  }
363  }
364  }
365 
366  return nextEdgeI;
367 }
368 
369 
370 // Finds connected feature edges by walking from prevEdgeI in direction of
371 // prevPointi. Marks feature edges visited in featVisited by assigning them
372 // the current feature line number. Returns cumulative length of edges walked.
373 // Works in one of two modes:
374 // - mark : step to edges with featVisited = -1.
375 // Mark edges visited with currentFeatI.
376 // - clear : step to edges with featVisited = currentFeatI
377 // Mark edges visited with -2 and erase from feature edges.
378 Foam::surfaceFeatures::labelScalar Foam::surfaceFeatures::walkSegment
379 (
380  const bool mark,
381  const List<edgeStatus>& edgeStat,
382  const label startEdgeI,
383  const label startPointi,
384  const label currentFeatI,
385  labelList& featVisited
386 )
387 {
388  label edgeI = startEdgeI;
389 
390  label vertI = startPointi;
391 
392  scalar visitedLength = 0.0;
393 
394  label nVisited = 0;
395 
396  if (featurePoints_.found(startPointi))
397  {
398  // Do not walk across feature points
399 
400  return labelScalar(nVisited, visitedLength);
401  }
402 
403  //
404  // Now we have:
405  // edgeI : first edge on this segment
406  // vertI : one of the endpoints of this segment
407  //
408  // Start walking, marking off edges (in featVisited)
409  // as we go along.
410  //
411 
412  label unsetVal;
413 
414  if (mark)
415  {
416  unsetVal = -1;
417  }
418  else
419  {
420  unsetVal = currentFeatI;
421  }
422 
423  do
424  {
425  // Step to next feature edge with value unsetVal
426  edgeI = nextFeatEdge(edgeStat, featVisited, unsetVal, edgeI, vertI);
427 
428  if (edgeI == -1 || edgeI == startEdgeI)
429  {
430  break;
431  }
432 
433  // Mark with current value. If in clean mode also remove feature edge.
434 
435  if (mark)
436  {
437  featVisited[edgeI] = currentFeatI;
438  }
439  else
440  {
441  featVisited[edgeI] = -2;
442  }
443 
444 
445  // Step to next vertex on edge
446  const edge& e = surf_.edges()[edgeI];
447 
448  vertI = e.otherVertex(vertI);
449 
450 
451  // Update cumulative length
452  visitedLength += e.mag(surf_.localPoints());
453 
454  nVisited++;
455 
456  if (nVisited > surf_.nEdges())
457  {
458  Warning<< "walkSegment : reached iteration limit in walking "
459  << "feature edges on surface from edge:" << startEdgeI
460  << " vertex:" << startPointi << nl
461  << "Returning with large length" << endl;
462 
463  return labelScalar(nVisited, GREAT);
464  }
465  }
466  while (true);
467 
468  return labelScalar(nVisited, visitedLength);
469 }
470 
471 
472 //- Divide into multiple normal bins
473 // - return REGION if != 2 normals
474 // - return REGION if 2 normals that make feature angle
475 // - otherwise return NONE and set normals,bins
476 //
477 // This has been relocated from surfaceFeatureExtract and could be cleaned
478 // up some more.
479 //
481 Foam::surfaceFeatures::surfaceFeatures::checkFlatRegionEdge
482 (
483  const scalar tol,
484  const scalar includedAngle,
485  const label edgeI
486 ) const
487 {
488  const triSurface& surf = surf_;
489 
490  const edge& e = surf.edges()[edgeI];
491  const labelList& eFaces = surf.edgeFaces()[edgeI];
492 
493  // Bin according to normal
494 
495  DynamicList<vector> normals(2);
496  DynamicList<labelList> bins(2);
497 
498  forAll(eFaces, eFacei)
499  {
500  const vector& n = surf.faceNormals()[eFaces[eFacei]];
501 
502  // Find the normal in normals
503  label index = -1;
504  forAll(normals, normalI)
505  {
506  if (mag(n & normals[normalI]) > (1-tol))
507  {
508  index = normalI;
509  break;
510  }
511  }
512 
513  if (index != -1)
514  {
515  bins[index].append(eFacei);
516  }
517  else if (normals.size() >= 2)
518  {
519  // Would be third normal. Mark as feature.
520  //Pout<< "** at edge:" << surf.localPoints()[e[0]]
521  // << surf.localPoints()[e[1]]
522  // << " have normals:" << normals
523  // << " and " << n << endl;
525  }
526  else
527  {
528  normals.append(n);
529  bins.append(labelList(1, eFacei));
530  }
531  }
532 
533 
534  // Check resulting number of bins
535  if (bins.size() == 1)
536  {
537  // Note: should check here whether they are two sets of faces
538  // that are planar or indeed 4 faces al coming together at an edge.
539  //Pout<< "** at edge:"
540  // << surf.localPoints()[e[0]]
541  // << surf.localPoints()[e[1]]
542  // << " have single normal:" << normals[0]
543  // << endl;
544  return surfaceFeatures::NONE;
545  }
546  else
547  {
548  // Two bins. Check if normals make an angle
549 
550  //Pout<< "** at edge:"
551  // << surf.localPoints()[e[0]]
552  // << surf.localPoints()[e[1]] << nl
553  // << " normals:" << normals << nl
554  // << " bins :" << bins << nl
555  // << endl;
556 
557  if (includedAngle >= 0)
558  {
559  scalar minCos = Foam::cos(degToRad(180.0 - includedAngle));
560 
561  forAll(eFaces, i)
562  {
563  const vector& ni = surf.faceNormals()[eFaces[i]];
564  for (label j=i+1; j<eFaces.size(); j++)
565  {
566  const vector& nj = surf.faceNormals()[eFaces[j]];
567  if (mag(ni & nj) < minCos)
568  {
569  //Pout<< "have sharp feature between normal:" << ni
570  // << " and " << nj << endl;
571 
572  // Is feature. Keep as region or convert to
573  // feature angle? For now keep as region.
575  }
576  }
577  }
578  }
579 
580  // So now we have two normals bins but need to make sure both
581  // bins have the same regions in it.
582 
583  // 1. store + or - region number depending
584  // on orientation of triangle in bins[0]
585  const labelList& bin0 = bins[0];
586  labelList regionAndNormal(bin0.size());
587  forAll(bin0, i)
588  {
589  const labelledTri& t = surf.localFaces()[eFaces[bin0[i]]];
590  const auto dir = t.edgeDirection(e);
591 
592  if (dir > 0)
593  {
594  regionAndNormal[i] = t.region()+1;
595  }
596  else if (dir == 0)
597  {
599  << exit(FatalError);
600  }
601  else
602  {
603  regionAndNormal[i] = -(t.region()+1);
604  }
605  }
606 
607  // 2. check against bin1
608  const labelList& bin1 = bins[1];
609  labelList regionAndNormal1(bin1.size());
610  forAll(bin1, i)
611  {
612  const labelledTri& t = surf.localFaces()[eFaces[bin1[i]]];
613  const auto dir = t.edgeDirection(e);
614 
615  label myRegionAndNormal;
616  if (dir > 0)
617  {
618  myRegionAndNormal = t.region()+1;
619  }
620  else
621  {
622  myRegionAndNormal = -(t.region()+1);
623  }
624 
625  regionAndNormal1[i] = myRegionAndNormal;
626 
627  label index = regionAndNormal.find(-myRegionAndNormal);
628  if (index == -1)
629  {
630  // Not found.
631  //Pout<< "cannot find region " << myRegionAndNormal
632  // << " in regions " << regionAndNormal << endl;
633 
635  }
636  }
637 
638  // Passed all checks, two normal bins with the same contents.
639  //Pout<< "regionAndNormal:" << regionAndNormal << endl;
640  //Pout<< "myRegionAndNormal:" << regionAndNormal1 << endl;
641  }
642 
643  return surfaceFeatures::NONE;
644 }
646 
647 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
648 
650 :
651  surf_(surf),
652  featurePoints_(0),
653  featureEdges_(0),
654  externalStart_(0),
655  internalStart_(0)
656 {}
658 
659 // Construct from components.
661 (
662  const triSurface& surf,
663  const labelList& featurePoints,
664  const labelList& featureEdges,
665  const label externalStart,
666  const label internalStart
667 )
668 :
669  surf_(surf),
670  featurePoints_(featurePoints),
671  featureEdges_(featureEdges),
672  externalStart_(externalStart),
673  internalStart_(externalStart)
674 {}
676 
677 // Construct from surface, angle and min length
679 (
680  const triSurface& surf,
681  const scalar includedAngle,
682  const scalar minLen,
683  const label minElems,
684  const bool geometricTestOnly
685 )
686 :
687  surf_(surf),
688  featurePoints_(0),
689  featureEdges_(0),
690  externalStart_(0),
691  internalStart_(0)
692 {
693  findFeatures(includedAngle, geometricTestOnly);
694 
695  if (minLen > 0 || minElems > 0)
696  {
697  trimFeatures(minLen, minElems, includedAngle);
698  }
699 }
700 
701 
703 (
704  const triSurface& surf,
705  const dictionary& featInfoDict
706 )
707 :
708  surf_(surf),
709  featurePoints_(featInfoDict.lookup("featurePoints")),
710  featureEdges_(featInfoDict.lookup("featureEdges")),
711  externalStart_(featInfoDict.get<label>("externalStart")),
712  internalStart_(featInfoDict.get<label>("internalStart"))
713 {}
714 
715 
717 (
718  const triSurface& surf,
719  const fileName& fName
720 )
721 :
722  surf_(surf),
723  featurePoints_(0),
724  featureEdges_(0),
725  externalStart_(0),
726  internalStart_(0)
727 {
728  IFstream str(fName);
729 
730  dictionary featInfoDict(str);
731 
732  featInfoDict.readEntry("featureEdges", featureEdges_);
733  featInfoDict.readEntry("featurePoints", featurePoints_);
734  featInfoDict.readEntry("externalStart", externalStart_);
735  featInfoDict.readEntry("internalStart", internalStart_);
736 }
737 
738 
740 (
741  const triSurface& surf,
742  const pointField& points,
743  const edgeList& edges,
744  const scalar mergeTol,
745  const bool geometricTestOnly
746 )
747 :
748  surf_(surf),
749  featurePoints_(0),
750  featureEdges_(0),
751  externalStart_(0),
752  internalStart_(0)
753 {
754  // Match edge mesh edges with the triSurface edges
755 
756  const labelListList& surfEdgeFaces = surf_.edgeFaces();
757  const edgeList& surfEdges = surf_.edges();
758 
759  scalar mergeTolSqr = sqr(mergeTol);
760 
761  EdgeMap<label> dynFeatEdges(2*edges.size());
762  DynamicList<labelList> dynFeatureEdgeFaces(edges.size());
763 
764  labelList edgeLabel;
765 
767  (
768  edges,
769  points,
770  mergeTolSqr,
771  edgeLabel // label of surface edge or -1
772  );
773 
774  label count = 0;
775  forAll(edgeLabel, sEdgeI)
776  {
777  const label sEdge = edgeLabel[sEdgeI];
778 
779  if (sEdge == -1)
780  {
781  continue;
782  }
783 
784  dynFeatEdges.insert(surfEdges[sEdge], count++);
785  dynFeatureEdgeFaces.append(surfEdgeFaces[sEdge]);
786  }
787 
788  // Find whether an edge is external or internal
789  List<edgeStatus> edgeStat(dynFeatEdges.size(), NONE);
790 
791  classifyFeatureAngles
792  (
793  dynFeatureEdgeFaces,
794  edgeStat,
795  GREAT,
796  geometricTestOnly
797  );
798 
799  // Transfer the edge status to a list encompassing all edges in the surface
800  // so that calcFeatPoints can be used.
801  List<edgeStatus> allEdgeStat(surf_.nEdges(), NONE);
802 
803  forAll(allEdgeStat, eI)
804  {
805  const auto iter = dynFeatEdges.cfind(surfEdges[eI]);
806 
807  if (iter.good())
808  {
809  allEdgeStat[eI] = edgeStat[iter.val()];
810  }
811  }
812 
813  edgeStat.clear();
814  dynFeatEdges.clear();
815 
816  setFromStatus(allEdgeStat, GREAT);
817 }
818 
819 
821 :
822  surf_(sf.surface()),
823  featurePoints_(sf.featurePoints()),
824  featureEdges_(sf.featureEdges()),
825  externalStart_(sf.externalStart()),
826  internalStart_(sf.internalStart())
827 {}
828 
830 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
831 
833 (
834  const bool regionEdges,
835  const bool externalEdges,
836  const bool internalEdges
837 ) const
838 {
839  DynamicList<label> selectedEdges;
840 
841  if (regionEdges)
842  {
843  selectedEdges.setCapacity(selectedEdges.size() + nRegionEdges());
844 
845  for (label i = 0; i < externalStart_; i++)
846  {
847  selectedEdges.append(featureEdges_[i]);
848  }
849  }
850 
851  if (externalEdges)
852  {
853  selectedEdges.setCapacity(selectedEdges.size() + nExternalEdges());
854 
855  for (label i = externalStart_; i < internalStart_; i++)
856  {
857  selectedEdges.append(featureEdges_[i]);
858  }
859  }
860 
861  if (internalEdges)
862  {
863  selectedEdges.setCapacity(selectedEdges.size() + nInternalEdges());
864 
865  for (label i = internalStart_; i < featureEdges_.size(); i++)
866  {
867  selectedEdges.append(featureEdges_[i]);
868  }
869  }
870 
871  return selectedEdges.shrink();
872 }
873 
874 
876 (
877  const scalar includedAngle,
878  const bool geometricTestOnly
879 )
880 {
881  scalar minCos = Foam::cos(degToRad(180.0 - includedAngle));
882 
883  // Per edge whether is feature edge.
884  List<edgeStatus> edgeStat(surf_.nEdges(), NONE);
885 
886  classifyFeatureAngles
887  (
888  surf_.edgeFaces(),
889  edgeStat,
890  minCos,
891  geometricTestOnly
892  );
893 
894  setFromStatus(edgeStat, includedAngle);
895 }
896 
898 // Remove small strings of edges. First assigns string number to
899 // every edge and then works out the length of them.
901 (
902  const scalar minLen,
903  const label minElems,
904  const scalar includedAngle
905 )
906 {
907  // Convert feature edge list to status per edge.
908  List<edgeStatus> edgeStat(toStatus());
909 
910 
911  // Mark feature edges according to connected lines.
912  // -1: unassigned
913  // -2: part of too small a feature line
914  // >0: feature line number
915  labelList featLines(surf_.nEdges(), -1);
916 
917  // Current featureline number.
918  label featI = 0;
919 
920  // Current starting edge
921  label startEdgeI = 0;
922 
923  do
924  {
925  // Find unset featureline
926  for (; startEdgeI < edgeStat.size(); startEdgeI++)
927  {
928  if
929  (
930  edgeStat[startEdgeI] != NONE
931  && featLines[startEdgeI] == -1
932  )
933  {
934  // Found unassigned feature edge.
935  break;
936  }
937  }
938 
939  if (startEdgeI == edgeStat.size())
940  {
941  // No unset feature edge found. All feature edges have line number
942  // assigned.
943  break;
944  }
945 
946  featLines[startEdgeI] = featI;
947 
948  const edge& startEdge = surf_.edges()[startEdgeI];
949 
950  // Walk 'left' and 'right' from startEdge.
951  labelScalar leftPath =
952  walkSegment
953  (
954  true, // 'mark' mode
955  edgeStat,
956  startEdgeI, // edge, not yet assigned to a featureLine
957  startEdge[0], // start of edge
958  featI, // Mark value
959  featLines // Mark for all edges
960  );
961 
962  labelScalar rightPath =
963  walkSegment
964  (
965  true,
966  edgeStat,
967  startEdgeI,
968  startEdge[1], // end of edge
969  featI,
970  featLines
971  );
972 
973  if
974  (
975  (
976  leftPath.len_
977  + rightPath.len_
978  + startEdge.mag(surf_.localPoints())
979  < minLen
980  )
981  || (leftPath.n_ + rightPath.n_ + 1 < minElems)
982  )
983  {
984  // Rewalk same route (recognizable by featLines == featI)
985  // to reset featLines.
986 
987  featLines[startEdgeI] = -2;
988 
989  walkSegment
990  (
991  false, // 'clean' mode
992  edgeStat,
993  startEdgeI, // edge, not yet assigned to a featureLine
994  startEdge[0], // start of edge
995  featI, // Unset value
996  featLines // Mark for all edges
997  );
998 
999  walkSegment
1000  (
1001  false,
1002  edgeStat,
1003  startEdgeI,
1004  startEdge[1], // end of edge
1005  featI,
1006  featLines
1007  );
1008  }
1009  else
1010  {
1011  featI++;
1012  }
1013  }
1014  while (true);
1015 
1016  // Unmark all feature lines that have featLines=-2
1017  forAll(featureEdges_, i)
1018  {
1019  label edgeI = featureEdges_[i];
1020 
1021  if (featLines[edgeI] == -2)
1022  {
1023  edgeStat[edgeI] = NONE;
1024  }
1025  }
1026 
1027  // Convert back to edge labels
1028  setFromStatus(edgeStat, includedAngle);
1029 
1030  return featLines;
1032 
1033 
1035 (
1036  List<edgeStatus>& edgeStat,
1037  const treeBoundBox& bb
1038 ) const
1039 {
1040  deleteBox(edgeStat, bb, true);
1042 
1043 
1045 (
1046  List<edgeStatus>& edgeStat,
1047  const treeBoundBox& bb
1048 ) const
1049 {
1050  deleteBox(edgeStat, bb, false);
1052 
1053 
1055 (
1056  List<edgeStatus>& edgeStat,
1057  const treeBoundBox& bb,
1058  const bool removeInside
1059 ) const
1060 {
1061  const edgeList& surfEdges = surf_.edges();
1062  const pointField& surfLocalPoints = surf_.localPoints();
1063 
1064  forAll(edgeStat, edgei)
1065  {
1066  const point eMid = surfEdges[edgei].centre(surfLocalPoints);
1067 
1068  if (removeInside ? bb.contains(eMid) : !bb.contains(eMid))
1069  {
1070  edgeStat[edgei] = surfaceFeatures::NONE;
1071  }
1072  }
1074 
1075 
1077 (
1078  List<edgeStatus>& edgeStat,
1079  const plane& cutPlane
1080 ) const
1081 {
1082  const edgeList& surfEdges = surf_.edges();
1083  const pointField& pts = surf_.points();
1084  const labelList& meshPoints = surf_.meshPoints();
1085 
1086  forAll(edgeStat, edgei)
1087  {
1088  const edge& e = surfEdges[edgei];
1089 
1090  const point& p0 = pts[meshPoints[e.start()]];
1091  const point& p1 = pts[meshPoints[e.end()]];
1092  const linePointRef line(p0, p1);
1093 
1094  // If edge does not intersect the plane, delete.
1095  scalar intersect = cutPlane.lineIntersect(line);
1096 
1097  point featPoint = intersect * (p1 - p0) + p0;
1098 
1099  if (!onLine(featPoint, line))
1100  {
1101  edgeStat[edgei] = surfaceFeatures::NONE;
1102  }
1103  }
1105 
1106 
1108 (
1109  List<edgeStatus>& edgeStat
1110 ) const
1111 {
1112  forAll(edgeStat, edgei)
1113  {
1114  if (surf_.edgeFaces()[edgei].size() == 1)
1115  {
1116  edgeStat[edgei] = surfaceFeatures::NONE;
1117  }
1118  }
1119 }
1120 
1121 
1122 //- Divide into multiple normal bins
1123 // - return REGION if != 2 normals
1124 // - return REGION if 2 normals that make feature angle
1125 // - otherwise return NONE and set normals,bins
1126 void Foam::surfaceFeatures::checkFlatRegionEdge
1127 (
1128  List<edgeStatus>& edgeStat,
1129  const scalar tol,
1130  const scalar includedAngle
1131 ) const
1132 {
1133  forAll(edgeStat, edgei)
1134  {
1135  if (edgeStat[edgei] == surfaceFeatures::REGION)
1136  {
1137  const labelList& eFaces = surf_.edgeFaces()[edgei];
1138 
1139  if (eFaces.size() > 2 && (eFaces.size() % 2) == 0)
1140  {
1141  edgeStat[edgei] = checkFlatRegionEdge
1142  (
1143  tol,
1144  includedAngle,
1145  edgei
1146  );
1147  }
1148  }
1149  }
1150 }
1151 
1154 {
1155  dictionary featInfoDict;
1156  featInfoDict.add("externalStart", externalStart_);
1157  featInfoDict.add("internalStart", internalStart_);
1158  featInfoDict.add("featureEdges", featureEdges_);
1159  featInfoDict.add("featurePoints", featurePoints_);
1160 
1161  featInfoDict.write(os);
1162 }
1163 
1165 void Foam::surfaceFeatures::write(const fileName& fName) const
1166 {
1167  OFstream os(fName);
1168  writeDict(os);
1169 }
1170 
1172 void Foam::surfaceFeatures::writeObj(const fileName& prefix) const
1173 {
1174  OFstream regionStr(prefix + "_regionEdges.obj");
1175  Pout<< "Writing region edges to " << regionStr.name() << endl;
1176 
1177  label verti = 0;
1178  for (label i = 0; i < externalStart_; i++)
1179  {
1180  const edge& e = surf_.edges()[featureEdges_[i]];
1181 
1182  meshTools::writeOBJ(regionStr, surf_.localPoints()[e[0]]); verti++;
1183  meshTools::writeOBJ(regionStr, surf_.localPoints()[e[1]]); verti++;
1184  regionStr << "l " << verti-1 << ' ' << verti << endl;
1185  }
1186 
1187 
1188  OFstream externalStr(prefix + "_externalEdges.obj");
1189  Pout<< "Writing external edges to " << externalStr.name() << endl;
1190 
1191  verti = 0;
1192  for (label i = externalStart_; i < internalStart_; i++)
1193  {
1194  const edge& e = surf_.edges()[featureEdges_[i]];
1195 
1196  meshTools::writeOBJ(externalStr, surf_.localPoints()[e[0]]); verti++;
1197  meshTools::writeOBJ(externalStr, surf_.localPoints()[e[1]]); verti++;
1198  externalStr << "l " << verti-1 << ' ' << verti << endl;
1199  }
1200 
1201  OFstream internalStr(prefix + "_internalEdges.obj");
1202  Pout<< "Writing internal edges to " << internalStr.name() << endl;
1203 
1204  verti = 0;
1205  for (label i = internalStart_; i < featureEdges_.size(); i++)
1206  {
1207  const edge& e = surf_.edges()[featureEdges_[i]];
1208 
1209  meshTools::writeOBJ(internalStr, surf_.localPoints()[e[0]]); verti++;
1210  meshTools::writeOBJ(internalStr, surf_.localPoints()[e[1]]); verti++;
1211  internalStr << "l " << verti-1 << ' ' << verti << endl;
1212  }
1213 
1214  OFstream pointStr(prefix + "_points.obj");
1215  Pout<< "Writing feature points to " << pointStr.name() << endl;
1216 
1217  for (const label pointi : featurePoints_)
1218  {
1219  meshTools::writeOBJ(pointStr, surf_.localPoints()[pointi]);
1220  }
1221 }
1222 
1225 {
1226  os << "Feature set:" << nl
1227  << " points : " << this->featurePoints().size() << nl
1228  << " edges : " << this->featureEdges().size() << nl
1229  << " of which" << nl
1230  << " region edges : " << this->nRegionEdges() << nl
1231  << " external edges : " << this->nExternalEdges() << nl
1232  << " internal edges : " << this->nInternalEdges() << endl;
1233 }
1234 
1235 
1236 // Get nearest vertex on patch for every point of surf in pointSet.
1238 (
1239  const labelList& pointLabels,
1240  const pointField& samples,
1241  const scalarField& maxDistSqr
1242 ) const
1243 {
1244  // Build tree out of all samples.
1245 
1246  // Define bound box here (gcc-4.8.5)
1247  const treeBoundBox overallBb(samples);
1248 
1250  (
1252  overallBb,
1253  8, // maxLevel
1254  10, // leafsize
1255  3.0 // duplicity
1256  );
1257  const auto& treeData = ppTree.shapes();
1258 
1259  // From patch point to surface point
1260  Map<label> nearest(2*pointLabels.size());
1261 
1262  const pointField& surfPoints = surf_.localPoints();
1263 
1264  forAll(pointLabels, i)
1265  {
1266  const label surfPointi = pointLabels[i];
1267 
1268  const point& surfPt = surfPoints[surfPointi];
1269 
1270  pointIndexHit info = ppTree.findNearest
1271  (
1272  surfPt,
1273  maxDistSqr[i]
1274  );
1275 
1276  if (!info.hit())
1277  {
1279  << "Problem for point "
1280  << surfPointi << " in tree " << ppTree.bb()
1281  << abort(FatalError);
1282  }
1283 
1284  label sampleI = info.index();
1285 
1286  if (treeData.centre(sampleI).distSqr(surfPt) < maxDistSqr[sampleI])
1287  {
1288  nearest.insert(sampleI, surfPointi);
1289  }
1290  }
1291 
1292 
1293  if (debug)
1294  {
1295  //
1296  // Dump to obj file
1297  //
1298 
1299  Pout<< "Dumping nearest surface feature points to nearestSamples.obj"
1300  << endl;
1301 
1302  OFstream objStream("nearestSamples.obj");
1303 
1304  label vertI = 0;
1305  forAllConstIters(nearest, iter)
1306  {
1307  meshTools::writeOBJ(objStream, samples[iter.key()]); vertI++;
1308  meshTools::writeOBJ(objStream, surfPoints[iter.val()]); vertI++;
1309  objStream<< "l " << vertI-1 << ' ' << vertI << endl;
1310  }
1311  }
1312 
1313  return nearest;
1314 }
1315 
1316 
1317 // Get nearest sample point for regularly sampled points along
1318 // selected edges. Return map from sample to edge label
1320 (
1321  const labelList& selectedEdges,
1322  const pointField& samples,
1323  const scalarField& sampleDist,
1324  const scalarField& maxDistSqr,
1325  const scalar minSampleDist
1326 ) const
1327 {
1328  const pointField& surfPoints = surf_.localPoints();
1329  const edgeList& surfEdges = surf_.edges();
1330 
1331  scalar maxSearchSqr = max(maxDistSqr);
1332 
1333  // Define bound box here (gcc-4.8.5)
1334  const treeBoundBox overallBb(samples);
1335 
1337  (
1339  overallBb,
1340  8, // maxLevel
1341  10, // leafsize
1342  3.0 // duplicity
1343  );
1344 
1345  // From patch point to surface edge.
1346  Map<label> nearest(2*selectedEdges.size());
1347 
1348  forAll(selectedEdges, i)
1349  {
1350  label surfEdgeI = selectedEdges[i];
1351 
1352  const edge& e = surfEdges[surfEdgeI];
1353 
1354  if (debug && (i % 1000) == 0)
1355  {
1356  Pout<< "looking at surface feature edge " << surfEdgeI
1357  << " verts:" << e << " points:" << surfPoints[e[0]]
1358  << ' ' << surfPoints[e[1]] << endl;
1359  }
1360 
1361  // Normalized edge vector
1362  vector eVec = e.vec(surfPoints);
1363  scalar eMag = mag(eVec);
1364  eVec /= eMag;
1365 
1366 
1367  //
1368  // Sample along edge
1369  //
1370 
1371  bool exit = false;
1372 
1373  // Coordinate along edge (0 .. eMag)
1374  scalar s = 0.0;
1375 
1376  while (true)
1377  {
1378  point edgePoint(surfPoints[e.start()] + s*eVec);
1379 
1380  pointIndexHit info = ppTree.findNearest
1381  (
1382  edgePoint,
1383  maxSearchSqr
1384  );
1385 
1386  if (!info.hit())
1387  {
1388  // No point close enough to surface edge.
1389  break;
1390  }
1391 
1392  label sampleI = info.index();
1393 
1394  if (info.point().distSqr(edgePoint) < maxDistSqr[sampleI])
1395  {
1396  nearest.insert(sampleI, surfEdgeI);
1397  }
1398 
1399  if (exit)
1400  {
1401  break;
1402  }
1403 
1404  // Step to next sample point using local distance.
1405  // Truncate to max 1/minSampleDist samples per feature edge.
1406  s += max(minSampleDist*eMag, sampleDist[sampleI]);
1407 
1408  if (s >= (1-minSampleDist)*eMag)
1409  {
1410  // Do one more sample, at edge endpoint
1411  s = eMag;
1412  exit = true;
1413  }
1414  }
1415  }
1416 
1417 
1418 
1419  if (debug)
1420  {
1421  // Dump to obj file
1422 
1423  Pout<< "Dumping nearest surface edges to nearestEdges.obj"
1424  << endl;
1425 
1426  OFstream objStream("nearestEdges.obj");
1427 
1428  label vertI = 0;
1429  forAllConstIters(nearest, iter)
1430  {
1431  const label sampleI = iter.key();
1432 
1433  const edge& e = surfEdges[iter.val()];
1434 
1435  meshTools::writeOBJ(objStream, samples[sampleI]); vertI++;
1436 
1437  point nearPt =
1438  e.line(surfPoints).nearestDist(samples[sampleI]).point();
1439 
1440  meshTools::writeOBJ(objStream, nearPt); vertI++;
1441 
1442  objStream<< "l " << vertI-1 << ' ' << vertI << endl;
1443  }
1444  }
1445 
1446  return nearest;
1447 }
1448 
1449 
1450 // Get nearest edge on patch for regularly sampled points along the feature
1451 // edges. Return map from patch edge to feature edge.
1452 //
1453 // Q: using point-based sampleDist and maxDist (distance to look around
1454 // each point). Should they be edge-based e.g. by averaging or max()?
1456 (
1457  const labelList& selectedEdges,
1458  const edgeList& sampleEdges,
1459  const labelList& selectedSampleEdges,
1460  const pointField& samplePoints,
1461  const scalarField& sampleDist,
1462  const scalarField& maxDistSqr,
1463  const scalar minSampleDist
1464 ) const
1465 {
1466  // Build tree out of selected sample edges.
1468  (
1469  treeDataEdge
1470  (
1471  sampleEdges,
1472  samplePoints,
1473  selectedSampleEdges
1474  ),
1475  treeBoundBox(samplePoints), // overall search domain
1476  8, // maxLevel
1477  10, // leafsize
1478  3.0 // duplicity
1479  );
1480 
1481  const pointField& surfPoints = surf_.localPoints();
1482  const edgeList& surfEdges = surf_.edges();
1483 
1484  scalar maxSearchSqr = max(maxDistSqr);
1485 
1486  Map<pointIndexHit> nearest(2*sampleEdges.size());
1487 
1488  //
1489  // Loop over all selected edges. Sample at regular intervals. Find nearest
1490  // sampleEdges (using octree)
1491  //
1492 
1493  forAll(selectedEdges, i)
1494  {
1495  label surfEdgeI = selectedEdges[i];
1496 
1497  const edge& e = surfEdges[surfEdgeI];
1498 
1499  if (debug && (i % 1000) == 0)
1500  {
1501  Pout<< "looking at surface feature edge " << surfEdgeI
1502  << " verts:" << e << " points:" << surfPoints[e[0]]
1503  << ' ' << surfPoints[e[1]] << endl;
1504  }
1505 
1506  // Normalized edge vector
1507  vector eVec = e.vec(surfPoints);
1508  scalar eMag = mag(eVec);
1509  eVec /= eMag;
1510 
1511 
1512  //
1513  // Sample along edge
1514  //
1515 
1516  bool exit = false;
1517 
1518  // Coordinate along edge (0 .. eMag)
1519  scalar s = 0.0;
1520 
1521  while (true)
1522  {
1523  point edgePoint(surfPoints[e.start()] + s*eVec);
1524 
1525  pointIndexHit info = ppTree.findNearest
1526  (
1527  edgePoint,
1528  maxSearchSqr
1529  );
1530 
1531  if (!info.hit())
1532  {
1533  // No edge close enough to surface edge.
1534  break;
1535  }
1536 
1537  label index = info.index();
1538 
1539  label sampleEdgeI = ppTree.shapes().objectIndex(index);
1540 
1541  const edge& e = sampleEdges[sampleEdgeI];
1542 
1543  if (info.point().distSqr(edgePoint) < maxDistSqr[e.start()])
1544  {
1545  nearest.insert
1546  (
1547  sampleEdgeI,
1548  pointIndexHit(true, edgePoint, surfEdgeI)
1549  );
1550  }
1551 
1552  if (exit)
1553  {
1554  break;
1555  }
1556 
1557  // Step to next sample point using local distance.
1558  // Truncate to max 1/minSampleDist samples per feature edge.
1559  // s += max(minSampleDist*eMag, sampleDist[e.start()]);
1560  s += 0.01*eMag;
1561 
1562  if (s >= (1-minSampleDist)*eMag)
1563  {
1564  // Do one more sample, at edge endpoint
1565  s = eMag;
1566  exit = true;
1567  }
1568  }
1569  }
1570 
1571 
1572  if (debug)
1573  {
1574  // Dump to obj file
1575 
1576  Pout<< "Dumping nearest surface feature edges to nearestEdges.obj"
1577  << endl;
1578 
1579  OFstream objStream("nearestEdges.obj");
1580 
1581  label vertI = 0;
1582  forAllConstIters(nearest, iter)
1583  {
1584  const label sampleEdgeI = iter.key();
1585 
1586  const edge& sampleEdge = sampleEdges[sampleEdgeI];
1587 
1588  // Write line from edgeMid to point on feature edge
1589  meshTools::writeOBJ(objStream, sampleEdge.centre(samplePoints));
1590  vertI++;
1591 
1592  meshTools::writeOBJ(objStream, iter.val().point());
1593  vertI++;
1594 
1595  objStream<< "l " << vertI-1 << ' ' << vertI << endl;
1596  }
1597  }
1598 
1599  return nearest;
1600 }
1601 
1602 
1603 // Get nearest surface edge for every sample. Return in form of
1604 // labelLists giving surfaceEdge label&intersectionpoint.
1606 (
1607  const labelList& selectedEdges,
1608  const pointField& samples,
1609  scalar searchSpanSqr, // Search span
1610  labelList& edgeLabel,
1611  labelList& edgeEndPoint,
1612  pointField& edgePoint
1613 ) const
1614 {
1615  edgeLabel.setSize(samples.size());
1616  edgeEndPoint.setSize(samples.size());
1617  edgePoint.setSize(samples.size());
1618 
1619  const pointField& localPoints = surf_.localPoints();
1620 
1621  treeBoundBox searchDomain(localPoints);
1622  searchDomain.inflate(0.1);
1623 
1625  (
1626  treeDataEdge
1627  (
1628  surf_.edges(),
1629  localPoints,
1630  selectedEdges
1631  ),
1632  searchDomain, // overall search domain
1633  8, // maxLevel
1634  10, // leafsize
1635  3.0 // duplicity
1636  );
1637  const auto& treeData = ppTree.shapes();
1638 
1639  forAll(samples, i)
1640  {
1641  const point& sample = samples[i];
1642 
1643  pointIndexHit info = ppTree.findNearest
1644  (
1645  sample,
1646  searchSpanSqr
1647  );
1648 
1649  if (!info.hit())
1650  {
1651  edgeLabel[i] = -1;
1652  }
1653  else
1654  {
1655  // Need to recalculate to classify edgeEndPoint
1656  pointIndexHit pHit = edgeNearest
1657  (
1658  treeData.line(info.index()),
1659  sample
1660  );
1661 
1662  edgeLabel[i] = treeData.objectIndex(info.index());
1663  edgePoint[i] = pHit.point();
1664  edgeEndPoint[i] = pHit.index();
1665  }
1666  }
1667 }
1668 
1669 
1670 // Get nearest point on nearest feature edge for every sample (is edge)
1672 (
1673  const labelList& selectedEdges,
1674 
1675  const edgeList& sampleEdges,
1676  const labelList& selectedSampleEdges,
1677  const pointField& samplePoints,
1678 
1679  const vector& searchSpan, // Search span
1680  labelList& edgeLabel, // label of surface edge or -1
1681  pointField& pointOnEdge, // point on above edge
1682  pointField& pointOnFeature // point on sample edge
1683 ) const
1684 {
1685  edgeLabel.setSize(selectedSampleEdges.size());
1686  pointOnEdge.setSize(selectedSampleEdges.size());
1687  pointOnFeature.setSize(selectedSampleEdges.size());
1688 
1689  treeBoundBox searchDomain(surf_.localPoints());
1690 
1692  (
1693  treeDataEdge
1694  (
1695  surf_.edges(),
1696  surf_.localPoints(),
1697  selectedEdges
1698  ),
1699  searchDomain, // overall search domain
1700  8, // maxLevel
1701  10, // leafsize
1702  3.0 // duplicity
1703  );
1704  const auto& treeData = ppTree.shapes();
1705 
1706  forAll(selectedSampleEdges, i)
1707  {
1708  const edge& e = sampleEdges[selectedSampleEdges[i]];
1709 
1710  linePointRef edgeLine = e.line(samplePoints);
1711 
1712  point eMid(edgeLine.centre());
1713 
1714  treeBoundBox tightest(eMid - searchSpan, eMid + searchSpan);
1715 
1716  pointIndexHit info = ppTree.findNearest
1717  (
1718  edgeLine,
1719  tightest,
1720  pointOnEdge[i]
1721  );
1722 
1723  if (!info.hit())
1724  {
1725  edgeLabel[i] = -1;
1726  }
1727  else
1728  {
1729  edgeLabel[i] = treeData.objectIndex(info.index());
1730  pointOnFeature[i] = info.point();
1731  }
1732  }
1733 }
1734 
1735 
1737 (
1738  const edgeList& edges,
1739  const pointField& points,
1740  scalar searchSpanSqr, // Search span
1741  labelList& edgeLabel
1742 ) const
1743 {
1744  edgeLabel = labelList(surf_.nEdges(), -1);
1745 
1746  treeBoundBox searchDomain(points);
1747  searchDomain.inflate(0.1);
1748 
1750  (
1751  treeDataEdge(edges, points), // All edges
1752 
1753  searchDomain, // overall search domain
1754  8, // maxLevel
1755  10, // leafsize
1756  3.0 // duplicity
1757  );
1758  const auto& treeData = ppTree.shapes();
1759 
1760  const edgeList& surfEdges = surf_.edges();
1761  const pointField& surfLocalPoints = surf_.localPoints();
1762 
1763  forAll(surfEdges, edgeI)
1764  {
1765  const edge& sample = surfEdges[edgeI];
1766 
1767  const point& startPoint = surfLocalPoints[sample.start()];
1768  const point& midPoint = sample.centre(surfLocalPoints);
1769 
1770  pointIndexHit infoMid = ppTree.findNearest
1771  (
1772  midPoint,
1773  searchSpanSqr
1774  );
1775 
1776  if (infoMid.hit())
1777  {
1778  const vector surfEdgeDir = midPoint - startPoint;
1779 
1780  const vector featEdgeDir = treeData.line(infoMid.index()).vec();
1781 
1782  // Check that the edges are nearly parallel
1783  if (mag(surfEdgeDir ^ featEdgeDir) < parallelTolerance)
1784  {
1785  edgeLabel[edgeI] = edgeI;
1786  }
1787  }
1788  }
1789 }
1790 
1791 
1792 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1795 {
1796  if (this == &rhs)
1797  {
1798  return; // Self-assignment is a no-op
1799  }
1800 
1801  if (&surf_ != &rhs.surface())
1802  {
1804  << "Operating on different surfaces"
1805  << abort(FatalError);
1806  }
1807 
1808  featurePoints_ = rhs.featurePoints();
1809  featureEdges_ = rhs.featureEdges();
1810  externalStart_ = rhs.externalStart();
1811  internalStart_ = rhs.internalStart();
1812 }
1813 
1814 
1815 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
A line primitive.
Definition: line.H:52
A class for handling file names.
Definition: fileName.H:72
Holds (reference to) pointField. Encapsulation of data needed for octree searches. Used for searching for nearest point. No bounding boxes around points. Only overlaps and calcNearest are implemented, rest makes little sense.
Definition: treeDataPoint.H:58
void deleteBox(List< edgeStatus > &edgeStat, const treeBoundBox &bb, const bool removeInside) const
Mark edge status as &#39;NONE&#39; for edges inside/outside box.
scalar mag(const UList< point > &pts) const
The length (L2-norm) of the edge vector.
Definition: edgeI.H:419
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
labelList pointLabels(nPoints, -1)
void write(const fileName &fName) const
Write as dictionary to file.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#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
const Type & shapes() const noexcept
Reference to shape.
List< edge > edgeList
List of edge.
Definition: edgeList.H:32
virtual const fileName & name() const override
Read/write access to the name of the stream.
Definition: OSstream.H:128
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
void inflate(const scalar factor)
Expand box by factor*mag(span) in all dimensions.
Definition: boundBoxI.H:381
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Unit conversion functions.
Output to file stream, using an OSstream.
Definition: OFstream.H:49
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
scalarField samples(nIntervals, Zero)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
Definition: pointIndexHit.H:58
Holds data for octree to work on an edges subset.
Definition: treeDataEdge.H:53
label externalStart() const
Start of external edges.
Map< label > nearestSamples(const labelList &selectedPoints, const pointField &samples, const scalarField &maxDistSqr) const
Find nearest sample for selected surface points.
void subsetBox(List< edgeStatus > &edgeStat, const treeBoundBox &bb) const
Mark edge status outside box as &#39;NONE&#39;.
Map< pointIndexHit > nearestEdges(const labelList &selectedEdges, const edgeList &sampleEdges, const labelList &selectedSampleEdges, const pointField &samplePoints, const scalarField &sampleDist, const scalarField &maxDistSqr, const scalar minSampleDist=0.1) const
Like nearestSamples but now gets nearest point on.
void setCapacity(const label len)
Alter the size of the underlying storage.
Definition: DynamicListI.H:303
This class describes the interaction of an object (often a face) and a point. It carries the info of ...
Definition: pointIndexHit.H:44
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:625
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:196
const labelList & featureEdges() const
Return feature edge list.
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition: plane.H:90
label objectIndex(const label index) const
Map from shape index to original (non-subset) edge label.
Definition: treeDataEdge.H:388
Lookup type of boundary radiation properties.
Definition: lookup.H:57
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
void writeObj(const fileName &prefix) const
Write to separate OBJ files (region, external, internal edges,.
scalar lineIntersect(const line< PointType, PointRef > &l) const
Return the cutting point between the plane and a line passing through the supplied points...
Definition: plane.H:317
labelList selectFeatureEdges(const bool regionEdges, const bool externalEdges, const bool internalEdges) const
Helper function: select a subset of featureEdges_.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
List< edgeStatus > toStatus() const
From member feature edges to status per edge.
bool insert(const edge &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:152
label start() const noexcept
The start (first) vertex label.
Definition: edge.H:147
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
Not a classified feature edge.
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of &#39;true&#39; entries.
Definition: BitOps.H:73
const point_type & point() const noexcept
Return point, no checks.
labelList trimFeatures(const scalar minLen, const label minElems, const scalar includedAngle)
Delete small sets of edges. Edges are stringed up and any.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dimensionedScalar cos(const dimensionedScalar &ds)
line< point, const point & > linePointRef
A line using referred points.
Definition: line.H:66
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
Mid-point interpolation (weighting factors = 0.5) scheme class.
Definition: midPoint.H:51
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition: edge.H:59
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
void setFromStatus(const List< edgeStatus > &edgeStat, const scalar includedAngle)
Set from status per edge.
void operator=(const surfaceFeatures &rhs)
void writeStats(Ostream &os) const
Write some information.
const triSurface & surface() const
void excludeBox(List< edgeStatus > &edgeStat, const treeBoundBox &bb) const
Mark edge status inside box as &#39;NONE&#39;.
const labelListList & edgeFaces() const
Return edge-face addressing.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
Vector< scalar > vector
Definition: vector.H:57
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicList.H:584
const treeBoundBox & bb() const
Top bounding box.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
label index() const noexcept
Return the hit index.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
const wordList surface
Standard surface field types (scalar, vector, tensor, etc)
bool contains(const vector &dir, const point &) const
Contains point (inside or on edge) and moving in direction.
Definition: treeBoundBox.C:413
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
dimensionedScalar sin(const dimensionedScalar &ds)
label nEdges() const
Number of edges in patch.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
int debug
Static debugging option.
const labelList & featurePoints() const
Return feature point list.
DynamicList< T, SizeMin > & shrink()
Calls shrink_to_fit() and returns a reference to the DynamicList.
Definition: DynamicListI.H:447
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
Input from file stream, using an ISstream.
Definition: IFstream.H:49
const Vector< Cmpt > & centre(const Foam::UList< Vector< Cmpt >> &) const noexcept
Return this (for point which is a typedef to Vector<scalar>)
Definition: VectorI.H:67
messageStream Warning
Warning stream (stdout output on master, null elsewhere), with additional &#39;FOAM Warning&#39; header text...
Point centre() const
Return centre (centroid)
Definition: lineI.H:83
void subsetPlane(List< edgeStatus > &edgeStat, const plane &cutPlane) const
If edge does not intersect the plane, mark as &#39;NONE&#39;.
void write(Ostream &os, const bool subDict=true) const
Write dictionary, normally with sub-dictionary formatting.
Definition: dictionaryIO.C:204
bool hit() const noexcept
Is there a hit?
vector point
Point is a vector.
Definition: point.H:37
Non-pointer based hierarchical recursive searching.
void findFeatures(const scalar includedAngle, const bool geometricTestOnly)
Find feature edges using provided included angle.
void nearestSurfEdge(const labelList &selectedEdges, const pointField &samples, scalar searchSpanSqr, labelList &edgeLabel, labelList &edgeEndPoint, pointField &edgePoint) const
Find nearest surface edge (out of selectedEdges) for.
surfaceFeatures(const triSurface &surf)
Construct from surface.
void excludeOpen(List< edgeStatus > &edgeStat) const
Mark edges with a single connected face as &#39;NONE&#39;.
void nearestFeatEdge(const edgeList &edges, const pointField &points, scalar searchSpanSqr, labelList &edgeLabel) const
Find nearest feature edge to each surface edge. Uses the.
point centre(const UList< point > &pts) const
Return centre point (centroid) of the edge.
Definition: edgeI.H:372
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:90
void writeDict(Ostream &os) const
Write as dictionary.
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
label internalStart() const
Start of internal edges.
Triangulated surface description with patch information.
Definition: triSurface.H:71
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))
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
const volScalarField & p0
Definition: EEqn.H:36
Holds feature edges/points of surface.
PointHit< point > pointHit
A PointHit with a 3D point.
Definition: pointHit.H:43
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
A HashTable to objects of type <T> with a label key.
const pointField & pts