1 // Copyright (c) 2013 GeometryFactory (France).
2 // All rights reserved.
3 //
4 // This file is part of CGAL (www.cgal.org).
5 //
6 // $URL: https://github.com/CGAL/cgal/blob/v5.3/Polygon_mesh_processing/include/CGAL/Polygon_mesh_processing/internal/Side_of_triangle_mesh/Ray_3_Triangle_3_traversal_traits.h $
7 // $Id: Ray_3_Triangle_3_traversal_traits.h d64faf3 2020-07-01T21:03:55+02:00 Sébastien Loriot
8 // SPDX-License-Identifier: GPL-3.0-or-later OR LicenseRef-Commercial
9 //
10 //
11 // Author(s)     : Sebastien Loriot
12 
13 
14 #ifndef CGAL_POINT_INSIDE_POLYHEDRON_RAY_3_TRIANGLE_3_TRAVERSAL_TRAITS_H
15 #define CGAL_POINT_INSIDE_POLYHEDRON_RAY_3_TRIANGLE_3_TRAVERSAL_TRAITS_H
16 
17 #include <CGAL/license/Polygon_mesh_processing/miscellaneous.h>
18 
19 
20 #include <boost/logic/tribool.hpp>
21 #include <CGAL/tags.h>
22 #include <CGAL/Bbox_3.h>
23 #include <CGAL/Intersections_3/Ray_3_Triangle_3.h>
24 #include <CGAL/internal/AABB_tree/Primitive_helper.h>
25 #include <CGAL/internal/AABB_tree/AABB_node.h>
26 
27 namespace CGAL {
28 namespace internal {
29 
30 template<typename AABBTraits, class Kernel, class Helper, class Tag_ray_is_vertical=Tag_false>
31 class Ray_3_Triangle_3_traversal_traits
32 {
33 protected:
34   //the status indicates whether the query point is strictly inside the polyhedron, and the number of intersected triangles if yes
35   std::pair<boost::logic::tribool,std::size_t>& m_status;
36   bool m_stop;
37   const AABBTraits& m_aabb_traits;
38   typedef typename AABBTraits::Primitive Primitive;
39   typedef CGAL::AABB_node<AABBTraits> Node;
40   Helper m_helper;
41 
42 public:
Ray_3_Triangle_3_traversal_traits(std::pair<boost::logic::tribool,std::size_t> & status,const AABBTraits & aabb_traits,const Helper & h)43   Ray_3_Triangle_3_traversal_traits(std::pair<boost::logic::tribool,std::size_t>& status,
44                                     const AABBTraits& aabb_traits,
45                                     const Helper& h)
46     :m_status(status),m_stop(false), m_aabb_traits(aabb_traits), m_helper(h)
47   {m_status.first=true;}
48 
go_further()49   bool go_further() const { return !m_stop; }
50 
51   template<class Query>
intersection(const Query & query,const Primitive & primitive)52   void intersection(const Query& query, const Primitive& primitive)
53   {
54     Intersections::internal::r3t3_do_intersect_endpoint_position_visitor visitor;
55     std::pair<bool,Intersections::internal::R3T3_intersection::type> res=
56       Intersections::internal::do_intersect(m_helper.get_primitive_datum(primitive, m_aabb_traits), query,Kernel(),visitor);
57 
58     if (res.first){
59       switch (res.second){
60         case Intersections::internal::R3T3_intersection::CROSS_FACET:
61           ++m_status.second;
62         break;
63         case Intersections::internal::R3T3_intersection::ENDPOINT_IN_TRIANGLE:
64           m_status.first=false;
65           m_stop=true;
66         break;
67         default:
68           m_status.first=boost::logic::indeterminate;
69           m_stop=true;
70       }
71     }
72   }
73 
74   template<class Query>
do_intersect(const Query & query,const Node & node)75   bool do_intersect(const Query& query, const Node& node) const
76   {
77     return m_aabb_traits.do_intersect_object()(query, m_helper.get_node_bbox(node));
78   }
79 };
80 
81 
82 //specialization for vertical ray
83 template<typename AABBTraits, class Kernel, class TraversalTraits>
84 class Ray_3_Triangle_3_traversal_traits<AABBTraits,Kernel,TraversalTraits,Tag_true>:
85   public Ray_3_Triangle_3_traversal_traits<AABBTraits,Kernel,TraversalTraits,Tag_false>
86 {
87   typedef Ray_3_Triangle_3_traversal_traits<AABBTraits,Kernel,TraversalTraits,Tag_false> Base;
88   typedef typename Kernel::Point_3 Point;
89   typedef typename Base::Primitive Primitive;
90   typedef CGAL::AABB_node<AABBTraits> Node;
91 
92 public:
Ray_3_Triangle_3_traversal_traits(std::pair<boost::logic::tribool,std::size_t> & status,const AABBTraits & aabb_traits,const TraversalTraits & tt)93   Ray_3_Triangle_3_traversal_traits(std::pair<boost::logic::tribool,std::size_t>& status, const AABBTraits& aabb_traits, const TraversalTraits& tt)
94     :Base(status, aabb_traits, tt){}
95 
96   template <class Query>
do_intersect(const Query & query,const Bbox_3 & bbox)97   bool do_intersect(const Query& query, const Bbox_3& bbox) const
98   {
99     const Point& source=query.point(0);
100     const Point& target=query.point(1);
101 
102     bool inc_z=target.z()>source.z();
103 
104     //the ray does not intersect the z-slab
105     if ( ( inc_z && source.z()>bbox.zmax() )|| (!inc_z && source.z()<bbox.zmin()) ) return false;
106 
107     //the source is not in the x-slab
108     if (source.x() > bbox.xmax() || source.x()<bbox.xmin()) return false;
109     //check if the source is not in the y-slab
110     return source.y() <= bbox.ymax() && source.y()>=bbox.ymin();
111   }
112 
113   template <class Query>
do_intersect(const Query & query,const Node & node)114   bool do_intersect(const Query& query, const Node& node) const
115   {
116     return do_intersect(query,this->m_helper.get_node_bbox(node));
117   }
118 
119 private:
x_project(const typename Kernel::Point_3 & p)120   typename Kernel::Point_2 x_project(const typename Kernel::Point_3& p) const{
121     return typename Kernel::Point_2(p.y(),p.z());
122   }
y_project(const typename Kernel::Point_3 & p)123   typename Kernel::Point_2 y_project(const typename Kernel::Point_3& p) const{
124     return typename Kernel::Point_2(p.x(),p.z());
125   }
z_project(const typename Kernel::Point_3 & p)126   typename Kernel::Point_2 z_project(const typename Kernel::Point_3& p) const{
127     return typename Kernel::Point_2(p.x(),p.y());
128   }
129 public:
130   template<class Query>
intersection(const Query & query,const Primitive & primitive)131   void intersection(const Query& query, const Primitive& primitive)
132   {
133     typename Kernel::Triangle_3 t = this->m_helper.get_primitive_datum(primitive, this->m_aabb_traits);
134     if ( !do_intersect(query,t.bbox()) ) return;
135 
136     typename Kernel::Point_2 p0=z_project(t[0]);
137     typename Kernel::Point_2 p1=z_project(t[1]);
138     typename Kernel::Point_2 p2=z_project(t[2]);
139     int indices[3]={0,1,2}; //to track whether triangle points have been swapt
140     typename Kernel::Point_2 q=z_project( query.source() );
141 
142     Orientation orient_2=orientation(p0,p1,p2);
143 
144     //check whether the face has a normal vector in the xy-plane
145     if (orient_2==COLLINEAR){
146       //in that case the projection of the triangle along the z-axis is a segment.
147       const typename Kernel::Point_2& other_point = p0!=p1?p1:p2;
148       if ( orientation(p0,other_point,q) != COLLINEAR ) return;///no intersection
149 
150       //check if the ray source is above or below the triangle and compare it
151       //with the direction of the ray
152       //TODO and if yes return
153       //this is just an optimisation, the current code is valid
154 
155       this->m_status.first=boost::logic::indeterminate;
156       this->m_stop=true;
157       return;
158     }
159 
160     //regular case
161     if (orient_2==NEGATIVE){
162       std::swap(p1,p2);
163       std::swap(indices[1],indices[2]);
164     }
165 
166     //check whether the ray intersect the supporting plane
167     Orientation orient_3 = orientation(t[indices[0]],t[indices[1]],t[indices[2]],query.source());
168     if ( orient_3!=COPLANAR &&
169           (
170             //indicates whether the ray is oriented toward the positive side of the plane
171             ( POSITIVE == CGAL::sign( query.to_vector().z() )  )
172               ==
173             //indicates whether the source of the ray is in the positive side of the plane
174             (orient_3==POSITIVE)
175           )
176     ) return; //no intersection
177 
178     //position against first segment
179     switch( orientation(p0,p1,q) ){
180       case COLLINEAR:
181         this->m_status.first=boost::logic::indeterminate;
182         this->m_stop=true;
183       case NEGATIVE:
184         return;
185       default:
186       {}
187     }
188     //position against second segment
189     switch( orientation(p1,p2,q) ){
190       case COLLINEAR:
191         this->m_status.first=boost::logic::indeterminate;
192         this->m_stop=true;
193       case NEGATIVE:
194         return;
195       default:
196       {}
197     }
198     //position against third segment
199     switch( orientation(p2,p0,q) ){
200       case COLLINEAR:
201         this->m_status.first=boost::logic::indeterminate;
202         this->m_stop=true;
203       case NEGATIVE:
204         return;
205       default:
206       {}
207     }
208 
209     if (orient_3==COPLANAR){
210       //the endpoint is inside the triangle
211       this->m_status.first=false;
212       this->m_stop=true;
213     }
214     else
215       ++(this->m_status.second);
216   }
217 };
218 
219 //special case when ray query is from another Kernel K1 is the kernel compatible with the AABB-tree
220 template<typename AABBTraits, class K1, class K2, class Helper>
221 class K2_Ray_3_K1_Triangle_3_traversal_traits
222 {
223   //the status indicates whether the query point is strictly inside the polyhedron, and the number of intersected triangles if yes
224   std::pair<boost::logic::tribool,std::size_t>& m_status;
225   bool m_stop;
226   const AABBTraits& m_aabb_traits;
227   typedef typename AABBTraits::Primitive Primitive;
228   typedef CGAL::AABB_node<AABBTraits> Node;
229   Helper m_helper;
230   CGAL::Cartesian_converter<K1,K2> to_K2;
231 
232 public:
K2_Ray_3_K1_Triangle_3_traversal_traits(std::pair<boost::logic::tribool,std::size_t> & status,const AABBTraits & aabb_traits,const Helper & h)233   K2_Ray_3_K1_Triangle_3_traversal_traits(std::pair<boost::logic::tribool,std::size_t>& status,
234                                           const AABBTraits& aabb_traits,
235                                           const Helper& h)
236     :m_status(status), m_stop(false), m_aabb_traits(aabb_traits), m_helper(h)
237   {m_status.first=true;}
238 
go_further()239   bool go_further() const { return !m_stop; }
240 
241   template<class Query>
intersection(const Query & query,const Primitive & primitive)242   void intersection(const Query& query, const Primitive& primitive)
243   {
244     Intersections::internal::r3t3_do_intersect_endpoint_position_visitor visitor;
245     std::pair<bool,Intersections::internal::R3T3_intersection::type> res=
246       Intersections::internal::do_intersect(to_K2(m_helper.get_primitive_datum(primitive, m_aabb_traits)),
247                                             query, K2(), visitor);
248 
249     if (res.first){
250       switch (res.second){
251         case Intersections::internal::R3T3_intersection::CROSS_FACET:
252           ++m_status.second;
253         break;
254         case Intersections::internal::R3T3_intersection::ENDPOINT_IN_TRIANGLE:
255           m_status.first=false;
256           m_stop=true;
257         break;
258         default:
259           m_status.first=boost::logic::indeterminate;
260           m_stop=true;
261       }
262     }
263   }
264 
265   template<class Query>
do_intersect(const Query & query,const Node & node)266   bool do_intersect(const Query& query, const Node& node) const
267   {
268     return CGAL::do_intersect(query, m_helper.get_node_bbox(node));
269   }
270 };
271 
272 
273 }// namespace internal
274 }// namespace CGAL
275 
276 #endif
277