treeDataFace.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-2013 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 "treeDataFace.H"
27 #include "polyMesh.H"
28 #include "triangleFuncs.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 defineTypeNameAndDebug(treeDataFace, 0);
35 
36 scalar treeDataFace::tolSqr = sqr(1e-6);
37 }
38 
39 
40 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41 
42 Foam::treeBoundBox Foam::treeDataFace::calcBb(const label faceI) const
43 {
44  const pointField& points = mesh_.points();
45 
46  const face& f = mesh_.faces()[faceI];
47 
48  treeBoundBox bb(points[f[0]], points[f[0]]);
49 
50  for (label fp = 1; fp < f.size(); fp++)
51  {
52  const point& p = points[f[fp]];
53 
54  bb.min() = min(bb.min(), p);
55  bb.max() = max(bb.max(), p);
56  }
57  return bb;
58 }
59 
60 
61 void Foam::treeDataFace::update()
62 {
63  forAll(faceLabels_, i)
64  {
65  isTreeFace_.set(faceLabels_[i], 1);
66  }
67 
68  if (cacheBb_)
69  {
70  bbs_.setSize(faceLabels_.size());
71 
72  forAll(faceLabels_, i)
73  {
74  bbs_[i] = calcBb(faceLabels_[i]);
75  }
76  }
77 }
78 
79 
80 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
81 
83 (
84  const bool cacheBb,
85  const primitiveMesh& mesh,
86  const labelUList& faceLabels
87 )
88 :
89  mesh_(mesh),
90  faceLabels_(faceLabels),
91  isTreeFace_(mesh.nFaces(), 0),
92  cacheBb_(cacheBb)
93 {
94  update();
95 }
96 
97 
99 (
100  const bool cacheBb,
101  const primitiveMesh& mesh,
102  const Xfer<labelList>& faceLabels
103 )
104 :
105  mesh_(mesh),
106  faceLabels_(faceLabels),
107  isTreeFace_(mesh.nFaces(), 0),
108  cacheBb_(cacheBb)
109 {
110  update();
111 }
112 
113 
115 (
116  const bool cacheBb,
117  const primitiveMesh& mesh
118 )
119 :
120  mesh_(mesh),
121  faceLabels_(identity(mesh_.nFaces())),
122  isTreeFace_(mesh.nFaces(), 0),
123  cacheBb_(cacheBb)
124 {
125  update();
126 }
127 
128 
130 (
131  const bool cacheBb,
132  const polyPatch& patch
133 )
134 :
135  mesh_(patch.boundaryMesh().mesh()),
136  faceLabels_
137  (
138  identity(patch.size())
139  + patch.start()
140  ),
141  isTreeFace_(mesh_.nFaces(), 0),
142  cacheBb_(cacheBb)
143 {
144  update();
145 }
146 
147 
149 (
150  const indexedOctree<treeDataFace>& tree
151 )
152 :
153  tree_(tree)
154 {}
155 
156 
158 (
159  const indexedOctree<treeDataFace>& tree
160 )
161 :
162  tree_(tree)
163 {}
164 
165 
166 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
167 
169 {
170  pointField cc(faceLabels_.size());
171 
172  forAll(faceLabels_, i)
173  {
174  cc[i] = mesh_.faceCentres()[faceLabels_[i]];
175  }
176 
177  return cc;
178 }
179 
180 
181 //- Get type (inside,outside,mixed,unknown) of point w.r.t. surface.
182 // Only makes sense for closed surfaces.
184 (
185  const indexedOctree<treeDataFace>& oc,
186  const point& sample
187 ) const
188 {
189  // Need to determine whether sample is 'inside' or 'outside'
190  // Done by finding nearest face. This gives back a face which is
191  // guaranteed to contain nearest point. This point can be
192  // - in interior of face: compare to face normal
193  // - on edge of face: compare to edge normal
194  // - on point of face: compare to point normal
195  // Unfortunately the octree does not give us back the intersection point
196  // or where on the face it has hit so we have to recreate all that
197  // information.
198 
199 
200  // Find nearest face to sample
201  pointIndexHit info = oc.findNearest(sample, sqr(GREAT));
202 
203  if (info.index() == -1)
204  {
206  (
207  "treeDataFace::getSampleType"
208  "(indexedOctree<treeDataFace>&, const point&)"
209  ) << "Could not find " << sample << " in octree."
210  << abort(FatalError);
211  }
212 
213 
214  // Get actual intersection point on face
215  label faceI = faceLabels_[info.index()];
216 
217  if (debug & 2)
218  {
219  Pout<< "getSampleType : sample:" << sample
220  << " nearest face:" << faceI;
221  }
222 
223  const pointField& points = mesh_.points();
224 
225  // Retest to classify where on face info is. Note: could be improved. We
226  // already have point.
227 
228  const face& f = mesh_.faces()[faceI];
229  const vector& area = mesh_.faceAreas()[faceI];
230  const point& fc = mesh_.faceCentres()[faceI];
231 
232  pointHit curHit = f.nearestPoint(sample, points);
233  const point& curPt = curHit.rawPoint();
234 
235  //
236  // 1] Check whether sample is above face
237  //
238 
239  if (curHit.hit())
240  {
241  // Nearest point inside face. Compare to face normal.
242 
243  if (debug & 2)
244  {
245  Pout<< " -> face hit:" << curPt
246  << " comparing to face normal " << area << endl;
247  }
248  return indexedOctree<treeDataFace>::getSide(area, sample - curPt);
249  }
250 
251  if (debug & 2)
252  {
253  Pout<< " -> face miss:" << curPt;
254  }
255 
256  //
257  // 2] Check whether intersection is on one of the face vertices or
258  // face centre
259  //
260 
261  const scalar typDimSqr = mag(area) + VSMALL;
262 
263  forAll(f, fp)
264  {
265  if ((magSqr(points[f[fp]] - curPt)/typDimSqr) < tolSqr)
266  {
267  // Face intersection point equals face vertex fp
268 
269  // Calculate point normal (wrong: uses face normals instead of
270  // triangle normals)
271  const labelList& pFaces = mesh_.pointFaces()[f[fp]];
272 
273  vector pointNormal(vector::zero);
274 
275  forAll(pFaces, i)
276  {
277  if (isTreeFace_.get(pFaces[i]) == 1)
278  {
279  vector n = mesh_.faceAreas()[pFaces[i]];
280  n /= mag(n) + VSMALL;
281 
282  pointNormal += n;
283  }
284  }
285 
286  if (debug & 2)
287  {
288  Pout<< " -> face point hit :" << points[f[fp]]
289  << " point normal:" << pointNormal
290  << " distance:"
291  << magSqr(points[f[fp]] - curPt)/typDimSqr << endl;
292  }
294  (
295  pointNormal,
296  sample - curPt
297  );
298  }
299  }
300  if ((magSqr(fc - curPt)/typDimSqr) < tolSqr)
301  {
302  // Face intersection point equals face centre. Normal at face centre
303  // is already average of face normals
304 
305  if (debug & 2)
306  {
307  Pout<< " -> centre hit:" << fc
308  << " distance:" << magSqr(fc - curPt)/typDimSqr << endl;
309  }
310 
311  return indexedOctree<treeDataFace>::getSide(area, sample - curPt);
312  }
313 
314 
315 
316  //
317  // 3] Get the 'real' edge the face intersection is on
318  //
319 
320  const labelList& myEdges = mesh_.faceEdges()[faceI];
321 
322  forAll(myEdges, myEdgeI)
323  {
324  const edge& e = mesh_.edges()[myEdges[myEdgeI]];
325 
326  pointHit edgeHit =
328  (
329  points[e.start()],
330  points[e.end()]
331  ).nearestDist(sample);
332 
333 
334  if ((magSqr(edgeHit.rawPoint() - curPt)/typDimSqr) < tolSqr)
335  {
336  // Face intersection point lies on edge e
337 
338  // Calculate edge normal (wrong: uses face normals instead of
339  // triangle normals)
340  const labelList& eFaces = mesh_.edgeFaces()[myEdges[myEdgeI]];
341 
342  vector edgeNormal(vector::zero);
343 
344  forAll(eFaces, i)
345  {
346  if (isTreeFace_.get(eFaces[i]) == 1)
347  {
348  vector n = mesh_.faceAreas()[eFaces[i]];
349  n /= mag(n) + VSMALL;
350 
351  edgeNormal += n;
352  }
353  }
354 
355  if (debug & 2)
356  {
357  Pout<< " -> real edge hit point:" << edgeHit.rawPoint()
358  << " comparing to edge normal:" << edgeNormal
359  << endl;
360  }
361 
362  // Found face intersection point on this edge. Compare to edge
363  // normal
365  (
366  edgeNormal,
367  sample - curPt
368  );
369  }
370  }
371 
372 
373  //
374  // 4] Get the internal edge the face intersection is on
375  //
376 
377  forAll(f, fp)
378  {
380  (
381  points[f[fp]],
382  fc
383  ).nearestDist(sample);
384 
385  if ((magSqr(edgeHit.rawPoint() - curPt)/typDimSqr) < tolSqr)
386  {
387  // Face intersection point lies on edge between two face triangles
388 
389  // Calculate edge normal as average of the two triangle normals
390  vector e = points[f[fp]] - fc;
391  vector ePrev = points[f[f.rcIndex(fp)]] - fc;
392  vector eNext = points[f[f.fcIndex(fp)]] - fc;
393 
394  vector nLeft = ePrev ^ e;
395  nLeft /= mag(nLeft) + VSMALL;
396 
397  vector nRight = e ^ eNext;
398  nRight /= mag(nRight) + VSMALL;
399 
400  if (debug & 2)
401  {
402  Pout<< " -> internal edge hit point:" << edgeHit.rawPoint()
403  << " comparing to edge normal "
404  << 0.5*(nLeft + nRight)
405  << endl;
406  }
407 
408  // Found face intersection point on this edge. Compare to edge
409  // normal
411  (
412  0.5*(nLeft + nRight),
413  sample - curPt
414  );
415  }
416  }
417 
418  if (debug & 2)
419  {
420  Pout<< "Did not find sample " << sample
421  << " anywhere related to nearest face " << faceI << endl
422  << "Face:";
423 
424  forAll(f, fp)
425  {
426  Pout<< " vertex:" << f[fp] << " coord:" << points[f[fp]]
427  << endl;
428  }
429  }
430 
431  // Can't determine status of sample with respect to nearest face.
432  // Either
433  // - tolerances are wrong. (if e.g. face has zero area)
434  // - or (more likely) surface is not closed.
435 
436  return volumeType::UNKNOWN;
437 }
438 
439 
440 // Check if any point on shape is inside cubeBb.
442 (
443  const label index,
444  const treeBoundBox& cubeBb
445 ) const
446 {
447  // 1. Quick rejection: bb does not intersect face bb at all
448  if (cacheBb_)
449  {
450  if (!cubeBb.overlaps(bbs_[index]))
451  {
452  return false;
453  }
454  }
455  else
456  {
457  if (!cubeBb.overlaps(calcBb(faceLabels_[index])))
458  {
459  return false;
460  }
461  }
462 
463  const pointField& points = mesh_.points();
464 
465 
466  // 2. Check if one or more face points inside
467  label faceI = faceLabels_[index];
468 
469  const face& f = mesh_.faces()[faceI];
470  if (cubeBb.containsAny(points, f))
471  {
472  return true;
473  }
474 
475  // 3. Difficult case: all points are outside but connecting edges might
476  // go through cube. Use triangle-bounding box intersection.
477  const point& fc = mesh_.faceCentres()[faceI];
478 
479  forAll(f, fp)
480  {
481  bool triIntersects = triangleFuncs::intersectBb
482  (
483  points[f[fp]],
484  points[f[f.fcIndex(fp)]],
485  fc,
486  cubeBb
487  );
488 
489  if (triIntersects)
490  {
491  return true;
492  }
493  }
494  return false;
495 }
496 
497 
498 void Foam::treeDataFace::findNearestOp::operator()
499 (
500  const labelUList& indices,
501  const point& sample,
502 
503  scalar& nearestDistSqr,
504  label& minIndex,
505  point& nearestPoint
506 ) const
507 {
508  const treeDataFace& shape = tree_.shapes();
509 
510  forAll(indices, i)
511  {
512  const label index = indices[i];
513 
514  const face& f = shape.mesh().faces()[shape.faceLabels()[index]];
515 
516  pointHit nearHit = f.nearestPoint(sample, shape.mesh().points());
517  scalar distSqr = sqr(nearHit.distance());
518 
519  if (distSqr < nearestDistSqr)
520  {
521  nearestDistSqr = distSqr;
522  minIndex = index;
523  nearestPoint = nearHit.rawPoint();
524  }
525  }
526 }
527 
528 
529 void Foam::treeDataFace::findNearestOp::operator()
530 (
531  const labelUList& indices,
532  const linePointRef& ln,
533 
534  treeBoundBox& tightest,
535  label& minIndex,
536  point& linePoint,
537  point& nearestPoint
538 ) const
539 {
541  (
542  "treeDataFace::findNearestOp::operator()"
543  "("
544  " const labelUList&,"
545  " const linePointRef&,"
546  " treeBoundBox&,"
547  " label&,"
548  " point&,"
549  " point&"
550  ") const"
551  );
552 }
553 
554 
555 bool Foam::treeDataFace::findIntersectOp::operator()
556 (
557  const label index,
558  const point& start,
559  const point& end,
560  point& intersectionPoint
561 ) const
562 {
563  const treeDataFace& shape = tree_.shapes();
564 
565  // Do quick rejection test
566  if (shape.cacheBb_)
567  {
568  const treeBoundBox& faceBb = shape.bbs_[index];
569 
570  if ((faceBb.posBits(start) & faceBb.posBits(end)) != 0)
571  {
572  // start and end in same block outside of faceBb.
573  return false;
574  }
575  }
576 
577  const label faceI = shape.faceLabels_[index];
578 
579  const vector dir(end - start);
580 
581  pointHit inter = shape.mesh_.faces()[faceI].intersection
582  (
583  start,
584  dir,
585  shape.mesh_.faceCentres()[faceI],
586  shape.mesh_.points(),
588  );
589 
590  if (inter.hit() && inter.distance() <= 1)
591  {
592  // Note: no extra test on whether intersection is in front of us
593  // since using half_ray
594  intersectionPoint = inter.hitPoint();
595  return true;
596  }
597  else
598  {
599  return false;
600  }
601 }
602 
603 
604 // ************************************************************************* //
const Point & rawPoint() const
Return point with no checking.
Definition: PointHit.H:158
virtual const pointField & points() const =0
Return mesh points.
treeDataFace(const bool cacheBb, const primitiveMesh &, const labelUList &)
Construct from mesh and subset of faces.
Definition: treeDataFace.C:83
label nFaces() const
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
vector point
Point is a vector.
Definition: point.H:41
dimensioned< scalar > mag(const dimensioned< Type > &)
findIntersectOp(const indexedOctree< treeDataFace > &tree)
Definition: treeDataFace.C:158
dimensioned< scalar > magSqr(const dimensioned< Type > &)
pointField shapePoints() const
Get representative point cloud for all shapes inside.
Definition: treeDataFace.C:168
const Point & hitPoint() const
Return hit point.
Definition: PointHit.H:126
A simple container for copying or transferring objects of type <T>.
Definition: Xfer.H:85
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
bool containsAny(const UList< point > &) const
Contains any of the points? (inside or on edge)
Definition: boundBox.C:261
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:75
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: UListI.H:58
bool overlaps(const label index, const treeBoundBox &sampleBb) const
Does (bb of) shape at index overlap bb.
Definition: treeDataFace.C:442
virtual const faceList & faces() const =0
Return faces.
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:278
Encapsulation of data needed to search for faces.
Definition: treeDataFace.H:58
volumeType getVolumeType(const indexedOctree< treeDataFace > &, const point &) const
Get type (inside,outside,mixed,unknown) of point w.r.t. surface.
Definition: treeDataFace.C:184
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
label index() const
Return index.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Namespace for OpenFOAM.
A line primitive.
Definition: line.H:56
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
label n
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:74
const double e
Elementary charge.
Definition: doubleFloat.H:78
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
volScalarField & p
Definition: createFields.H:51
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: POSIX.C:855
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const labelList & faceLabels() const
Definition: treeDataFace.H:179
findNearestOp(const indexedOctree< treeDataFace > &tree)
Definition: treeDataFace.C:149
#define forAll(list, i)
Definition: UList.H:421
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
errorManip< error > abort(error &err)
Definition: errorManip.H:131
direction posBits(const point &) const
Position of point relative to bounding box.
Definition: treeBoundBox.C:476
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
Non-pointer based hierarchical recursive searching.
Definition: treeDataEdge.H:47
scalar distance() const
Return distance to hit.
Definition: PointHit.H:139
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:300
static volumeType getSide(const vector &outsideNormal, const vector &vec)
Helper function to return the side. Returns outside if.
static const Vector min
Definition: Vector.H:83
error FatalError
pointHit nearestPoint(const point &p, const pointField &) const
Return nearest point to face.
const polyMesh & mesh() const
Return the mesh reference.
label end() const
Return end vertex label.
Definition: edgeI.H:92
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
static const Vector zero
Definition: Vector.H:80
bool hit() const
Is there a hit.
Definition: PointHit.H:120
static bool intersectBb(const point &p0, const point &p1, const point &p2, const treeBoundBox &cubeBb)
Does triangle intersect bounding box.
label start() const
Return start vertex label.
Definition: edgeI.H:81
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
#define notImplemented(functionName)
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:356
bool overlaps(const boundBox &) const
Overlaps/touches boundingBox?
Definition: boundBoxI.H:120
defineTypeNameAndDebug(combustionModel, 0)
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
const primitiveMesh & mesh() const
Definition: treeDataFace.H:184
const vectorField & faceCentres() const