BWAPI
|
00001 // Copyright (c) 2000,2001 Utrecht University (The Netherlands), 00002 // ETH Zurich (Switzerland), Freie Universitaet Berlin (Germany), 00003 // INRIA Sophia-Antipolis (France), Martin-Luther-University Halle-Wittenberg 00004 // (Germany), Max-Planck-Institute Saarbruecken (Germany), RISC Linz (Austria), 00005 // and Tel-Aviv University (Israel). All rights reserved. 00006 // 00007 // This file is part of CGAL (www.cgal.org); you can redistribute it and/or 00008 // modify it under the terms of the GNU Lesser General Public License as 00009 // published by the Free Software Foundation; version 2.1 of the License. 00010 // See the file LICENSE.LGPL distributed with CGAL. 00011 // 00012 // Licensees holding a valid commercial license may use this file in 00013 // accordance with the commercial license agreement provided with the software. 00014 // 00015 // This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE 00016 // WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE. 00017 // 00018 // $URL: svn+ssh://scm.gforge.inria.fr/svn/cgal/branches/CGAL-3.5-branch/Conic_2/include/CGAL/Conic_misc.h $ 00019 // $Id: Conic_misc.h 28567 2006-02-16 14:30:13Z lsaboret $ 00020 // 00021 // 00022 // Author(s) : Bernd Gaertner, Sven Schoenherr <sven@inf.ethz.ch> 00023 00024 #ifndef CGAL_CONIC_MISC_H 00025 #define CGAL_CONIC_MISC_H 00026 00027 #include <cmath> 00028 #include <CGAL/kernel_assertions.h> 00029 00030 CGAL_BEGIN_NAMESPACE 00031 00032 template < class R> 00033 class Conic_2; 00034 00035 00036 enum Conic_type 00037 { 00038 HYPERBOLA = -1, 00039 PARABOLA, 00040 ELLIPSE 00041 }; 00042 00043 00044 typedef CGAL::Bounded_side Convex_side; 00045 const Convex_side ON_CONVEX_SIDE = CGAL::ON_BOUNDED_SIDE; 00046 const Convex_side ON_NONCONVEX_SIDE = CGAL::ON_UNBOUNDED_SIDE; 00047 00048 00049 00050 00051 template < class NT > 00052 NT best_value (NT *values, int nr_values, 00053 NT a2, NT a1, NT a0, 00054 NT b3, NT b2, NT b1, NT b0) 00055 { 00056 bool det_positive = false; 00057 NT d, q, max_det = 0, det, best = -1; 00058 for (int i=0; i<nr_values; ++i) { 00059 NT x = values[i]; 00060 d = (a2*x+a1)*x+a0; 00061 q = ((b3*x+b2)*x+b1)*x+b0; 00062 det = d*d*d/(q*q); 00063 // if q==0, this root value doesn't qualify for the 00064 // best value. Under roundoff errors, q might be very 00065 // small but nonzero, so that the value is erroneously 00066 // being considered; however, d should be very small 00067 // in this case as well, so that det won't compete 00068 // for max_det below. 00069 if (CGAL_NTS is_positive(det) && !CGAL_NTS is_zero(q)) 00070 if (!det_positive || (det > max_det)) { 00071 max_det = det; 00072 best = x; 00073 det_positive = true; 00074 } 00075 } 00076 CGAL_kernel_precondition (det_positive); 00077 return best; 00078 } 00079 00080 00081 template < class NT > 00082 int solve_cubic (NT c3, NT c2, NT c1, NT c0, 00083 NT& r1, NT& r2, NT& r3) 00084 { 00085 if (c3 == 0.0) { 00086 // quadratic equation 00087 if (c2 == 0) { 00088 // linear equation 00089 CGAL_kernel_precondition (c1 != 0); 00090 r1 = -c0/c1; 00091 return 1; 00092 } 00093 NT D = c1*c1-4*c2*c0; 00094 if (D < 0.0) 00095 // only complex roots 00096 return 0; 00097 if (D == 0.0) { 00098 // one real root 00099 r1 = -c1/(2.0*c2); 00100 return 1; 00101 } 00102 // two real roots 00103 r1 = (-c1 + CGAL_NTS sqrt(D))/(2.0*c2); 00104 r2 = (-c1 - CGAL_NTS sqrt(D))/(2.0*c2); 00105 return 2; 00106 } 00107 00108 // cubic equation 00109 // define the gamma_i 00110 NT g2 = c2/c3, 00111 g1 = c1/c3, 00112 g0 = c0/c3; 00113 00114 // define a, b 00115 NT a = g1 - g2*g2/3.0, 00116 b = 2.0*g2*g2*g2/27.0 - g1*g2/3.0 + g0; 00117 00118 if (a == 0) { 00119 // one real root 00120 /***** r1 = cbrt(-b) - g2/3.0; *****/ 00121 r1 = exp(log(-b)/3.0) - g2/3.0; 00122 return 1; 00123 } 00124 00125 // define D 00126 NT D = a*a*a/27.0 + b*b/4.0; 00127 if (D >= 0.0) { 00128 // real case 00129 /***** NT u = cbrt(-b/2.0 + CGAL_NTS sqrt(D)), *****/ 00130 NT u = exp(log(-b/2.0 + CGAL_NTS sqrt(D))), 00131 alpha = 1.0 - a/(3.0*u*u); 00132 if (D == 0) { 00133 // two distinct real roots 00134 r1 = u*alpha - g2/3.0; 00135 r2 = -0.5*alpha*u - g2/3.0; 00136 return 2; 00137 } 00138 // one real root 00139 r1 = u*alpha - g2/3.0; 00140 return 1; 00141 } 00142 // complex case 00143 NT r_prime = CGAL_NTS sqrt(-a/3), 00144 phi_prime = acos (-b/(2.0*r_prime*r_prime*r_prime))/3.0, 00145 u_R = r_prime * cos (phi_prime), 00146 u_I = r_prime * sin (phi_prime); 00147 // three distinct real roots 00148 r1 = 2.0*u_R - g2/3.0; 00149 r2 = -u_R + u_I*std::sqrt(3.0) - g2/3.0; 00150 r3 = -u_R - u_I*std::sqrt(3.0) - g2/3.0; 00151 return 3; 00152 } 00153 00154 00155 00156 CGAL_END_NAMESPACE 00157 00158 #endif // CGAL_CONIC_MISC_H 00159 00160 // ===== EOF ==================================================================