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