meshSurfaceMapperCornersAndEdges.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | cfMesh: A library for mesh generation
4  \\ / O peration |
5  \\ / A nd | Author: Franjo Juretic (franjo.juretic@c-fields.com)
6  \\/ M anipulation | Copyright (C) Creative Fields, Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of cfMesh.
10 
11  cfMesh is free software; you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by the
13  Free Software Foundation; either version 3 of the License, or (at your
14  option) any later version.
15 
16  cfMesh 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 cfMesh. If not, see <http://www.gnu.org/licenses/>.
23 
24 Description
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "meshOctree.H"
29 #include "triSurf.H"
30 #include "triSurfacePartitioner.H"
31 #include "meshSurfaceMapper.H"
32 #include "meshSurfaceEngine.H"
34 #include "meshSurfacePartitioner.H"
35 #include "labelledScalar.H"
36 
37 #include "helperFunctions.H"
38 
39 # ifdef USE_OMP
40 #include <omp.h>
41 # endif
42 
43 //#define DEBUGMapping
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
53 (
54  const labelLongList& nodesToMap,
55  scalarList& mappingDistance
56 ) const
57 {
58  const vectorField& faceCentres = surfaceEngine_.faceCentres();
59  const VRWGraph& pFaces = surfaceEngine_.pointFaces();
60  const labelList& bPoints = surfaceEngine_.boundaryPoints();
61  const pointFieldPMG& points = surfaceEngine_.points();
62 
63  //- generate search distance for corner nodes
64  mappingDistance.setSize(nodesToMap.size());
65  # ifdef USE_OMP
66  # pragma omp parallel for schedule(dynamic, 50)
67  # endif
68  forAll(nodesToMap, i)
69  {
70  const label bpI = nodesToMap[i];
71 
72  mappingDistance[i] = 0.0;
73 
74  const point& p = points[bPoints[bpI]];
75  forAllRow(pFaces, bpI, pfI)
76  {
77  const scalar d = magSqr(faceCentres[pFaces(bpI, pfI)] - p);
78  mappingDistance[i] = Foam::max(mappingDistance[i], d);
79  }
80 
81  //- safety factor
82  mappingDistance[i] *= 4.0;
83  }
84 
85  if( Pstream::parRun() )
86  {
87  //- make sure that corner nodesd at parallel boundaries
88  //- have the same range in which they accept the corners
89  const VRWGraph& bpAtProcs = surfaceEngine_.bpAtProcs();
90  const labelList& globalPointLabel =
91  surfaceEngine_.globalBoundaryPointLabel();
92 
93  //- create the map for exchanging data
94  std::map<label, DynList<labelledScalar> > exchangeData;
95  const DynList<label>& neiProcs = surfaceEngine_.bpNeiProcs();
96  forAll(neiProcs, i)
97  exchangeData.insert
98  (
99  std::make_pair(neiProcs[i], DynList<labelledScalar>())
100  );
101 
102  Map<label> globalToLocal;
103 
104  forAll(nodesToMap, nI)
105  {
106  const label bpI = nodesToMap[nI];
107 
108  if( bpAtProcs.sizeOfRow(bpI) != 0 )
109  globalToLocal.insert(globalPointLabel[bpI], nI);
110 
111  forAllRow(bpAtProcs, bpI, i)
112  {
113  const label neiProc = bpAtProcs(bpI, i);
114  if( neiProc == Pstream::myProcNo() )
115  continue;
116 
117  exchangeData[neiProc].append
118  (
119  labelledScalar(globalPointLabel[bpI], mappingDistance[nI])
120  );
121  }
122  }
123 
124  //- exchange data between processors
125  LongList<labelledScalar> receivedData;
126  help::exchangeMap(exchangeData, receivedData);
127 
128  //- select the maximum mapping distance for processor points
129  forAll(receivedData, i)
130  {
131  const labelledScalar& ls = receivedData[i];
132 
133  const label nI = globalToLocal[ls.scalarLabel()];
134 
135  //- choose the maximum value for the mapping distance
136  mappingDistance[nI] = Foam::max(mappingDistance[nI], ls.value());
137  }
138  }
139 }
140 
142 (
143  const label bfI,
144  const label patchI
145 ) const
146 {
147  const face& bf = surfaceEngine_.boundaryFaces()[bfI];
148 
149  const pointFieldPMG& points = surfaceEngine_.points();
150 
151  const point centre = bf.centre(points);
152  const vector area = bf.normal(points);
153 
154  point projCentre;
155  scalar dSq;
156  label nt;
157 
158  meshOctree_.findNearestSurfacePointInRegion
159  (
160  projCentre,
161  dSq,
162  nt,
163  patchI,
164  centre
165  );
166 
167  DynList<point> projPoints(bf.size());
168  forAll(bf, pI)
169  {
170  meshOctree_.findNearestSurfacePointInRegion
171  (
172  projPoints[pI],
173  dSq,
174  nt,
175  patchI,
176  points[bf[pI]]
177  );
178  }
179 
180  vector projArea(vector::zero);
181  forAll(bf, pI)
182  {
183  projArea +=
185  (
186  projPoints[pI],
187  projPoints[bf.fcIndex(pI)],
188  projCentre
189  ).normal();
190  }
191 
192  return magSqr(centre - projCentre) + mag(mag(projArea) - mag(area));
193 }
194 
196 {
197  const triSurfacePartitioner& sPartitioner = surfacePartitioner();
198  const labelList& surfCorners = sPartitioner.corners();
199  const List<DynList<label> >& cornerPatches = sPartitioner.cornerPatches();
200 
201  const meshSurfacePartitioner& mPart = meshPartitioner();
202  const labelHashSet& corners = mPart.corners();
203  const VRWGraph& pPatches = mPart.pointPatches();
204 
206  const labelList& bPoints = surfaceEngine_.boundaryPoints();
207 
208  const triSurf& surf = meshOctree_.surface();
209  const pointField& sPoints = surf.points();
210 
211  //std::map<label, scalar> mappingDistance;
212  scalarList mappingDistance;
213  findMappingDistance(nodesToMap, mappingDistance);
214 
215  //- for every corner in the mesh surface find the nearest corner in the
216  //- triSurface
218 
219  # ifdef USE_OMP
220  # pragma omp parallel for schedule(dynamic, 50)
221  # endif
222  forAll(nodesToMap, cornerI)
223  {
224  const label bpI = nodesToMap[cornerI];
225  if( !corners.found(bpI) )
227  (
228  "meshSurfaceMapper::mapCorners(const labelLongList&)"
229  ) << "Trying to map a point that is not a corner"
230  << abort(FatalError);
231 
232  const point& p = points[bPoints[bpI]];
233  const scalar maxDist = mappingDistance[cornerI];
234 
235  //- find the nearest position to the given point patches
236  const DynList<label> patches = pPatches[bpI];
237 
238  point mapPointApprox(p);
239  scalar distSqApprox;
240 
241  label iter(0);
242  while( iter++ < 20 )
243  {
244  point newP(vector::zero);
245  forAll(patches, patchI)
246  {
247  point np;
248  label nt;
250  (
251  np,
252  distSqApprox,
253  nt,
254  patches[patchI],
255  mapPointApprox
256  );
257 
258  newP += np;
259  }
260 
261  newP /= patches.size();
262 
263  if( magSqr(newP - mapPointApprox) < 1e-8 * maxDist )
264  break;
265 
266  mapPointApprox = newP;
267  }
268  distSqApprox = magSqr(mapPointApprox - p);
269 
270  //- find the nearest triSurface corner for the given corner
271  scalar distSq(mappingDistance[cornerI]);
272  point mapPoint(p);
273  forAll(surfCorners, scI)
274  {
275  const label cornerID = surfCorners[scI];
276  const point& sp = sPoints[cornerID];
277 
278  if( Foam::magSqr(sp - p) < distSq )
279  {
280  bool store(true);
281  const DynList<label>& cPatches = cornerPatches[scI];
282  forAll(patches, i)
283  {
284  if( !cPatches.contains(patches[i]) )
285  {
286  store = false;
287  break;
288  }
289  }
290 
291  if( store )
292  {
293  mapPoint = sp;
294  distSq = Foam::magSqr(sp - p);
295  }
296  }
297  }
298 
299  if( distSq > 1.2 * distSqApprox )
300  {
301  mapPoint = mapPointApprox;
302  }
303 
304  //- move the point to the nearest corner
305  sMod.moveBoundaryVertexNoUpdate(bpI, mapPoint);
306  }
307 
308  sMod.updateGeometry(nodesToMap);
309 }
310 
312 {
314  const labelList& bPoints = surfaceEngine_.boundaryPoints();
315 
316  const meshSurfacePartitioner& mPart = meshPartitioner();
317  const VRWGraph& pPatches = mPart.pointPatches();
318 
319  //const triSurf& surf = meshOctree_.surface();
320  //const pointField& sPoints = surf.points();
321 
322  //- find mapping distance for selected vertices
323  scalarList mappingDistance;
324  findMappingDistance(nodesToMap, mappingDistance);
325 
326  const VRWGraph* bpAtProcsPtr(NULL);
327  if( Pstream::parRun() )
328  bpAtProcsPtr = &surfaceEngine_.bpAtProcs();
329 
330  LongList<parMapperHelper> parallelBndNodes;
331 
333 
334  //- map point to the nearest vertex on the surface mesh
335  # ifdef USE_OMP
336  # pragma omp parallel for schedule(dynamic, 50)
337  # endif
338  forAll(nodesToMap, i)
339  {
340  const label bpI = nodesToMap[i];
341  const point& p = points[bPoints[bpI]];
342 
343  //- find patches at this edge point
344  const DynList<label> patches = pPatches[bpI];
345 
346  const scalar maxDist = mappingDistance[i];
347 
348  //- find approximate position of the vertex on the edge
349  point mapPointApprox(p);
350  scalar distSqApprox;
351  label iter(0);
352  while( iter++ < 20 )
353  {
354  point newP(vector::zero);
355 
356  forAll(patches, patchI)
357  {
358  point np;
359  label nt;
361  (
362  np,
363  distSqApprox,
364  nt,
365  patches[patchI],
366  mapPointApprox
367  );
368 
369  newP += np;
370  }
371 
372  newP /= patches.size();
373 
374  if( magSqr(newP - mapPointApprox) < 1e-8 * maxDist )
375  break;
376 
377  mapPointApprox = newP;
378  }
379  distSqApprox = magSqr(mapPointApprox - p);
380 
381  //- find the nearest vertex on the triSurface feature edge
382  point mapPoint(mapPointApprox);
383  scalar distSq(distSqApprox);
384  label nse;
385  meshOctree_.findNearestEdgePoint(mapPoint, distSq, nse, p, patches);
386 
387  //- use the vertex with the smallest mapping distance
388  if( distSq > 1.2 * distSqApprox )
389  {
390  mapPoint = mapPointApprox;
391  distSq = distSqApprox;
392  }
393 
394  //- check if the mapping distance is within the given tolerances
395  if( distSq > maxDist )
396  {
397  //- this indicates possible problems
398  //- reduce the mapping distance
399  const scalar f = Foam::sqrt(maxDist / distSq);
400  distSq = mappingDistance[i];
401  mapPoint = f * (mapPoint - p) + p;
402  }
403 
404  //- move the point to the nearest edge vertex
405  sMod.moveBoundaryVertexNoUpdate(bpI, mapPoint);
406 
407  if( bpAtProcsPtr && bpAtProcsPtr->sizeOfRow(bpI) )
408  {
409  # ifdef USE_OMP
410  # pragma omp critical
411  # endif
412  {
413  parallelBndNodes.append
414  (
416  (
417  mapPoint,
418  distSq,
419  bpI,
420  -1
421  )
422  );
423  }
424  }
425  }
426 
427  sMod.updateGeometry(nodesToMap);
428 
429  mapToSmallestDistance(parallelBndNodes);
430 }
431 
433 {
434  const meshSurfacePartitioner& mPart = meshPartitioner();
435  const labelHashSet& cornerPoints = mPart.corners();
436  labelLongList selectedPoints;
437  forAllConstIter(labelHashSet, cornerPoints, it)
438  selectedPoints.append(it.key());
439 
440  mapCorners(selectedPoints);
441 
442  selectedPoints.clear();
443  const labelHashSet& edgePoints = mPart.edgePoints();
444  forAllConstIter(labelHashSet, edgePoints, it)
445  selectedPoints.append(it.key());
446 
447  mapEdgeNodes(selectedPoints);
448 }
449 
450 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
451 
452 } // End namespace Foam
453 
454 // ************************************************************************* //
Foam::meshSurfaceEngine::bpAtProcs
const VRWGraph & bpAtProcs() const
processors which contain the vertex
Definition: meshSurfaceEngineI.H:451
Foam::face::normal
vector normal(const pointField &) const
Vector normal; magnitude is equal to area of face.
Definition: face.C:552
Foam::triSurfacePartitioner
Definition: triSurfacePartitioner.H:53
Foam::LongList::append
void append(const T &e)
Append an element at the end of the list.
Definition: LongListI.H:265
Foam::meshSurfaceMapper::surfaceEngine_
const meshSurfaceEngine & surfaceEngine_
mesh surface
Definition: meshSurfaceMapper.H:63
triSurf.H
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
Foam::parMapperHelper
Definition: parMapperHelper.H:51
p
p
Definition: pEqn.H:62
meshSurfaceMapper.H
Foam::meshSurfaceMapper::faceMetricInPatch
scalar faceMetricInPatch(const label bfI, const label patchI) const
calculate face metric
Definition: meshSurfaceMapperCornersAndEdges.C:142
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
triSurfacePartitioner.H
Foam::triSurfPoints::points
const pointField & points() const
access to points
Definition: triSurfPointsI.H:44
Foam::LongList::clear
void clear()
Clear the list, i.e. set next free to zero.
Definition: LongListI.H:230
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:377
Foam::face::centre
point centre(const pointField &) const
Centre point of face.
Definition: face.C:488
Foam::Map< label >
Foam::meshSurfaceMapper::surfacePartitioner
const triSurfacePartitioner & surfacePartitioner() const
Definition: meshSurfaceMapper.H:88
Foam::LongList::size
label size() const
Size of the active part of the list.
Definition: LongListI.H:203
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::HashSet< label, Hash< label > >
Foam::DynList::contains
bool contains(const T &e) const
check if the element is in the list (takes linear time)
Definition: DynListI.H:339
Foam::meshSurfaceMapper::findMappingDistance
void findMappingDistance(const labelLongList &nodesToMap, scalarList &mappingDistance) const
find mapping distance for selected points
Definition: meshSurfaceMapperCornersAndEdges.C:53
Foam::labelledScalar::scalarLabel
label scalarLabel() const
return scalar label
Definition: labelledScalar.H:82
Foam::meshSurfacePartitioner::pointPatches
const VRWGraph & pointPatches() const
Definition: meshSurfacePartitioner.H:123
meshOctree.H
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
pFaces
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=*(cellModeller::lookup("hex"));labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells]=cellShape(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:235
Foam::labelledScalar::value
const scalar & value() const
return the value
Definition: labelledScalar.H:88
Foam::LongList< label >
Foam::meshSurfaceEngineModifier
Definition: meshSurfaceEngineModifier.H:48
Foam::meshSurfaceMapper::mapToSmallestDistance
void mapToSmallestDistance(LongList< parMapperHelper > &)
Definition: meshSurfaceMapperMapVertices.C:104
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:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::meshOctree::findNearestEdgePoint
bool findNearestEdgePoint(point &edgePoint, scalar &distSq, label &nearestEdge, const point &p, const DynList< label > &regions) const
find nearest feature-edges vertex to a given vertex
Definition: meshOctreeFindNearestSurfacePoint.C:215
Foam::meshSurfaceEngineModifier::moveBoundaryVertexNoUpdate
void moveBoundaryVertexNoUpdate(const label bpI, const point &newP)
relocate the selected boundary vertex
Definition: meshSurfaceEngineModifier.C:68
forAllRow
#define forAllRow(graph, rowI, index)
Definition: VRWGraph.H:277
Foam::meshSurfaceEngine::points
const pointFieldPMG & points() const
Definition: meshSurfaceEngineI.H:49
Foam::VRWGraph::sizeOfRow
label sizeOfRow(const label rowI) const
Returns the number of elements in the given row.
Definition: VRWGraphI.H:127
Foam::FatalError
error FatalError
labelledScalar.H
Foam::HashTable::found
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::meshSurfaceEngineModifier::updateGeometry
void updateGeometry(const labelLongList &)
Definition: meshSurfaceEngineModifier.C:234
meshSurfaceEngine.H
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::DynList< label >
Foam::meshSurfacePartitioner
Definition: meshSurfacePartitioner.H:52
Foam::meshSurfaceMapper::mapCornersAndEdges
void mapCornersAndEdges()
Definition: meshSurfaceMapperCornersAndEdges.C:432
Foam::labelledScalar
Definition: labelledScalar.H:50
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:405
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::meshSurfaceMapper::mapEdgeNodes
void mapEdgeNodes(const labelLongList &nodesToMap)
Definition: meshSurfaceMapperCornersAndEdges.C:311
Foam::meshSurfacePartitioner::edgePoints
const labelHashSet & edgePoints() const
return labels of edge points (from the list of boundary points)
Definition: meshSurfacePartitioner.H:135
helperFunctions.H
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::help::exchangeMap
void exchangeMap(const std::map< label, ListType > &m, LongList< T > &data, const Pstream::commsTypes commsType)
Definition: helperFunctionsPar.C:129
Foam::triSurfacePartitioner::corners
const labelList & corners() const
return corner nodes
Definition: triSurfacePartitioner.C:64
Foam::meshSurfaceMapper::meshPartitioner
const meshSurfacePartitioner & meshPartitioner() const
Definition: meshSurfaceMapper.H:78
meshSurfaceEngineModifier.H
Foam::meshSurfaceEngine::boundaryPoints
const labelList & boundaryPoints() const
Definition: meshSurfaceEngineI.H:84
patches
patches[0]
Definition: createSingleCellMesh.H:36
FatalErrorIn
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:313
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::meshOctree::findNearestSurfacePointInRegion
void findNearestSurfacePointInRegion(point &nearest, scalar &distSq, label &nearestTriangle, const label region, const point &p) const
find nearest surface point for vertex in a given region
Definition: meshOctreeFindNearestSurfacePoint.C:131
Foam::meshOctree::surface
const triSurf & surface() const
return a reference to the surface
Definition: meshOctreeI.H:130
Foam::meshSurfacePartitioner::corners
const labelHashSet & corners() const
return labels of corner points (from the list of boundary points)
Definition: meshSurfacePartitioner.H:129
Foam::pointFieldPMG
Definition: pointFieldPMG.H:50
Foam::VRWGraph
Definition: VRWGraph.H:101
Foam::meshSurfaceMapper::mapCorners
void mapCorners(const labelLongList &nodesToMap)
map corner nodes to the boundary
Definition: meshSurfaceMapperCornersAndEdges.C:195
Foam::triSurf
Definition: triSurf.H:59
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Foam::meshSurfaceMapper::meshOctree_
const meshOctree & meshOctree_
reference to the octree
Definition: meshSurfaceMapper.H:66
Foam::triSurfacePartitioner::cornerPatches
const List< DynList< label > > & cornerPatches() const
return corner patches
Definition: triSurfacePartitioner.C:69
Foam::triPointRef
triangle< point, const point & > triPointRef
Definition: triangle.H:75
meshSurfacePartitioner.H