Loading...
Searching...
No Matches
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-2024 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "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"
37#include "edgeHashes.H"
38
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
41namespace Foam
42{
44
45 const scalar surfaceFeatures::parallelTolerance = sin(degToRad(1.0));
46
47
49// Check if the point is on the line
50static 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
74Foam::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 }
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;
137
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
206void 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
267void 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.
330Foam::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.
378Foam::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//
481Foam::surfaceFeatures::surfaceFeatures::checkFlatRegionEdge
482(
483 const scalar tol,
484 const scalar includedAngle,
485 const label edgeI,
486 const point& leftPoint
487) const
488{
489 const triSurface& surf = surf_;
490
491 const edge& e = surf.edges()[edgeI];
492 const labelList& eFaces = surf.edgeFaces()[edgeI];
493
494 // Bin according to normal
495
496 DynamicList<vector> normals(2);
498
499 forAll(eFaces, eFacei)
500 {
501 const vector& n = surf.faceNormals()[eFaces[eFacei]];
502
503 // Find the normal in normals
504 label index = -1;
505 forAll(normals, normalI)
506 {
507 if (mag(n & normals[normalI]) > (1-tol))
508 {
509 index = normalI;
510 break;
511 }
512 }
513
514 if (index != -1)
515 {
516 bins[index].append(eFacei);
517 }
518 else if (normals.size() >= 2)
519 {
520 // Would be third normal. Mark as feature.
521 //Pout<< "** at edge:" << surf.localPoints()[e[0]]
522 // << surf.localPoints()[e[1]]
523 // << " have normals:" << normals
524 // << " and " << n << endl;
526 }
527 else
528 {
529 normals.append(n);
530 bins.append(labelList(1, eFacei));
531 }
532 }
533
534
535 // Check resulting number of bins
536 if (bins.size() == 1)
537 {
538 // Note: should check here whether they are two sets of faces
539 // that are planar or indeed 4 faces al coming together at an edge.
540 //Pout<< "** at edge:"
541 // << surf.localPoints()[e[0]]
542 // << surf.localPoints()[e[1]]
543 // << " have single normal:" << normals[0]
544 // << endl;
546 }
547 else
548 {
549 // Two bins. Check if normals make an angle
550
551 //Pout<< "** at edge:"
552 // << surf.localPoints()[e[0]]
553 // << surf.localPoints()[e[1]] << nl
554 // << " normals:" << normals << nl
555 // << " bins :" << bins << nl
556 // << endl;
557
558 if (includedAngle >= 0)
559 {
560 scalar minCos = Foam::cos(degToRad(180.0 - includedAngle));
561
562 forAll(eFaces, i)
563 {
564 const vector& ni = surf.faceNormals()[eFaces[i]];
565 for (label j=i+1; j<eFaces.size(); j++)
566 {
567 const vector& nj = surf.faceNormals()[eFaces[j]];
568 if (mag(ni & nj) < minCos)
569 {
570 //Pout<< "have sharp feature between normal:" << ni
571 // << " and " << nj << endl;
572
573 // Is feature. Keep as region or convert to
574 // feature angle? For now keep as region.
576 }
577 }
578 }
579 }
580
581 // So now we have two normals bins but need to make sure both
582 // bins have the same regions in it.
583
584 // 1. store + or - region number depending
585 // on orientation of triangle in bins[0]
586 const labelList& bin0 = bins[0];
587 labelList regionAndNormal(bin0.size());
588 forAll(bin0, i)
589 {
590 const labelledTri& t = surf.localFaces()[eFaces[bin0[i]]];
591 const auto dir = t.edgeDirection(e);
592
593 if (dir > 0)
594 {
595 regionAndNormal[i] = t.region()+1;
596 }
597 else if (dir == 0)
598 {
600 << exit(FatalError);
601 }
602 else
603 {
604 regionAndNormal[i] = -(t.region()+1);
605 }
606 }
607
608 // 2. check against bin1
609 const labelList& bin1 = bins[1];
610 labelList regionAndNormal1(bin1.size());
611 forAll(bin1, i)
612 {
613 const labelledTri& t = surf.localFaces()[eFaces[bin1[i]]];
614 const auto dir = t.edgeDirection(e);
615
616 label myRegionAndNormal;
617 if (dir > 0)
618 {
619 myRegionAndNormal = t.region()+1;
620 }
621 else
622 {
623 myRegionAndNormal = -(t.region()+1);
624 }
625
626 regionAndNormal1[i] = myRegionAndNormal;
627
628 label index = regionAndNormal.find(-myRegionAndNormal);
629 if (index == -1)
630 {
631 // Not found.
632 //Pout<< "cannot find region " << myRegionAndNormal
633 // << " in regions " << regionAndNormal << endl;
634
636 }
637 }
638
639 // Passed all checks, two normal bins with the same contents.
640 //Pout<< "regionAndNormal:" << regionAndNormal << endl;
641 //Pout<< "myRegionAndNormal:" << regionAndNormal1 << endl;
642 }
643
645}
646
647/*
648
649// TBD. Big problem for duplicate triangles with opposing normals: we
650// don't know which one of the duplicates gets found on either side so
651// the normal might be + or -. Hence on e.g. motorBike.obj we see feature
652// lines where there shouldn't be
653Foam::surfaceFeatures::edgeStatus
654Foam::surfaceFeatures::surfaceFeatures::checkBinRegion
655(
656 const label edgei,
657 const labelList& bin0,
658 const labelList& bin1
659) const
660{
661 const triSurface& surf = surf_;
662 const labelList& eFaces = surf.edgeFaces()[edgei];
663 const edge& e = surf.edges()[edgei];
664
665 // 1. store + or - region number depending
666 // on orientation of triangle in bins[0]
667 labelList regionAndNormal(bin0.size());
668 forAll(bin0, i)
669 {
670 const labelledTri& t = surf.localFaces()[eFaces[bin0[i]]];
671 const auto dir = t.edgeDirection(e);
672
673 if (dir > 0)
674 {
675 regionAndNormal[i] = t.region()+1;
676 }
677 else if (dir == 0)
678 {
679 FatalErrorInFunction
680 << exit(FatalError);
681 }
682 else
683 {
684 regionAndNormal[i] = -(t.region()+1);
685 }
686 }
687
688 // 2. check against bin1
689 labelList regionAndNormal1(bin1.size());
690 forAll(bin1, i)
691 {
692 const labelledTri& t = surf.localFaces()[eFaces[bin1[i]]];
693 const auto dir = t.edgeDirection(e);
694
695 label myRegionAndNormal;
696 if (dir > 0)
697 {
698 myRegionAndNormal = t.region()+1;
699 }
700 else if (dir == 0)
701 {
702 myRegionAndNormal = 0;
703 FatalErrorInFunction
704 << exit(FatalError);
705 }
706 else
707 {
708 myRegionAndNormal = -(t.region()+1);
709 }
710
711 regionAndNormal1[i] = myRegionAndNormal;
712
713 label index = regionAndNormal.find(-myRegionAndNormal);
714 if (index == -1)
715 {
716 // Not found.
717 //Pout<< "cannot find region " << myRegionAndNormal
718 // << " in regions " << regionAndNormal << endl;
719
720 return surfaceFeatures::REGION;
721 }
722 }
723 return surfaceFeatures::NONE;
724}
725
726
727Foam::surfaceFeatures::edgeStatus
728Foam::surfaceFeatures::surfaceFeatures::checkFlatRegionEdge
729(
730 const scalar tol,
731 const scalar includedAngle,
732 const label edgei,
733 const point& leftPoint
734) const
735{
736 const triSurface& surf = surf_;
737 const edge& e = surf.edges()[edgei];
738 const auto& mp = surf.meshPoints();
739 const point eMid(edge(mp[e[0]], mp[e[1]]).centre(surf.points()));
740 const labelList& eFaces = surf.edgeFaces()[edgei];
741
742 vector leftVec(leftPoint-eMid);
743 leftVec.normalise();
744
745 // Bin according to normal and location w.r.t. first face
746 plane pl(eMid, leftVec);
747
748
749 DynamicList<labelList> leftBins(2);
750 DynamicList<labelList> rightBins(2);
751 DynamicList<vector> leftNormals(2);
752 DynamicList<vector> rightNormals(2);
753
754 // Append first face since this is what leftPoint was created from in the
755 // first place
756 leftNormals.append(surf.faceNormals()[eFaces[0]]);
757 leftBins.append(labelList(1, 0));
758
759 for (label eFacei = 1; eFacei < eFaces.size(); ++eFacei)
760 {
761 const label facei = eFaces[eFacei];
762
763 const bool isLeft(pl.signedDistance(surf.faceCentres()[facei]) > 0);
764
765 DynamicList<labelList>& bins = (isLeft ? leftBins : rightBins);
766 DynamicList<vector>& normals = (isLeft ? leftNormals : rightNormals);
767
768 const vector& n = surf.faceNormals()[facei];
769
770 // Find the normal in normals
771 label index = -1;
772 forAll(normals, normalI)
773 {
774 if (mag(n & normals[normalI]) > (1-tol))
775 {
776 index = normalI;
777 break;
778 }
779 }
780
781 if (index != -1)
782 {
783 bins[index].append(eFacei);
784// Pout<< "edge:" << edgei << " verts:" << e
785// << " found existing normal bin:" << index
786// << " after:" << flatOutput(bins[index])
787// << endl;
788 }
789 else if ((leftNormals.size()+rightNormals.size()) >= 2)
790 {
791 // Would be third normal. Mark as feature.
792 //Pout<< "** at edge:" << surf.localPoints()[e[0]]
793 // << surf.localPoints()[e[1]]
794 // << " have normals:" << normals
795 // << " and " << n << endl;
796 return surfaceFeatures::REGION;
797 }
798 else
799 {
800 normals.append(n);
801 bins.append(labelList(1, eFacei));
802 }
803 }
804
805 // Check resulting number of bins
806 if ((leftBins.size()+rightBins.size()) == 1)
807 {
808 // Note: should check here whether they are two sets of faces
809 // that are planar or indeed 4 faces al coming together at an edge.
810 //Pout<< "** at edge:"
811 // << surf.localPoints()[e[0]]
812 // << surf.localPoints()[e[1]]
813 // << " have single normal:" << normals[0]
814 // << endl;
815 return surfaceFeatures::NONE;
816 }
817 else
818 {
819 // Two bins. Check if normals make an angle
820
821 //Pout<< "** at edge:"
822 // << surf.localPoints()[e[0]]
823 // << surf.localPoints()[e[1]] << nl
824 // << " normals:" << normals << nl
825 // << " bins :" << bins << nl
826 // << endl;
827
828 if (includedAngle >= 0)
829 {
830 scalar minCos = Foam::cos(degToRad(180.0 - includedAngle));
831
832 forAll(eFaces, i)
833 {
834 const vector& ni = surf.faceNormals()[eFaces[i]];
835 for (label j=i+1; j<eFaces.size(); j++)
836 {
837 const vector& nj = surf.faceNormals()[eFaces[j]];
838 if (mag(ni & nj) < minCos)
839 {
840 //Pout<< "have sharp feature between normal:" << ni
841 // << " and " << nj << endl;
842
843 // Is feature. Keep as region or convert to
844 // feature angle? For now keep as region.
845 return surfaceFeatures::REGION;
846 }
847 }
848 }
849 }
850
851 // So now we have two normals bins but need to make sure both
852 // bins have the same regions in it.
853
854 if (leftBins.size() == 2)
855 {
856 return checkBinRegion
857 (
858 edgei,
859 leftBins[0],
860 leftBins[1]
861 );
862 }
863 else if (leftBins.size() == 1 && rightBins.size() == 1)
864 {
865 return checkBinRegion
866 (
867 edgei,
868 leftBins[0],
869 rightBins[0]
870 );
871 }
872 else if (rightBins.size() == 2)
873 {
874 return checkBinRegion
875 (
876 edgei,
877 rightBins[0],
878 rightBins[1]
879 );
880 }
881 else
882 {
883 FatalErrorInFunction << "leftBins:" << leftBins
884 << " rightBins:" << rightBins << exit(FatalError);
885 }
886 }
887
888 return surfaceFeatures::NONE;
889}
891
892// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
893
895:
896 surf_(surf),
897 featurePoints_(0),
898 featureEdges_(0),
899 externalStart_(0),
900 internalStart_(0)
902
903
904// Construct from components.
906(
907 const triSurface& surf,
908 const labelList& featurePoints,
909 const labelList& featureEdges,
910 const label externalStart,
911 const label internalStart
912)
913:
914 surf_(surf),
915 featurePoints_(featurePoints),
916 featureEdges_(featureEdges),
917 externalStart_(externalStart),
918 internalStart_(externalStart)
920
921
922// Construct from surface, angle and min length
924(
925 const triSurface& surf,
926 const scalar includedAngle,
927 const scalar minLen,
928 const label minElems,
929 const bool geometricTestOnly
930)
931:
932 surf_(surf),
933 featurePoints_(0),
934 featureEdges_(0),
935 externalStart_(0),
936 internalStart_(0)
937{
938 findFeatures(includedAngle, geometricTestOnly);
939
940 if (minLen > 0 || minElems > 0)
941 {
942 trimFeatures(minLen, minElems, includedAngle);
944}
945
946
948(
949 const triSurface& surf,
950 const dictionary& featInfoDict
951)
952:
953 surf_(surf),
954 featurePoints_(featInfoDict.lookup("featurePoints")),
955 featureEdges_(featInfoDict.lookup("featureEdges")),
956 externalStart_(featInfoDict.get<label>("externalStart")),
957 internalStart_(featInfoDict.get<label>("internalStart"))
958{}
959
960
962(
963 const triSurface& surf,
964 const fileName& fName
965)
966:
967 surf_(surf),
968 featurePoints_(0),
969 featureEdges_(0),
970 externalStart_(0),
971 internalStart_(0)
972{
973 IFstream str(fName);
974
975 dictionary featInfoDict(str);
976
977 featInfoDict.readEntry("featureEdges", featureEdges_);
978 featInfoDict.readEntry("featurePoints", featurePoints_);
979 featInfoDict.readEntry("externalStart", externalStart_);
980 featInfoDict.readEntry("internalStart", internalStart_);
981}
982
983
985(
986 const triSurface& surf,
987 const pointField& points,
988 const edgeList& edges,
989 const scalar mergeTol,
990 const bool geometricTestOnly
991)
992:
993 surf_(surf),
994 featurePoints_(0),
995 featureEdges_(0),
996 externalStart_(0),
997 internalStart_(0)
998{
999 // Match edge mesh edges with the triSurface edges
1000
1001 const labelListList& surfEdgeFaces = surf_.edgeFaces();
1002 const edgeList& surfEdges = surf_.edges();
1003
1004 scalar mergeTolSqr = sqr(mergeTol);
1005
1006 EdgeMap<label> dynFeatEdges(2*edges.size());
1007 DynamicList<labelList> dynFeatureEdgeFaces(edges.size());
1008
1009 labelList edgeLabel;
1010
1012 (
1013 edges,
1014 points,
1015 mergeTolSqr,
1016 edgeLabel // label of surface edge or -1
1017 );
1018
1019 label count = 0;
1020 forAll(edgeLabel, sEdgeI)
1021 {
1022 const label sEdge = edgeLabel[sEdgeI];
1023
1024 if (sEdge == -1)
1025 {
1026 continue;
1027 }
1028
1029 dynFeatEdges.insert(surfEdges[sEdge], count++);
1030 dynFeatureEdgeFaces.append(surfEdgeFaces[sEdge]);
1031 }
1032
1033 // Find whether an edge is external or internal
1034 List<edgeStatus> edgeStat(dynFeatEdges.size(), NONE);
1035
1036 classifyFeatureAngles
1037 (
1038 dynFeatureEdgeFaces,
1039 edgeStat,
1040 GREAT,
1041 geometricTestOnly
1042 );
1043
1044 // Transfer the edge status to a list encompassing all edges in the surface
1045 // so that calcFeatPoints can be used.
1046 List<edgeStatus> allEdgeStat(surf_.nEdges(), NONE);
1047
1048 forAll(allEdgeStat, eI)
1049 {
1050 const auto iter = dynFeatEdges.cfind(surfEdges[eI]);
1051
1052 if (iter.good())
1053 {
1054 allEdgeStat[eI] = edgeStat[iter.val()];
1055 }
1056 }
1057
1058 edgeStat.clear();
1059 dynFeatEdges.clear();
1060
1061 setFromStatus(allEdgeStat, GREAT);
1062}
1063
1064
1066:
1067 surf_(sf.surface()),
1068 featurePoints_(sf.featurePoints()),
1069 featureEdges_(sf.featureEdges()),
1070 externalStart_(sf.externalStart()),
1071 internalStart_(sf.internalStart())
1072{}
1074
1075// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1076
1078(
1079 const bool regionEdges,
1080 const bool externalEdges,
1081 const bool internalEdges
1082) const
1083{
1084 DynamicList<label> selectedEdges;
1085
1086 // Presizing
1087 {
1088 label count = 0;
1089
1090 if (regionEdges) count += nRegionEdges();
1091 if (externalEdges) count += nExternalEdges();
1092 if (internalEdges) count += nInternalEdges();
1093
1094 selectedEdges.reserve_exact(count);
1095 }
1096
1097 if (regionEdges)
1098 {
1099 for (label i = 0; i < externalStart_; i++)
1100 {
1101 selectedEdges.append(featureEdges_[i]);
1102 }
1103 }
1104
1105 if (externalEdges)
1106 {
1107 for (label i = externalStart_; i < internalStart_; i++)
1108 {
1109 selectedEdges.append(featureEdges_[i]);
1110 }
1111 }
1112
1113 if (internalEdges)
1114 {
1115 for (label i = internalStart_; i < featureEdges_.size(); i++)
1116 {
1117 selectedEdges.append(featureEdges_[i]);
1118 }
1119 }
1120
1121 return labelList(std::move(selectedEdges));
1122}
1123
1124
1126(
1127 const scalar includedAngle,
1128 const bool geometricTestOnly
1129)
1130{
1131 scalar minCos = Foam::cos(degToRad(180.0 - includedAngle));
1132
1133 // Per edge whether is feature edge.
1134 List<edgeStatus> edgeStat(surf_.nEdges(), NONE);
1135
1136 classifyFeatureAngles
1137 (
1138 surf_.edgeFaces(),
1139 edgeStat,
1140 minCos,
1141 geometricTestOnly
1142 );
1143
1144 setFromStatus(edgeStat, includedAngle);
1145}
1147
1148// Remove small strings of edges. First assigns string number to
1149// every edge and then works out the length of them.
1151(
1152 const scalar minLen,
1153 const label minElems,
1154 const scalar includedAngle
1155)
1156{
1157 // Convert feature edge list to status per edge.
1158 List<edgeStatus> edgeStat(toStatus());
1159
1160
1161 // Mark feature edges according to connected lines.
1162 // -1: unassigned
1163 // -2: part of too small a feature line
1164 // >0: feature line number
1165 labelList featLines(surf_.nEdges(), -1);
1166
1167 // Current featureline number.
1168 label featI = 0;
1169
1170 // Current starting edge
1171 label startEdgeI = 0;
1172
1173 do
1174 {
1175 // Find unset featureline
1176 for (; startEdgeI < edgeStat.size(); startEdgeI++)
1177 {
1178 if
1179 (
1180 edgeStat[startEdgeI] != NONE
1181 && featLines[startEdgeI] == -1
1182 )
1183 {
1184 // Found unassigned feature edge.
1185 break;
1186 }
1187 }
1188
1189 if (startEdgeI == edgeStat.size())
1190 {
1191 // No unset feature edge found. All feature edges have line number
1192 // assigned.
1193 break;
1194 }
1195
1196 featLines[startEdgeI] = featI;
1197
1198 const edge& startEdge = surf_.edges()[startEdgeI];
1199
1200 // Walk 'left' and 'right' from startEdge.
1201 labelScalar leftPath =
1202 walkSegment
1203 (
1204 true, // 'mark' mode
1205 edgeStat,
1206 startEdgeI, // edge, not yet assigned to a featureLine
1207 startEdge[0], // start of edge
1208 featI, // Mark value
1209 featLines // Mark for all edges
1210 );
1211
1212 labelScalar rightPath =
1213 walkSegment
1214 (
1215 true,
1216 edgeStat,
1217 startEdgeI,
1218 startEdge[1], // end of edge
1219 featI,
1220 featLines
1221 );
1222
1223 if
1224 (
1225 (
1226 leftPath.len_
1227 + rightPath.len_
1228 + startEdge.mag(surf_.localPoints())
1229 < minLen
1230 )
1231 || (leftPath.n_ + rightPath.n_ + 1 < minElems)
1232 )
1233 {
1234 // Rewalk same route (recognizable by featLines == featI)
1235 // to reset featLines.
1236
1237 featLines[startEdgeI] = -2;
1238
1239 walkSegment
1240 (
1241 false, // 'clean' mode
1242 edgeStat,
1243 startEdgeI, // edge, not yet assigned to a featureLine
1244 startEdge[0], // start of edge
1245 featI, // Unset value
1246 featLines // Mark for all edges
1247 );
1248
1249 walkSegment
1250 (
1251 false,
1252 edgeStat,
1253 startEdgeI,
1254 startEdge[1], // end of edge
1255 featI,
1256 featLines
1257 );
1258 }
1259 else
1260 {
1261 featI++;
1262 }
1263 }
1264 while (true);
1265
1266 // Unmark all feature lines that have featLines=-2
1267 forAll(featureEdges_, i)
1268 {
1269 label edgeI = featureEdges_[i];
1270
1271 if (featLines[edgeI] == -2)
1272 {
1273 edgeStat[edgeI] = NONE;
1274 }
1275 }
1276
1277 // Convert back to edge labels
1278 setFromStatus(edgeStat, includedAngle);
1279
1280 return featLines;
1281}
1282
1283
1285(
1286 List<edgeStatus>& edgeStat,
1287 const treeBoundBox& bb
1288) const
1289{
1290 deleteBox(edgeStat, bb, true);
1291}
1292
1293
1295(
1296 List<edgeStatus>& edgeStat,
1297 const treeBoundBox& bb
1298) const
1299{
1300 deleteBox(edgeStat, bb, false);
1301}
1302
1303
1305(
1306 List<edgeStatus>& edgeStat,
1307 const treeBoundBox& bb,
1308 const bool removeInside
1309) const
1310{
1311 const edgeList& surfEdges = surf_.edges();
1312 const pointField& surfLocalPoints = surf_.localPoints();
1313
1314 forAll(edgeStat, edgei)
1315 {
1316 const point eMid = surfEdges[edgei].centre(surfLocalPoints);
1317
1318 if (removeInside ? bb.contains(eMid) : !bb.contains(eMid))
1319 {
1320 edgeStat[edgei] = surfaceFeatures::NONE;
1321 }
1323}
1324
1325
1327(
1328 List<edgeStatus>& edgeStat,
1329 const plane& cutPlane
1330) const
1331{
1332 const edgeList& surfEdges = surf_.edges();
1333 const pointField& pts = surf_.points();
1334 const labelList& meshPoints = surf_.meshPoints();
1335
1336 forAll(edgeStat, edgei)
1337 {
1338 const edge& e = surfEdges[edgei];
1339
1340 const point& p0 = pts[meshPoints[e.start()]];
1341 const point& p1 = pts[meshPoints[e.end()]];
1342 const linePointRef line(p0, p1);
1343
1344 // If edge does not intersect the plane, delete.
1345 scalar intersect = cutPlane.lineIntersect(line);
1346
1347 point featPoint = intersect * (p1 - p0) + p0;
1348
1349 if (!onLine(featPoint, line))
1350 {
1351 edgeStat[edgei] = surfaceFeatures::NONE;
1352 }
1354}
1355
1356
1358(
1359 List<edgeStatus>& edgeStat
1360) const
1361{
1362 forAll(edgeStat, edgei)
1363 {
1364 if (surf_.edgeFaces()[edgei].size() == 1)
1365 {
1366 edgeStat[edgei] = surfaceFeatures::NONE;
1367 }
1369}
1370
1371
1373(
1374 List<edgeStatus>& edgeStat
1375) const
1376{
1377 forAll(edgeStat, edgei)
1378 {
1379 if (surf_.edgeFaces()[edgei].size() > 2)
1380 {
1381 edgeStat[edgei] = surfaceFeatures::NONE;
1382 }
1383 }
1384}
1385
1386
1387//- Divide into multiple normal bins
1388// - return REGION if != 2 normals
1389// - return REGION if 2 normals that make feature angle
1390// - otherwise return NONE and set normals,bins
1391void Foam::surfaceFeatures::checkFlatRegionEdge
1392(
1393 List<edgeStatus>& edgeStat,
1394 const scalar tol,
1395 const scalar includedAngle
1396) const
1397{
1398 forAll(edgeStat, edgei)
1399 {
1400 if (edgeStat[edgei] == surfaceFeatures::REGION)
1401 {
1402 const labelList& eFaces = surf_.edgeFaces()[edgei];
1403
1404 if (eFaces.size() > 2 && (eFaces.size() % 2) == 0)
1405 {
1406 const point& leftPoint = surf_.faceCentres()[eFaces[0]];
1407
1408 edgeStat[edgei] = checkFlatRegionEdge
1409 (
1410 tol,
1411 includedAngle,
1412 edgei,
1413 leftPoint
1414 );
1415 }
1416 }
1417 }
1418}
1419
1422{
1423 dictionary featInfoDict;
1424 featInfoDict.add("externalStart", externalStart_);
1425 featInfoDict.add("internalStart", internalStart_);
1426 featInfoDict.add("featureEdges", featureEdges_);
1427 featInfoDict.add("featurePoints", featurePoints_);
1428
1429 featInfoDict.write(os);
1430}
1431
1433void Foam::surfaceFeatures::write(const fileName& fName) const
1434{
1435 OFstream os(fName);
1436 writeDict(os);
1437}
1438
1440void Foam::surfaceFeatures::writeObj(const fileName& prefix) const
1441{
1442 OFstream regionStr(prefix + "_regionEdges.obj");
1443 Pout<< "Writing region edges to " << regionStr.name() << endl;
1444
1445 label verti = 0;
1446 for (label i = 0; i < externalStart_; i++)
1447 {
1448 const edge& e = surf_.edges()[featureEdges_[i]];
1449
1450 meshTools::writeOBJ(regionStr, surf_.localPoints()[e[0]]); verti++;
1451 meshTools::writeOBJ(regionStr, surf_.localPoints()[e[1]]); verti++;
1452 regionStr << "l " << verti-1 << ' ' << verti << endl;
1453 }
1454
1455
1456 OFstream externalStr(prefix + "_externalEdges.obj");
1457 Pout<< "Writing external edges to " << externalStr.name() << endl;
1458
1459 verti = 0;
1460 for (label i = externalStart_; i < internalStart_; i++)
1461 {
1462 const edge& e = surf_.edges()[featureEdges_[i]];
1463
1464 meshTools::writeOBJ(externalStr, surf_.localPoints()[e[0]]); verti++;
1465 meshTools::writeOBJ(externalStr, surf_.localPoints()[e[1]]); verti++;
1466 externalStr << "l " << verti-1 << ' ' << verti << endl;
1467 }
1468
1469 OFstream internalStr(prefix + "_internalEdges.obj");
1470 Pout<< "Writing internal edges to " << internalStr.name() << endl;
1471
1472 verti = 0;
1473 for (label i = internalStart_; i < featureEdges_.size(); i++)
1474 {
1475 const edge& e = surf_.edges()[featureEdges_[i]];
1476
1477 meshTools::writeOBJ(internalStr, surf_.localPoints()[e[0]]); verti++;
1478 meshTools::writeOBJ(internalStr, surf_.localPoints()[e[1]]); verti++;
1479 internalStr << "l " << verti-1 << ' ' << verti << endl;
1480 }
1481
1482 OFstream pointStr(prefix + "_points.obj");
1483 Pout<< "Writing feature points to " << pointStr.name() << endl;
1484
1485 for (const label pointi : featurePoints_)
1486 {
1487 meshTools::writeOBJ(pointStr, surf_.localPoints()[pointi]);
1488 }
1489}
1490
1493{
1494 os << "Feature set:" << nl
1495 << " points : " << this->featurePoints().size() << nl
1496 << " edges : " << this->featureEdges().size() << nl
1497 << " of which" << nl
1498 << " region edges : " << this->nRegionEdges() << nl
1499 << " external edges : " << this->nExternalEdges() << nl
1500 << " internal edges : " << this->nInternalEdges() << endl;
1501}
1502
1503
1504// Get nearest vertex on patch for every point of surf in pointSet.
1506(
1507 const labelList& pointLabels,
1508 const pointField& samples,
1509 const scalarField& maxDistSqr
1510) const
1511{
1512 // Build tree out of all samples.
1513
1514 // Define bound box here (gcc-4.8.5)
1515 const treeBoundBox overallBb(samples);
1516
1518 (
1520 overallBb,
1521 8, // maxLevel
1522 10, // leafsize
1523 3.0 // duplicity
1524 );
1525 const auto& treeData = ppTree.shapes();
1526
1527 // From patch point to surface point
1528 Map<label> nearest(2*pointLabels.size());
1529
1530 const pointField& surfPoints = surf_.localPoints();
1531
1533 {
1534 const label surfPointi = pointLabels[i];
1535
1536 const point& surfPt = surfPoints[surfPointi];
1537
1538 pointIndexHit info = ppTree.findNearest
1539 (
1540 surfPt,
1541 maxDistSqr[i]
1542 );
1543
1544 if (!info.hit())
1545 {
1547 << "Problem for point "
1548 << surfPointi << " in tree " << ppTree.bb()
1549 << abort(FatalError);
1550 }
1551
1552 label sampleI = info.index();
1553
1554 if (treeData.centre(sampleI).distSqr(surfPt) < maxDistSqr[sampleI])
1555 {
1556 nearest.insert(sampleI, surfPointi);
1557 }
1558 }
1559
1560
1561 if (debug)
1562 {
1563 //
1564 // Dump to obj file
1565 //
1566
1567 Pout<< "Dumping nearest surface feature points to nearestSamples.obj"
1568 << endl;
1569
1570 OFstream objStream("nearestSamples.obj");
1571
1572 label vertI = 0;
1573 forAllConstIters(nearest, iter)
1574 {
1575 meshTools::writeOBJ(objStream, samples[iter.key()]); vertI++;
1576 meshTools::writeOBJ(objStream, surfPoints[iter.val()]); vertI++;
1577 objStream<< "l " << vertI-1 << ' ' << vertI << endl;
1578 }
1579 }
1580
1581 return nearest;
1582}
1583
1584
1585// Get nearest sample point for regularly sampled points along
1586// selected edges. Return map from sample to edge label
1588(
1589 const labelList& selectedEdges,
1590 const pointField& samples,
1591 const scalarField& sampleDist,
1592 const scalarField& maxDistSqr,
1593 const scalar minSampleDist
1594) const
1595{
1596 const pointField& surfPoints = surf_.localPoints();
1597 const edgeList& surfEdges = surf_.edges();
1598
1599 scalar maxSearchSqr = max(maxDistSqr);
1600
1601 // Define bound box here (gcc-4.8.5)
1602 const treeBoundBox overallBb(samples);
1603
1605 (
1607 overallBb,
1608 8, // maxLevel
1609 10, // leafsize
1610 3.0 // duplicity
1611 );
1612
1613 // From patch point to surface edge.
1614 Map<label> nearest(2*selectedEdges.size());
1615
1616 forAll(selectedEdges, i)
1617 {
1618 label surfEdgeI = selectedEdges[i];
1619
1620 const edge& e = surfEdges[surfEdgeI];
1621
1622 if (debug && (i % 1000) == 0)
1623 {
1624 Pout<< "looking at surface feature edge " << surfEdgeI
1625 << " verts:" << e << " points:" << surfPoints[e[0]]
1626 << ' ' << surfPoints[e[1]] << endl;
1627 }
1628
1629 // Normalized edge vector
1630 vector eVec = e.vec(surfPoints);
1631 scalar eMag = mag(eVec);
1632 eVec /= eMag;
1633
1634
1635 //
1636 // Sample along edge
1637 //
1638
1639 bool exit = false;
1640
1641 // Coordinate along edge (0 .. eMag)
1642 scalar s = 0.0;
1643
1644 while (true)
1645 {
1646 point edgePoint(surfPoints[e.start()] + s*eVec);
1647
1648 pointIndexHit info = ppTree.findNearest
1649 (
1650 edgePoint,
1651 maxSearchSqr
1652 );
1653
1654 if (!info.hit())
1655 {
1656 // No point close enough to surface edge.
1657 break;
1658 }
1659
1660 label sampleI = info.index();
1661
1662 if (info.point().distSqr(edgePoint) < maxDistSqr[sampleI])
1663 {
1664 nearest.insert(sampleI, surfEdgeI);
1665 }
1666
1667 if (exit)
1668 {
1669 break;
1670 }
1671
1672 // Step to next sample point using local distance.
1673 // Truncate to max 1/minSampleDist samples per feature edge.
1674 s += max(minSampleDist*eMag, sampleDist[sampleI]);
1675
1676 if (s >= (1-minSampleDist)*eMag)
1677 {
1678 // Do one more sample, at edge endpoint
1679 s = eMag;
1680 exit = true;
1681 }
1682 }
1683 }
1684
1685
1686
1687 if (debug)
1688 {
1689 // Dump to obj file
1690
1691 Pout<< "Dumping nearest surface edges to nearestEdges.obj"
1692 << endl;
1693
1694 OFstream objStream("nearestEdges.obj");
1695
1696 label vertI = 0;
1697 forAllConstIters(nearest, iter)
1698 {
1699 const label sampleI = iter.key();
1700
1701 const edge& e = surfEdges[iter.val()];
1702
1703 meshTools::writeOBJ(objStream, samples[sampleI]); vertI++;
1704
1705 point nearPt =
1706 e.line(surfPoints).nearestDist(samples[sampleI]).point();
1707
1708 meshTools::writeOBJ(objStream, nearPt); vertI++;
1709
1710 objStream<< "l " << vertI-1 << ' ' << vertI << endl;
1711 }
1712 }
1713
1714 return nearest;
1715}
1716
1717
1718// Get nearest edge on patch for regularly sampled points along the feature
1719// edges. Return map from patch edge to feature edge.
1720//
1721// Q: using point-based sampleDist and maxDist (distance to look around
1722// each point). Should they be edge-based e.g. by averaging or max()?
1724(
1725 const labelList& selectedEdges,
1726 const edgeList& sampleEdges,
1727 const labelList& selectedSampleEdges,
1728 const pointField& samplePoints,
1729 const scalarField& sampleDist,
1730 const scalarField& maxDistSqr,
1731 const scalar minSampleDist
1732) const
1733{
1734 // Build tree out of selected sample edges.
1736 (
1738 (
1739 sampleEdges,
1740 samplePoints,
1741 selectedSampleEdges
1742 ),
1743 treeBoundBox(samplePoints), // overall search domain
1744 8, // maxLevel
1745 10, // leafsize
1746 3.0 // duplicity
1747 );
1748
1749 const pointField& surfPoints = surf_.localPoints();
1750 const edgeList& surfEdges = surf_.edges();
1751
1752 scalar maxSearchSqr = max(maxDistSqr);
1753
1754 Map<pointIndexHit> nearest(2*sampleEdges.size());
1755
1756 //
1757 // Loop over all selected edges. Sample at regular intervals. Find nearest
1758 // sampleEdges (using octree)
1759 //
1760
1761 forAll(selectedEdges, i)
1762 {
1763 label surfEdgeI = selectedEdges[i];
1764
1765 const edge& e = surfEdges[surfEdgeI];
1766
1767 if (debug && (i % 1000) == 0)
1768 {
1769 Pout<< "looking at surface feature edge " << surfEdgeI
1770 << " verts:" << e << " points:" << surfPoints[e[0]]
1771 << ' ' << surfPoints[e[1]] << endl;
1772 }
1773
1774 // Normalized edge vector
1775 vector eVec = e.vec(surfPoints);
1776 scalar eMag = mag(eVec);
1777 eVec /= eMag;
1778
1779
1780 //
1781 // Sample along edge
1782 //
1783
1784 bool exit = false;
1785
1786 // Coordinate along edge (0 .. eMag)
1787 scalar s = 0.0;
1788
1789 while (true)
1790 {
1791 point edgePoint(surfPoints[e.start()] + s*eVec);
1792
1793 pointIndexHit info = ppTree.findNearest
1794 (
1795 edgePoint,
1796 maxSearchSqr
1797 );
1798
1799 if (!info.hit())
1800 {
1801 // No edge close enough to surface edge.
1802 break;
1803 }
1804
1805 label index = info.index();
1806
1807 label sampleEdgeI = ppTree.shapes().objectIndex(index);
1808
1809 const edge& e = sampleEdges[sampleEdgeI];
1810
1811 if (info.point().distSqr(edgePoint) < maxDistSqr[e.start()])
1812 {
1813 nearest.insert
1814 (
1815 sampleEdgeI,
1816 pointIndexHit(true, edgePoint, surfEdgeI)
1817 );
1818 }
1819
1820 if (exit)
1821 {
1822 break;
1823 }
1824
1825 // Step to next sample point using local distance.
1826 // Truncate to max 1/minSampleDist samples per feature edge.
1827 // s += max(minSampleDist*eMag, sampleDist[e.start()]);
1828 s += 0.01*eMag;
1829
1830 if (s >= (1-minSampleDist)*eMag)
1831 {
1832 // Do one more sample, at edge endpoint
1833 s = eMag;
1834 exit = true;
1835 }
1836 }
1837 }
1838
1839
1840 if (debug)
1841 {
1842 // Dump to obj file
1843
1844 Pout<< "Dumping nearest surface feature edges to nearestEdges.obj"
1845 << endl;
1846
1847 OFstream objStream("nearestEdges.obj");
1848
1849 label vertI = 0;
1850 forAllConstIters(nearest, iter)
1851 {
1852 const label sampleEdgeI = iter.key();
1853
1854 const edge& sampleEdge = sampleEdges[sampleEdgeI];
1855
1856 // Write line from edgeMid to point on feature edge
1857 meshTools::writeOBJ(objStream, sampleEdge.centre(samplePoints));
1858 vertI++;
1859
1860 meshTools::writeOBJ(objStream, iter.val().point());
1861 vertI++;
1862
1863 objStream<< "l " << vertI-1 << ' ' << vertI << endl;
1864 }
1865 }
1866
1867 return nearest;
1868}
1869
1870
1871// Get nearest surface edge for every sample. Return in form of
1872// labelLists giving surfaceEdge label&intersectionpoint.
1874(
1875 const labelList& selectedEdges,
1876 const pointField& samples,
1877 scalar searchSpanSqr, // Search span
1878 labelList& edgeLabel,
1879 labelList& edgeEndPoint,
1880 pointField& edgePoint
1881) const
1882{
1883 edgeLabel.setSize(samples.size());
1884 edgeEndPoint.setSize(samples.size());
1885 edgePoint.setSize(samples.size());
1886
1887 const pointField& localPoints = surf_.localPoints();
1888
1889 treeBoundBox searchDomain(localPoints);
1890 searchDomain.inflate(0.1);
1891
1893 (
1895 (
1896 surf_.edges(),
1897 localPoints,
1898 selectedEdges
1899 ),
1900 searchDomain, // overall search domain
1901 8, // maxLevel
1902 10, // leafsize
1903 3.0 // duplicity
1904 );
1905 const auto& treeData = ppTree.shapes();
1906
1907 forAll(samples, i)
1908 {
1909 const point& sample = samples[i];
1910
1911 pointIndexHit info = ppTree.findNearest
1912 (
1913 sample,
1914 searchSpanSqr
1915 );
1916
1917 if (!info.hit())
1918 {
1919 edgeLabel[i] = -1;
1920 }
1921 else
1922 {
1923 // Need to recalculate to classify edgeEndPoint
1924 pointIndexHit pHit = edgeNearest
1925 (
1926 treeData.line(info.index()),
1927 sample
1928 );
1929
1930 edgeLabel[i] = treeData.objectIndex(info.index());
1931 edgePoint[i] = pHit.point();
1932 edgeEndPoint[i] = pHit.index();
1933 }
1934 }
1935}
1936
1937
1938// Get nearest point on nearest feature edge for every sample (is edge)
1940(
1941 const labelList& selectedEdges,
1942
1943 const edgeList& sampleEdges,
1944 const labelList& selectedSampleEdges,
1945 const pointField& samplePoints,
1946
1947 const vector& searchSpan, // Search span
1948 labelList& edgeLabel, // label of surface edge or -1
1949 pointField& pointOnEdge, // point on above edge
1950 pointField& pointOnFeature // point on sample edge
1951) const
1952{
1953 edgeLabel.setSize(selectedSampleEdges.size());
1954 pointOnEdge.setSize(selectedSampleEdges.size());
1955 pointOnFeature.setSize(selectedSampleEdges.size());
1956
1957 treeBoundBox searchDomain(surf_.localPoints());
1958
1960 (
1962 (
1963 surf_.edges(),
1964 surf_.localPoints(),
1965 selectedEdges
1966 ),
1967 searchDomain, // overall search domain
1968 8, // maxLevel
1969 10, // leafsize
1970 3.0 // duplicity
1971 );
1972 const auto& treeData = ppTree.shapes();
1973
1974 forAll(selectedSampleEdges, i)
1975 {
1976 const edge& e = sampleEdges[selectedSampleEdges[i]];
1977
1978 linePointRef edgeLine = e.line(samplePoints);
1979
1980 point eMid(edgeLine.centre());
1981
1982 treeBoundBox tightest(eMid - searchSpan, eMid + searchSpan);
1983
1984 pointIndexHit info = ppTree.findNearest
1985 (
1986 edgeLine,
1987 tightest,
1988 pointOnEdge[i]
1989 );
1990
1991 if (!info.hit())
1992 {
1993 edgeLabel[i] = -1;
1994 }
1995 else
1996 {
1997 edgeLabel[i] = treeData.objectIndex(info.index());
1998 pointOnFeature[i] = info.point();
1999 }
2000 }
2001}
2002
2005(
2006 const edgeList& edges,
2007 const pointField& points,
2008 scalar searchSpanSqr, // Search span
2009 labelList& edgeLabel
2010) const
2011{
2012 edgeLabel = labelList(surf_.nEdges(), -1);
2013
2014 treeBoundBox searchDomain(points);
2015 searchDomain.inflate(0.1);
2016
2018 (
2019 treeDataEdge(edges, points), // All edges
2020
2021 searchDomain, // overall search domain
2022 8, // maxLevel
2023 10, // leafsize
2024 3.0 // duplicity
2025 );
2026 const auto& treeData = ppTree.shapes();
2027
2028 const edgeList& surfEdges = surf_.edges();
2029 const pointField& surfLocalPoints = surf_.localPoints();
2030
2031 forAll(surfEdges, edgeI)
2032 {
2033 const edge& sample = surfEdges[edgeI];
2034
2035 const point& startPoint = surfLocalPoints[sample.start()];
2036 const point& midPoint = sample.centre(surfLocalPoints);
2037
2038 pointIndexHit infoMid = ppTree.findNearest
2039 (
2040 midPoint,
2041 searchSpanSqr
2042 );
2043
2044 if (infoMid.hit())
2045 {
2046 const vector surfEdgeDir = midPoint - startPoint;
2047
2048 const vector featEdgeDir = treeData.line(infoMid.index()).vec();
2049
2050 // Check that the edges are nearly parallel
2051 if (mag(surfEdgeDir ^ featEdgeDir) < parallelTolerance)
2052 {
2053 edgeLabel[edgeI] = edgeI;
2054 }
2055 }
2056 }
2057}
2058
2059
2060// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
2063{
2064 if (this == &rhs)
2065 {
2066 return; // Self-assignment is a no-op
2067 }
2068
2069 if (&surf_ != &rhs.surface())
2070 {
2072 << "Operating on different surfaces"
2073 << abort(FatalError);
2074 }
2075
2076 featurePoints_ = rhs.featurePoints();
2077 featureEdges_ = rhs.featureEdges();
2078 externalStart_ = rhs.externalStart();
2079 internalStart_ = rhs.internalStart();
2080}
2081
2082
2083// ************************************************************************* //
label n
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void reserve_exact(const label len)
Reserve allocation space for at least this size, allocating new space if required and retaining old c...
void append(const T &val)
Copy append an element to the end of this list.
Map from edge (expressed as its endpoints) to value. Hashing (and ==) on an edge is symmetric.
Definition edgeHashes.H:59
Minimal example by using system/controlDict.functions:
const_iterator cfind(const Key &key) const
Find and return an const_iterator set at the hashed entry.
Definition HashTableI.H:113
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition HashTableI.H:152
label size() const noexcept
The number of elements in table.
Definition HashTable.H:358
void clear()
Remove all entries from table.
Definition HashTable.C:742
Input from file stream as an ISstream, normally using std::ifstream for the actual input.
Definition IFstream.H:55
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition List.C:347
void setSize(label n)
Alias for resize().
Definition List.H:536
void clear()
Clear the list, i.e. set size to zero.
Definition ListI.H:133
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
virtual const fileName & name() const override
Read/write access to the name of the stream.
Definition OSstream.H:134
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
label index() const noexcept
Return the hit index.
bool hit() const noexcept
Is there a hit?
const point_type & point() const noexcept
Return point, no checks.
label nEdges() const
Number of edges in patch.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
const labelListList & edgeFaces() const
Return edge-face addressing.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
const Cmpt & x() const noexcept
Access to the vector x component.
Definition Vector.H:135
const Cmpt & z() const noexcept
Access to the vector z component.
Definition Vector.H:145
scalar distSqr(const Vector< Cmpt > &v2) const
The L2-norm distance squared from another vector. The magSqr() of the difference.
Definition VectorI.H:95
const Cmpt & y() const noexcept
Access to the vector y component.
Definition Vector.H:140
void inflate(const scalar factor)
Expand box by factor*mag(span) in all dimensions.
Definition boundBoxI.H:381
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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,...
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition dictionary.C:625
void write(Ostream &os, const bool subDict=true) const
Write dictionary, normally with sub-dictionary formatting.
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
point centre(const UList< point > &pts) const
Return centre point (centroid) of the edge.
Definition edgeI.H:388
scalar mag(const UList< point > &pts) const
The length (L2-norm) of the edge vector.
Definition edgeI.H:435
A class for handling file names.
Definition fileName.H:75
Non-pointer based hierarchical recursive searching.
const Type & shapes() const noexcept
Reference to shape.
const treeBoundBox & bb() const
Top bounding box.
A triFace with additional (region) index.
Definition labelledTri.H:56
A line primitive.
Definition line.H:180
PointRef end() const noexcept
The end (second) point.
Definition line.H:252
Point centre() const
Return centre (centroid).
Definition lineI.H:83
PointRef start() const noexcept
The start (first) point.
Definition line.H:247
Mid-point interpolation (weighting factors = 0.5) scheme class.
Definition midPoint.H:54
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition plane.H:91
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
Lookup type of boundary radiation properties.
Definition lookup.H:60
Holds feature edges/points of surface.
void findFeatures(const scalar includedAngle, const bool geometricTestOnly)
Find feature edges using provided included angle.
labelList selectFeatureEdges(const bool regionEdges, const bool externalEdges, const bool internalEdges) const
Helper function: select a subset of featureEdges_.
void excludeNonManifold(List< edgeStatus > &edgeStat) const
Mark edges with >2 connected faces as 'NONE'.
void setFromStatus(const List< edgeStatus > &edgeStat, const scalar includedAngle)
Set from status per edge.
void write(const fileName &fName) const
Write as dictionary to file.
const triSurface & surface() const
void excludeOpen(List< edgeStatus > &edgeStat) const
Mark edges with a single connected face as 'NONE'.
void subsetPlane(List< edgeStatus > &edgeStat, const plane &cutPlane) const
If edge does not intersect the plane, mark as 'NONE'.
label nRegionEdges() const
Return number of region edges.
label nExternalEdges() const
Return number of external edges.
void operator=(const surfaceFeatures &rhs)
void writeObj(const fileName &prefix) const
Write to separate OBJ files (region, external, internal edges,.
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.
List< edgeStatus > toStatus() const
From member feature edges to status per edge.
void subsetBox(List< edgeStatus > &edgeStat, const treeBoundBox &bb) const
Mark edge status outside box as 'NONE'.
@ EXTERNAL
"Convex" edge
@ NONE
Not a classified feature edge.
@ INTERNAL
"Concave" edge
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.
void deleteBox(List< edgeStatus > &edgeStat, const treeBoundBox &bb, const bool removeInside) const
Mark edge status as 'NONE' for edges inside/outside box.
void writeDict(Ostream &os) const
Write as dictionary.
label internalStart() const
Start of internal edges.
label externalStart() const
Start of external edges.
void writeStats(Ostream &os) const
Write some information.
Map< label > nearestSamples(const labelList &selectedPoints, const pointField &samples, const scalarField &maxDistSqr) const
Find nearest sample for selected surface points.
void nearestFeatEdge(const edgeList &edges, const pointField &points, scalar searchSpanSqr, labelList &edgeLabel) const
Find nearest feature edge to each surface edge. Uses the.
label nInternalEdges() const
Return number of internal edges.
const labelList & featureEdges() const
Return feature edge list.
surfaceFeatures(const triSurface &surf)
Construct from surface.
const labelList & featurePoints() const
Return feature point list.
labelList trimFeatures(const scalar minLen, const label minElems, const scalar includedAngle)
Delete small sets of edges. Edges are stringed up and any.
void excludeBox(List< edgeStatus > &edgeStat, const treeBoundBox &bb) const
Mark edge status inside box as 'NONE'.
Standard boundBox with extra functionality for use in octree.
bool contains(const vector &dir, const point &) const
Contains point (inside or on edge) and moving in direction.
Holds data for octree to work on an edges subset.
Holds (reference to) pointField. Encapsulation of data needed for octree searches....
int edgeDirection(const Foam::edge &e) const
Test the edge direction on the face.
Definition triFaceI.H:444
Triangulated surface description with patch information.
Definition triSurface.H:74
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
const volScalarField & p0
Definition EEqn.H:36
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Namespace for handling debugging switches.
Definition debug.C:45
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition meshTools.C:196
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
PointHit< point > pointHit
A PointHit with a 3D point.
Definition pointHit.H:51
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sin(const dimensionedScalar &ds)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
line< point, const point & > linePointRef
A line using referred points.
Definition line.H:66
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
void rhs(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensionedScalar cos(const dimensionedScalar &ds)
messageStream Warning
Warning stream (stdout output on master, null elsewhere), with additional 'FOAM Warning' header text.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
labelList pointLabels(nPoints, -1)
const pointField & pts
volScalarField & b
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235
Unit conversion functions.
scalarField samples(nIntervals, Zero)