1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153
|
/* poly/zsolve_cubic.c
*
* Copyright (C) 1996, 1997, 1998, 1999, 2000, 2007, 2009 Brian Gough
*
* This program is free software; you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation; either version 3 of the License, or (at
* your option) any later version.
*
* This program is distributed in the hope that it will be useful, but
* WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
* General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with this program; if not, write to the Free Software
* Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
*/
/* zsolve_cubic.c - finds the complex roots of x^3 + a x^2 + b x + c = 0 */
#include <config.h>
#include <math.h>
#include <gsl/gsl_math.h>
#include <gsl/gsl_complex.h>
#include <gsl/gsl_poly.h>
#define SWAP(a,b) do { double tmp = b ; b = a ; a = tmp ; } while(0)
int
gsl_poly_complex_solve_cubic (double a, double b, double c,
gsl_complex *z0, gsl_complex *z1,
gsl_complex *z2)
{
double q = (a * a - 3 * b);
double r = (2 * a * a * a - 9 * a * b + 27 * c);
double Q = q / 9;
double R = r / 54;
double Q3 = Q * Q * Q;
double R2 = R * R;
double CR2 = 729 * r * r;
double CQ3 = 2916 * q * q * q;
if (R == 0 && Q == 0)
{
GSL_REAL (*z0) = -a / 3;
GSL_IMAG (*z0) = 0;
GSL_REAL (*z1) = -a / 3;
GSL_IMAG (*z1) = 0;
GSL_REAL (*z2) = -a / 3;
GSL_IMAG (*z2) = 0;
return 3;
}
else if (CR2 == CQ3)
{
/* this test is actually R2 == Q3, written in a form suitable
for exact computation with integers */
/* Due to finite precision some double roots may be missed, and
will be considered to be a pair of complex roots z = x +/-
epsilon i close to the real axis. */
double sqrtQ = sqrt (Q);
if (R > 0)
{
GSL_REAL (*z0) = -2 * sqrtQ - a / 3;
GSL_IMAG (*z0) = 0;
GSL_REAL (*z1) = sqrtQ - a / 3;
GSL_IMAG (*z1) = 0;
GSL_REAL (*z2) = sqrtQ - a / 3;
GSL_IMAG (*z2) = 0;
}
else
{
GSL_REAL (*z0) = -sqrtQ - a / 3;
GSL_IMAG (*z0) = 0;
GSL_REAL (*z1) = -sqrtQ - a / 3;
GSL_IMAG (*z1) = 0;
GSL_REAL (*z2) = 2 * sqrtQ - a / 3;
GSL_IMAG (*z2) = 0;
}
return 3;
}
else if (R2 < Q3)
{
double sgnR = (R >= 0 ? 1 : -1);
double ratio = sgnR * sqrt (R2 / Q3);
double theta = acos (ratio);
double norm = -2 * sqrt (Q);
double r0 = norm * cos (theta / 3) - a / 3;
double r1 = norm * cos ((theta + 2.0 * M_PI) / 3) - a / 3;
double r2 = norm * cos ((theta - 2.0 * M_PI) / 3) - a / 3;
/* Sort r0, r1, r2 into increasing order */
if (r0 > r1)
SWAP (r0, r1);
if (r1 > r2)
{
SWAP (r1, r2);
if (r0 > r1)
SWAP (r0, r1);
}
GSL_REAL (*z0) = r0;
GSL_IMAG (*z0) = 0;
GSL_REAL (*z1) = r1;
GSL_IMAG (*z1) = 0;
GSL_REAL (*z2) = r2;
GSL_IMAG (*z2) = 0;
return 3;
}
else
{
double sgnR = (R >= 0 ? 1 : -1);
double A = -sgnR * pow (fabs (R) + sqrt (R2 - Q3), 1.0 / 3.0);
double B = Q / A;
if (A + B < 0)
{
GSL_REAL (*z0) = A + B - a / 3;
GSL_IMAG (*z0) = 0;
GSL_REAL (*z1) = -0.5 * (A + B) - a / 3;
GSL_IMAG (*z1) = -(sqrt (3.0) / 2.0) * fabs(A - B);
GSL_REAL (*z2) = -0.5 * (A + B) - a / 3;
GSL_IMAG (*z2) = (sqrt (3.0) / 2.0) * fabs(A - B);
}
else
{
GSL_REAL (*z0) = -0.5 * (A + B) - a / 3;
GSL_IMAG (*z0) = -(sqrt (3.0) / 2.0) * fabs(A - B);
GSL_REAL (*z1) = -0.5 * (A + B) - a / 3;
GSL_IMAG (*z1) = (sqrt (3.0) / 2.0) * fabs(A - B);
GSL_REAL (*z2) = A + B - a / 3;
GSL_IMAG (*z2) = 0;
}
return 3;
}
}
|