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-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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  << "Could not find " << sample << " in octree."
176  << abort(FatalError);
177  }
178 
179  // Get actual intersection point on face
180  label faceI = info.index();
181 
182  if (debug & 2)
183  {
184  Pout<< "getSampleType : sample:" << sample
185  << " nearest face:" << faceI;
186  }
187 
188  const typename PatchType::FaceType& localF = patch_.localFaces()[faceI];
189  const typename PatchType::FaceType& f = patch_[faceI];
190  const pointField& points = patch_.points();
191  const labelList& mp = patch_.meshPoints();
192 
193  // Retest to classify where on face info is. Note: could be improved. We
194  // already have point.
195 
196  pointHit curHit = f.nearestPoint(sample, points);
197  const vector area = f.normal(points);
198  const point& curPt = curHit.rawPoint();
199 
200  //
201  // 1] Check whether sample is above face
202  //
203 
204  if (curHit.hit())
205  {
206  // Nearest point inside face. Compare to face normal.
207 
208  if (debug & 2)
209  {
210  Pout<< " -> face hit:" << curPt
211  << " comparing to face normal " << area << endl;
212  }
214  (
215  area,
216  sample - curPt
217  );
218  }
219 
220  if (debug & 2)
221  {
222  Pout<< " -> face miss:" << curPt;
223  }
224 
225  //
226  // 2] Check whether intersection is on one of the face vertices or
227  // face centre
228  //
229 
230  const scalar typDimSqr = mag(area) + VSMALL;
231 
232 
233  forAll(f, fp)
234  {
235  if ((magSqr(points[f[fp]] - curPt)/typDimSqr) < planarTol_)
236  {
237  // Face intersection point equals face vertex fp
238 
239  // Calculate point normal (wrong: uses face normals instead of
240  // triangle normals)
241 
243  (
244  patch_.pointNormals()[localF[fp]],
245  sample - curPt
246  );
247  }
248  }
249 
250  const point fc(f.centre(points));
251 
252  if ((magSqr(fc - curPt)/typDimSqr) < planarTol_)
253  {
254  // Face intersection point equals face centre. Normal at face centre
255  // is already average of face normals
256 
257  if (debug & 2)
258  {
259  Pout<< " -> centre hit:" << fc
260  << " distance:" << magSqr(fc - curPt)/typDimSqr << endl;
261  }
262 
264  (
265  area,
266  sample - curPt
267  );
268  }
269 
270 
271 
272  //
273  // 3] Get the 'real' edge the face intersection is on
274  //
275 
276  const labelList& fEdges = patch_.faceEdges()[faceI];
277 
278  forAll(fEdges, fEdgeI)
279  {
280  label edgeI = fEdges[fEdgeI];
281  const edge& e = patch_.edges()[edgeI];
282  const linePointRef ln(points[mp[e.start()]], points[mp[e.end()]]);
283  pointHit edgeHit = ln.nearestDist(sample);
284 
285  if ((magSqr(edgeHit.rawPoint() - curPt)/typDimSqr) < planarTol_)
286  {
287  // Face intersection point lies on edge e
288 
289  // Calculate edge normal (wrong: uses face normals instead of
290  // triangle normals)
291  const labelList& eFaces = patch_.edgeFaces()[edgeI];
292 
293  vector edgeNormal(vector::zero);
294 
295  forAll(eFaces, i)
296  {
297  edgeNormal += patch_.faceNormals()[eFaces[i]];
298  }
299 
300  if (debug & 2)
301  {
302  Pout<< " -> real edge hit point:" << edgeHit.rawPoint()
303  << " comparing to edge normal:" << edgeNormal
304  << endl;
305  }
306 
307  // Found face intersection point on this edge. Compare to edge
308  // normal
310  (
311  edgeNormal,
312  sample - curPt
313  );
314  }
315  }
316 
317 
318  //
319  // 4] Get the internal edge the face intersection is on
320  //
321 
322  forAll(f, fp)
323  {
324  pointHit edgeHit = linePointRef(points[f[fp]], fc).nearestDist(sample);
325 
326  if ((magSqr(edgeHit.rawPoint() - curPt)/typDimSqr) < planarTol_)
327  {
328  // Face intersection point lies on edge between two face triangles
329 
330  // Calculate edge normal as average of the two triangle normals
331  vector e = points[f[fp]] - fc;
332  vector ePrev = points[f[f.rcIndex(fp)]] - fc;
333  vector eNext = points[f[f.fcIndex(fp)]] - fc;
334 
335  vector nLeft = ePrev ^ e;
336  nLeft /= mag(nLeft) + VSMALL;
337 
338  vector nRight = e ^ eNext;
339  nRight /= mag(nRight) + VSMALL;
340 
341  if (debug & 2)
342  {
343  Pout<< " -> internal edge hit point:" << edgeHit.rawPoint()
344  << " comparing to edge normal "
345  << 0.5*(nLeft + nRight)
346  << endl;
347  }
348 
349  // Found face intersection point on this edge. Compare to edge
350  // normal
352  (
353  0.5*(nLeft + nRight),
354  sample - curPt
355  );
356  }
357  }
358 
359  if (debug & 2)
360  {
361  Pout<< "Did not find sample " << sample
362  << " anywhere related to nearest face " << faceI << endl
363  << "Face:";
364 
365  forAll(f, fp)
366  {
367  Pout<< " vertex:" << f[fp]
368  << " coord:" << points[f[fp]]
369  << endl;
370  }
371  }
372 
373  // Can't determine status of sample with respect to nearest face.
374  // Either
375  // - tolerances are wrong. (if e.g. face has zero area)
376  // - or (more likely) surface is not closed.
377 
378  return volumeType::UNKNOWN;
379 }
380 
381 
382 // Check if any point on shape is inside cubeBb.
383 template<class PatchType>
385 (
386  const label index,
387  const treeBoundBox& cubeBb
388 ) const
389 {
390  // 1. Quick rejection: bb does not intersect face bb at all
391  if (cacheBb_)
392  {
393  if (!cubeBb.overlaps(bbs_[index]))
394  {
395  return false;
396  }
397  }
398  else
399  {
400  if (!cubeBb.overlaps(calcBb(patch_.points(), patch_[index])))
401  {
402  return false;
403  }
404  }
405 
406 
407  // 2. Check if one or more face points inside
408 
409  const pointField& points = patch_.points();
410  const typename PatchType::FaceType& f = patch_[index];
411 
412  if (cubeBb.containsAny(points, f))
413  {
414  return true;
415  }
416 
417  // 3. Difficult case: all points are outside but connecting edges might
418  // go through cube. Use triangle-bounding box intersection.
419  const point fc = f.centre(points);
420 
421  if (f.size() == 3)
422  {
423  return triangleFuncs::intersectBb
424  (
425  points[f[0]],
426  points[f[1]],
427  points[f[2]],
428  cubeBb
429  );
430  }
431  else
432  {
433  forAll(f, fp)
434  {
435  bool triIntersects = triangleFuncs::intersectBb
436  (
437  points[f[fp]],
438  points[f[f.fcIndex(fp)]],
439  fc,
440  cubeBb
441  );
442 
443  if (triIntersects)
444  {
445  return true;
446  }
447  }
448  }
449 
450  return false;
451 }
452 
453 
454 // Check if any point on shape is inside sphere.
455 template<class PatchType>
457 (
458  const label index,
459  const point& centre,
460  const scalar radiusSqr
461 ) const
462 {
463  // 1. Quick rejection: sphere does not intersect face bb at all
464  if (cacheBb_)
465  {
466  if (!bbs_[index].overlaps(centre, radiusSqr))
467  {
468  return false;
469  }
470  }
471  else
472  {
473  if (!calcBb(patch_.points(), patch_[index]).overlaps(centre, radiusSqr))
474  {
475  return false;
476  }
477  }
478 
479  const pointField& points = patch_.points();
480  const face& f = patch_[index];
481 
482  pointHit nearHit = f.nearestPoint(centre, points);
483 
484  // If the distance to the nearest point on the face from the sphere centres
485  // is within the radius, then the sphere touches the face.
486  if (sqr(nearHit.distance()) < radiusSqr)
487  {
488  return true;
489  }
490 
491  return false;
492 }
493 
494 
495 template<class PatchType>
497 (
498  const labelUList& indices,
499  const point& sample,
500 
501  scalar& nearestDistSqr,
502  label& minIndex,
503  point& nearestPoint
504 ) const
505 {
506  const treeDataPrimitivePatch<PatchType>& shape = tree_.shapes();
507  const PatchType& patch = shape.patch();
508 
509  const pointField& points = patch.points();
510 
511  forAll(indices, i)
512  {
513  const label index = indices[i];
514  const typename PatchType::FaceType& f = patch[index];
515 
516  pointHit nearHit = f.nearestPoint(sample, 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 template<class PatchType>
531 (
532  const labelUList& indices,
533  const linePointRef& ln,
534 
535  treeBoundBox& tightest,
536  label& minIndex,
537  point& linePoint,
538  point& nearestPoint
539 ) const
540 {
542 }
543 
544 
545 template<class PatchType>
547 (
548  const label index,
549  const point& start,
550  const point& end,
551  point& intersectionPoint
552 ) const
553 {
554  return findIntersection(tree_, index, start, end, intersectionPoint);
555 }
556 
557 
558 template<class PatchType>
560 (
561  const label index,
562  const point& start,
563  const point& end,
564  point& intersectionPoint
565 ) const
566 {
567  if (!shapeMask_.empty() && findIndex(shapeMask_, index) != -1)
568  {
569  return false;
570  }
571 
572  return findIntersection(tree_, index, start, end, intersectionPoint);
573 }
574 
575 
576 template<class PatchType>
578 (
579  const label index,
580  const point& start,
581  const point& end,
582  point& intersectionPoint
583 ) const
584 {
585  if (edgeID_ == -1)
586  {
588  << "EdgeID not set. Please set edgeID to the index of"
589  << " the edge you are testing"
590  << exit(FatalError);
591  }
592 
593  const treeDataPrimitivePatch<PatchType>& shape = tree_.shapes();
594  const PatchType& patch = shape.patch();
595 
596  const typename PatchType::FaceType& f = patch.localFaces()[index];
597  const edge& e = patch.edges()[edgeID_];
598 
599  if (findIndex(f, e[0]) == -1 && findIndex(f, e[1]) == -1)
600  {
601  return findIntersection(tree_, index, start, end, intersectionPoint);
602  }
603  else
604  {
605  return false;
606  }
607 }
608 
609 
610 template<class PatchType>
612 (
614  const label index,
615  const point& start,
616  const point& end,
617  point& intersectionPoint
618 )
619 {
620  const treeDataPrimitivePatch<PatchType>& shape = tree.shapes();
621  const PatchType& patch = shape.patch();
622 
623  const pointField& points = patch.points();
624  const typename PatchType::FaceType& f = patch[index];
625 
626  // Do quick rejection test
627  if (shape.cacheBb_)
628  {
629  const treeBoundBox& faceBb = shape.bbs_[index];
630 
631  if ((faceBb.posBits(start) & faceBb.posBits(end)) != 0)
632  {
633  // start and end in same block outside of faceBb.
634  return false;
635  }
636  }
637 
638  const vector dir(end - start);
639  pointHit inter;
640 
641  if (f.size() == 3)
642  {
643  inter = triPointRef
644  (
645  points[f[0]],
646  points[f[1]],
647  points[f[2]]
648  ).intersection(start, dir, intersection::HALF_RAY, shape.planarTol_);
649  }
650  else
651  {
652  const pointField& faceCentres = patch.faceCentres();
653 
654  inter = f.intersection
655  (
656  start,
657  dir,
658  faceCentres[index],
659  points,
660  intersection::HALF_RAY,
661  shape.planarTol_
662  );
663  }
664 
665  if (inter.hit() && inter.distance() <= 1)
666  {
667  // Note: no extra test on whether intersection is in front of us
668  // since using half_ray
669  intersectionPoint = inter.hitPoint();
670 
671  return true;
672  }
673  else
674  {
675  return false;
676  }
677 }
678 
679 
680 // ************************************************************************* //
Foam::treeDataPrimitivePatch::calcBb
static treeBoundBox calcBb(const pointField &, const face &)
Calculate face bounding box.
Definition: treeDataPrimitivePatch.C:36
Foam::constant::atomic::mp
const dimensionedScalar mp
Proton mass.
Foam::PointIndexHit::index
label index() const
Return index.
Definition: PointIndexHit.H:111
Foam::treeDataPrimitivePatch::findIntersectOp::findIntersectOp
findIntersectOp(const indexedOctree< treeDataPrimitivePatch > &tree)
Definition: treeDataPrimitivePatch.C:100
Foam::PointHit::hit
bool hit() const
Is there a hit.
Definition: PointHit.H:120
p
p
Definition: pEqn.H:62
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
Foam::treeDataPrimitivePatch::patch
const PatchType & patch() const
Return access to the underlying patch.
Definition: treeDataPrimitivePatch.H:222
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::PointHit
This class describes the interaction of a face and a point. It carries the info of a successful hit a...
Definition: PointHit.H:51
Foam::findIndex
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
Foam::DynamicList< label >
Foam::PointHit::rawPoint
const Point & rawPoint() const
Return point with no checking.
Definition: PointHit.H:158
Foam::treeBoundBox
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:75
indexedOctree.H
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Foam::treeDataPrimitivePatch::overlaps
bool overlaps(const label index, const treeBoundBox &sampleBb) const
Does shape at index overlap bb.
Definition: treeDataPrimitivePatch.C:385
Foam::treeBoundBox::posBits
direction posBits(const point &) const
Position of point relative to bounding box.
Definition: treeBoundBox.C:469
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
triFace.H
Foam::treeDataPrimitivePatch::update
void update()
Initialise all member data.
Definition: treeDataPrimitivePatch.C:55
Foam::treeDataPrimitivePatch::getVolumeType
volumeType getVolumeType(const indexedOctree< treeDataPrimitivePatch< PatchType > > &, const point &) const
Get type (inside,outside,mixed,unknown) of point w.r.t. surface.
Definition: treeDataPrimitivePatch.C:153
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:365
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
Foam::treeDataPrimitivePatch::cacheBb_
const bool cacheBb_
Whether to precalculate and store face bounding box.
Definition: treeDataPrimitivePatch.H:71
Foam::treeDataPrimitivePatch::findAllIntersectOp::findAllIntersectOp
findAllIntersectOp(const indexedOctree< treeDataPrimitivePatch > &tree, DynamicList< label > &shapeMask)
Definition: treeDataPrimitivePatch.C:110
Foam::PointHit::distance
scalar distance() const
Return distance to hit.
Definition: PointHit.H:139
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::volumeType
Definition: volumeType.H:54
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::indexedOctree
Non-pointer based hierarchical recursive searching.
Definition: treeDataTriSurface.H:47
treeDataPrimitivePatch.H
Foam::boundBox::containsAny
bool containsAny(const UList< point > &) const
Contains any of the points? (inside or on edge)
Definition: boundBox.C:261
Foam::treeDataPrimitivePatch::shapePoints
pointField shapePoints() const
Get representative point cloud for all shapes inside.
Definition: treeDataPrimitivePatch.C:136
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
Foam::FatalError
error FatalError
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
triangleFuncs.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::treeDataPrimitivePatch::treeDataPrimitivePatch
treeDataPrimitivePatch(const bool cacheBb, const PatchType &, const scalar planarTol)
Construct from patch.
Definition: treeDataPrimitivePatch.C:74
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::treeDataPrimitivePatch::findNearestOp::findNearestOp
findNearestOp(const indexedOctree< treeDataPrimitivePatch > &tree)
Definition: treeDataPrimitivePatch.C:90
Foam::linePointRef
line< point, const point & > linePointRef
Line using referred points.
Definition: linePointRef.H:45
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::treeDataPrimitivePatch::findIntersection
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.
Definition: treeDataPrimitivePatch.C:612
Foam::treeDataPrimitivePatch
Encapsulation of data needed to search on PrimitivePatches.
Definition: treeDataPrimitivePatch.H:61
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::treeDataPrimitivePatch::bbs_
treeBoundBoxList bbs_
Face bounding boxes (valid only if cacheBb_)
Definition: treeDataPrimitivePatch.H:77
Foam::UList
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
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:94
Foam::treeDataPrimitivePatch::findSelfIntersectOp::findSelfIntersectOp
findSelfIntersectOp(const indexedOctree< treeDataPrimitivePatch > &tree, const label edgeID)
Definition: treeDataPrimitivePatch.C:123
Foam::line
A line primitive.
Definition: line.H:56
Foam::ln
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: POSIX.C:854
Foam::treeDataPrimitivePatch::planarTol_
const scalar planarTol_
Tolerance to use for intersection tests.
Definition: treeDataPrimitivePatch.H:74
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::treeBoundBox::overlaps
bool overlaps(const boundBox &) const
Overlaps other bounding box?
Definition: boundBoxI.H:120
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
triSurfaceTools.H
Foam::PointHit::hitPoint
const Point & hitPoint() const
Return hit point.
Definition: PointHit.H:126
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Foam::triPointRef
triangle< point, const point & > triPointRef
Definition: triangle.H:75