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