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