Loading...
Searching...
No Matches
refinementFeatures.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-2015 OpenFOAM Foundation
9 Copyright (C) 2020-2022 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 "refinementFeatures.H"
30#include "Time.H"
31#include "Tuple2.H"
32#include "DynamicField.H"
33#include "featureEdgeMesh.H"
34#include "meshRefinement.H"
35
36// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
37
38void Foam::refinementFeatures::read
39(
40 const objectRegistry& io,
41 const PtrList<dictionary>& featDicts
42)
43{
44 forAll(featDicts, featI)
45 {
46 const dictionary& dict = featDicts[featI];
47
48 fileName featFileName
49 (
51 (
52 dict,
53 "file",
54 dryRun_,
57 )
58 );
59
60
61 // Try reading extendedEdgeMesh first
62
63 IOobject extFeatObj
64 (
65 featFileName, // name
66 io.time().constant(), // instance
67 "extendedFeatureEdgeMesh", // local
68 io.time(), // registry
72 );
73
74 const fileName fName
75 (
76 extFeatObj.typeFilePath<extendedFeatureEdgeMesh>()
77 );
78
79 if (!fName.empty() && extendedEdgeMesh::canRead(fName))
80 {
81 autoPtr<extendedEdgeMesh> eMeshPtr = extendedEdgeMesh::New
82 (
83 fName
84 );
85
86 if (!dryRun_)
87 {
88 Info<< "Read extendedFeatureEdgeMesh " << extFeatObj.name()
89 << nl << incrIndent;
90 eMeshPtr().writeStats(Info);
91 Info<< decrIndent << endl;
92 }
93
94 set(featI, new extendedFeatureEdgeMesh(extFeatObj, eMeshPtr()));
95 }
96 else
97 {
98 // Try reading edgeMesh
99
100 IOobject featObj
101 (
102 featFileName, // name
103 io.time().constant(), // instance
104 "triSurface", // local
105 io.time(), // registry
109 );
110
111 const fileName fName
112 (
113 featObj.typeFilePath<featureEdgeMesh>()
114 );
115
116 if (fName.empty())
117 {
119 << "Could not open " << featObj.objectPath()
120 << exit(FatalIOError);
121 }
122
123 // Read as edgeMesh
124 autoPtr<edgeMesh> eMeshPtr = edgeMesh::New(fName);
125 const edgeMesh& eMesh = eMeshPtr();
126
127 if (!dryRun_)
128 {
129 Info<< "Read edgeMesh " << featObj.name() << nl
130 << incrIndent;
131 eMesh.writeStats(Info);
132 Info<< decrIndent << endl;
133 }
134
135 // Analyse for feature points. These are all classified as mixed
136 // points for lack of anything better
137 const labelListList& pointEdges = eMesh.pointEdges();
138
139 labelList oldToNew(eMesh.points().size(), -1);
140 DynamicField<point> newPoints(eMesh.points().size());
141 forAll(pointEdges, pointi)
142 {
143 if (pointEdges[pointi].size() > 2)
144 {
145 oldToNew[pointi] = newPoints.size();
146 newPoints.append(eMesh.points()[pointi]);
147 }
148 //else if (pointEdges[pointi].size() == 2)
149 //MEJ: do something based on a feature angle?
150 }
151 label nFeatures = newPoints.size();
152 forAll(oldToNew, pointi)
153 {
154 if (oldToNew[pointi] == -1)
155 {
156 oldToNew[pointi] = newPoints.size();
157 newPoints.append(eMesh.points()[pointi]);
158 }
159 }
160
161
162 const edgeList& edges = eMesh.edges();
163 edgeList newEdges(edges.size());
164 forAll(edges, edgeI)
165 {
166 const edge& e = edges[edgeI];
167 newEdges[edgeI] = edge
168 (
169 oldToNew[e[0]],
170 oldToNew[e[1]]
171 );
172 }
173
174 // Construct an extendedEdgeMesh with
175 // - all points on more than 2 edges : mixed feature points
176 // - all edges : external edges
177
178 extendedEdgeMesh eeMesh
179 (
180 newPoints, // pts
181 newEdges, // eds
182 0, // (point) concaveStart
183 0, // (point) mixedStart
184 nFeatures, // (point) nonFeatureStart
185 edges.size(), // (edge) internalStart
186 edges.size(), // (edge) flatStart
187 edges.size(), // (edge) openStart
188 edges.size(), // (edge) multipleStart
189 vectorField(0), // normals
190 List<extendedEdgeMesh::sideVolumeType>(0),// normalVolumeTypes
191 vectorField(0), // edgeDirections
192 labelListList(0), // normalDirections
193 labelListList(0), // edgeNormals
194 labelListList(0), // featurePointNormals
195 labelListList(0), // featurePointEdges
196 identity(newEdges.size()) // regionEdges
197 );
198
199 //Info<< "Constructed extendedFeatureEdgeMesh " << featObj.name()
200 // << nl << incrIndent;
201 //eeMesh.writeStats(Info);
202 //Info<< decrIndent << endl;
203
204 set(featI, new extendedFeatureEdgeMesh(featObj, eeMesh));
205 }
206
207 const extendedEdgeMesh& eMesh = operator[](featI);
208
209 if (dict.found("levels"))
210 {
211 List<Tuple2<scalar, label>> distLevels(dict.lookup("levels"));
212
213 if (dict.size() < 1)
214 {
216 << " : levels should be at least size 1" << endl
217 << "levels : " << dict.lookup("levels")
218 << exit(FatalError);
219 }
220
221 distances_[featI].setSize(distLevels.size());
222 levels_[featI].setSize(distLevels.size());
223
224 forAll(distLevels, j)
225 {
226 distances_[featI][j] = distLevels[j].first();
227 levels_[featI][j] = distLevels[j].second();
228
229 if (levels_[featI][j] < 0)
230 {
232 << "Feature " << featFileName
233 << " has illegal refinement level " << levels_[featI][j]
234 << exit(FatalError);
235 }
236
237 // Check in incremental order
238 if (j > 0)
239 {
240 if
241 (
242 (distances_[featI][j] <= distances_[featI][j-1])
243 || (levels_[featI][j] > levels_[featI][j-1])
244 )
245 {
247 << " : Refinement should be specified in order"
248 << " of increasing distance"
249 << " (and decreasing refinement level)." << endl
250 << "Distance:" << distances_[featI][j]
251 << " refinementLevel:" << levels_[featI][j]
252 << exit(FatalError);
253 }
254 }
255 }
256 }
257 else
258 {
259 // Look up 'level' for single level
260 levels_[featI] =
262 (
263 1,
265 (
266 dict,
267 "level",
268 dryRun_,
270 0
271 )
272 );
273 distances_[featI] = scalarField(1, Zero);
274 }
275
276 if (!dryRun_)
277 {
278 Info<< "Refinement level according to distance to "
279 << featFileName << " (" << eMesh.points().size() << " points, "
280 << eMesh.edges().size() << " edges)." << endl;
281 forAll(levels_[featI], j)
282 {
283 Info<< " level " << levels_[featI][j]
284 << " for all cells within " << distances_[featI][j]
285 << " metre." << endl;
286 }
287 }
288 }
289}
290
291
292void Foam::refinementFeatures::buildTrees(const label featI)
293{
294 const extendedEdgeMesh& eMesh = operator[](featI);
295 const pointField& points = eMesh.points();
296 const edgeList& edges = eMesh.edges();
297
298 // Calculate bb of all points
300
301 // Random number generator. Bit dodgy since not exactly random ;-)
302 Random rndGen(65431);
303
304 // Slightly extended bb. Slightly off-centred just so on symmetric
305 // geometry there are less face/edge aligned items.
306 bb.inflate(rndGen, 1e-4, ROOTVSMALL);
307
308 edgeTrees_.set
309 (
310 featI,
312 (
313 treeDataEdge(edges, points), // All edges
314
315 bb, // overall search domain
316 8, // maxLevel
317 10, // leafsize
318 3.0 // duplicity
319 )
320 );
321
322
323 labelList featurePoints(identity(eMesh.nonFeatureStart()));
324
325 pointTrees_.set
326 (
327 featI,
329 (
330 treeDataPoint(points, featurePoints),
331
332 bb, // overall search domain
333 8, // maxLevel
334 10, // leafsize
335 3.0 // duplicity
336 )
337 );
338}
339
340
341// Find maximum level of a shell.
342void Foam::refinementFeatures::findHigherLevel
343(
344 const pointField& pt,
345 const label featI,
346 labelList& maxLevel
347) const
348{
349 const labelList& levels = levels_[featI];
350
351 const scalarField& distances = distances_[featI];
352
353 // Collect all those points that have a current maxLevel less than
354 // (any of) the shell. Also collect the furthest distance allowable
355 // to any shell with a higher level.
356
357 pointField candidates(pt.size());
358 labelList candidateMap(pt.size());
359 scalarField candidateDistSqr(pt.size());
360 label candidateI = 0;
361
362 forAll(maxLevel, pointi)
363 {
364 forAllReverse(levels, levelI)
365 {
366 if (levels[levelI] > maxLevel[pointi])
367 {
368 candidates[candidateI] = pt[pointi];
369 candidateMap[candidateI] = pointi;
370 candidateDistSqr[candidateI] = sqr(distances[levelI]);
371 candidateI++;
372 break;
373 }
374 }
375 }
376 candidates.setSize(candidateI);
377 candidateMap.setSize(candidateI);
378 candidateDistSqr.setSize(candidateI);
379
380 // Do the expensive nearest test only for the candidate points.
381 const indexedOctree<treeDataEdge>& tree = edgeTrees_[featI];
382
383 List<pointIndexHit> nearInfo(candidates.size());
384 forAll(candidates, candidateI)
385 {
386 nearInfo[candidateI] = tree.findNearest
387 (
388 candidates[candidateI],
389 candidateDistSqr[candidateI]
390 );
391 }
392
393 // Update maxLevel
394 forAll(nearInfo, candidateI)
395 {
396 if (nearInfo[candidateI].hit())
397 {
398 // Check which level it actually is in.
399 label minDistI = findLower
400 (
401 distances,
402 nearInfo[candidateI].point().dist(candidates[candidateI])
403 );
404
405 label pointi = candidateMap[candidateI];
406
407 // pt is inbetween shell[minDistI] and shell[minDistI+1]
408 maxLevel[pointi] = levels[minDistI+1];
409 }
411}
412
413
414// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
415
418{
419 if (!regionEdgeTreesPtr_)
420 {
421 regionEdgeTreesPtr_.reset
422 (
423 new PtrList<indexedOctree<treeDataEdge>>(size())
424 );
425 PtrList<indexedOctree<treeDataEdge>>& trees = regionEdgeTreesPtr_();
426
427 forAll(*this, featI)
428 {
429 const extendedEdgeMesh& eMesh = operator[](featI);
430 const pointField& points = eMesh.points();
431 const edgeList& edges = eMesh.edges();
432
433 // Calculate bb of all points
434 treeBoundBox bb(points);
435
436 // Random number generator. Bit dodgy since not exactly random ;-)
437 Random rndGen(65431);
438
439 // Slightly extended bb. Slightly off-centred just so on symmetric
440 // geometry there are less face/edge aligned items.
441 bb.inflate(rndGen, 1e-4, ROOTVSMALL);
442
443 trees.set
444 (
445 featI,
446 new indexedOctree<treeDataEdge>
447 (
448 treeDataEdge(edges, points, eMesh.regionEdges()),
449
450 bb, // overall search domain
451 8, // maxLevel
452 10, // leafsize
453 3.0 // duplicity
454 )
455 );
456 }
457 }
459 return *regionEdgeTreesPtr_;
460}
461
462
463// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
464
466(
467 const objectRegistry& io,
468 const PtrList<dictionary>& featDicts,
469 const bool dryRun
470)
471:
472 PtrList<extendedFeatureEdgeMesh>(featDicts.size()),
473 distances_(featDicts.size()),
474 levels_(featDicts.size()),
475 edgeTrees_(featDicts.size()),
476 pointTrees_(featDicts.size()),
477 dryRun_(dryRun)
478{
479 // Read features
480 read(io, featDicts);
481
482 // Search engines
483 forAll(*this, i)
484 {
485 buildTrees(i);
486 }
487}
488
489
490//Foam::refinementFeatures::refinementFeatures
491//(
492// const objectRegistry& io,
493// const PtrList<dictionary>& featDicts,
494// const scalar minCos
495//)
496//:
497// PtrList<extendedFeatureEdgeMesh>(featDicts.size()),
498// distances_(featDicts.size()),
499// levels_(featDicts.size()),
500// edgeTrees_(featDicts.size()),
501// pointTrees_(featDicts.size())
502//{
503// // Read features
504// read(io, featDicts);
505//
506// // Search engines
507// forAll(*this, i)
508// {
509// const edgeMesh& eMesh = operator[](i);
510// const pointField& points = eMesh.points();
511// const edgeList& edges = eMesh.edges();
512// const labelListList& pointEdges = eMesh.pointEdges();
513//
514// DynamicList<label> featurePoints;
515// forAll(pointEdges, pointi)
516// {
517// const labelList& pEdges = pointEdges[pointi];
518// if (pEdges.size() > 2)
519// {
520// featurePoints.append(pointi);
521// }
522// else if (pEdges.size() == 2)
523// {
524// // Check the angle
525// const edge& e0 = edges[pEdges[0]];
526// const edge& e1 = edges[pEdges[1]];
527//
528// const point& p = points[pointi];
529// const point& p0 = points[e0.otherVertex(pointi)];
530// const point& p1 = points[e1.otherVertex(pointi)];
531//
532// vector v0 = p-p0;
533// scalar v0Mag = mag(v0);
534//
535// vector v1 = p1-p;
536// scalar v1Mag = mag(v1);
537//
538// if
539// (
540// v0Mag > SMALL
541// && v1Mag > SMALL
542// && ((v0/v0Mag & v1/v1Mag) < minCos)
543// )
544// {
545// featurePoints.append(pointi);
546// }
547// }
548// }
549//
550// Info<< "Detected " << featurePoints.size()
551// << " featurePoints out of " << points.size()
552// << " points on feature " << i //eMesh.name()
553// << " when using feature cos " << minCos << endl;
554//
555// buildTrees(i, featurePoints);
556// }
557//}
558
559
560// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
561
563(
564 const scalar maxRatio,
565 const boundBox& meshBb,
566 const bool report,
567 Ostream& os
568) const
569{
570 if (report)
571 {
572 os<< "Checking for size." << endl;
573 }
574
575 bool hasError = false;
576
577 forAll(*this, i)
578 {
579 const extendedFeatureEdgeMesh& em = operator[](i);
580 const boundBox bb(em.points(), true);
581
582 for (label j = i+1; j < size(); j++)
583 {
584 const extendedFeatureEdgeMesh& em2 = operator[](j);
585 const boundBox bb2(em2.points(), true);
586
587 scalar ratio = bb.mag()/bb2.mag();
588
589 if (ratio > maxRatio || ratio < 1.0/maxRatio)
590 {
591 hasError = true;
592
593 if (report)
594 {
595 os << " " << em.name()
596 << " bounds differ from " << em2.name()
597 << " by more than a factor 100:" << nl
598 << " bounding box : " << bb << nl
599 << " bounding box : " << bb2
600 << endl;
601 }
602 }
603 }
604 }
605
606 forAll(*this, i)
607 {
608 const extendedFeatureEdgeMesh& em = operator[](i);
609 const boundBox bb(em.points(), true);
610 if (!meshBb.contains(bb))
611 {
612 if (report)
613 {
614 os << " " << em.name()
615 << " bounds not fully contained in mesh"<< nl
616 << " bounding box : " << bb << nl
617 << " mesh bounding box : " << meshBb
618 << endl;
619 }
620 }
621 }
622
623 if (report)
624 {
626 }
627
628 return returnReduceOr(hasError);
629}
630
631
633(
634 const pointField& samples,
635 const scalarField& nearestDistSqr,
636 labelList& nearFeature,
637 List<pointIndexHit>& nearInfo,
638 vectorField& nearNormal
639) const
640{
641 nearFeature.setSize(samples.size());
642 nearFeature = -1;
643 nearInfo.setSize(samples.size());
644 nearInfo = pointIndexHit();
645 nearNormal.setSize(samples.size());
646 nearNormal = Zero;
647
648 forAll(edgeTrees_, featI)
649 {
650 const indexedOctree<treeDataEdge>& tree = edgeTrees_[featI];
651 const treeDataEdge& treeData = tree.shapes();
652
653 if (!treeData.empty())
654 {
655 forAll(samples, sampleI)
656 {
657 const point& sample = samples[sampleI];
658
659 scalar distSqr;
660 if (nearInfo[sampleI].hit())
661 {
662 distSqr = nearInfo[sampleI].point().distSqr(sample);
663 }
664 else
665 {
666 distSqr = nearestDistSqr[sampleI];
667 }
668
669 pointIndexHit info = tree.findNearest(sample, distSqr);
670
671 if (info.hit())
672 {
673 nearFeature[sampleI] = featI;
674 nearInfo[sampleI] = pointIndexHit
675 (
676 info.hit(),
677 info.point(),
678 treeData.objectIndex(info.index())
679 );
680 nearNormal[sampleI] = treeData.line(info.index()).unitVec();
682 }
683 }
684 }
685}
686
687
689(
690 const pointField& samples,
691 const scalarField& nearestDistSqr,
692 labelList& nearFeature,
693 List<pointIndexHit>& nearInfo,
694 vectorField& nearNormal
695) const
696{
697 nearFeature.setSize(samples.size());
698 nearFeature = -1;
699 nearInfo.setSize(samples.size());
700 nearInfo = pointIndexHit();
701 nearNormal.setSize(samples.size());
702 nearNormal = Zero;
703
704
705 const PtrList<indexedOctree<treeDataEdge>>& regionTrees =
706 regionEdgeTrees();
707
708 forAll(regionTrees, featI)
709 {
710 const indexedOctree<treeDataEdge>& regionTree = regionTrees[featI];
711 const treeDataEdge& treeData = regionTree.shapes();
712
713 forAll(samples, sampleI)
714 {
715 const point& sample = samples[sampleI];
716
717 scalar distSqr;
718 if (nearInfo[sampleI].hit())
719 {
720 distSqr = nearInfo[sampleI].point().distSqr(sample);
721 }
722 else
723 {
724 distSqr = nearestDistSqr[sampleI];
725 }
726
727 // Find anything closer than current best
728 pointIndexHit info = regionTree.findNearest(sample, distSqr);
729
730 if (info.hit())
731 {
732 nearFeature[sampleI] = featI;
733 nearInfo[sampleI] = pointIndexHit
734 (
735 info.hit(),
736 info.point(),
737 treeData.objectIndex(info.index())
738 );
739 nearNormal[sampleI] = treeData.line(info.index()).unitVec();
740 }
741 }
742 }
743}
744
745
746//void Foam::refinementFeatures::findNearestPoint
747//(
748// const pointField& samples,
749// const scalarField& nearestDistSqr,
750// labelList& nearFeature,
751// labelList& nearIndex
752//) const
753//{
754// nearFeature.setSize(samples.size());
755// nearFeature = -1;
756// nearIndex.setSize(samples.size());
757// nearIndex = -1;
758//
759// forAll(pointTrees_, featI)
760// {
761// const indexedOctree<treeDataPoint>& tree = pointTrees_[featI];
762//
763// if (!tree.shapes().empty())
764// {
765// forAll(samples, sampleI)
766// {
767// const point& sample = samples[sampleI];
768//
769// scalar distSqr;
770// if (nearFeature[sampleI] != -1)
771// {
772// const nearTree = pointTrees_[nearFeature[sampleI]];
773// distSqr = sample.distSqr
774// (
775// nearTree.shapes()[nearIndex[sampleI]]
776// );
777// }
778// else
779// {
780// distSqr = nearestDistSqr[sampleI];
781// }
782//
783// pointIndexHit info = tree.findNearest(sample, distSqr);
784//
785// if (info.hit())
786// {
787// nearFeature[sampleI] = featI;
788// nearIndex[sampleI] = info.index();
789// }
790// }
791// }
792// }
793//}
794
795
797(
798 const pointField& samples,
799 const scalarField& nearestDistSqr,
800 labelList& nearFeature,
801 List<pointIndexHit>& nearInfo
802) const
803{
804 nearFeature.setSize(samples.size());
805 nearFeature = -1;
806 nearInfo.setSize(samples.size());
807 nearInfo = pointIndexHit();
808
809 forAll(pointTrees_, featI)
810 {
811 const indexedOctree<treeDataPoint>& tree = pointTrees_[featI];
812 const auto& treeData = tree.shapes();
813
814 if (!treeData.empty())
815 {
816 forAll(samples, sampleI)
817 {
818 const point& sample = samples[sampleI];
819
820 scalar distSqr;
821 if (nearFeature[sampleI] != -1)
822 {
823 distSqr = nearInfo[sampleI].point().distSqr(sample);
824 }
825 else
826 {
827 distSqr = nearestDistSqr[sampleI];
828 }
829
830 pointIndexHit info = tree.findNearest(sample, distSqr);
831
832 if (info.hit())
833 {
834 nearFeature[sampleI] = featI;
835 nearInfo[sampleI] = pointIndexHit
836 (
837 info.hit(),
838 info.point(),
839 treeData.objectIndex(info.index())
840 );
842 }
843 }
844 }
845}
846
847
848void Foam::refinementFeatures::findHigherLevel
849(
850 const pointField& pt,
851 const labelList& ptLevel,
852 labelList& maxLevel
853) const
854{
855 // Maximum level of any shell. Start off with level of point.
856 maxLevel = ptLevel;
857
858 forAll(*this, featI)
859 {
860 findHigherLevel(pt, featI, maxLevel);
861 }
862}
863
864
865Foam::scalar Foam::refinementFeatures::maxDistance() const
866{
867 scalar overallMax = -GREAT;
868 forAll(distances_, featI)
869 {
870 overallMax = max(overallMax, max(distances_[featI]));
871 }
872 return overallMax;
873}
874
875
876// ************************************************************************* //
Minimal example by using system/controlDict.functions:
@ NO_REGISTER
Do not request registration (bool: false).
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
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 setSize(label n)
Alias for resize().
Definition List.H:536
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.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
const extendedFeatureEdgeMesh * set(const label i) const
Definition PtrList.H:171
Random number generator.
Definition Random.H:56
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
const extendedFeatureEdgeMesh & operator[](const label i) const
Definition UPtrListI.H:289
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
scalar mag() const
The magnitude/length of the bounding box diagonal.
Definition boundBoxI.H:198
void inflate(const scalar factor)
Expand box by factor*mag(span) in all dimensions.
Definition boundBoxI.H:381
const pointField & points() const noexcept
Return points.
Definition edgeMeshI.H:92
const edgeList & edges() const noexcept
Return edges.
Definition edgeMeshI.H:98
static autoPtr< edgeMesh > New(const fileName &name, const word &fileType)
Read construct from filename with given format.
Definition edgeMeshNew.C:27
Description of feature edges and points.
const labelList & regionEdges() const
Return the feature edges which are on the boundary between.
static autoPtr< extendedEdgeMesh > New(const fileName &name, const word &fileType)
Select constructed from filename with given file format.
static bool canRead(const fileName &name, bool verbose=false)
Can we read this file format?
static const fileName null
An empty fileName.
Definition fileName.H:111
Non-pointer based hierarchical recursive searching.
const Type & shapes() const noexcept
Reference to shape.
@ REGEX
Regular expression.
Definition keyType.H:83
Point unitVec() const
Return the unit vector (start-to-end).
Definition lineI.H:135
static Type get(const dictionary &dict, const word &keyword, const bool noExit, enum keyType::option matchOpt=keyType::REGEX, const Type &deflt=Zero)
Wrapper around dictionary::get which does not exit.
Registry of regIOobjects.
void findNearestPoint(const pointField &samples, const scalarField &nearestDistSqr, labelList &nearFeature, List< pointIndexHit > &nearInfo) const
Find nearest feature point. Sets.
bool checkSizes(const scalar maxRatio, const boundBox &meshBb, const bool report, Ostream &os) const
Check sizes - return true if error and stream to os.
scalar maxDistance() const
Highest distance of all features.
void findNearestRegionEdge(const pointField &samples, const scalarField &nearestDistSqr, labelList &nearFeature, List< pointIndexHit > &nearInfo, vectorField &nearNormal) const
Find nearest point on nearest region edge. Sets.
void findNearestEdge(const pointField &samples, const scalarField &nearestDistSqr, labelList &nearFeature, List< pointIndexHit > &nearInfo, vectorField &nearNormal) const
Find nearest point on nearest feature edge. Sets.
const PtrList< indexedOctree< treeDataEdge > > & regionEdgeTrees() const
refinementFeatures(const objectRegistry &io, const PtrList< dictionary > &featDicts, const bool dryRun=false)
Construct from description.
Standard boundBox with extra functionality for use in octree.
Holds data for octree to work on an edges subset.
bool empty() const noexcept
Is the effective edge selection empty?
const linePointRef line(const label index) const
Geometric line for edge at specified shape index. Frequently used.
label objectIndex(const label index) const
Map from shape index to original (non-subset) edge label.
Holds (reference to) pointField. Encapsulation of data needed for octree searches....
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const auto & io
const pointField & points
List< edge > edgeList
List of edge.
Definition edgeList.H:32
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
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
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
messageStream Info
Information stream (stdout output on master, null elsewhere).
label findLower(const ListType &input, const T &val, const label start, const ComparePredicate &comp)
Binary search to find the index of the last element in a sorted list that is less than value.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition Ostream.H:490
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
vectorField pointField
pointField is a vectorField.
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition Ostream.H:499
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
List< treeBoundBox > meshBb(1, treeBoundBox(coarseMesh.points()).extend(rndGen, 1e-3))
Tree tree(triangles.begin(), triangles.end())
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition stdFoam.H:315
scalarField samples(nIntervals, Zero)
Random rndGen