sampledTriSurfaceMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "sampledTriSurfaceMesh.H"
27 #include "meshSearch.H"
28 #include "Tuple2.H"
29 #include "globalIndex.H"
30 #include "treeDataCell.H"
31 #include "treeDataFace.H"
32 #include "meshTools.H"
33 
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(sampledTriSurfaceMesh, 0);
42  (
43  sampledSurface,
44  sampledTriSurfaceMesh,
45  word
46  );
47 
48  template<>
50  {
51  "cells",
52  "insideCells",
53  "boundaryFaces"
54  };
55 
57  sampledTriSurfaceMesh::samplingSourceNames_;
58 
59 
60  //- Private class for finding nearest
61  // Comprising:
62  // - global index
63  // - sqr(distance)
65 
67  {
68 
69  public:
70 
71  void operator()(nearInfo& x, const nearInfo& y) const
72  {
73  if (y.first() < x.first())
74  {
75  x = y;
76  }
77  }
78  };
79 }
80 
81 
82 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
83 
85 Foam::sampledTriSurfaceMesh::nonCoupledboundaryTree() const
86 {
87  // Variant of meshSearch::boundaryTree() that only does non-coupled
88  // boundary faces.
89 
90  if (!boundaryTreePtr_.valid())
91  {
92  // all non-coupled boundary faces (not just walls)
93  const polyBoundaryMesh& patches = mesh().boundaryMesh();
94 
95  labelList bndFaces(mesh().nFaces()-mesh().nInternalFaces());
96  label bndI = 0;
97  forAll(patches, patchi)
98  {
99  const polyPatch& pp = patches[patchi];
100  if (!pp.coupled())
101  {
102  forAll(pp, i)
103  {
104  bndFaces[bndI++] = pp.start()+i;
105  }
106  }
107  }
108  bndFaces.setSize(bndI);
109 
110 
111  treeBoundBox overallBb(mesh().points());
112  Random rndGen(123456);
113  overallBb = overallBb.extend(rndGen, 1e-4);
114  overallBb.min() -= point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
115  overallBb.max() += point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
116 
117  boundaryTreePtr_.reset
118  (
120  (
121  treeDataFace // all information needed to search faces
122  (
123  false, // do not cache bb
124  mesh(),
125  bndFaces // boundary faces only
126  ),
127  overallBb, // overall search domain
128  8, // maxLevel
129  10, // leafsize
130  3.0 // duplicity
131  )
132  );
133  }
134 
135  return boundaryTreePtr_();
136 }
137 
138 
139 bool Foam::sampledTriSurfaceMesh::update(const meshSearch& meshSearcher)
140 {
141  // Find the cells the triangles of the surface are in.
142  // Does approximation by looking at the face centres only
143  const pointField& fc = surface_.faceCentres();
144 
145  List<nearInfo> nearest(fc.size());
146 
147  // Global numbering for cells/faces - only used to uniquely identify local
148  // elements
149  globalIndex globalCells
150  (
151  (sampleSource_ == cells || sampleSource_ == insideCells)
152  ? mesh().nCells()
153  : mesh().nFaces()
154  );
155 
156  forAll(nearest, i)
157  {
158  nearest[i].first() = GREAT;
159  nearest[i].second() = labelMax;
160  }
161 
162  if (sampleSource_ == cells)
163  {
164  // Search for nearest cell
165 
166  const indexedOctree<treeDataCell>& cellTree = meshSearcher.cellTree();
167 
168  forAll(fc, triI)
169  {
170  pointIndexHit nearInfo = cellTree.findNearest
171  (
172  fc[triI],
173  sqr(GREAT)
174  );
175  if (nearInfo.hit())
176  {
177  nearest[triI].first() = magSqr(nearInfo.hitPoint()-fc[triI]);
178  nearest[triI].second() = globalCells.toGlobal(nearInfo.index());
179  }
180  }
181  }
182  else if (sampleSource_ == insideCells)
183  {
184  // Search for cell containing point
185 
186  const indexedOctree<treeDataCell>& cellTree = meshSearcher.cellTree();
187 
188  forAll(fc, triI)
189  {
190  if (cellTree.bb().contains(fc[triI]))
191  {
192  label index = cellTree.findInside(fc[triI]);
193  if (index != -1)
194  {
195  nearest[triI].first() = 0.0;
196  nearest[triI].second() = globalCells.toGlobal(index);
197  }
198  }
199  }
200  }
201  else
202  {
203  // Search for nearest boundaryFace
204 
206  //const indexedOctree<treeDataFace>& bTree = meshSearcher.boundaryTree()
207  //- Search on all non-coupled boundary faces
208  const indexedOctree<treeDataFace>& bTree = nonCoupledboundaryTree();
209 
210  forAll(fc, triI)
211  {
212  pointIndexHit nearInfo = bTree.findNearest
213  (
214  fc[triI],
215  sqr(GREAT)
216  );
217  if (nearInfo.hit())
218  {
219  nearest[triI].first() = magSqr(nearInfo.hitPoint()-fc[triI]);
220  nearest[triI].second() = globalCells.toGlobal
221  (
222  bTree.shapes().faceLabels()[nearInfo.index()]
223  );
224  }
225  }
226  }
227 
228 
229  // See which processor has the nearest. Mark and subset
230  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
231 
234 
235  labelList cellOrFaceLabels(fc.size(), -1);
236 
237  label nFound = 0;
238  forAll(nearest, triI)
239  {
240  if (nearest[triI].second() == labelMax)
241  {
242  // Not found on any processor. How to map?
243  }
244  else if (globalCells.isLocal(nearest[triI].second()))
245  {
246  cellOrFaceLabels[triI] = globalCells.toLocal
247  (
248  nearest[triI].second()
249  );
250  nFound++;
251  }
252  }
253 
254 
255  if (debug)
256  {
257  Pout<< "Local out of faces:" << cellOrFaceLabels.size()
258  << " keeping:" << nFound << endl;
259  }
260 
261  // Now subset the surface. Do not use triSurface::subsetMesh since requires
262  // original surface to be in compact numbering.
263 
264  const triSurface& s = surface_;
265 
266  // Compact to original triangle
267  labelList faceMap(s.size());
268  // Compact to original points
269  labelList pointMap(s.points().size());
270  // From original point to compact points
271  labelList reversePointMap(s.points().size(), -1);
272 
273  {
274  label newPointi = 0;
275  label newFacei = 0;
276 
277  forAll(s, facei)
278  {
279  if (cellOrFaceLabels[facei] != -1)
280  {
281  faceMap[newFacei++] = facei;
282 
283  const triSurface::FaceType& f = s[facei];
284  forAll(f, fp)
285  {
286  if (reversePointMap[f[fp]] == -1)
287  {
288  pointMap[newPointi] = f[fp];
289  reversePointMap[f[fp]] = newPointi++;
290  }
291  }
292  }
293  }
294  faceMap.setSize(newFacei);
295  pointMap.setSize(newPointi);
296  }
297 
298 
299  // Subset cellOrFaceLabels
300  cellOrFaceLabels = UIndirectList<label>(cellOrFaceLabels, faceMap)();
301 
302  // Store any face per point (without using pointFaces())
303  labelList pointToFace(pointMap.size());
304 
305  // Create faces and points for subsetted surface
306  faceList& faces = this->storedFaces();
307  faces.setSize(faceMap.size());
308  forAll(faceMap, i)
309  {
310  const triFace& f = s[faceMap[i]];
311  triFace newF
312  (
313  reversePointMap[f[0]],
314  reversePointMap[f[1]],
315  reversePointMap[f[2]]
316  );
317  faces[i] = newF.triFaceFace();
318 
319  forAll(newF, fp)
320  {
321  pointToFace[newF[fp]] = i;
322  }
323  }
324 
325  this->storedPoints() = pointField(s.points(), pointMap);
326 
327  if (debug)
328  {
329  print(Pout);
330  Pout<< endl;
331  }
332 
333 
334 
335  // Collect the samplePoints and sampleElements
336  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
337 
339  {
340  samplePoints_.setSize(pointMap.size());
341  sampleElements_.setSize(pointMap.size(), -1);
342 
343  if (sampleSource_ == cells)
344  {
345  // samplePoints_ : per surface point a location inside the cell
346  // sampleElements_ : per surface point the cell
347 
348  forAll(points(), pointi)
349  {
350  const point& pt = points()[pointi];
351  label celli = cellOrFaceLabels[pointToFace[pointi]];
352  sampleElements_[pointi] = celli;
353 
354  // Check if point inside cell
355  if
356  (
357  mesh().pointInCell
358  (
359  pt,
360  sampleElements_[pointi],
361  meshSearcher.decompMode()
362  )
363  )
364  {
365  samplePoints_[pointi] = pt;
366  }
367  else
368  {
369  // Find nearest point on faces of cell
370  const cell& cFaces = mesh().cells()[celli];
371 
372  scalar minDistSqr = VGREAT;
373 
374  forAll(cFaces, i)
375  {
376  const face& f = mesh().faces()[cFaces[i]];
377  pointHit info = f.nearestPoint(pt, mesh().points());
378  if (info.distance() < minDistSqr)
379  {
380  minDistSqr = info.distance();
381  samplePoints_[pointi] = info.rawPoint();
382  }
383  }
384  }
385  }
386  }
387  else if (sampleSource_ == insideCells)
388  {
389  // samplePoints_ : per surface point a location inside the cell
390  // sampleElements_ : per surface point the cell
391 
392  forAll(points(), pointi)
393  {
394  const point& pt = points()[pointi];
395  label celli = cellOrFaceLabels[pointToFace[pointi]];
396  sampleElements_[pointi] = celli;
397  samplePoints_[pointi] = pt;
398  }
399  }
400  else
401  {
402  // samplePoints_ : per surface point a location on the boundary
403  // sampleElements_ : per surface point the boundary face containing
404  // the location
405 
406  forAll(points(), pointi)
407  {
408  const point& pt = points()[pointi];
409  label facei = cellOrFaceLabels[pointToFace[pointi]];
410  sampleElements_[pointi] = facei;
411  samplePoints_[pointi] = mesh().faces()[facei].nearestPoint
412  (
413  pt,
414  mesh().points()
415  ).rawPoint();
416  }
417  }
418  }
419  else
420  {
421  // if sampleSource_ == cells:
422  // samplePoints_ : n/a
423  // sampleElements_ : per surface triangle the cell
424  // if sampleSource_ == insideCells:
425  // samplePoints_ : n/a
426  // sampleElements_ : -1 or per surface triangle the cell
427  // else:
428  // samplePoints_ : n/a
429  // sampleElements_ : per surface triangle the boundary face
430  samplePoints_.clear();
431  sampleElements_.transfer(cellOrFaceLabels);
432  }
433 
434 
435  if (debug)
436  {
437  this->clearOut();
438  OFstream str(mesh().time().path()/"surfaceToMesh.obj");
439  Info<< "Dumping correspondence from local surface (points or faces)"
440  << " to mesh (cells or faces) to " << str.name() << endl;
441  label vertI = 0;
442 
444  {
445  if (sampleSource_ == cells || sampleSource_ == insideCells)
446  {
447  forAll(samplePoints_, pointi)
448  {
449  meshTools::writeOBJ(str, points()[pointi]);
450  vertI++;
451 
452  meshTools::writeOBJ(str, samplePoints_[pointi]);
453  vertI++;
454 
455  label celli = sampleElements_[pointi];
456  meshTools::writeOBJ(str, mesh().cellCentres()[celli]);
457  vertI++;
458  str << "l " << vertI-2 << ' ' << vertI-1 << ' ' << vertI
459  << nl;
460  }
461  }
462  else
463  {
464  forAll(samplePoints_, pointi)
465  {
466  meshTools::writeOBJ(str, points()[pointi]);
467  vertI++;
468 
469  meshTools::writeOBJ(str, samplePoints_[pointi]);
470  vertI++;
471 
472  label facei = sampleElements_[pointi];
473  meshTools::writeOBJ(str, mesh().faceCentres()[facei]);
474  vertI++;
475  str << "l " << vertI-2 << ' ' << vertI-1 << ' ' << vertI
476  << nl;
477  }
478  }
479  }
480  else
481  {
482  if (sampleSource_ == cells || sampleSource_ == insideCells)
483  {
484  forAll(sampleElements_, triI)
485  {
486  meshTools::writeOBJ(str, faceCentres()[triI]);
487  vertI++;
488 
489  label celli = sampleElements_[triI];
490  meshTools::writeOBJ(str, mesh().cellCentres()[celli]);
491  vertI++;
492  str << "l " << vertI-1 << ' ' << vertI << nl;
493  }
494  }
495  else
496  {
497  forAll(sampleElements_, triI)
498  {
499  meshTools::writeOBJ(str, faceCentres()[triI]);
500  vertI++;
501 
502  label facei = sampleElements_[triI];
503  meshTools::writeOBJ(str, mesh().faceCentres()[facei]);
504  vertI++;
505  str << "l " << vertI-1 << ' ' << vertI << nl;
506  }
507  }
508  }
509  }
510 
511 
512  needsUpdate_ = false;
513  return true;
514 }
515 
516 
517 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
518 
520 (
521  const word& name,
522  const polyMesh& mesh,
523  const word& surfaceName,
524  const samplingSource sampleSource
525 )
526 :
527  sampledSurface(name, mesh),
528  surface_
529  (
530  IOobject
531  (
532  surfaceName,
533  mesh.time().constant(), // instance
534  "triSurface", // local
535  mesh, // registry
538  false
539  )
540  ),
541  sampleSource_(sampleSource),
542  needsUpdate_(true),
543  sampleElements_(0),
544  samplePoints_(0)
545 {}
546 
547 
549 (
550  const word& name,
551  const polyMesh& mesh,
552  const dictionary& dict
553 )
554 :
555  sampledSurface(name, mesh, dict),
556  surface_
557  (
558  IOobject
559  (
560  dict.lookup("surface"),
561  mesh.time().constant(), // instance
562  "triSurface", // local
563  mesh, // registry
566  false
567  )
568  ),
569  sampleSource_(samplingSourceNames_[dict.lookup("source")]),
570  needsUpdate_(true),
571  sampleElements_(0),
572  samplePoints_(0)
573 {}
574 
575 
577 (
578  const word& name,
579  const polyMesh& mesh,
580  const triSurface& surface,
581  const word& sampleSourceName
582 )
583 :
584  sampledSurface(name, mesh),
585  surface_
586  (
587  IOobject
588  (
589  name,
590  mesh.time().constant(), // instance
591  "triSurface", // local
592  mesh, // registry
595  false
596  ),
597  surface
598  ),
599  sampleSource_(samplingSourceNames_[sampleSourceName]),
600  needsUpdate_(true),
601  sampleElements_(0),
602  samplePoints_(0)
603 {}
604 
605 
606 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
607 
609 {}
610 
611 
612 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
613 
615 {
616  return needsUpdate_;
617 }
618 
619 
621 {
622  // already marked as expired
623  if (needsUpdate_)
624  {
625  return false;
626  }
627 
630 
631  boundaryTreePtr_.clear();
632  sampleElements_.clear();
633  samplePoints_.clear();
634 
635  needsUpdate_ = true;
636  return true;
637 }
638 
639 
641 {
642  if (!needsUpdate_)
643  {
644  return false;
645  }
646 
647  // Mesh search engine, no triangulation of faces.
648  meshSearch meshSearcher(mesh(), polyMesh::FACE_PLANES);
649 
650  return update(meshSearcher);
651 }
652 
653 
655 {
656  if (!needsUpdate_)
657  {
658  return false;
659  }
660 
661  // Mesh search engine on subset, no triangulation of faces.
662  meshSearch meshSearcher(mesh(), bb, polyMesh::FACE_PLANES);
663 
664  return update(meshSearcher);
665 }
666 
667 
669 (
670  const volScalarField& vField
671 ) const
672 {
673  return sampleField(vField);
674 }
675 
676 
678 (
679  const volVectorField& vField
680 ) const
681 {
682  return sampleField(vField);
683 }
684 
686 (
687  const volSphericalTensorField& vField
688 ) const
689 {
690  return sampleField(vField);
691 }
692 
693 
695 (
696  const volSymmTensorField& vField
697 ) const
698 {
699  return sampleField(vField);
700 }
701 
702 
704 (
705  const volTensorField& vField
706 ) const
707 {
708  return sampleField(vField);
709 }
710 
711 
713 (
714  const interpolation<scalar>& interpolator
715 ) const
716 {
717  return interpolateField(interpolator);
718 }
719 
720 
722 (
723  const interpolation<vector>& interpolator
724 ) const
725 {
726  return interpolateField(interpolator);
727 }
728 
730 (
731  const interpolation<sphericalTensor>& interpolator
732 ) const
733 {
734  return interpolateField(interpolator);
735 }
736 
737 
739 (
740  const interpolation<symmTensor>& interpolator
741 ) const
742 {
743  return interpolateField(interpolator);
744 }
745 
746 
748 (
749  const interpolation<tensor>& interpolator
750 ) const
751 {
752  return interpolateField(interpolator);
753 }
754 
755 
757 {
758  os << "sampledTriSurfaceMesh: " << name() << " :"
759  << " surface:" << surface_.objectRegistry::name()
760  << " faces:" << faces().size()
761  << " points:" << points().size();
762 }
763 
764 
765 // ************************************************************************* //
polyMesh::cellDecomposition decompMode() const
Definition: meshSearch.H:199
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search...
Definition: meshSearch.H:57
virtual ~sampledTriSurfaceMesh()
Destructor.
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
#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
virtual bool expire()
Mark the surface as needing an update.
const indexedOctree< treeDataCell > & cellTree() const
Get (demand driven) reference to octree holding all cells.
Definition: meshSearch.C:598
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
const double e
Elementary charge.
Definition: doubleFloat.H:78
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
An abstract class for surfaces with sampling.
A 2-tuple for storing two objects of different types.
Definition: HashTable.H:66
Encapsulation of data needed to search for faces.
Definition: treeDataFace.H:58
tUEqn clear()
const Type1 & first() const
Return first.
Definition: Tuple2.H:94
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Output to file stream.
Definition: OFstream.H:82
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
bool interpolate() const
Interpolation requested for surface.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
patches[0]
const fileName & name() const
Return the name of the stream.
Definition: OFstream.H:120
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:52
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Macros for easy insertion into run-time selection tables.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
virtual bool update()
Update the surface as required.
samplingSource
Types of communications.
virtual bool needsUpdate() const
Does the surface need an update?
scalar y
const Point & hitPoint() const
Return hit point.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
dynamicFvMesh & mesh
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:310
virtual void clearGeom() const
const cellShapeList & cells
const pointField & points
A class for handling words, derived from string.
Definition: word.H:59
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:68
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
cachedRandom rndGen(label(0), -1)
bool hit() const
Is there a hit.
const Field< PointType > & points() const
Return reference to global points.
const Type & shapes() const
Reference to shape.
static const label labelMax
Definition: label.H:62
Triangle with additional region number.
Definition: labelledTri.H:57
const treeBoundBox & bb() const
Top bounding box.
Simple random number generator.
Definition: Random.H:49
Foam::polyBoundaryMesh.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
const Point & rawPoint() const
Return point with no checking.
Definition: PointHit.H:158
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
static const char nl
Definition: Ostream.H:262
label findInside(const point &) const
Find shape containing point. Only implemented for certain.
const Time & time() const
Return time.
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
defineTypeNameAndDebug(combustionModel, 0)
void operator()(nearInfo &x, const nearInfo &y) const
labelList f(nPoints)
treeBoundBox extend(Random &, const scalar s) const
Return slightly wider bounding box.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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
label patchi
vector point
Point is a vector.
Definition: point.H:41
pointHit nearestPoint(const point &p, const pointField &) const
Return nearest point to face.
label newPointi
Definition: readKivaGrid.H:501
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
A List with indirect addressing.
Definition: fvMatrix.H:106
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:300
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:87
messageStream Info
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
virtual void print(Ostream &) const
Write.
scalar distance() const
Return distance to hit.
Definition: PointHit.H:139
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
A class for managing temporary objects.
Definition: PtrList.H:53
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Triangulated surface description with patch information.
Definition: triSurface.H:65
A topoSetSource to select faces based on use of points.
Definition: pointToFace.H:49
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
sampledTriSurfaceMesh(const word &name, const polyMesh &mesh, const word &surfaceName, const samplingSource sampleSource)
Construct from components.
label index() const
Return index.
Tuple2< scalar, label > nearInfo
Private class for finding nearest.
face triFaceFace() const
Return triangle as a face.
Definition: triFaceI.H:140
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576