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 | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "refinementFeatures.H"
27 #include "Time.H"
28 #include "Tuple2.H"
29 #include "DynamicField.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 void Foam::refinementFeatures::read
34 (
35  const objectRegistry& io,
36  const PtrList<dictionary>& featDicts
37 )
38 {
39  forAll(featDicts, featI)
40  {
41  const dictionary& dict = featDicts[featI];
42 
43  fileName featFileName(dict.lookup("file"));
44 
45 
46  // Try reading extendedEdgeMesh first
47 
48  IOobject extFeatObj
49  (
50  featFileName, // name
51  io.time().constant(), // instance
52  "extendedFeatureEdgeMesh", // local
53  io.time(), // registry
56  false
57  );
58 
59  const fileName fName(extFeatObj.filePath());
60 
61  if (!fName.empty() && extendedEdgeMesh::canRead(fName))
62  {
63  autoPtr<extendedEdgeMesh> eMeshPtr = extendedEdgeMesh::New
64  (
65  fName
66  );
67 
68  Info<< "Read extendedFeatureEdgeMesh " << extFeatObj.name()
69  << nl << incrIndent;
70  eMeshPtr().writeStats(Info);
71  Info<< decrIndent << endl;
72 
73  set(featI, new extendedFeatureEdgeMesh(extFeatObj, eMeshPtr()));
74  }
75  else
76  {
77  // Try reading edgeMesh
78 
79  IOobject featObj
80  (
81  featFileName, // name
82  io.time().constant(), // instance
83  "triSurface", // local
84  io.time(), // registry
87  false
88  );
89 
90  const fileName fName(featObj.filePath());
91 
92  if (fName.empty())
93  {
95  (
96  dict
97  ) << "Could not open " << featObj.objectPath()
98  << exit(FatalIOError);
99  }
100 
101 
102  // Read as edgeMesh
103  autoPtr<edgeMesh> eMeshPtr = edgeMesh::New(fName);
104  const edgeMesh& eMesh = eMeshPtr();
105 
106  Info<< "Read edgeMesh " << featObj.name() << nl
107  << incrIndent;
108  eMesh.writeStats(Info);
109  Info<< decrIndent << endl;
110 
111 
112  // Analyse for feature points. These are all classified as mixed
113  // points for lack of anything better
114  const labelListList& pointEdges = eMesh.pointEdges();
115 
116  labelList oldToNew(eMesh.points().size(), -1);
117  DynamicField<point> newPoints(eMesh.points().size());
118  forAll(pointEdges, pointi)
119  {
120  if (pointEdges[pointi].size() > 2)
121  {
122  oldToNew[pointi] = newPoints.size();
123  newPoints.append(eMesh.points()[pointi]);
124  }
125  //else if (pointEdges[pointi].size() == 2)
126  //MEJ: do something based on a feature angle?
127  }
128  label nFeatures = newPoints.size();
129  forAll(oldToNew, pointi)
130  {
131  if (oldToNew[pointi] == -1)
132  {
133  oldToNew[pointi] = newPoints.size();
134  newPoints.append(eMesh.points()[pointi]);
135  }
136  }
137 
138 
139  const edgeList& edges = eMesh.edges();
140  edgeList newEdges(edges.size());
141  forAll(edges, edgeI)
142  {
143  const edge& e = edges[edgeI];
144  newEdges[edgeI] = edge
145  (
146  oldToNew[e[0]],
147  oldToNew[e[1]]
148  );
149  }
150 
151  // Construct an extendedEdgeMesh with
152  // - all points on more than 2 edges : mixed feature points
153  // - all edges : external edges
154 
155  extendedEdgeMesh eeMesh
156  (
157  newPoints, // pts
158  newEdges, // eds
159  0, // (point) concaveStart
160  0, // (point) mixedStart
161  nFeatures, // (point) nonFeatureStart
162  edges.size(), // (edge) internalStart
163  edges.size(), // (edge) flatStart
164  edges.size(), // (edge) openStart
165  edges.size(), // (edge) multipleStart
166  vectorField(0), // normals
167  List<extendedEdgeMesh::sideVolumeType>(0),// normalVolumeTypes
168  vectorField(0), // edgeDirections
169  labelListList(0), // normalDirections
170  labelListList(0), // edgeNormals
171  labelListList(0), // featurePointNormals
172  labelListList(0), // featurePointEdges
173  identity(newEdges.size()) // regionEdges
174  );
175 
176  //Info<< "Constructed extendedFeatureEdgeMesh " << featObj.name()
177  // << nl << incrIndent;
178  //eeMesh.writeStats(Info);
179  //Info<< decrIndent << endl;
180 
181  set(featI, new extendedFeatureEdgeMesh(featObj, eeMesh));
182  }
183 
184  const extendedEdgeMesh& eMesh = operator[](featI);
185 
186  if (dict.found("levels"))
187  {
188  List<Tuple2<scalar, label>> distLevels(dict["levels"]);
189 
190  if (dict.size() < 1)
191  {
193  << " : levels should be at least size 1" << endl
194  << "levels : " << dict["levels"]
195  << exit(FatalError);
196  }
197 
198  distances_[featI].setSize(distLevels.size());
199  levels_[featI].setSize(distLevels.size());
200 
201  forAll(distLevels, j)
202  {
203  distances_[featI][j] = distLevels[j].first();
204  levels_[featI][j] = distLevels[j].second();
205 
206  // Check in incremental order
207  if (j > 0)
208  {
209  if
210  (
211  (distances_[featI][j] <= distances_[featI][j-1])
212  || (levels_[featI][j] > levels_[featI][j-1])
213  )
214  {
216  << " : Refinement should be specified in order"
217  << " of increasing distance"
218  << " (and decreasing refinement level)." << endl
219  << "Distance:" << distances_[featI][j]
220  << " refinementLevel:" << levels_[featI][j]
221  << exit(FatalError);
222  }
223  }
224  }
225  }
226  else
227  {
228  // Look up 'level' for single level
229  levels_[featI] = labelList(1, readLabel(dict.lookup("level")));
230  distances_[featI] = scalarField(1, 0.0);
231  }
232 
233  Info<< "Refinement level according to distance to "
234  << featFileName << " (" << eMesh.points().size() << " points, "
235  << eMesh.edges().size() << " edges)." << endl;
236  forAll(levels_[featI], j)
237  {
238  Info<< " level " << levels_[featI][j]
239  << " for all cells within " << distances_[featI][j]
240  << " metre." << endl;
241  }
242  }
243 }
244 
245 
246 void Foam::refinementFeatures::buildTrees(const label featI)
247 {
248  const extendedEdgeMesh& eMesh = operator[](featI);
249  const pointField& points = eMesh.points();
250  const edgeList& edges = eMesh.edges();
251 
252  // Calculate bb of all points
253  treeBoundBox bb(points);
254 
255  // Random number generator. Bit dodgy since not exactly random ;-)
256  Random rndGen(65431);
257 
258  // Slightly extended bb. Slightly off-centred just so on symmetric
259  // geometry there are less face/edge aligned items.
260  bb = bb.extend(rndGen, 1e-4);
261  bb.min() -= point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
262  bb.max() += point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
263 
264  edgeTrees_.set
265  (
266  featI,
267  new indexedOctree<treeDataEdge>
268  (
269  treeDataEdge
270  (
271  false, // do not cache bb
272  edges,
273  points,
274  identity(edges.size())
275  ),
276  bb, // overall search domain
277  8, // maxLevel
278  10, // leafsize
279  3.0 // duplicity
280  )
281  );
282 
283 
284  labelList featurePoints(identity(eMesh.nonFeatureStart()));
285 
286  pointTrees_.set
287  (
288  featI,
289  new indexedOctree<treeDataPoint>
290  (
291  treeDataPoint(points, featurePoints),
292  bb, // overall search domain
293  8, // maxLevel
294  10, // leafsize
295  3.0 // duplicity
296  )
297  );
298 }
299 
300 
303 {
304  if (!regionEdgeTreesPtr_.valid())
305  {
306  regionEdgeTreesPtr_.reset
307  (
309  );
310  PtrList<indexedOctree<treeDataEdge>>& trees = regionEdgeTreesPtr_();
311 
312  forAll(*this, featI)
313  {
314  const extendedEdgeMesh& eMesh = operator[](featI);
315  const pointField& points = eMesh.points();
316  const edgeList& edges = eMesh.edges();
317 
318  // Calculate bb of all points
319  treeBoundBox bb(points);
320 
321  // Random number generator. Bit dodgy since not exactly random ;-)
322  Random rndGen(65431);
323 
324  // Slightly extended bb. Slightly off-centred just so on symmetric
325  // geometry there are less face/edge aligned items.
326  bb = bb.extend(rndGen, 1e-4);
327  bb.min() -= point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
328  bb.max() += point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
329 
330  trees.set
331  (
332  featI,
334  (
336  (
337  false, // do not cache bb
338  edges,
339  points,
340  eMesh.regionEdges()
341  ),
342  bb, // overall search domain
343  8, // maxLevel
344  10, // leafsize
345  3.0 // duplicity
346  )
347  );
348  }
349  }
350  return regionEdgeTreesPtr_();
351 }
352 
353 
354 // Find maximum level of a shell.
355 void Foam::refinementFeatures::findHigherLevel
356 (
357  const pointField& pt,
358  const label featI,
359  labelList& maxLevel
360 ) const
361 {
362  const labelList& levels = levels_[featI];
363 
364  const scalarField& distances = distances_[featI];
365 
366  // Collect all those points that have a current maxLevel less than
367  // (any of) the shell. Also collect the furthest distance allowable
368  // to any shell with a higher level.
369 
370  pointField candidates(pt.size());
371  labelList candidateMap(pt.size());
372  scalarField candidateDistSqr(pt.size());
373  label candidateI = 0;
374 
375  forAll(maxLevel, pointi)
376  {
377  forAllReverse(levels, levelI)
378  {
379  if (levels[levelI] > maxLevel[pointi])
380  {
381  candidates[candidateI] = pt[pointi];
382  candidateMap[candidateI] = pointi;
383  candidateDistSqr[candidateI] = sqr(distances[levelI]);
384  candidateI++;
385  break;
386  }
387  }
388  }
389  candidates.setSize(candidateI);
390  candidateMap.setSize(candidateI);
391  candidateDistSqr.setSize(candidateI);
392 
393  // Do the expensive nearest test only for the candidate points.
394  const indexedOctree<treeDataEdge>& tree = edgeTrees_[featI];
395 
396  List<pointIndexHit> nearInfo(candidates.size());
397  forAll(candidates, candidateI)
398  {
399  nearInfo[candidateI] = tree.findNearest
400  (
401  candidates[candidateI],
402  candidateDistSqr[candidateI]
403  );
404  }
405 
406  // Update maxLevel
407  forAll(nearInfo, candidateI)
408  {
409  if (nearInfo[candidateI].hit())
410  {
411  // Check which level it actually is in.
412  label minDistI = findLower
413  (
414  distances,
415  mag(nearInfo[candidateI].hitPoint()-candidates[candidateI])
416  );
417 
418  label pointi = candidateMap[candidateI];
419 
420  // pt is inbetween shell[minDistI] and shell[minDistI+1]
421  maxLevel[pointi] = levels[minDistI+1];
422  }
423  }
424 }
425 
426 
427 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
428 
430 (
431  const objectRegistry& io,
432  const PtrList<dictionary>& featDicts
433 )
434 :
436  distances_(featDicts.size()),
437  levels_(featDicts.size()),
438  edgeTrees_(featDicts.size()),
439  pointTrees_(featDicts.size())
440 {
441  // Read features
442  read(io, featDicts);
443 
444  // Search engines
445  forAll(*this, i)
446  {
447  buildTrees(i);
448  }
449 }
450 
451 
452 //Foam::refinementFeatures::refinementFeatures
453 //(
454 // const objectRegistry& io,
455 // const PtrList<dictionary>& featDicts,
456 // const scalar minCos
457 //)
458 //:
459 // PtrList<extendedFeatureEdgeMesh>(featDicts.size()),
460 // distances_(featDicts.size()),
461 // levels_(featDicts.size()),
462 // edgeTrees_(featDicts.size()),
463 // pointTrees_(featDicts.size())
464 //{
465 // // Read features
466 // read(io, featDicts);
467 //
468 // // Search engines
469 // forAll(*this, i)
470 // {
471 // const edgeMesh& eMesh = operator[](i);
472 // const pointField& points = eMesh.points();
473 // const edgeList& edges = eMesh.edges();
474 // const labelListList& pointEdges = eMesh.pointEdges();
475 //
476 // DynamicList<label> featurePoints;
477 // forAll(pointEdges, pointi)
478 // {
479 // const labelList& pEdges = pointEdges[pointi];
480 // if (pEdges.size() > 2)
481 // {
482 // featurePoints.append(pointi);
483 // }
484 // else if (pEdges.size() == 2)
485 // {
486 // // Check the angle
487 // const edge& e0 = edges[pEdges[0]];
488 // const edge& e1 = edges[pEdges[1]];
489 //
490 // const point& p = points[pointi];
491 // const point& p0 = points[e0.otherVertex(pointi)];
492 // const point& p1 = points[e1.otherVertex(pointi)];
493 //
494 // vector v0 = p-p0;
495 // scalar v0Mag = mag(v0);
496 //
497 // vector v1 = p1-p;
498 // scalar v1Mag = mag(v1);
499 //
500 // if
501 // (
502 // v0Mag > SMALL
503 // && v1Mag > SMALL
504 // && ((v0/v0Mag & v1/v1Mag) < minCos)
505 // )
506 // {
507 // featurePoints.append(pointi);
508 // }
509 // }
510 // }
511 //
512 // Info<< "Detected " << featurePoints.size()
513 // << " featurePoints out of " << points.size()
514 // << " points on feature " << i //eMesh.name()
515 // << " when using feature cos " << minCos << endl;
516 //
517 // buildTrees(i, featurePoints);
518 // }
519 //}
520 
521 
522 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
523 
525 (
526  const pointField& samples,
527  const scalarField& nearestDistSqr,
528  labelList& nearFeature,
529  List<pointIndexHit>& nearInfo,
530  vectorField& nearNormal
531 ) const
532 {
533  nearFeature.setSize(samples.size());
534  nearFeature = -1;
535  nearInfo.setSize(samples.size());
536  nearInfo = pointIndexHit();
537  nearNormal.setSize(samples.size());
538  nearNormal = Zero;
539 
540  forAll(edgeTrees_, featI)
541  {
542  const indexedOctree<treeDataEdge>& tree = edgeTrees_[featI];
543 
544  if (tree.shapes().size() > 0)
545  {
546  forAll(samples, sampleI)
547  {
548  const point& sample = samples[sampleI];
549 
550  scalar distSqr;
551  if (nearInfo[sampleI].hit())
552  {
553  distSqr = magSqr(nearInfo[sampleI].hitPoint()-sample);
554  }
555  else
556  {
557  distSqr = nearestDistSqr[sampleI];
558  }
559 
560  pointIndexHit info = tree.findNearest(sample, distSqr);
561 
562  if (info.hit())
563  {
564  nearFeature[sampleI] = featI;
565  nearInfo[sampleI] = pointIndexHit
566  (
567  info.hit(),
568  info.hitPoint(),
569  tree.shapes().edgeLabels()[info.index()]
570  );
571 
572  const treeDataEdge& td = tree.shapes();
573  const edge& e = td.edges()[nearInfo[sampleI].index()];
574  nearNormal[sampleI] = e.vec(td.points());
575  nearNormal[sampleI] /= mag(nearNormal[sampleI])+VSMALL;
576  }
577  }
578  }
579  }
580 }
581 
582 
584 (
585  const pointField& samples,
586  const scalarField& nearestDistSqr,
587  labelList& nearFeature,
588  List<pointIndexHit>& nearInfo,
589  vectorField& nearNormal
590 ) const
591 {
592  nearFeature.setSize(samples.size());
593  nearFeature = -1;
594  nearInfo.setSize(samples.size());
595  nearInfo = pointIndexHit();
596  nearNormal.setSize(samples.size());
597  nearNormal = Zero;
598 
599 
600  const PtrList<indexedOctree<treeDataEdge>>& regionTrees =
601  regionEdgeTrees();
602 
603  forAll(regionTrees, featI)
604  {
605  const indexedOctree<treeDataEdge>& regionTree = regionTrees[featI];
606 
607  forAll(samples, sampleI)
608  {
609  const point& sample = samples[sampleI];
610 
611  scalar distSqr;
612  if (nearInfo[sampleI].hit())
613  {
614  distSqr = magSqr(nearInfo[sampleI].hitPoint()-sample);
615  }
616  else
617  {
618  distSqr = nearestDistSqr[sampleI];
619  }
620 
621  // Find anything closer than current best
622  pointIndexHit info = regionTree.findNearest(sample, distSqr);
623 
624  if (info.hit())
625  {
626  const treeDataEdge& td = regionTree.shapes();
627 
628  nearFeature[sampleI] = featI;
629  nearInfo[sampleI] = pointIndexHit
630  (
631  info.hit(),
632  info.hitPoint(),
633  regionTree.shapes().edgeLabels()[info.index()]
634  );
635 
636  const edge& e = td.edges()[nearInfo[sampleI].index()];
637  nearNormal[sampleI] = e.vec(td.points());
638  nearNormal[sampleI] /= mag(nearNormal[sampleI])+VSMALL;
639  }
640  }
641  }
642 }
643 
644 
645 //void Foam::refinementFeatures::findNearestPoint
646 //(
647 // const pointField& samples,
648 // const scalarField& nearestDistSqr,
649 // labelList& nearFeature,
650 // labelList& nearIndex
651 //) const
652 //{
653 // nearFeature.setSize(samples.size());
654 // nearFeature = -1;
655 // nearIndex.setSize(samples.size());
656 // nearIndex = -1;
657 //
658 // forAll(pointTrees_, featI)
659 // {
660 // const indexedOctree<treeDataPoint>& tree = pointTrees_[featI];
661 //
662 // if (tree.shapes().pointLabels().size() > 0)
663 // {
664 // forAll(samples, sampleI)
665 // {
666 // const point& sample = samples[sampleI];
667 //
668 // scalar distSqr;
669 // if (nearFeature[sampleI] != -1)
670 // {
671 // label nearFeatI = nearFeature[sampleI];
672 // const indexedOctree<treeDataPoint>& nearTree =
673 // pointTrees_[nearFeatI];
674 // label featPointi =
675 // nearTree.shapes().pointLabels()[nearIndex[sampleI]];
676 // const point& featPt =
677 // operator[](nearFeatI).points()[featPointi];
678 // distSqr = magSqr(featPt-sample);
679 // }
680 // else
681 // {
682 // distSqr = nearestDistSqr[sampleI];
683 // }
684 //
685 // pointIndexHit info = tree.findNearest(sample, distSqr);
686 //
687 // if (info.hit())
688 // {
689 // nearFeature[sampleI] = featI;
690 // nearIndex[sampleI] = info.index();
691 // }
692 // }
693 // }
694 // }
695 //}
696 
697 
699 (
700  const pointField& samples,
701  const scalarField& nearestDistSqr,
702  labelList& nearFeature,
703  List<pointIndexHit>& nearInfo
704 ) const
705 {
706  nearFeature.setSize(samples.size());
707  nearFeature = -1;
708  nearInfo.setSize(samples.size());
709  nearInfo = pointIndexHit();
710 
711  forAll(pointTrees_, featI)
712  {
713  const indexedOctree<treeDataPoint>& tree = pointTrees_[featI];
714 
715  if (tree.shapes().pointLabels().size() > 0)
716  {
717  forAll(samples, sampleI)
718  {
719  const point& sample = samples[sampleI];
720 
721  scalar distSqr;
722  if (nearFeature[sampleI] != -1)
723  {
724  distSqr = magSqr(nearInfo[sampleI].hitPoint()-sample);
725  }
726  else
727  {
728  distSqr = nearestDistSqr[sampleI];
729  }
730 
731  pointIndexHit info = tree.findNearest(sample, distSqr);
732 
733  if (info.hit())
734  {
735  nearFeature[sampleI] = featI;
736  nearInfo[sampleI] = pointIndexHit
737  (
738  info.hit(),
739  info.hitPoint(),
740  tree.shapes().pointLabels()[info.index()]
741  );
742  }
743  }
744  }
745  }
746 }
747 
748 
749 void Foam::refinementFeatures::findHigherLevel
750 (
751  const pointField& pt,
752  const labelList& ptLevel,
753  labelList& maxLevel
754 ) const
755 {
756  // Maximum level of any shell. Start off with level of point.
757  maxLevel = ptLevel;
758 
759  forAll(*this, featI)
760  {
761  findHigherLevel(pt, featI, maxLevel);
762  }
763 }
764 
765 
767 {
768  scalar overallMax = -GREAT;
769  forAll(distances_, featI)
770  {
771  overallMax = max(overallMax, max(distances_[featI]));
772  }
773  return overallMax;
774 }
775 
776 
777 // ************************************************************************* //
cachedRandom rndGen(label(0),-1)
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
static autoPtr< extendedEdgeMesh > New(const fileName &, const word &ext)
Select constructed from filename (explicit extension)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
refinementFeatures(const objectRegistry &io, const PtrList< dictionary > &featDicts)
Construct from description.
const double e
Elementary charge.
Definition: doubleFloat.H:78
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
static autoPtr< edgeMesh > New(const fileName &, const word &ext)
Select constructed from filename (explicit extension)
Definition: edgeMeshNew.C:31
bool hit() const
Is there a hit.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
Description of feature edges and points.
treeBoundBox extend(Random &, const scalar s) const
Return slightly wider bounding box.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:440
Holds data for octree to work on an edges subset.
Definition: treeDataEdge.H:53
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
T & first()
Return the first element of the list.
Definition: UListI.H:114
const edgeList & edges() const
Definition: treeDataEdge.H:169
const pointField & points() const
Return points.
Definition: edgeMeshI.H:39
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
List< edge > edgeList
Definition: edgeList.H:38
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
scalar maxDistance() const
Highest distance of all features.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void findNearestPoint(const pointField &samples, const scalarField &nearestDistSqr, labelList &nearFeature, List< pointIndexHit > &nearInfo) const
Find nearest feature point. Sets:
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:91
label readLabel(Istream &is)
Definition: label.H:64
Simple random number generator.
Definition: Random.H:49
const labelList & regionEdges() const
Return the feature edges which are on the boundary between.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
static const char nl
Definition: Ostream.H:262
const List< scalarField > & distances() const
Per featureEdgeMesh the list of ranges.
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:237
const T & operator[](const label) const
Return element const reference.
Definition: UPtrListI.H:103
void setSize(const label)
Reset size of List.
Definition: List.C:295
vector point
Point is a vector.
Definition: point.H:41
Non-pointer based hierarchical recursive searching.
Definition: treeDataEdge.H:47
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:62
const Type & shapes() const
Reference to shape.
label index() const
Return index.
static bool canRead(const fileName &, const bool verbose=false)
Can we read this file format?
const edgeList & edges() const
Return edges.
Definition: edgeMeshI.H:45
const Point & hitPoint() const
Return hit point.
vector vec(const pointField &) const
Return the vector (end - start)
Definition: edgeI.H:157
const pointField & points() const
Definition: treeDataEdge.H:174
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:87
messageStream Info
label findLower(const ListType &, typename ListType::const_reference, const label stary, const BinaryOp &bop)
Find last element < given value in sorted list and return index,.
dimensioned< scalar > mag(const dimensioned< Type > &)
const PtrList< indexedOctree< treeDataEdge > > & regionEdgeTrees() const
Field< vector > vectorField
Specialisation of Field<T> for vector.
Registry of regIOobjects.
const labelListList & levels() const
Per featureEdgeMesh the list of level.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:230
void findNearestEdge(const pointField &samples, const scalarField &nearestDistSqr, labelList &nearFeature, List< pointIndexHit > &nearInfo, vectorField &nearNormal) const
Find nearest point on nearest feature edge. Sets:
void findNearestRegionEdge(const pointField &samples, const scalarField &nearestDistSqr, labelList &nearFeature, List< pointIndexHit > &nearInfo, vectorField &nearNormal) const
Find nearest point on nearest region edge. Sets:
Tuple2< scalar, label > nearInfo
Private class for finding nearest.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
IOerror FatalIOError