Factorize code of Do_intersect(Bbox_3, <something>)

Uniform use of do_intersect_bbox_segment_aux(..) with various Boolean
template arguments.
This commit is contained in:
Laurent Rineau 2012-06-15 15:06:34 +00:00
parent 1c6ba1851b
commit 41ba29e19a
3 changed files with 38 additions and 205 deletions

View File

@ -28,8 +28,10 @@
#include <CGAL/Profile_counter.h> #include <CGAL/Profile_counter.h>
#include <CGAL/internal/Static_filters/Static_filter_error.h> #include <CGAL/internal/Static_filters/Static_filter_error.h>
#include <CGAL/internal/Static_filters/tools.h> #include <CGAL/internal/Static_filters/tools.h>
#include <CGAL/internal/Intersections_3/Bbox_3_Segment_3_do_intersect.h> #include <CGAL/internal/Intersections_3/Bbox_3_Segment_3_do_intersect.h>
#include <cmath> // for CGAL::internal::do_intersect_bbox_segment_aux
#include <iostream> #include <iostream>
// inspired from http://cag.csail.mit.edu/~amy/papers/box-jgt.pdf // inspired from http://cag.csail.mit.edu/~amy/papers/box-jgt.pdf
@ -123,7 +125,11 @@ public:
CGAL_BRANCH_PROFILER_BRANCH_1(tmp); CGAL_BRANCH_PROFILER_BRANCH_1(tmp);
const Uncertain<result_type> ub = const Uncertain<result_type> ub =
do_intersect_bbox_segment_aux<double, true, true, true> do_intersect_bbox_segment_aux
<double,
true, // bounded at t=0
true, // bounded at t=0
true> // do use static filters
(px, py, pz, (px, py, pz,
qx, qy, qz, qx, qy, qz,
b); b);
@ -165,7 +171,11 @@ public:
CGAL_BRANCH_PROFILER_BRANCH_1(tmp); CGAL_BRANCH_PROFILER_BRANCH_1(tmp);
const Uncertain<result_type> ub = const Uncertain<result_type> ub =
do_intersect_bbox_segment_aux<double, true, false, true> do_intersect_bbox_segment_aux
<double,
true, // bounded at t=0
false,// not bounded at t=0
true> // do use static filters
(px, py, pz, (px, py, pz,
qx, qy, qz, qx, qy, qz,
b); b);

View File

