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 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
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 (
50 meshRefinement::get<fileName>
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
71 false
72 );
73
74 const fileName fName(typeFilePath<extendedFeatureEdgeMesh>(extFeatObj));
75
76 if (!fName.empty() && extendedEdgeMesh::canRead(fName))
77 {
78 autoPtr<extendedEdgeMesh> eMeshPtr = extendedEdgeMesh::New
79 (
80 fName
81 );
82
83 if (!dryRun_)
84 {
85 Info<< "Read extendedFeatureEdgeMesh " << extFeatObj.name()
86 << nl << incrIndent;
87 eMeshPtr().writeStats(Info);
88 Info<< decrIndent << endl;
89 }
90
91 set(featI, new extendedFeatureEdgeMesh(extFeatObj, eMeshPtr()));
92 }
93 else
94 {
95 // Try reading edgeMesh
96
97 IOobject featObj
98 (
99 featFileName, // name
100 io.time().constant(), // instance
101 "triSurface", // local
102 io.time(), // registry
105 false
106 );
107
108 const fileName fName(typeFilePath<featureEdgeMesh>(featObj));
109
110 if (fName.empty())
111 {
113 << "Could not open " << featObj.objectPath()
114 << exit(FatalIOError);
115 }
116
117 // Read as edgeMesh
118 autoPtr<edgeMesh> eMeshPtr = edgeMesh::New(fName);
119 const edgeMesh& eMesh = eMeshPtr();
120
121 if (!dryRun_)
122 {
123 Info<< "Read edgeMesh " << featObj.name() << nl
124 << incrIndent;
125 eMesh.writeStats(Info);
126 Info<< decrIndent << endl;
127 }
128
129 // Analyse for feature points. These are all classified as mixed
130 // points for lack of anything better
131 const labelListList& pointEdges = eMesh.pointEdges();
132
133 labelList oldToNew(eMesh.points().size(), -1);
134 DynamicField<point> newPoints(eMesh.points().size());
135 forAll(pointEdges, pointi)
136 {
137 if (pointEdges[pointi].size() > 2)
138 {
139 oldToNew[pointi] = newPoints.size();
140 newPoints.append(eMesh.points()[pointi]);
141 }
142 //else if (pointEdges[pointi].size() == 2)
143 //MEJ: do something based on a feature angle?
144 }
145 label nFeatures = newPoints.size();
146 forAll(oldToNew, pointi)
147 {
148 if (oldToNew[pointi] == -1)
149 {
150 oldToNew[pointi] = newPoints.size();
151 newPoints.append(eMesh.points()[pointi]);
152 }
153 }
154
155
156 const edgeList& edges = eMesh.edges();
157 edgeList newEdges(edges.size());
158 forAll(edges, edgeI)
159 {
160 const edge& e = edges[edgeI];
161 newEdges[edgeI] = edge
162 (
163 oldToNew[e[0]],
164 oldToNew[e[1]]
165 );
166 }
167
168 // Construct an extendedEdgeMesh with
169 // - all points on more than 2 edges : mixed feature points
170 // - all edges : external edges
171
172 extendedEdgeMesh eeMesh
173 (
174 newPoints, // pts
175 newEdges, // eds
176 0, // (point) concaveStart
177 0, // (point) mixedStart
178 nFeatures, // (point) nonFeatureStart
179 edges.size(), // (edge) internalStart
180 edges.size(), // (edge) flatStart
181 edges.size(), // (edge) openStart
182 edges.size(), // (edge) multipleStart
183 vectorField(0), // normals
184 List<extendedEdgeMesh::sideVolumeType>(0),// normalVolumeTypes
185 vectorField(0), // edgeDirections
186 labelListList(0), // normalDirections
187 labelListList(0), // edgeNormals
188 labelListList(0), // featurePointNormals
189 labelListList(0), // featurePointEdges
190 identity(newEdges.size()) // regionEdges
191 );
192
193 //Info<< "Constructed extendedFeatureEdgeMesh " << featObj.name()
194 // << nl << incrIndent;
195 //eeMesh.writeStats(Info);
196 //Info<< decrIndent << endl;
197
198 set(featI, new extendedFeatureEdgeMesh(featObj, eeMesh));
199 }
200
201 const extendedEdgeMesh& eMesh = operator[](featI);
202
203 if (dict.found("levels"))
204 {
205 List<Tuple2<scalar, label>> distLevels(dict.lookup("levels"));
206
207 if (dict.size() < 1)
208 {
210 << " : levels should be at least size 1" << endl
211 << "levels : " << dict.lookup("levels")
212 << exit(FatalError);
213 }
214
215 distances_[featI].setSize(distLevels.size());
216 levels_[featI].setSize(distLevels.size());
217
218 forAll(distLevels, j)
219 {
220 distances_[featI][j] = distLevels[j].first();
221 levels_[featI][j] = distLevels[j].second();
222
223 if (levels_[featI][j] < 0)
224 {
226 << "Feature " << featFileName
227 << " has illegal refinement level " << levels_[featI][j]
228 << exit(FatalError);
229 }
230
231 // Check in incremental order
232 if (j > 0)
233 {
234 if
235 (
236 (distances_[featI][j] <= distances_[featI][j-1])
237 || (levels_[featI][j] > levels_[featI][j-1])
238 )
239 {
241 << " : Refinement should be specified in order"
242 << " of increasing distance"
243 << " (and decreasing refinement level)." << endl
244 << "Distance:" << distances_[featI][j]
245 << " refinementLevel:" << levels_[featI][j]
246 << exit(FatalError);
247 }
248 }
249 }
250 }
251 else
252 {
253 // Look up 'level' for single level
254 levels_[featI] =
256 (
257 1,
258 meshRefinement::get<label>
259 (
260 dict,
261 "level",
262 dryRun_,
264 0
265 )
266 );
267 distances_[featI] = scalarField(1, Zero);
268 }
269
270 if (!dryRun_)
271 {
272 Info<< "Refinement level according to distance to "
273 << featFileName << " (" << eMesh.points().size() << " points, "
274 << eMesh.edges().size() << " edges)." << endl;
275 forAll(levels_[featI], j)
276 {
277 Info<< " level " << levels_[featI][j]
278 << " for all cells within " << distances_[featI][j]
279 << " metre." << endl;
280 }
281 }
282 }
283}
284
285
286void Foam::refinementFeatures::buildTrees(const label featI)
287{
288 const extendedEdgeMesh& eMesh = operator[](featI);
289 const pointField& points = eMesh.points();
290 const edgeList& edges = eMesh.edges();
291
292 // Calculate bb of all points
293 treeBoundBox bb(points);
294
295 // Random number generator. Bit dodgy since not exactly random ;-)
296 Random rndGen(65431);
297
298 // Slightly extended bb. Slightly off-centred just so on symmetric
299 // geometry there are less face/edge aligned items.
300 bb = bb.extend(rndGen, 1e-4);
301 bb.min() -= point::uniform(ROOTVSMALL);
302 bb.max() += point::uniform(ROOTVSMALL);
303
304 edgeTrees_.set
305 (
306 featI,
307 new indexedOctree<treeDataEdge>
308 (
309 treeDataEdge
310 (
311 false, // do not cache bb
312 edges,
313 points,
314 identity(edges.size())
315 ),
316 bb, // overall search domain
317 8, // maxLevel
318 10, // leafsize
319 3.0 // duplicity
320 )
321 );
322
323
324 labelList featurePoints(identity(eMesh.nonFeatureStart()));
325
326 pointTrees_.set
327 (
328 featI,
329 new indexedOctree<treeDataPoint>
330 (
331 treeDataPoint(points, featurePoints),
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 mag(nearInfo[candidateI].hitPoint()-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 }
410 }
411}
412
413
414// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
415
418{
419 if (!regionEdgeTreesPtr_)
420 {
421 regionEdgeTreesPtr_.reset
422 (
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
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 = bb.extend(rndGen, 1e-4);
442 bb.min() -= point::uniform(ROOTVSMALL);
443 bb.max() += point::uniform(ROOTVSMALL);
444
445 trees.set
446 (
447 featI,
449 (
451 (
452 false, // do not cache bb
453 edges,
454 points,
455 eMesh.regionEdges()
456 ),
457 bb, // overall search domain
458 8, // maxLevel
459 10, // leafsize
460 3.0 // duplicity
461 )
462 );
463 }
464 }
465
466 return *regionEdgeTreesPtr_;
467}
468
469
470// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
471
473(
474 const objectRegistry& io,
475 const PtrList<dictionary>& featDicts,
476 const bool dryRun
477)
478:
479 PtrList<extendedFeatureEdgeMesh>(featDicts.size()),
480 distances_(featDicts.size()),
481 levels_(featDicts.size()),
482 edgeTrees_(featDicts.size()),
483 pointTrees_(featDicts.size()),
484 dryRun_(dryRun)
485{
486 // Read features
487 read(io, featDicts);
488
489 // Search engines
490 forAll(*this, i)
491 {
492 buildTrees(i);
493 }
494}
495
496
497//Foam::refinementFeatures::refinementFeatures
498//(
499// const objectRegistry& io,
500// const PtrList<dictionary>& featDicts,
501// const scalar minCos
502//)
503//:
504// PtrList<extendedFeatureEdgeMesh>(featDicts.size()),
505// distances_(featDicts.size()),
506// levels_(featDicts.size()),
507// edgeTrees_(featDicts.size()),
508// pointTrees_(featDicts.size())
509//{
510// // Read features
511// read(io, featDicts);
512//
513// // Search engines
514// forAll(*this, i)
515// {
516// const edgeMesh& eMesh = operator[](i);
517// const pointField& points = eMesh.points();
518// const edgeList& edges = eMesh.edges();
519// const labelListList& pointEdges = eMesh.pointEdges();
520//
521// DynamicList<label> featurePoints;
522// forAll(pointEdges, pointi)
523// {
524// const labelList& pEdges = pointEdges[pointi];
525// if (pEdges.size() > 2)
526// {
527// featurePoints.append(pointi);
528// }
529// else if (pEdges.size() == 2)
530// {
531// // Check the angle
532// const edge& e0 = edges[pEdges[0]];
533// const edge& e1 = edges[pEdges[1]];
534//
535// const point& p = points[pointi];
536// const point& p0 = points[e0.otherVertex(pointi)];
537// const point& p1 = points[e1.otherVertex(pointi)];
538//
539// vector v0 = p-p0;
540// scalar v0Mag = mag(v0);
541//
542// vector v1 = p1-p;
543// scalar v1Mag = mag(v1);
544//
545// if
546// (
547// v0Mag > SMALL
548// && v1Mag > SMALL
549// && ((v0/v0Mag & v1/v1Mag) < minCos)
550// )
551// {
552// featurePoints.append(pointi);
553// }
554// }
555// }
556//
557// Info<< "Detected " << featurePoints.size()
558// << " featurePoints out of " << points.size()
559// << " points on feature " << i //eMesh.name()
560// << " when using feature cos " << minCos << endl;
561//
562// buildTrees(i, featurePoints);
563// }
564//}
565
566
567// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
568
570(
571 const scalar maxRatio,
572 const boundBox& meshBb,
573 const bool report,
574 Ostream& os
575) const
576{
577 if (report)
578 {
579 os<< "Checking for size." << endl;
580 }
581
582 bool hasError = false;
583
584 forAll(*this, i)
585 {
586 const extendedFeatureEdgeMesh& em = operator[](i);
587 const boundBox bb(em.points(), true);
588
589 for (label j = i+1; j < size(); j++)
590 {
591 const extendedFeatureEdgeMesh& em2 = operator[](j);
592 const boundBox bb2(em2.points(), true);
593
594 scalar ratio = bb.mag()/bb2.mag();
595
596 if (ratio > maxRatio || ratio < 1.0/maxRatio)
597 {
598 hasError = true;
599
600 if (report)
601 {
602 os << " " << em.name()
603 << " bounds differ from " << em2.name()
604 << " by more than a factor 100:" << nl
605 << " bounding box : " << bb << nl
606 << " bounding box : " << bb2
607 << endl;
608 }
609 }
610 }
611 }
612
613 forAll(*this, i)
614 {
615 const extendedFeatureEdgeMesh& em = operator[](i);
616 const boundBox bb(em.points(), true);
617 if (!meshBb.contains(bb))
618 {
619 if (report)
620 {
621 os << " " << em.name()
622 << " bounds not fully contained in mesh"<< nl
623 << " bounding box : " << bb << nl
624 << " mesh bounding box : " << meshBb
625 << endl;
626 }
627 }
628 }
629
630 if (report)
631 {
632 os<< endl;
633 }
634
635 return returnReduce(hasError, orOp<bool>());
636}
637
638
640(
641 const pointField& samples,
642 const scalarField& nearestDistSqr,
643 labelList& nearFeature,
644 List<pointIndexHit>& nearInfo,
645 vectorField& nearNormal
646) const
647{
648 nearFeature.setSize(samples.size());
649 nearFeature = -1;
650 nearInfo.setSize(samples.size());
651 nearInfo = pointIndexHit();
652 nearNormal.setSize(samples.size());
653 nearNormal = Zero;
654
655 forAll(edgeTrees_, featI)
656 {
657 const indexedOctree<treeDataEdge>& tree = edgeTrees_[featI];
658
659 if (tree.shapes().size() > 0)
660 {
661 forAll(samples, sampleI)
662 {
663 const point& sample = samples[sampleI];
664
665 scalar distSqr;
666 if (nearInfo[sampleI].hit())
667 {
668 distSqr = magSqr(nearInfo[sampleI].hitPoint()-sample);
669 }
670 else
671 {
672 distSqr = nearestDistSqr[sampleI];
673 }
674
675 pointIndexHit info = tree.findNearest(sample, distSqr);
676
677 if (info.hit())
678 {
679 nearFeature[sampleI] = featI;
680 nearInfo[sampleI] = pointIndexHit
681 (
682 info.hit(),
683 info.hitPoint(),
684 tree.shapes().edgeLabels()[info.index()]
685 );
686
687 const treeDataEdge& td = tree.shapes();
688 const edge& e = td.edges()[nearInfo[sampleI].index()];
689
690 nearNormal[sampleI] = e.unitVec(td.points());
691 }
692 }
693 }
694 }
695}
696
697
699(
700 const pointField& samples,
701 const scalarField& nearestDistSqr,
702 labelList& nearFeature,
703 List<pointIndexHit>& nearInfo,
704 vectorField& nearNormal
705) const
706{
707 nearFeature.setSize(samples.size());
708 nearFeature = -1;
709 nearInfo.setSize(samples.size());
710 nearInfo = pointIndexHit();
711 nearNormal.setSize(samples.size());
712 nearNormal = Zero;
713
714
715 const PtrList<indexedOctree<treeDataEdge>>& regionTrees =
716 regionEdgeTrees();
717
718 forAll(regionTrees, featI)
719 {
720 const indexedOctree<treeDataEdge>& regionTree = regionTrees[featI];
721
722 forAll(samples, sampleI)
723 {
724 const point& sample = samples[sampleI];
725
726 scalar distSqr;
727 if (nearInfo[sampleI].hit())
728 {
729 distSqr = magSqr(nearInfo[sampleI].hitPoint()-sample);
730 }
731 else
732 {
733 distSqr = nearestDistSqr[sampleI];
734 }
735
736 // Find anything closer than current best
737 pointIndexHit info = regionTree.findNearest(sample, distSqr);
738
739 if (info.hit())
740 {
741 const treeDataEdge& td = regionTree.shapes();
742
743 nearFeature[sampleI] = featI;
744 nearInfo[sampleI] = pointIndexHit
745 (
746 info.hit(),
747 info.hitPoint(),
748 regionTree.shapes().edgeLabels()[info.index()]
749 );
750
751 const edge& e = td.edges()[nearInfo[sampleI].index()];
752
753 nearNormal[sampleI] = e.unitVec(td.points());
754 }
755 }
756 }
757}
758
759
760//void Foam::refinementFeatures::findNearestPoint
761//(
762// const pointField& samples,
763// const scalarField& nearestDistSqr,
764// labelList& nearFeature,
765// labelList& nearIndex
766//) const
767//{
768// nearFeature.setSize(samples.size());
769// nearFeature = -1;
770// nearIndex.setSize(samples.size());
771// nearIndex = -1;
772//
773// forAll(pointTrees_, featI)
774// {
775// const indexedOctree<treeDataPoint>& tree = pointTrees_[featI];
776//
777// if (tree.shapes().pointLabels().size() > 0)
778// {
779// forAll(samples, sampleI)
780// {
781// const point& sample = samples[sampleI];
782//
783// scalar distSqr;
784// if (nearFeature[sampleI] != -1)
785// {
786// label nearFeatI = nearFeature[sampleI];
787// const indexedOctree<treeDataPoint>& nearTree =
788// pointTrees_[nearFeatI];
789// label featPointi =
790// nearTree.shapes().pointLabels()[nearIndex[sampleI]];
791// const point& featPt =
792// operator[](nearFeatI).points()[featPointi];
793// distSqr = magSqr(featPt-sample);
794// }
795// else
796// {
797// distSqr = nearestDistSqr[sampleI];
798// }
799//
800// pointIndexHit info = tree.findNearest(sample, distSqr);
801//
802// if (info.hit())
803// {
804// nearFeature[sampleI] = featI;
805// nearIndex[sampleI] = info.index();
806// }
807// }
808// }
809// }
810//}
811
812
814(
815 const pointField& samples,
816 const scalarField& nearestDistSqr,
817 labelList& nearFeature,
818 List<pointIndexHit>& nearInfo
819) const
820{
821 nearFeature.setSize(samples.size());
822 nearFeature = -1;
823 nearInfo.setSize(samples.size());
824 nearInfo = pointIndexHit();
825
826 forAll(pointTrees_, featI)
827 {
828 const indexedOctree<treeDataPoint>& tree = pointTrees_[featI];
829
830 if (tree.shapes().pointLabels().size() > 0)
831 {
832 forAll(samples, sampleI)
833 {
834 const point& sample = samples[sampleI];
835
836 scalar distSqr;
837 if (nearFeature[sampleI] != -1)
838 {
839 distSqr = magSqr(nearInfo[sampleI].hitPoint()-sample);
840 }
841 else
842 {
843 distSqr = nearestDistSqr[sampleI];
844 }
845
846 pointIndexHit info = tree.findNearest(sample, distSqr);
847
848 if (info.hit())
849 {
850 nearFeature[sampleI] = featI;
851 nearInfo[sampleI] = pointIndexHit
852 (
853 info.hit(),
854 info.hitPoint(),
855 tree.shapes().pointLabels()[info.index()]
856 );
857 }
858 }
859 }
860 }
861}
862
863
864void Foam::refinementFeatures::findHigherLevel
865(
866 const pointField& pt,
867 const labelList& ptLevel,
868 labelList& maxLevel
869) const
870{
871 // Maximum level of any shell. Start off with level of point.
872 maxLevel = ptLevel;
873
874 forAll(*this, featI)
875 {
876 findHigherLevel(pt, featI, maxLevel);
877 }
878}
879
880
882{
883 scalar overallMax = -GREAT;
884 forAll(distances_, featI)
885 {
886 overallMax = max(overallMax, max(distances_[featI]));
887 }
888 return overallMax;
889}
890
891
892// ************************************************************************* //
Minimal example by using system/controlDict.functions:
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
const Time & time() const
Return Time associated with the objectRegistry.
Definition: IOobject.C:506
void setSize(const label n)
Alias for resize()
Definition: List.H:218
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
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.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
const extendedFeatureEdgeMesh * set(const label i) const
Definition: PtrList.H:138
virtual bool read()
Re-read model coefficients if they have changed.
Random number generator.
Definition: Random.H:60
const word & constant() const
Return constant name.
Definition: TimePathsI.H:96
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
T & first()
Return the first element of the list.
Definition: UListI.H:202
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
const T & operator[](const label i) const
Return const reference to the element.
Definition: UPtrListI.H:234
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:64
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:133
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:386
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
Particle-size distribution model wherein random samples are drawn from the doubly-truncated uniform p...
Definition: uniform.H:164
const pointField & points() const noexcept
Return points.
Definition: edgeMeshI.H:99
const edgeList & edges() const noexcept
Return edges.
Definition: edgeMeshI.H:105
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
Description of feature edges and points.
const labelList & regionEdges() const
Return the feature edges which are on the boundary between.
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:102
void checkSizes() const
Check that all components of sizes() are non-negative.
Non-pointer based hierarchical recursive searching.
Definition: indexedOctree.H:74
const Type & shapes() const
Reference to shape.
@ REGEX
Regular expression.
Definition: keyType.H:82
Registry of regIOobjects.
Encapsulates queries for features.
void findNearestPoint(const pointField &samples, const scalarField &nearestDistSqr, labelList &nearFeature, List< pointIndexHit > &nearInfo) const
Find nearest feature point. Sets.
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
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:89
treeBoundBox extend(Random &rndGen, const scalar s) const
Return slightly wider bounding box.
Holds data for octree to work on an edges subset.
Definition: treeDataEdge.H:57
const edgeList & edges() const
Definition: treeDataEdge.H:173
const pointField & points() const
Definition: treeDataEdge.H:178
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
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
messageStream Info
Information stream (stdout output on master, null elsewhere)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
label findLower(const ListType &input, const T &val, const label start, const ComparePredicate &comp)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:349
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
List< edge > edgeList
A List of edges.
Definition: edgeList.H:63
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:356
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
List< treeBoundBox > meshBb(1, treeBoundBox(boundBox(coarseMesh.points(), false)).extend(rndGen, 1e-3))
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: stdFoam.H:346
scalarField samples(nIntervals, Zero)
Random rndGen
Definition: createFields.H:23