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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2015-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "treeDataPrimitivePatch.H"
30 #include "indexedOctree.H"
31 #include "triangle.H"
33 #include "triFace.H"
34 
35 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
36 
37 template<class PatchType>
40 {
41  const auto& points = pp.points();
42 
43  treeBoundBoxList bbs(pp.size());
44 
45  // Like std::transform with [&](const auto& f)
46  // which may not work with older C++ versions
47 
48  {
49  auto iter = bbs.begin();
50 
51  for (const auto& f : pp)
52  {
53  *iter = treeBoundBox(points, f);
54  ++iter;
55  }
56  }
57 
58  return bbs;
59 }
60 
61 
62 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
63 
64 template<class PatchType>
65 inline Foam::treeBoundBox
66 Foam::treeDataPrimitivePatch<PatchType>::getBounds(const label patchFacei) const
67 {
68  return treeBoundBox(patch_.points(), patch_[patchFacei]);
69 }
70 
71 
72 template<class PatchType>
74 {
75  if (cacheBb_)
76  {
77  bbs_ = treeDataPrimitivePatch<PatchType>::boxes(patch_);
78  }
79 }
80 
81 
82 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
83 
84 template<class PatchType>
86 (
87  const bool cacheBb,
88  const PatchType& patch,
89  const scalar planarTol
90 )
91 :
92  patch_(patch),
93  cacheBb_(cacheBb),
94  planarTol_(planarTol)
95 {
96  update();
97 }
98 
99 
100 template<class PatchType>
102 (
104 )
105 :
106  tree_(tree)
107 {}
108 
109 
110 template<class PatchType>
112 (
114 )
115 :
116  tree_(tree)
117 {}
118 
119 
120 template<class PatchType>
122 (
124  DynamicList<label>& shapeMask
125 )
126 :
127  tree_(tree),
128  shapeMask_(shapeMask)
129 {}
130 
131 
132 template<class PatchType>
135 (
137  const label edgeID
138 )
139 :
140  tree_(tree),
141  edgeID_(edgeID)
142 {}
144 
145 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
146 
147 template<class PatchType>
150 (
151  const labelUList& indices
152 ) const
153 {
154  treeBoundBox bb;
155 
156  for (const label patchFacei : indices)
157  {
158  bb.add(patch_.points(), patch_[patchFacei]);
159  }
161  return bb;
162 }
163 
164 
165 template<class PatchType>
167 (
169  const point& sample
170 ) const
171 {
172  // Need to determine whether sample is 'inside' or 'outside'
173  // Done by finding nearest face. This gives back a face which is
174  // guaranteed to contain nearest point. This point can be
175  // - in interior of face: compare to face normal
176  // - on edge of face: compare to edge normal
177  // - on point of face: compare to point normal
178  // Unfortunately the octree does not give us back the intersection point
179  // or where on the face it has hit so we have to recreate all that
180  // information.
181 
182 
183  // Find nearest face to sample
184  pointIndexHit info = oc.findNearest(sample, sqr(GREAT));
185 
186  if (info.index() == -1)
187  {
189  << "Could not find " << sample << " in octree."
190  << abort(FatalError);
191  }
192 
193  // Get actual intersection point on face
194  label facei = info.index();
195 
196  if (debug & 2)
197  {
198  Pout<< "getSampleType : sample:" << sample
199  << " nearest face:" << facei;
200  }
201 
202  const auto& localF = patch_.localFaces()[facei];
203  const auto& f = patch_[facei];
204  const auto& points = patch_.points();
205 
206  // Retest to classify where on face info is. Note: could be improved. We
207  // already have point.
208 
209  pointHit curHit = f.nearestPoint(sample, points);
210  const vector area = f.areaNormal(points);
211  const point& curPt = curHit.point();
212 
213  //
214  // 1] Check whether sample is above face
215  //
216 
217  if (curHit.hit())
218  {
219  // Nearest point inside face. Compare to face normal.
220 
221  if (debug & 2)
222  {
223  Pout<< " -> face hit:" << curPt
224  << " comparing to face normal " << area << endl;
225  }
227  (
228  area,
229  sample - curPt
230  );
231  }
232 
233  if (debug & 2)
234  {
235  Pout<< " -> face miss:" << curPt;
236  }
237 
238  //
239  // 2] Check whether intersection is on one of the face vertices or
240  // face centre
241  //
242 
243  const scalar typDimSqr = mag(area) + VSMALL;
244 
245 
246  forAll(f, fp)
247  {
248  const scalar relDistSqr = (magSqr(points[f[fp]] - curPt)/typDimSqr);
249 
250  if (relDistSqr < planarTol_)
251  {
252  // Face intersection point equals face vertex fp
253 
254  // Calculate point normal (wrong: uses face normals instead of
255  // triangle normals)
256 
258  (
259  patch_.pointNormals()[localF[fp]],
260  sample - curPt
261  );
262  }
263  }
264 
265  const point fc(f.centre(points));
266 
267  const scalar relDistSqr = (magSqr(fc - curPt)/typDimSqr);
268 
269  if (relDistSqr < planarTol_)
270  {
271  // Face intersection point equals face centre. Normal at face centre
272  // is already average of face normals
273 
274  if (debug & 2)
275  {
276  Pout<< " -> centre hit:" << fc
277  << " distance:" << relDistSqr << endl;
278  }
279 
281  (
282  area,
283  sample - curPt
284  );
285  }
286 
287 
288  //
289  // 3] Get the 'real' edge the face intersection is on
290  //
291 
292  for (const label edgei : patch_.faceEdges()[facei])
293  {
294  pointHit edgeHit =
295  patch_.meshEdge(edgei).line(points).nearestDist(sample);
296 
297  const scalar relDistSqr = (edgeHit.point().distSqr(curPt)/typDimSqr);
298 
299  if (relDistSqr < planarTol_)
300  {
301  // Face intersection point lies on edge e
302 
303  // Calculate edge normal (wrong: uses face normals instead of
304  // triangle normals)
305  vector edgeNormal(Zero);
306 
307  for (const label eFacei : patch_.edgeFaces()[edgei])
308  {
309  edgeNormal += patch_.faceNormals()[eFacei];
310  }
311 
312  if (debug & 2)
313  {
314  Pout<< " -> real edge hit point:" << edgeHit.point()
315  << " comparing to edge normal:" << edgeNormal
316  << endl;
317  }
318 
319  // Found face intersection point on this edge. Compare to edge
320  // normal
322  (
323  edgeNormal,
324  sample - curPt
325  );
326  }
327  }
328 
329 
330  //
331  // 4] Get the internal edge the face intersection is on
332  //
333 
334  forAll(f, fp)
335  {
336  pointHit edgeHit = linePointRef(points[f[fp]], fc).nearestDist(sample);
337 
338  const scalar relDistSqr = (edgeHit.point().distSqr(curPt)/typDimSqr);
339 
340  if (relDistSqr < planarTol_)
341  {
342  // Face intersection point lies on edge between two face triangles
343 
344  // Calculate edge normal as average of the two triangle normals
345  vector e = points[f[fp]] - fc;
346  vector ePrev = points[f[f.rcIndex(fp)]] - fc;
347  vector eNext = points[f[f.fcIndex(fp)]] - fc;
348 
349  vector nLeft = normalised(ePrev ^ e);
350  vector nRight = normalised(e ^ eNext);
351 
352  if (debug & 2)
353  {
354  Pout<< " -> internal edge hit point:" << edgeHit.point()
355  << " comparing to edge normal "
356  << 0.5*(nLeft + nRight)
357  << endl;
358  }
359 
360  // Found face intersection point on this edge. Compare to edge
361  // normal
363  (
364  0.5*(nLeft + nRight),
365  sample - curPt
366  );
367  }
368  }
369 
370  if (debug & 2)
371  {
372  Pout<< "Did not find sample " << sample
373  << " anywhere related to nearest face " << facei << endl
374  << "Face:";
375 
376  forAll(f, fp)
377  {
378  Pout<< " vertex:" << f[fp]
379  << " coord:" << points[f[fp]]
380  << endl;
381  }
382  }
383 
384  // Can't determine status of sample with respect to nearest face.
385  // Either
386  // - tolerances are wrong. (if e.g. face has zero area)
387  // - or (more likely) surface is not closed.
388 
390 }
391 
392 
393 // Check if any point on shape is inside searchBox.
394 template<class PatchType>
396 (
397  const label index,
398  const treeBoundBox& searchBox
399 ) const
400 {
401  // 1. Quick rejection: bb does not intersect face bb at all
402  if
403  (
404  cacheBb_
405  ? !searchBox.overlaps(bbs_[index])
406  : !searchBox.overlaps(getBounds(index))
407  )
408  {
409  return false;
410  }
411 
412 
413  // 2. Check if one or more face points inside
414 
415  const auto& points = patch_.points();
416  const auto& f = patch_[index];
417 
418  if (f.size() == 3)
419  {
420  const triPointRef tri(points[f[0]], points[f[1]], points[f[2]]);
421 
422  return searchBox.intersects(tri);
423  }
424 
425  if (searchBox.containsAny(points, f))
426  {
427  return true;
428  }
429 
430  // 3. Difficult case: all points are outside but connecting edges might
431  // go through cube. Use triangle-bounding box intersection.
432 
433  const point fc = f.centre(points);
434 
435  forAll(f, fp)
436  {
437  const triPointRef tri
438  (
439  points[f.thisLabel(fp)], points[f.nextLabel(fp)], fc
440  );
441 
442  if (searchBox.intersects(tri))
443  {
444  return true;
445  }
446  }
447 
448  return false;
449 }
450 
451 
452 // Check if any point on shape is inside sphere.
453 template<class PatchType>
455 (
456  const label index,
457  const point& centre,
458  const scalar radiusSqr
459 ) const
460 {
461  // 1. Quick rejection: sphere does not intersect face bb at all
462  if
463  (
464  cacheBb_
465  ? !bbs_[index].overlaps(centre, radiusSqr)
466  : !getBounds(index).overlaps(centre, radiusSqr)
467  )
468  {
469  return false;
470  }
471 
472  const auto& points = patch_.points();
473  const auto& f = patch_[index];
474 
475  pointHit nearHit = f.nearestPoint(centre, points);
476 
477  // If the distance to the nearest point on the face from the sphere centres
478  // is within the radius, then the sphere touches the face.
479  if (sqr(nearHit.distance()) < radiusSqr)
480  {
481  return true;
482  }
483 
484  return false;
485 }
486 
487 
488 // * * * * * * * * * * * * * * * * Searching * * * * * * * * * * * * * * * * //
489 
490 template<class PatchType>
492 (
493  const labelUList& indices,
494  const point& sample,
495 
496  scalar& nearestDistSqr,
497  label& minIndex,
498  point& nearestPoint
499 ) const
500 {
501  const auto& points = patch_.points();
502 
503  for (const label index : indices)
504  {
505  const auto& f = patch_[index];
506 
507  const pointHit nearHit = f.nearestPoint(sample, points);
508  const scalar distSqr = sqr(nearHit.distance());
509 
510  if (distSqr < nearestDistSqr)
511  {
512  nearestDistSqr = distSqr;
513  minIndex = index;
514  nearestPoint = nearHit.point();
515  }
516  }
517 }
518 
519 
520 template<class PatchType>
522 (
523  const labelUList& indices,
524  const point& sample,
525 
526  scalar& nearestDistSqr,
527  label& minIndex,
528  point& nearestPoint
529 ) const
530 {
531  tree_.shapes().findNearest
532  (
533  indices,
534  sample,
535  nearestDistSqr,
536  minIndex,
537  nearestPoint
538  );
539 }
540 
541 
542 template<class PatchType>
544 (
545  const labelUList& indices,
546  const linePointRef& ln,
547 
548  treeBoundBox& tightest,
549  label& minIndex,
550  point& linePoint,
551  point& nearestPoint
552 ) const
553 {
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  return findIntersection(tree_, index, start, end, intersectionPoint);
568 }
569 
570 
571 template<class PatchType>
573 (
574  const label index,
575  const point& start,
576  const point& end,
577  point& intersectionPoint
578 ) const
579 {
580  if (shapeMask_.found(index))
581  {
582  return false;
583  }
585  return findIntersection(tree_, index, start, end, intersectionPoint);
586 }
587 
588 
589 template<class PatchType>
591 (
592  const label index,
593  const point& start,
594  const point& end,
595  point& intersectionPoint
596 ) const
597 {
598  if (edgeID_ == -1)
599  {
601  << "EdgeID not set. Please set edgeID to the index of"
602  << " the edge you are testing"
603  << exit(FatalError);
604  }
605 
606  const treeDataPrimitivePatch<PatchType>& shape = tree_.shapes();
607  const PatchType& patch = shape.patch();
608 
609  const auto& f = patch.localFaces()[index];
610  const edge& e = patch.edges()[edgeID_];
611 
612  if (!f.found(e[0]) && !f.found(e[1]))
613  {
614  return findIntersection(tree_, index, start, end, intersectionPoint);
615  }
617  return false;
618 }
619 
620 
621 template<class PatchType>
623 (
625  const label index,
626  const point& start,
627  const point& end,
628  point& intersectionPoint
629 )
630 {
631  const treeDataPrimitivePatch<PatchType>& shape = tree.shapes();
632  const PatchType& patch = shape.patch();
633 
634  const auto& points = patch.points();
635  const auto& f = patch[index];
636 
637  // Do quick rejection test
638  if (shape.cacheBb_)
639  {
640  const treeBoundBox& faceBb = shape.bbs_[index];
641 
642  if ((faceBb.posBits(start) & faceBb.posBits(end)) != 0)
643  {
644  // start and end in same block outside of faceBb.
645  return false;
646  }
647  }
648 
649  const vector dir(end - start);
650  pointHit inter;
651 
652  if (f.size() == 3)
653  {
654  inter = triPointRef
655  (
656  points[f[0]],
657  points[f[1]],
658  points[f[2]]
659  ).intersection(start, dir, intersection::HALF_RAY, shape.planarTol_);
660  }
661  else
662  {
663  const pointField& faceCentres = patch.faceCentres();
664 
665  inter = f.intersection
666  (
667  start,
668  dir,
669  faceCentres[index],
670  points,
672  shape.planarTol_
673  );
674  }
675 
676  if (inter.hit() && inter.distance() <= 1)
677  {
678  // Note: no extra test on whether intersection is in front of us
679  // since using half_ray
680  intersectionPoint = inter.point();
681 
682  return true;
683  }
684 
685  return false;
686 }
687 
688 
689 // ************************************************************************* //
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.
bool overlaps(const label index, const treeBoundBox &searchBox) const
Does shape at index overlap searchBox.
A line primitive.
Definition: line.H:52
const point & centre(const label index) const
Representative point (face centre) at shape index.
const PatchType & patch() const noexcept
The underlying patch.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
volumeType getVolumeType(const indexedOctree< treeDataPrimitivePatch< PatchType >> &, const point &) const
Get type (inside,outside,mixed,unknown) of point w.r.t. surface.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
bool overlaps(const boundBox &bb) const
Overlaps with other bounding box, sphere etc?
Definition: boundBoxI.H:439
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void findNearest(const labelUList &indices, const point &sample, scalar &nearestDistSqr, label &nearestIndex, point &nearestPoint) const
Calculates nearest (to sample) point in shape.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
This class describes the interaction of an object (often a face) and a point. It carries the info of ...
Definition: pointIndexHit.H:44
static volumeType getSide(const vector &outsideNormal, const vector &vec)
Helper function to return the side. Returns outside if.
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:674
An enumeration wrapper for classification of a location as being inside/outside of a volume...
Definition: volumeType.H:55
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
Unknown state.
Definition: volumeType.H:64
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
direction posBits(const point &pt) const
Position of point relative to bounding box.
Definition: treeBoundBox.C:486
void add(const boundBox &bb)
Extend to include the second box.
Definition: boundBoxI.H:323
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
bool intersects(const point &overallStart, const vector &overallVec, const point &start, const point &end, point &pt, direction &ptBits) const
Intersects segment; set point to intersection position and face,.
Definition: treeBoundBox.C:255
line< point, const point & > linePointRef
A line using referred points.
Definition: line.H:66
findNearestOp(const indexedOctree< treeDataPrimitivePatch > &tree)
const pointField & points
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition: edge.H:59
const wordList area
Standard area field types (scalar, vector, tensor, etc)
Encapsulation of data needed to search on PrimitivePatches.
treeBoundBox bounds(const labelUList &indices) const
Return bounding box for the specified face indices.
Tree tree(triangles.begin(), triangles.end())
Vector< scalar > vector
Definition: vector.H:57
Describes the interaction of a object and a (templated) point. It carries the info of a successful hi...
Definition: pointHit.H:43
label index() const noexcept
Return the hit index.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition: UListI.H:391
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: POSIX.C:1237
findSelfIntersectOp(const indexedOctree< treeDataPrimitivePatch > &tree, const label edgeID)
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
int debug
Static debugging option.
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:201
mesh update()
labelList f(nPoints)
bool containsAny(const UList< point > &points) const
Contains any of the points? (inside or on edge)
Definition: boundBox.C:385
static treeBoundBoxList boxes(const PatchType &patch)
Calculate and return bounding boxes for each patch face.
vector point
Point is a vector.
Definition: point.H:37
Non-pointer based hierarchical recursive searching.
findAllIntersectOp(const indexedOctree< treeDataPrimitivePatch > &tree, DynamicList< label > &shapeMask)
scalar distance() const noexcept
Return distance to hit.
Definition: pointHit.H:169
const std::string patch
OpenFOAM patch number as a std::string.
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition: triangleFwd.H:39
const point_type & point() const noexcept
Return the point, no checks.
Definition: pointHit.H:161
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:90
treeDataPrimitivePatch(const bool cacheBb, const PatchType &patch, const scalar planarTol)
Construct from patch.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:686
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
PointHit< point > pointHit
A PointHit with a 3D point.
Definition: pointHit.H:43
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
findIntersectOp(const indexedOctree< treeDataPrimitivePatch > &tree)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127