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-2018 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  // Slightly extended bb. Slightly off-centred just so on symmetric
257  // geometry there are less face/edge aligned items.
258  bb = bb.extend(1e-4);
259 
260  edgeTrees_.set
261  (
262  featI,
263  new indexedOctree<treeDataEdge>
264  (
265  treeDataEdge
266  (
267  false, // do not cache bb
268  edges,
269  points,
270  identity(edges.size())
271  ),
272  bb, // overall search domain
273  8, // maxLevel
274  10, // leafsize
275  3.0 // duplicity
276  )
277  );
278 
279 
280  labelList featurePoints(identity(eMesh.nonFeatureStart()));
281 
282  pointTrees_.set
283  (
284  featI,
285  new indexedOctree<treeDataPoint>
286  (
287  treeDataPoint(points, featurePoints),
288  bb, // overall search domain
289  8, // maxLevel
290  10, // leafsize
291  3.0 // duplicity
292  )
293  );
294 }
295 
296 
297 // Find maximum level of a shell.
298 void Foam::refinementFeatures::findHigherLevel
299 (
300  const pointField& pt,
301  const label featI,
302  labelList& maxLevel
303 ) const
304 {
305  const labelList& levels = levels_[featI];
306 
307  const scalarField& distances = distances_[featI];
308 
309  // Collect all those points that have a current maxLevel less than
310  // (any of) the shell. Also collect the furthest distance allowable
311  // to any shell with a higher level.
312 
313  pointField candidates(pt.size());
314  labelList candidateMap(pt.size());
315  scalarField candidateDistSqr(pt.size());
316  label candidateI = 0;
317 
318  forAll(maxLevel, pointi)
319  {
320  forAllReverse(levels, levelI)
321  {
322  if (levels[levelI] > maxLevel[pointi])
323  {
324  candidates[candidateI] = pt[pointi];
325  candidateMap[candidateI] = pointi;
326  candidateDistSqr[candidateI] = sqr(distances[levelI]);
327  candidateI++;
328  break;
329  }
330  }
331  }
332  candidates.setSize(candidateI);
333  candidateMap.setSize(candidateI);
334  candidateDistSqr.setSize(candidateI);
335 
336  // Do the expensive nearest test only for the candidate points.
337  const indexedOctree<treeDataEdge>& tree = edgeTrees_[featI];
338 
339  List<pointIndexHit> nearInfo(candidates.size());
340  forAll(candidates, candidateI)
341  {
342  nearInfo[candidateI] = tree.findNearest
343  (
344  candidates[candidateI],
345  candidateDistSqr[candidateI]
346  );
347  }
348 
349  // Update maxLevel
350  forAll(nearInfo, candidateI)
351  {
352  if (nearInfo[candidateI].hit())
353  {
354  // Check which level it actually is in.
355  label minDistI = findLower
356  (
357  distances,
358  mag(nearInfo[candidateI].hitPoint()-candidates[candidateI])
359  );
360 
361  label pointi = candidateMap[candidateI];
362 
363  // pt is in between shell[minDistI] and shell[minDistI+1]
364  maxLevel[pointi] = levels[minDistI+1];
365  }
366  }
367 }
368 
369 
370 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
371 
374 {
375  if (!regionEdgeTreesPtr_.valid())
376  {
377  regionEdgeTreesPtr_.reset
378  (
380  );
381  PtrList<indexedOctree<treeDataEdge>>& trees = regionEdgeTreesPtr_();
382 
383  forAll(*this, featI)
384  {
385  const extendedEdgeMesh& eMesh = operator[](featI);
386  const pointField& points = eMesh.points();
387  const edgeList& edges = eMesh.edges();
388 
389  // Calculate bb of all points
390  treeBoundBox bb(points);
391 
392  // Slightly extended bb. Slightly off-centred just so on symmetric
393  // geometry there are less face/edge aligned items.
394  bb = bb.extend(1e-4);
395 
396  trees.set
397  (
398  featI,
400  (
402  (
403  false, // do not cache bb
404  edges,
405  points,
406  eMesh.regionEdges()
407  ),
408  bb, // overall search domain
409  8, // maxLevel
410  10, // leafsize
411  3.0 // duplicity
412  )
413  );
414  }
415  }
416  return regionEdgeTreesPtr_();
417 }
418 
419 
420 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
421 
423 (
424  const objectRegistry& io,
425  const PtrList<dictionary>& featDicts
426 )
427 :
429  distances_(featDicts.size()),
430  levels_(featDicts.size()),
431  edgeTrees_(featDicts.size()),
432  pointTrees_(featDicts.size())
433 {
434  // Read features
435  read(io, featDicts);
436 
437  // Search engines
438  forAll(*this, i)
439  {
440  buildTrees(i);
441  }
442 }
443 
444 
445 //Foam::refinementFeatures::refinementFeatures
446 //(
447 // const objectRegistry& io,
448 // const PtrList<dictionary>& featDicts,
449 // const scalar minCos
450 //)
451 //:
452 // PtrList<extendedFeatureEdgeMesh>(featDicts.size()),
453 // distances_(featDicts.size()),
454 // levels_(featDicts.size()),
455 // edgeTrees_(featDicts.size()),
456 // pointTrees_(featDicts.size())
457 //{
458 // // Read features
459 // read(io, featDicts);
460 //
461 // // Search engines
462 // forAll(*this, i)
463 // {
464 // const edgeMesh& eMesh = operator[](i);
465 // const pointField& points = eMesh.points();
466 // const edgeList& edges = eMesh.edges();
467 // const labelListList& pointEdges = eMesh.pointEdges();
468 //
469 // DynamicList<label> featurePoints;
470 // forAll(pointEdges, pointi)
471 // {
472 // const labelList& pEdges = pointEdges[pointi];
473 // if (pEdges.size() > 2)
474 // {
475 // featurePoints.append(pointi);
476 // }
477 // else if (pEdges.size() == 2)
478 // {
479 // // Check the angle
480 // const edge& e0 = edges[pEdges[0]];
481 // const edge& e1 = edges[pEdges[1]];
482 //
483 // const point& p = points[pointi];
484 // const point& p0 = points[e0.otherVertex(pointi)];
485 // const point& p1 = points[e1.otherVertex(pointi)];
486 //
487 // vector v0 = p-p0;
488 // scalar v0Mag = mag(v0);
489 //
490 // vector v1 = p1-p;
491 // scalar v1Mag = mag(v1);
492 //
493 // if
494 // (
495 // v0Mag > small
496 // && v1Mag > small
497 // && ((v0/v0Mag & v1/v1Mag) < minCos)
498 // )
499 // {
500 // featurePoints.append(pointi);
501 // }
502 // }
503 // }
504 //
505 // Info<< "Detected " << featurePoints.size()
506 // << " featurePoints out of " << points.size()
507 // << " points on feature " << i // eMesh.name()
508 // << " when using feature cos " << minCos << endl;
509 //
510 // buildTrees(i, featurePoints);
511 // }
512 //}
513 
514 
515 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
516 
518 (
519  const pointField& samples,
520  const scalarField& nearestDistSqr,
521  labelList& nearFeature,
522  List<pointIndexHit>& nearInfo,
523  vectorField& nearNormal
524 ) const
525 {
526  nearFeature.setSize(samples.size());
527  nearFeature = -1;
528  nearInfo.setSize(samples.size());
529  nearInfo = pointIndexHit();
530  nearNormal.setSize(samples.size());
531  nearNormal = Zero;
532 
533  forAll(edgeTrees_, featI)
534  {
535  const indexedOctree<treeDataEdge>& tree = edgeTrees_[featI];
536 
537  if (tree.shapes().size() > 0)
538  {
539  forAll(samples, sampleI)
540  {
541  const point& sample = samples[sampleI];
542 
543  scalar distSqr;
544  if (nearInfo[sampleI].hit())
545  {
546  distSqr = magSqr(nearInfo[sampleI].hitPoint()-sample);
547  }
548  else
549  {
550  distSqr = nearestDistSqr[sampleI];
551  }
552 
553  pointIndexHit info = tree.findNearest(sample, distSqr);
554 
555  if (info.hit())
556  {
557  nearFeature[sampleI] = featI;
558  nearInfo[sampleI] = pointIndexHit
559  (
560  info.hit(),
561  info.hitPoint(),
562  tree.shapes().edgeLabels()[info.index()]
563  );
564 
565  const treeDataEdge& td = tree.shapes();
566  const edge& e = td.edges()[nearInfo[sampleI].index()];
567  nearNormal[sampleI] = e.vec(td.points());
568  nearNormal[sampleI] /= mag(nearNormal[sampleI])+vSmall;
569  }
570  }
571  }
572  }
573 }
574 
575 
577 (
578  const pointField& samples,
579  const scalarField& nearestDistSqr,
580  labelList& nearFeature,
581  List<pointIndexHit>& nearInfo,
582  vectorField& nearNormal
583 ) const
584 {
585  nearFeature.setSize(samples.size());
586  nearFeature = -1;
587  nearInfo.setSize(samples.size());
588  nearInfo = pointIndexHit();
589  nearNormal.setSize(samples.size());
590  nearNormal = Zero;
591 
592 
593  const PtrList<indexedOctree<treeDataEdge>>& regionTrees =
594  regionEdgeTrees();
595 
596  forAll(regionTrees, featI)
597  {
598  const indexedOctree<treeDataEdge>& regionTree = regionTrees[featI];
599 
600  forAll(samples, sampleI)
601  {
602  const point& sample = samples[sampleI];
603 
604  scalar distSqr;
605  if (nearInfo[sampleI].hit())
606  {
607  distSqr = magSqr(nearInfo[sampleI].hitPoint()-sample);
608  }
609  else
610  {
611  distSqr = nearestDistSqr[sampleI];
612  }
613 
614  // Find anything closer than current best
615  pointIndexHit info = regionTree.findNearest(sample, distSqr);
616 
617  if (info.hit())
618  {
619  const treeDataEdge& td = regionTree.shapes();
620 
621  nearFeature[sampleI] = featI;
622  nearInfo[sampleI] = pointIndexHit
623  (
624  info.hit(),
625  info.hitPoint(),
626  regionTree.shapes().edgeLabels()[info.index()]
627  );
628 
629  const edge& e = td.edges()[nearInfo[sampleI].index()];
630  nearNormal[sampleI] = e.vec(td.points());
631  nearNormal[sampleI] /= mag(nearNormal[sampleI])+vSmall;
632  }
633  }
634  }
635 }
636 
637 
638 //void Foam::refinementFeatures::findNearestPoint
639 //(
640 // const pointField& samples,
641 // const scalarField& nearestDistSqr,
642 // labelList& nearFeature,
643 // labelList& nearIndex
644 //) const
645 //{
646 // nearFeature.setSize(samples.size());
647 // nearFeature = -1;
648 // nearIndex.setSize(samples.size());
649 // nearIndex = -1;
650 //
651 // forAll(pointTrees_, featI)
652 // {
653 // const indexedOctree<treeDataPoint>& tree = pointTrees_[featI];
654 //
655 // if (tree.shapes().pointLabels().size() > 0)
656 // {
657 // forAll(samples, sampleI)
658 // {
659 // const point& sample = samples[sampleI];
660 //
661 // scalar distSqr;
662 // if (nearFeature[sampleI] != -1)
663 // {
664 // label nearFeatI = nearFeature[sampleI];
665 // const indexedOctree<treeDataPoint>& nearTree =
666 // pointTrees_[nearFeatI];
667 // label featPointi =
668 // nearTree.shapes().pointLabels()[nearIndex[sampleI]];
669 // const point& featPt =
670 // operator[](nearFeatI).points()[featPointi];
671 // distSqr = magSqr(featPt-sample);
672 // }
673 // else
674 // {
675 // distSqr = nearestDistSqr[sampleI];
676 // }
677 //
678 // pointIndexHit info = tree.findNearest(sample, distSqr);
679 //
680 // if (info.hit())
681 // {
682 // nearFeature[sampleI] = featI;
683 // nearIndex[sampleI] = info.index();
684 // }
685 // }
686 // }
687 // }
688 //}
689 
690 
692 (
693  const pointField& samples,
694  const scalarField& nearestDistSqr,
695  labelList& nearFeature,
696  List<pointIndexHit>& nearInfo
697 ) const
698 {
699  nearFeature.setSize(samples.size());
700  nearFeature = -1;
701  nearInfo.setSize(samples.size());
702  nearInfo = pointIndexHit();
703 
704  forAll(pointTrees_, featI)
705  {
706  const indexedOctree<treeDataPoint>& tree = pointTrees_[featI];
707 
708  if (tree.shapes().pointLabels().size() > 0)
709  {
710  forAll(samples, sampleI)
711  {
712  const point& sample = samples[sampleI];
713 
714  scalar distSqr;
715  if (nearFeature[sampleI] != -1)
716  {
717  distSqr = magSqr(nearInfo[sampleI].hitPoint()-sample);
718  }
719  else
720  {
721  distSqr = nearestDistSqr[sampleI];
722  }
723 
724  pointIndexHit info = tree.findNearest(sample, distSqr);
725 
726  if (info.hit())
727  {
728  nearFeature[sampleI] = featI;
729  nearInfo[sampleI] = pointIndexHit
730  (
731  info.hit(),
732  info.hitPoint(),
733  tree.shapes().pointLabels()[info.index()]
734  );
735  }
736  }
737  }
738  }
739 }
740 
741 
742 void Foam::refinementFeatures::findHigherLevel
743 (
744  const pointField& pt,
745  const labelList& ptLevel,
746  labelList& maxLevel
747 ) const
748 {
749  // Maximum level of any shell. Start off with level of point.
750  maxLevel = ptLevel;
751 
752  forAll(*this, featI)
753  {
754  findHigherLevel(pt, featI, maxLevel);
755  }
756 }
757 
758 
760 {
761  scalar overallMax = -great;
762  forAll(distances_, featI)
763  {
764  overallMax = max(overallMax, max(distances_[featI]));
765  }
766  return overallMax;
767 }
768 
769 
770 // ************************************************************************* //
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.
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:256
#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
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
label readLabel(Istream &is)
Definition: label.H:64
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:265
const pointField & points() const
Return points.
Definition: edgeMeshI.H:53
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:240
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: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 doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
Registry of regIOobjects.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:233
label index() const
Return index.
Tuple2< scalar, label > nearInfo
Private class for finding nearest.
treeBoundBox extend(const scalar s) const
Return asymetrically extended bounding box, with guaranteed.
IOerror FatalIOError