Factorize code in CGAL/internal/Static_filters/Do_intersect_3.h

The static filters of Do_intersect(BBox_3, Segment_3) and
Do_intersect(BBox_3, Segment_3) now use do_intersect_bbox_segment_aux(..)
from CGAL/internal/Intersections_3/Bbox_3_Segment_3_do_intersect.h
This commit is contained in:
Laurent Rineau 2012-06-15 14:57:36 +00:00
parent a72bd80380
commit 1c6ba1851b
1 changed files with 7 additions and 321 deletions

View File

@ -130,179 +130,6 @@ public:
if(is_indeterminate(ub)) return Base::operator()(s,b);
else return ub.sup();
// -----------------------------------
// treat x coord
// -----------------------------------
double dmin, dmax, tmin, tmax;
if ( qx >= px ) // this is input and needs no epsilon
{
if(px > bxmax) return false; // segment on the right of bbox
if(qx < bxmin) return false; // segment on the left of bbox
tmax = bxmax - px;
dmax = qx - px;
if ( bxmin < px ) // tmin < 0 means px is in the x-range of bbox
{
tmin = 0;
dmin = 1;
} else {
tmin = bxmin - px;
dmin = dmax;
}
}
else
{
if(qx > bxmax) return false; // segment on the right of bbox
if(px < bxmin) return false; // segment on the left of bbox
tmax = px - bxmin;
dmax = px - qx;
if ( px < bxmax ) // tmin < 0 means px is in the x-range of bbox
{
tmin = 0;
dmin = 1;
} else {
tmin = px - bxmax;
dmin = dmax;
}
}
double m = CGAL::abs(tmin), m2;
m2 = CGAL::abs(tmax); if(m2 > m) { m = m2; }
m2 = CGAL::abs(dmin); if(m2 > m) { m = m2; }
if(m < 7e-294) {
// underflow in the computation of 'error'
return Base::operator()(s,b);
}
const double EPS_1 = 3.55618e-15;
double error = EPS_1 * m;
switch(sign_with_error( tmax - dmax, error)) {
case POSITIVE:
tmax = 1;
dmax = 1;
break;
case NEGATIVE:
break;
default:
// ambiguity of the comparison tmax > dmin
// let's call the exact predicate
return Base::operator()(s,b);
}
// -----------------------------------
// treat y coord
// -----------------------------------
double d_, tmin_, tmax_;
if ( qy >= py ) // this is input and needs no epsilon
{
tmin_ = bymin - py;
tmax_ = bymax - py;
d_ = qy - py;
}
else
{
tmin_ = py - bymax;
tmax_ = py - bymin;
d_ = py - qy;
}
m2 = CGAL::abs(tmin_); if(m2 > m) { m = m2; }
m2 = CGAL::abs(tmax_); if(m2 > m) { m = m2; }
m2 = CGAL::abs(d_); if(m2 > m) { m = m2; }
if(m < 3e-147) {
// underflow in the computation of 'error'
return Base::operator()(s,b);
}
error = EPS_1 * m * m;
if(m > 1e153) { /* sqrt(max_double [hadamard]/2) */
// potential overflow on the computation of 'sign1' and 'sign2'
return Base::operator()(s,b);
}
Sign sign1 = sign_with_error( (d_*tmin) - (dmin*tmax_) , error);
Sign sign2 = sign_with_error( (dmax*tmin_) - (d_*tmax) , error);
if(sign1 == POSITIVE || sign2 == POSITIVE)
return false; // We are *sure* the segment is outside the box, on one
// side or the other.
if(sign1 == ZERO || sign2 == ZERO) {
return Base::operator()(s,b); // We are *unsure*: one *may be*
// positive.
}
// Here we are sure the two signs are negative. We can continue with
// the rest of the function...
// epsilons needed
switch(sign_with_error((dmin*tmin_) - (d_*tmin) , error)) {
case POSITIVE:
tmin = tmin_;
dmin = d_;
break;
case NEGATIVE:
break;
default: // uncertainty
return Base::operator()(s,b);
}
// epsilons needed
switch(sign_with_error((d_*tmax) - (dmax*tmax_) , error)) {
case POSITIVE:
tmax = tmax_;
dmax = d_;
case NEGATIVE:
break;
default: // uncertainty
return Base::operator()(s,b);
}
// -----------------------------------
// treat z coord
// -----------------------------------
if ( qz >= pz ) // this is input and needs no epsilon
{
tmin_ = bzmin - pz;
tmax_ = bzmax - pz;
d_ = qz - pz;
}
else
{
tmin_ = pz - bzmax;
tmax_ = pz - bzmin;
d_ = pz - qz;
}
m2 = CGAL::abs(tmin_); if(m2 > m) { m = m2; }
m2 = CGAL::abs(tmax_); if(m2 > m) { m = m2; }
m2 = CGAL::abs(d_); if(m2 > m) { m = m2; }
// m may have changed
error = EPS_1 * m * m;
if(m > 1e153) { /* sqrt(max_double [hadamard]/2) */
// potential overflow on the computation of 'sign1' and 'sign2'
return Base::operator()(s,b);
}
sign1 = sign_with_error( (dmin*tmax_) - (d_*tmin) , error);
sign2 = sign_with_error( (d_*tmax) - (dmax*tmin_) , error);
if(sign1 == NEGATIVE || sign2 == NEGATIVE) {
return false; // We are *sure* the segment is outside the box, on one
// side or the other.
}
if(sign1 == ZERO || sign2 == ZERO) {
return Base::operator()(s,b); // We are *unsure*: one *may be*
// negative.
}
return true; // We are *sure* the two signs are positive.
}
return Base::operator()(s,b);
}
@ -337,155 +164,14 @@ public:
{
CGAL_BRANCH_PROFILER_BRANCH_1(tmp);
// -----------------------------------
// treat x coord
// -----------------------------------
double dmin, dmax, tmin, tmax;
if ( qx >= px ) // this is input and needs no epsilon
{
if(px > bxmax) return false; // if tmax < 0, ray on the right of bbox
tmax = bxmax - px;
dmax = qx - px;
if ( bxmin < px ) // tmin < 0 means px is in the x-range of bbox
{
tmin = 0;
dmin = 1;
} else {
if( px == qx ) return false; // if dmin == 0
tmin = bxmin - px;
dmin = dmax;
}
}
else
{
if(px < bxmin) return false; // if tmax < 0, ray on the left of bbox
tmax = px - bxmin;
const Uncertain<result_type> ub =
do_intersect_bbox_segment_aux<double, true, false, true>
(px, py, pz,
qx, qy, qz,
b);
dmax = px - qx;
if ( px < bxmax ) // tmin < 0 means px is in the x-range of bbox
{
tmin = 0;
dmin = 1;
} else {
if( px == qx ) return false; // if dmin == 0
tmin = px - bxmax;
dmin = dmax;
}
}
// -----------------------------------
// treat y coord
// -----------------------------------
double d_, tmin_, tmax_;
if ( qy >= py ) // this is input and needs no epsilon
{
tmin_ = bymin - py;
tmax_ = bymax - py;
d_ = qy - py;
}
else
{
tmin_ = py - bymax;
tmax_ = py - bymin;
d_ = py - qy;
}
double m = CGAL::abs(tmin), m2;
m2 = CGAL::abs(tmax); if(m2 > m) { m = m2; }
m2 = CGAL::abs(dmin); if(m2 > m) { m = m2; }
m2 = CGAL::abs(tmin_); if(m2 > m) { m = m2; }
m2 = CGAL::abs(tmax_); if(m2 > m) { m = m2; }
m2 = CGAL::abs(d_); if(m2 > m) { m = m2; }
if(m < 3e-147) {
// underflow in the computation of 'error'
return Base::operator()(r,b);
}
const double EPS_1 = 3.55618e-15;
double error = EPS_1 * m * m;
if(m > 1e153) { /* sqrt(max_double [hadamard]/2) */
// potential overflow on the computation of 'sign1' and 'sign2'
return Base::operator()(r,b);
}
Sign sign1 = sign_with_error( (d_*tmin) - (dmin*tmax_) , error);
Sign sign2 = sign_with_error( (dmax*tmin_) - (d_*tmax) , error);
if(sign1 == POSITIVE || sign2 == POSITIVE)
return false; // We are *sure* the ray is outside the box, on one
// side or the other.
if(sign1 == ZERO || sign2 == ZERO) {
return Base::operator()(r,b); // We are *unsure*: one *may be*
// positive.
}
// Here we are sure the two signs are negative. We can continue with
// the rest of the function...
// epsilons needed
switch(sign_with_error((dmin*tmin_) - (d_*tmin) , error)) {
case POSITIVE:
tmin = tmin_;
dmin = d_;
break;
case NEGATIVE:
break;
default: // uncertainty
return Base::operator()(r,b);
}
// epsilons needed
switch(sign_with_error((d_*tmax) - (dmax*tmax_) , error)) {
case POSITIVE:
tmax = tmax_;
dmax = d_;
case NEGATIVE:
break;
default: // uncertainty
return Base::operator()(r,b);
}
// -----------------------------------
// treat z coord
// -----------------------------------
if ( qz >= pz ) // this is input and needs no epsilon
{
tmin_ = bzmin - pz;
tmax_ = bzmax - pz;
d_ = qz - pz;
}
else
{
tmin_ = pz - bzmax;
tmax_ = pz - bzmin;
d_ = pz - qz;
}
m2 = CGAL::abs(tmin_); if(m2 > m) { m = m2; }
m2 = CGAL::abs(tmax_); if(m2 > m) { m = m2; }
m2 = CGAL::abs(d_); if(m2 > m) { m = m2; }
// m may have changed
error = EPS_1 * m * m;
if(m > 1e153) { /* sqrt(max_double [hadamard]/2) */
// potential overflow on the computation of 'sign1' and 'sign2'
return Base::operator()(r,b);
}
sign1 = sign_with_error( (dmin*tmax_) - (d_*tmin) , error);
sign2 = sign_with_error( (d_*tmax) - (dmax*tmin_) , error);
if(sign1 == NEGATIVE || sign2 == NEGATIVE) {
return false; // We are *sure* the ray is outside the box, on one
// side or the other.
}
if(sign1 == ZERO || sign2 == ZERO) {
return Base::operator()(r,b); // We are *unsure*: one *may be*
// negative.
}
return true; // We are *sure* the two signs are positive.
if(is_indeterminate(ub)) return Base::operator()(r,b);
else return ub.sup();
}
return Base::operator()(r,b);
}