sampledTriSurfaceMesh.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 "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  overallBb = overallBb.extend(1e-4);
113 
114  boundaryTreePtr_.reset
115  (
117  (
118  treeDataFace // all information needed to search faces
119  (
120  false, // do not cache bb
121  mesh(),
122  bndFaces // boundary faces only
123  ),
124  overallBb, // overall search domain
125  8, // maxLevel
126  10, // leafsize
127  3.0 // duplicity
128  )
129  );
130  }
131 
132  return boundaryTreePtr_();
133 }
134 
135 
136 bool Foam::sampledTriSurfaceMesh::update(const meshSearch& meshSearcher)
137 {
138  // Find the cells the triangles of the surface are in.
139  // Does approximation by looking at the face centres only
140  const pointField& fc = surface_.faceCentres();
141 
142  List<nearInfo> nearest(fc.size());
143 
144  // Global numbering for cells/faces - only used to uniquely identify local
145  // elements
146  globalIndex globalCells
147  (
148  (sampleSource_ == cells || sampleSource_ == insideCells)
149  ? mesh().nCells()
150  : mesh().nFaces()
151  );
152 
153  forAll(nearest, i)
154  {
155  nearest[i].first() = great;
156  nearest[i].second() = labelMax;
157  }
158 
159  if (sampleSource_ == cells)
160  {
161  // Search for nearest cell
162 
163  const indexedOctree<treeDataCell>& cellTree = meshSearcher.cellTree();
164 
165  forAll(fc, triI)
166  {
167  pointIndexHit nearInfo = cellTree.findNearest
168  (
169  fc[triI],
170  sqr(great)
171  );
172  if (nearInfo.hit())
173  {
174  nearest[triI].first() = magSqr(nearInfo.hitPoint()-fc[triI]);
175  nearest[triI].second() = globalCells.toGlobal(nearInfo.index());
176  }
177  }
178  }
179  else if (sampleSource_ == insideCells)
180  {
181  // Search for cell containing point
182 
183  const indexedOctree<treeDataCell>& cellTree = meshSearcher.cellTree();
184 
185  forAll(fc, triI)
186  {
187  if (cellTree.bb().contains(fc[triI]))
188  {
189  label index = cellTree.findInside(fc[triI]);
190  if (index != -1)
191  {
192  nearest[triI].first() = 0.0;
193  nearest[triI].second() = globalCells.toGlobal(index);
194  }
195  }
196  }
197  }
198  else
199  {
200  // Search on all non-coupled boundary faces
201 
202  const indexedOctree<treeDataFace>& bTree = nonCoupledboundaryTree();
203 
204  forAll(fc, triI)
205  {
206  pointIndexHit nearInfo = bTree.findNearest
207  (
208  fc[triI],
209  sqr(great)
210  );
211  if (nearInfo.hit())
212  {
213  nearest[triI].first() = magSqr(nearInfo.hitPoint()-fc[triI]);
214  nearest[triI].second() = globalCells.toGlobal
215  (
216  bTree.shapes().faceLabels()[nearInfo.index()]
217  );
218  }
219  }
220  }
221 
222 
223  // See which processor has the nearest. Mark and subset
224  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
225 
228 
229  labelList cellOrFaceLabels(fc.size(), -1);
230 
231  label nFound = 0;
232  forAll(nearest, triI)
233  {
234  if (nearest[triI].second() == labelMax)
235  {
236  // Not found on any processor. How to map?
237  }
238  else if (globalCells.isLocal(nearest[triI].second()))
239  {
240  cellOrFaceLabels[triI] = globalCells.toLocal
241  (
242  nearest[triI].second()
243  );
244  nFound++;
245  }
246  }
247 
248 
249  if (debug)
250  {
251  Pout<< "Local out of faces:" << cellOrFaceLabels.size()
252  << " keeping:" << nFound << endl;
253  }
254 
255  // Now subset the surface. Do not use triSurface::subsetMesh since requires
256  // original surface to be in compact numbering.
257 
258  const triSurface& s = surface_;
259 
260  // Compact to original triangle
261  labelList faceMap(s.size());
262  // Compact to original points
263  labelList pointMap(s.points().size());
264  // From original point to compact points
265  labelList reversePointMap(s.points().size(), -1);
266 
267  {
268  label newPointi = 0;
269  label newFacei = 0;
270 
271  forAll(s, facei)
272  {
273  if (cellOrFaceLabels[facei] != -1)
274  {
275  faceMap[newFacei++] = facei;
276 
277  const triSurface::FaceType& f = s[facei];
278  forAll(f, fp)
279  {
280  if (reversePointMap[f[fp]] == -1)
281  {
282  pointMap[newPointi] = f[fp];
283  reversePointMap[f[fp]] = newPointi++;
284  }
285  }
286  }
287  }
288  faceMap.setSize(newFacei);
289  pointMap.setSize(newPointi);
290  }
291 
292 
293  // Subset cellOrFaceLabels
294  cellOrFaceLabels = UIndirectList<label>(cellOrFaceLabels, faceMap)();
295 
296  // Store any face per point (without using pointFaces())
297  labelList pointToFace(pointMap.size());
298 
299  // Create faces and points for subsetted surface
300  faceList& faces = this->storedFaces();
301  faces.setSize(faceMap.size());
302  forAll(faceMap, i)
303  {
304  const triFace& f = s[faceMap[i]];
305  triFace newF
306  (
307  reversePointMap[f[0]],
308  reversePointMap[f[1]],
309  reversePointMap[f[2]]
310  );
311  faces[i] = newF.triFaceFace();
312 
313  forAll(newF, fp)
314  {
315  pointToFace[newF[fp]] = i;
316  }
317  }
318 
319  this->storedPoints() = pointField(s.points(), pointMap);
320 
321  if (debug)
322  {
323  print(Pout);
324  Pout<< endl;
325  }
326 
327 
328 
329  // Collect the samplePoints and sampleElements
330  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
331 
333  {
334  samplePoints_.setSize(pointMap.size());
335  sampleElements_.setSize(pointMap.size(), -1);
336 
337  if (sampleSource_ == cells)
338  {
339  // samplePoints_ : per surface point a location inside the cell
340  // sampleElements_ : per surface point the cell
341 
342  forAll(points(), pointi)
343  {
344  const point& pt = points()[pointi];
345  label celli = cellOrFaceLabels[pointToFace[pointi]];
346  sampleElements_[pointi] = celli;
347 
348  // Check if point inside cell
349  if
350  (
351  mesh().pointInCell
352  (
353  pt,
354  sampleElements_[pointi],
355  meshSearcher.decompMode()
356  )
357  )
358  {
359  samplePoints_[pointi] = pt;
360  }
361  else
362  {
363  // Find nearest point on faces of cell
364  const cell& cFaces = mesh().cells()[celli];
365 
366  scalar minDistSqr = vGreat;
367 
368  forAll(cFaces, i)
369  {
370  const face& f = mesh().faces()[cFaces[i]];
371  pointHit info = f.nearestPoint(pt, mesh().points());
372  if (info.distance() < minDistSqr)
373  {
374  minDistSqr = info.distance();
375  samplePoints_[pointi] = info.rawPoint();
376  }
377  }
378  }
379  }
380  }
381  else if (sampleSource_ == insideCells)
382  {
383  // samplePoints_ : per surface point a location inside the cell
384  // sampleElements_ : per surface point the cell
385 
386  forAll(points(), pointi)
387  {
388  const point& pt = points()[pointi];
389  label celli = cellOrFaceLabels[pointToFace[pointi]];
390  sampleElements_[pointi] = celli;
391  samplePoints_[pointi] = pt;
392  }
393  }
394  else
395  {
396  // samplePoints_ : per surface point a location on the boundary
397  // sampleElements_ : per surface point the boundary face containing
398  // the location
399 
400  forAll(points(), pointi)
401  {
402  const point& pt = points()[pointi];
403  label facei = cellOrFaceLabels[pointToFace[pointi]];
404  sampleElements_[pointi] = facei;
405  samplePoints_[pointi] = mesh().faces()[facei].nearestPoint
406  (
407  pt,
408  mesh().points()
409  ).rawPoint();
410  }
411  }
412  }
413  else
414  {
415  // if sampleSource_ == cells:
416  // samplePoints_ : n/a
417  // sampleElements_ : per surface triangle the cell
418  // if sampleSource_ == insideCells:
419  // samplePoints_ : n/a
420  // sampleElements_ : -1 or per surface triangle the cell
421  // else:
422  // samplePoints_ : n/a
423  // sampleElements_ : per surface triangle the boundary face
424  samplePoints_.clear();
425  sampleElements_.transfer(cellOrFaceLabels);
426  }
427 
428 
429  if (debug)
430  {
431  this->clearOut();
432  OFstream str(mesh().time().path()/"surfaceToMesh.obj");
433  Info<< "Dumping correspondence from local surface (points or faces)"
434  << " to mesh (cells or faces) to " << str.name() << endl;
435  label vertI = 0;
436 
438  {
439  if (sampleSource_ == cells || sampleSource_ == insideCells)
440  {
441  forAll(samplePoints_, pointi)
442  {
443  meshTools::writeOBJ(str, points()[pointi]);
444  vertI++;
445 
446  meshTools::writeOBJ(str, samplePoints_[pointi]);
447  vertI++;
448 
449  label celli = sampleElements_[pointi];
450  meshTools::writeOBJ(str, mesh().cellCentres()[celli]);
451  vertI++;
452  str << "l " << vertI-2 << ' ' << vertI-1 << ' ' << vertI
453  << nl;
454  }
455  }
456  else
457  {
458  forAll(samplePoints_, pointi)
459  {
460  meshTools::writeOBJ(str, points()[pointi]);
461  vertI++;
462 
463  meshTools::writeOBJ(str, samplePoints_[pointi]);
464  vertI++;
465 
466  label facei = sampleElements_[pointi];
467  meshTools::writeOBJ(str, mesh().faceCentres()[facei]);
468  vertI++;
469  str << "l " << vertI-2 << ' ' << vertI-1 << ' ' << vertI
470  << nl;
471  }
472  }
473  }
474  else
475  {
476  if (sampleSource_ == cells || sampleSource_ == insideCells)
477  {
478  forAll(sampleElements_, triI)
479  {
480  meshTools::writeOBJ(str, faceCentres()[triI]);
481  vertI++;
482 
483  label celli = sampleElements_[triI];
484  meshTools::writeOBJ(str, mesh().cellCentres()[celli]);
485  vertI++;
486  str << "l " << vertI-1 << ' ' << vertI << nl;
487  }
488  }
489  else
490  {
491  forAll(sampleElements_, triI)
492  {
493  meshTools::writeOBJ(str, faceCentres()[triI]);
494  vertI++;
495 
496  label facei = sampleElements_[triI];
497  meshTools::writeOBJ(str, mesh().faceCentres()[facei]);
498  vertI++;
499  str << "l " << vertI-1 << ' ' << vertI << nl;
500  }
501  }
502  }
503  }
504 
505 
506  needsUpdate_ = false;
507  return true;
508 }
509 
510 
511 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
512 
514 (
515  const word& name,
516  const polyMesh& mesh,
517  const word& surfaceName,
518  const samplingSource sampleSource
519 )
520 :
521  sampledSurface(name, mesh),
522  surface_
523  (
524  IOobject
525  (
526  surfaceName,
527  mesh.time().constant(), // instance
528  "triSurface", // local
529  mesh, // registry
532  false
533  )
534  ),
535  sampleSource_(sampleSource),
536  needsUpdate_(true),
537  sampleElements_(0),
538  samplePoints_(0)
539 {}
540 
541 
543 (
544  const word& name,
545  const polyMesh& mesh,
546  const dictionary& dict
547 )
548 :
549  sampledSurface(name, mesh, dict),
550  surface_
551  (
552  IOobject
553  (
554  dict.lookup("surface"),
555  mesh.time().constant(), // instance
556  "triSurface", // local
557  mesh, // registry
560  false
561  )
562  ),
563  sampleSource_(samplingSourceNames_[dict.lookup("source")]),
564  needsUpdate_(true),
565  sampleElements_(0),
566  samplePoints_(0)
567 {}
568 
569 
571 (
572  const word& name,
573  const polyMesh& mesh,
574  const triSurface& surface,
575  const word& sampleSourceName
576 )
577 :
578  sampledSurface(name, mesh),
579  surface_
580  (
581  IOobject
582  (
583  name,
584  mesh.time().constant(), // instance
585  "triSurface", // local
586  mesh, // registry
589  false
590  ),
591  surface
592  ),
593  sampleSource_(samplingSourceNames_[sampleSourceName]),
594  needsUpdate_(true),
595  sampleElements_(0),
596  samplePoints_(0)
597 {}
598 
599 
600 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
601 
603 {}
604 
605 
606 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
607 
609 {
610  return needsUpdate_;
611 }
612 
613 
615 {
616  // already marked as expired
617  if (needsUpdate_)
618  {
619  return false;
620  }
621 
624 
625  boundaryTreePtr_.clear();
626  sampleElements_.clear();
627  samplePoints_.clear();
628 
629  needsUpdate_ = true;
630  return true;
631 }
632 
633 
635 {
636  if (!needsUpdate_)
637  {
638  return false;
639  }
640 
641  // Mesh search engine, no triangulation of faces.
642  meshSearch meshSearcher(mesh(), polyMesh::FACE_PLANES);
643 
644  return update(meshSearcher);
645 }
646 
647 
649 {
650  if (!needsUpdate_)
651  {
652  return false;
653  }
654 
655  // Mesh search engine on subset, no triangulation of faces.
656  meshSearch meshSearcher(mesh(), bb, polyMesh::FACE_PLANES);
657 
658  return update(meshSearcher);
659 }
660 
661 
663 (
664  const volScalarField& vField
665 ) const
666 {
667  return sampleField(vField);
668 }
669 
670 
672 (
673  const volVectorField& vField
674 ) const
675 {
676  return sampleField(vField);
677 }
678 
680 (
681  const volSphericalTensorField& vField
682 ) const
683 {
684  return sampleField(vField);
685 }
686 
687 
689 (
690  const volSymmTensorField& vField
691 ) const
692 {
693  return sampleField(vField);
694 }
695 
696 
698 (
699  const volTensorField& vField
700 ) const
701 {
702  return sampleField(vField);
703 }
704 
705 
707 (
708  const interpolation<scalar>& interpolator
709 ) const
710 {
711  return interpolateField(interpolator);
712 }
713 
714 
716 (
717  const interpolation<vector>& interpolator
718 ) const
719 {
720  return interpolateField(interpolator);
721 }
722 
724 (
725  const interpolation<sphericalTensor>& interpolator
726 ) const
727 {
728  return interpolateField(interpolator);
729 }
730 
731 
733 (
734  const interpolation<symmTensor>& interpolator
735 ) const
736 {
737  return interpolateField(interpolator);
738 }
739 
740 
742 (
743  const interpolation<tensor>& interpolator
744 ) const
745 {
746  return interpolateField(interpolator);
747 }
748 
749 
751 {
752  os << "sampledTriSurfaceMesh: " << name() << " :"
753  << " surface:" << surface_.objectRegistry::name()
754  << " faces:" << faces().size()
755  << " points:" << points().size();
756 }
757 
758 
759 // ************************************************************************* //
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:596
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
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:256
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:51
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
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.
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:265
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)
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
label patchi
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
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
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
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:66
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.
treeBoundBox extend(const scalar s) const
Return asymetrically extended bounding box, with guaranteed.
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