@ -26,107 +26,15 @@
#include <CGAL/Line_3.h> #include <CGAL/Line_3.h>
#include <CGAL/Bbox_3.h> #include <CGAL/Bbox_3.h>
#include <CGAL/internal/Intersections_3/Bbox_3_Segment_3_do_intersect.h>
// for CGAL::internal::do_intersect_bbox_segment_aux
// inspired from http://cag.csail.mit.edu/~amy/papers/box-jgt.pdf // inspired from http://cag.csail.mit.edu/~amy/papers/box-jgt.pdf
namespace CGAL { namespace CGAL {
namespace internal { namespace internal {
template <typename FT>
inline
bool
bbox_line_do_intersect_aux(const FT& px, const FT& py, const FT& pz,
const FT& vx, const FT& vy, const FT& vz,
const FT& bxmin, const FT& bymin, const FT& bzmin,
const FT& bxmax, const FT& bymax, const FT& bzmax)
{
// -----------------------------------
// treat x coord
// -----------------------------------
FT dmin, tmin, tmax;
if ( vx >= 0 )
{
tmin = bxmin - px;
tmax = bxmax - px;
dmin = vx;
}
else
{
tmin = px - bxmax;
tmax = px - bxmin;
dmin = -vx;
}
//if px is not in the x-slab
if ( dmin == FT(0) && (tmin > FT(0) || tmax < FT(0)) ) return false;
FT dmax = dmin;
// -----------------------------------
// treat y coord
// -----------------------------------
FT d_, tmin_, tmax_;
if ( vy >= 0 )
{
tmin_ = bymin - py;
tmax_ = bymax - py;
d_ = vy;
}
else
{
tmin_ = py - bymax;
tmax_ = py - bymin;
d_ = -vy;
}
if ( d_ == FT(0) ){
//if py is not in the y-slab
if( (tmin_ > FT(0) || tmax_ < FT(0)) ) return false;
}
else
if ( (dmin*tmax_) < (d_*tmin) || (dmax*tmin_) > (d_*tmax) )
return false;
if( (dmin*tmin_) > (d_*tmin) )
{
tmin = tmin_;
dmin = d_;
}
if( (dmax*tmax_) < (d_*tmax) )
{
tmax = tmax_;
dmax = d_;
}
// -----------------------------------
// treat z coord
// -----------------------------------
if ( vz >= 0 )
{
tmin_ = bzmin - pz;
tmax_ = bzmax - pz;
d_ = vz;
}
else
{
tmin_ = pz - bzmax;
tmax_ = pz - bzmin;
d_ = -vz;
}
//if pz is not in the z-slab
//if ( d_ == FT(0) && (tmin_ > FT(0) || tmax_ < FT(0)) ) return false;
//The previous line is not needed as either dmin or d_ are not 0
//(otherwise the direction of the line would be null).
// The following is equivalent to the in z-slab test if d_=0.
return ( (dmin*tmax_) >= (d_*tmin) && (dmax*tmin_) <= (d_*tmax) );
}
template <class K> template <class K>
bool do_intersect(const typename K::Line_3& line, bool do_intersect(const typename K::Line_3& line,
const CGAL::Bbox_3& bbox, const CGAL::Bbox_3& bbox,
@ -139,11 +47,16 @@ namespace internal {
const Point_3& point = line.point(); const Point_3& point = line.point();
const Vector_3& v = line.to_vector(); const Vector_3& v = line.to_vector();
return bbox_line_do_intersect_aux( return do_intersect_bbox_segment_aux
<FT,
false, // not bounded at t=0
false, // not bounded at t=1
false> // do not use static filters
(
point.x(), point.y(), point.z(), point.x(), point.y(), point.z(),
v.x(), v.y(), v.z(), v.x(), v.y(), v.z(),
FT(bbox.xmin()), FT(bbox.ymin()), FT(bbox.zmin()), bbox
FT(bbox.xmax()), FT(bbox.ymax()), FT(bbox.zmax()) ); );
} }
} // namespace internal } // namespace internal

View File

@ -24,7 +24,9 @@
#include <CGAL/Ray_3.h> #include <CGAL/Ray_3.h>
#include <CGAL/Bbox_3.h> #include <CGAL/Bbox_3.h>
#include <CGAL/internal/Intersections_3/Bbox_3_Segment_3_do_intersect.h> #include <CGAL/internal/Intersections_3/Bbox_3_Segment_3_do_intersect.h>
// for CGAL::internal::do_intersect_bbox_segment_aux
// inspired from http://cag.csail.mit.edu/~amy/papers/box-jgt.pdf // inspired from http://cag.csail.mit.edu/~amy/papers/box-jgt.pdf
@ -32,99 +34,6 @@ namespace CGAL {
namespace internal { namespace internal {
template <typename FT>
inline
bool
bbox_ray_do_intersect_aux(const FT& px, const FT& py, const FT& pz,
const FT& qx, const FT& qy, const FT& qz,
const FT& bxmin, const FT& bymin, const FT& bzmin,
const FT& bxmax, const FT& bymax, const FT& bzmax)
{
// (px, py, pz) is the source
// -----------------------------------
// treat x coord
// -----------------------------------
FT dmin, tmin, tmax;
if ( qx >= px )
{
tmin = bxmin - px;
tmax = bxmax - px;
dmin = qx - px;
if ( tmax < FT(0) )
return false;
}
else
{
tmin = px - bxmax;
tmax = px - bxmin;
dmin = px - qx;
if ( tmax < FT(0) )
return false;
}
FT dmax = dmin;
if ( tmin > FT(0) )
{
if ( dmin == FT(0) )
return false;
}
else
{
tmin = FT(0);
dmin = FT(1);
}
// -----------------------------------
// treat y coord
// -----------------------------------
FT d_, tmin_, tmax_;
if ( qy >= py )
{
tmin_ = bymin - py;
tmax_ = bymax - py;
d_ = qy - py;
}
else
{
tmin_ = py - bymax;
tmax_ = py - bymin;
d_ = py - qy;
}
if ( (dmin*tmax_) < (d_*tmin) || (dmax*tmin_) > (d_*tmax) )
return false;
if( (dmin*tmin_) > (d_*tmin) )
{
tmin = tmin_;
dmin = d_;
}
if( (dmax*tmax_) < (d_*tmax) )
{
tmax = tmax_;
dmax = d_;
}
// -----------------------------------
// treat z coord
// -----------------------------------
if ( qz >= pz )
{
tmin_ = bzmin - pz;
tmax_ = bzmax - pz;
d_ = qz - pz;
}
else
{
tmin_ = pz - bzmax;
tmax_ = pz - bzmin;
d_ = pz - qz;
}
return ( (dmin*tmax_) >= (d_*tmin) && (dmax*tmin_) <= (d_*tmax) );
}
template <class K> template <class K>
bool do_intersect(const typename K::Ray_3& ray, bool do_intersect(const typename K::Ray_3& ray,
const CGAL::Bbox_3& bbox, const CGAL::Bbox_3& bbox,
@ -136,15 +45,16 @@ namespace internal {
const Point_3& source = ray.source(); const Point_3& source = ray.source();
const Point_3& point_on_ray = ray.second_point(); const Point_3& point_on_ray = ray.second_point();
return do_intersect_bbox_segment_aux<FT, true, false, false>( return do_intersect_bbox_segment_aux
<FT,
true, // bounded at t=0
false, // not bounded at t=1
false> // do not use static filters
(
source.x(), source.y(), source.z(), source.x(), source.y(), source.z(),
point_on_ray.x(), point_on_ray.y(), point_on_ray.z(), point_on_ray.x(), point_on_ray.y(), point_on_ray.z(),
bbox); bbox
// const CGAL::cpp0x::array<FT, 6> rray = {source.x(), source.y(), source.z(), );
// point_on_ray.x(), point_on_ray.y(), point_on_ray.z() };
// return do_intersect_bbox_segment_aux<FT, true, false>
// (rray,
// *reinterpret_cast<const CGAL::cpp0x::array<double, 6>*>(&*bbox.cartesian_begin()) );
} }
} // namespace internal } // namespace internal