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