treeDataPrimitivePatch.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 "treeDataPrimitivePatch.H"
27 #include "indexedOctree.H"
28 #include "triangleFuncs.H"
29 #include "triSurfaceTools.H"
30 #include "triFace.H"
31 
32 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
34 template<class PatchType>
36 (
37  const pointField& points,
38  const face& f
39 )
40 {
41  treeBoundBox bb(points[f[0]], points[f[0]]);
42 
43  for (label fp = 1; fp < f.size(); fp++)
44  {
45  const point& p = points[f[fp]];
46 
47  bb.min() = min(bb.min(), p);
48  bb.max() = max(bb.max(), p);
49  }
50  return bb;
51 }
52 
53 
54 template<class PatchType>
56 {
57  if (cacheBb_)
58  {
59  bbs_.setSize(patch_.size());
60 
61  forAll(patch_, i)
62  {
63  bbs_[i] = calcBb(patch_.points(), patch_[i]);
64  }
65  }
66 }
67 
68 
69 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
70 
71 // Construct from components
72 template<class PatchType>
74 (
75  const bool cacheBb,
76  const PatchType& patch,
77  const scalar planarTol
78 )
79 :
80  patch_(patch),
81  cacheBb_(cacheBb),
82  planarTol_(planarTol)
83 {
84  update();
85 }
86 
87 
88 template<class PatchType>
90 (
92 )
93 :
94  tree_(tree)
95 {}
96 
97 
98 template<class PatchType>
100 (
102 )
103 :
104  tree_(tree)
105 {}
106 
107 
108 template<class PatchType>
110 (
112  DynamicList<label>& shapeMask
113 )
114 :
115  tree_(tree),
116  shapeMask_(shapeMask)
117 {}
118 
119 
120 template<class PatchType>
123 (
125  const label edgeID
126 )
127 :
128  tree_(tree),
129  edgeID_(edgeID)
130 {}
131 
132 
133 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
134 
135 template<class PatchType>
137 {
138  pointField cc(patch_.size());
139 
140  forAll(patch_, i)
141  {
142  cc[i] = patch_[i].centre(patch_.points());
143  }
144 
145  return cc;
146 }
147 
148 
149 //- Get type (inside,outside,mixed,unknown) of point w.r.t. surface.
150 // Only makes sense for closed surfaces.
151 template<class PatchType>
153 (
155  const point& sample
156 ) const
157 {
158  // Need to determine whether sample is 'inside' or 'outside'
159  // Done by finding nearest face. This gives back a face which is
160  // guaranteed to contain nearest point. This point can be
161  // - in interior of face: compare to face normal
162  // - on edge of face: compare to edge normal
163  // - on point of face: compare to point normal
164  // Unfortunately the octree does not give us back the intersection point
165  // or where on the face it has hit so we have to recreate all that
166  // information.
167 
168 
169  // Find nearest face to sample
170  pointIndexHit info = oc.findNearest(sample, sqr(GREAT));
171 
172  if (info.index() == -1)
173  {
175  (
176  "treeDataPrimitivePatch::getSampleType"
177  "(indexedOctree<treeDataPrimitivePatch>&, const point&)"
178  ) << "Could not find " << sample << " in octree."
179  << abort(FatalError);
180  }
181 
182  // Get actual intersection point on face
183  label faceI = info.index();
184 
185  if (debug & 2)
186  {
187  Pout<< "getSampleType : sample:" << sample
188  << " nearest face:" << faceI;
189  }
190 
191  const pointField& points = patch_.localPoints();
192  const typename PatchType::FaceType& f = patch_.localFaces()[faceI];
193 
194  // Retest to classify where on face info is. Note: could be improved. We
195  // already have point.
196 
197  pointHit curHit = f.nearestPoint(sample, points);
198  const vector area = f.normal(points);
199  const point& curPt = curHit.rawPoint();
200 
201  //
202  // 1] Check whether sample is above face
203  //
204 
205  if (curHit.hit())
206  {
207  // Nearest point inside face. Compare to face normal.
208 
209  if (debug & 2)
210  {
211  Pout<< " -> face hit:" << curPt
212  << " comparing to face normal " << area << endl;
213  }
215  (
216  area,
217  sample - curPt
218  );
219  }
220 
221  if (debug & 2)
222  {
223  Pout<< " -> face miss:" << curPt;
224  }
225 
226  //
227  // 2] Check whether intersection is on one of the face vertices or
228  // face centre
229  //
230 
231  const scalar typDimSqr = mag(area) + VSMALL;
232 
233 
234  forAll(f, fp)
235  {
236  if ((magSqr(points[f[fp]] - curPt)/typDimSqr) < planarTol_)
237  {
238  // Face intersection point equals face vertex fp
239 
240  // Calculate point normal (wrong: uses face normals instead of
241  // triangle normals)
242 
244  (
245  patch_.pointNormals()[f[fp]],
246  sample - curPt
247  );
248  }
249  }
250 
251  const point fc(f.centre(points));
252 
253  if ((magSqr(fc - curPt)/typDimSqr) < planarTol_)
254  {
255  // Face intersection point equals face centre. Normal at face centre
256  // is already average of face normals
257 
258  if (debug & 2)
259  {
260  Pout<< " -> centre hit:" << fc
261  << " distance:" << magSqr(fc - curPt)/typDimSqr << endl;
262  }
263 
265  (
266  area,
267  sample - curPt
268  );
269  }
270 
271 
272 
273  //
274  // 3] Get the 'real' edge the face intersection is on
275  //
276 
277  const labelList& fEdges = patch_.faceEdges()[faceI];
278 
279  forAll(fEdges, fEdgeI)
280  {
281  label edgeI = fEdges[fEdgeI];
282  const edge& e = patch_.edges()[edgeI];
283 
284  pointHit edgeHit = e.line(points).nearestDist(sample);
285 
286  if ((magSqr(edgeHit.rawPoint() - curPt)/typDimSqr) < planarTol_)
287  {
288  // Face intersection point lies on edge e
289 
290  // Calculate edge normal (wrong: uses face normals instead of
291  // triangle normals)
292  const labelList& eFaces = patch_.edgeFaces()[edgeI];
293 
294  vector edgeNormal(vector::zero);
295 
296  forAll(eFaces, i)
297  {
298  edgeNormal += patch_.faceNormals()[eFaces[i]];
299  }
300 
301  if (debug & 2)
302  {
303  Pout<< " -> real edge hit point:" << edgeHit.rawPoint()
304  << " comparing to edge normal:" << edgeNormal
305  << endl;
306  }
307 
308  // Found face intersection point on this edge. Compare to edge
309  // normal
311  (
312  edgeNormal,
313  sample - curPt
314  );
315  }
316  }
317 
318 
319  //
320  // 4] Get the internal edge the face intersection is on
321  //
322 
323  forAll(f, fp)
324  {
325  pointHit edgeHit = linePointRef
326  (
327  points[f[fp]],
328  fc
329  ).nearestDist(sample);
330 
331  if ((magSqr(edgeHit.rawPoint() - curPt)/typDimSqr) < planarTol_)
332  {
333  // Face intersection point lies on edge between two face triangles
334 
335  // Calculate edge normal as average of the two triangle normals
336  vector e = points[f[fp]] - fc;
337  vector ePrev = points[f[f.rcIndex(fp)]] - fc;
338  vector eNext = points[f[f.fcIndex(fp)]] - fc;
339 
340  vector nLeft = ePrev ^ e;
341  nLeft /= mag(nLeft) + VSMALL;
342 
343  vector nRight = e ^ eNext;
344  nRight /= mag(nRight) + VSMALL;
345 
346  if (debug & 2)
347  {
348  Pout<< " -> internal edge hit point:" << edgeHit.rawPoint()
349  << " comparing to edge normal "
350  << 0.5*(nLeft + nRight)
351  << endl;
352  }
353 
354  // Found face intersection point on this edge. Compare to edge
355  // normal
357  (
358  0.5*(nLeft + nRight),
359  sample - curPt
360  );
361  }
362  }
363 
364  if (debug & 2)
365  {
366  Pout<< "Did not find sample " << sample
367  << " anywhere related to nearest face " << faceI << endl
368  << "Face:";
369 
370  forAll(f, fp)
371  {
372  Pout<< " vertex:" << f[fp] << " coord:" << points[f[fp]]
373  << endl;
374  }
375  }
376 
377  // Can't determine status of sample with respect to nearest face.
378  // Either
379  // - tolerances are wrong. (if e.g. face has zero area)
380  // - or (more likely) surface is not closed.
381 
382  return volumeType::UNKNOWN;
383 }
384 
385 
386 // Check if any point on shape is inside cubeBb.
387 template<class PatchType>
389 (
390  const label index,
391  const treeBoundBox& cubeBb
392 ) const
393 {
394  // 1. Quick rejection: bb does not intersect face bb at all
395  if (cacheBb_)
396  {
397  if (!cubeBb.overlaps(bbs_[index]))
398  {
399  return false;
400  }
401  }
402  else
403  {
404  if (!cubeBb.overlaps(calcBb(patch_.points(), patch_[index])))
405  {
406  return false;
407  }
408  }
409 
410 
411  // 2. Check if one or more face points inside
412 
413  const pointField& points = patch_.points();
414  const typename PatchType::FaceType& f = patch_[index];
415 
416  if (cubeBb.containsAny(points, f))
417  {
418  return true;
419  }
420 
421  // 3. Difficult case: all points are outside but connecting edges might
422  // go through cube. Use triangle-bounding box intersection.
423  const point fc = f.centre(points);
424 
425  if (f.size() == 3)
426  {
427  return triangleFuncs::intersectBb
428  (
429  points[f[0]],
430  points[f[1]],
431  points[f[2]],
432  cubeBb
433  );
434  }
435  else
436  {
437  forAll(f, fp)
438  {
439  bool triIntersects = triangleFuncs::intersectBb
440  (
441  points[f[fp]],
442  points[f[f.fcIndex(fp)]],
443  fc,
444  cubeBb
445  );
446 
447  if (triIntersects)
448  {
449  return true;
450  }
451  }
452  }
453 
454  return false;
455 }
456 
457 
458 // Check if any point on shape is inside sphere.
459 template<class PatchType>
461 (
462  const label index,
463  const point& centre,
464  const scalar radiusSqr
465 ) const
466 {
467  // 1. Quick rejection: sphere does not intersect face bb at all
468  if (cacheBb_)
469  {
470  if (!bbs_[index].overlaps(centre, radiusSqr))
471  {
472  return false;
473  }
474  }
475  else
476  {
477  if (!calcBb(patch_.points(), patch_[index]).overlaps(centre, radiusSqr))
478  {
479  return false;
480  }
481  }
482 
483  const pointField& points = patch_.points();
484  const face& f = patch_[index];
485 
486  pointHit nearHit = f.nearestPoint(centre, points);
487 
488  // If the distance to the nearest point on the face from the sphere centres
489  // is within the radius, then the sphere touches the face.
490  if (sqr(nearHit.distance()) < radiusSqr)
491  {
492  return true;
493  }
494 
495  return false;
496 }
497 
498 
499 template<class PatchType>
501 (
502  const labelUList& indices,
503  const point& sample,
504 
505  scalar& nearestDistSqr,
506  label& minIndex,
507  point& nearestPoint
508 ) const
509 {
510  const treeDataPrimitivePatch<PatchType>& shape = tree_.shapes();
511  const PatchType& patch = shape.patch();
512 
513  const pointField& points = patch.points();
514 
515  forAll(indices, i)
516  {
517  const label index = indices[i];
518  const typename PatchType::FaceType& f = patch[index];
519 
520  pointHit nearHit = f.nearestPoint(sample, points);
521  scalar distSqr = sqr(nearHit.distance());
522 
523  if (distSqr < nearestDistSqr)
524  {
525  nearestDistSqr = distSqr;
526  minIndex = index;
527  nearestPoint = nearHit.rawPoint();
528  }
529  }
530 }
531 
532 
533 template<class PatchType>
535 (
536  const labelUList& indices,
537  const linePointRef& ln,
538 
539  treeBoundBox& tightest,
540  label& minIndex,
541  point& linePoint,
542  point& nearestPoint
543 ) const
544 {
546  (
547  "treeDataPrimitivePatch<PatchType>::findNearestOp::operator()"
548  "("
549  " const labelUList&,"
550  " const linePointRef&,"
551  " treeBoundBox&,"
552  " label&,"
553  " point&,"
554  " point&"
555  ") const"
556  );
557 }
558 
559 
560 template<class PatchType>
562 (
563  const label index,
564  const point& start,
565  const point& end,
566  point& intersectionPoint
567 ) const
568 {
569  return findIntersection(tree_, index, start, end, intersectionPoint);
570 }
571 
572 
573 template<class PatchType>
575 (
576  const label index,
577  const point& start,
578  const point& end,
579  point& intersectionPoint
580 ) const
581 {
582  if (!shapeMask_.empty() && findIndex(shapeMask_, index) != -1)
583  {
584  return false;
585  }
586 
587  return findIntersection(tree_, index, start, end, intersectionPoint);
588 }
589 
590 
591 template<class PatchType>
593 (
594  const label index,
595  const point& start,
596  const point& end,
597  point& intersectionPoint
598 ) const
599 {
600  if (edgeID_ == -1)
601  {
603  (
604  "findSelfIntersectOp::operator()\n"
605  "(\n"
606  " const label index,\n"
607  " const point& start,\n"
608  " const point& end,\n"
609  " point& intersectionPoint\n"
610  ") const"
611  ) << "EdgeID not set. Please set edgeID to the index of"
612  << " the edge you are testing"
613  << exit(FatalError);
614  }
615 
616  const treeDataPrimitivePatch<PatchType>& shape = tree_.shapes();
617  const PatchType& patch = shape.patch();
618 
619  const typename PatchType::FaceType& f = patch.localFaces()[index];
620  const edge& e = patch.edges()[edgeID_];
621 
622  if (findIndex(f, e[0]) == -1 && findIndex(f, e[1]) == -1)
623  {
624  return findIntersection(tree_, index, start, end, intersectionPoint);
625  }
626  else
627  {
628  return false;
629  }
630 }
631 
632 
633 template<class PatchType>
635 (
637  const label index,
638  const point& start,
639  const point& end,
640  point& intersectionPoint
641 )
642 {
643  const treeDataPrimitivePatch<PatchType>& shape = tree.shapes();
644  const PatchType& patch = shape.patch();
645 
646  const pointField& points = patch.points();
647  const typename PatchType::FaceType& f = patch[index];
648 
649  // Do quick rejection test
650  if (shape.cacheBb_)
651  {
652  const treeBoundBox& faceBb = shape.bbs_[index];
653 
654  if ((faceBb.posBits(start) & faceBb.posBits(end)) != 0)
655  {
656  // start and end in same block outside of faceBb.
657  return false;
658  }
659  }
660 
661  const vector dir(end - start);
662  pointHit inter;
663 
664  if (f.size() == 3)
665  {
666  inter = triPointRef
667  (
668  points[f[0]],
669  points[f[1]],
670  points[f[2]]
671  ).intersection(start, dir, intersection::HALF_RAY, shape.planarTol_);
672  }
673  else
674  {
675  const pointField& faceCentres = patch.faceCentres();
676 
677  inter = f.intersection
678  (
679  start,
680  dir,
681  faceCentres[index],
682  points,
683  intersection::HALF_RAY,
684  shape.planarTol_
685  );
686  }
687 
688  if (inter.hit() && inter.distance() <= 1)
689  {
690  // Note: no extra test on whether intersection is in front of us
691  // since using half_ray
692  intersectionPoint = inter.hitPoint();
693 
694  return true;
695  }
696  else
697  {
698  return false;
699  }
700 }
701 
702 
703 // ************************************************************************* //
PointHit< Point > nearestDist(const Point &p) const
Return nearest distance to line from a given point.
Definition: lineI.H:95
const Point & rawPoint() const
Return point with no checking.
Definition: PointHit.H:158
findAllIntersectOp(const indexedOctree< treeDataPrimitivePatch > &tree, DynamicList< label > &shapeMask)
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
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
dimensioned< scalar > mag(const dimensioned< Type > &)
const PatchType & patch() const
Return access to the underlying patch.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
const Point & hitPoint() const
Return hit point.
Definition: PointHit.H:126
findNearestOp(const indexedOctree< treeDataPrimitivePatch > &tree)
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
A line primitive.
Definition: line.H:56
triangle< point, const point & > triPointRef
Definition: triPointRef.H:44
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
bool overlaps(const label index, const treeBoundBox &sampleBb) const
Does shape at index overlap bb.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
volScalarField & p
Definition: createFields.H:51
findSelfIntersectOp(const indexedOctree< treeDataPrimitivePatch > &tree, const label edgeID)
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: POSIX.C:855
linePointRef line(const pointField &) const
Return edge line.
Definition: edgeI.H:169
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:78
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define forAll(list, i)
Definition: UList.H:421
static bool findIntersection(const indexedOctree< treeDataPrimitivePatch< PatchType > > &tree, const label index, const point &start, const point &end, point &intersectionPoint)
Helper: find intersection of line with shapes.
findIntersectOp(const indexedOctree< treeDataPrimitivePatch > &tree)
Encapsulation of data needed to search on PrimitivePatches.
volumeType getVolumeType(const indexedOctree< treeDataPrimitivePatch< PatchType > > &, const point &) const
Get type (inside,outside,mixed,unknown) of point w.r.t. surface.
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
const Vector< Cmpt > & centre(const Foam::List< Vector< Cmpt > > &) const
Return *this (used for point which is a typedef to Vector<scalar>.
Definition: VectorI.H:106
treeDataPrimitivePatch(const bool cacheBb, const PatchType &, const scalar planarTol)
Construct from patch.
#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
static const Vector min
Definition: Vector.H:83
error FatalError
pointHit nearestPoint(const point &p, const pointField &) const
Return nearest point to face.
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
bool hit() const
Is there a hit.
Definition: PointHit.H:120
pointField shapePoints() const
Get representative point cloud for all shapes inside.
line< point, const point & > linePointRef
Line using referred points.
Definition: linePointRef.H:45
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
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53