/* rpoly.cpp -- Jenkins-Traub real polynomial root finder. * * (C) 2000, C. Bond. All rights reserved. * * Translation of TOMS493 from FORTRAN to C. This * implementation of Jenkins-Traub partially adapts * the original code to a C environment by restruction * many of the 'goto' controls to better fit a block * structured form. It also eliminates the global memory * allocation in favor of local, dynamic memory management. * * The calling conventions are slightly modified to return * the number of roots found as the function value. * * INPUT: * op - double precision vector of coefficients in order of * decreasing powers. * degree - integer degree of polynomial * * OUTPUT: * zeror,zeroi - output double precision vectors of the * real and imaginary parts of the zeros. * * RETURN: * returnval: -1 if leading coefficient is zero, otherwise * number of roots found. */ #include "math.h" void quad(double a,double b1,double c,double *sr,double *si, double *lr,double *li); void fxshfr(int l2, int *nz); void quadit(double *uu,double *vv,int *nz); void realit(double *sss, int *nz, int *iflag); void calcsc(int *type); void nextk(int *type); void newest(int type,double *uu,double *vv); void quadsd(int n,double *u,double *v,double *p,double *q, double *a,double *b); double *p,*qp,*k,*qk,*svk; double sr,si,u,v,a,b,c,d,a1,a2; double a3,a6,a7,e,f,g,h,szr,szi,lzr,lzi; double eta,are,mre; int n,nn,nmi,zerok; int rpoly(double *op, int degree, double *zeror, double *zeroi) { double t,aa,bb,cc,*temp,factor,rot; double *pt; double lo,max,min,xx,yy,cosr,sinr,xxx,x,sc,bnd; double xm,ff,df,dx,infin,smalno,base; int cnt,nz,i,j,jj,l,nm1,zerok; /* The following statements set machine constants. */ base = 2.0; eta = 2.22e-16; infin = 3.4e38; smalno = 1.2e-38; are = eta; mre = eta; lo = smalno/eta; /* Initialization of constants for shift rotation. */ xx = sqrt(0.5); yy = -xx; rot = 94.0; rot *= 0.017453293; cosr = cos(rot); sinr = sin(rot); n = degree; /* Algorithm fails of the leading coefficient is zero. */ if (op[0] == 0.0) return -1; /* Remove the zeros at the origin, if any. */ while (op[n] == 0.0) { j = degree - n; zeror[j] = 0.0; zeroi[j] = 0.0; n--; } if (n < 1) return -1; /* * Allocate memory here */ temp = new double [degree+1]; pt = new double [degree+1]; p = new double [degree+1]; qp = new double [degree+1]; k = new double [degree+1]; qk = new double [degree+1]; svk = new double [degree+1]; /* Make a copy of the coefficients. */ for (i=0;i<=n;i++) p[i] = op[i]; /* Start the algorithm for one zero. */ _40: if (n == 1) { zeror[degree-1] = -p[1]/p[0]; zeroi[degree-1] = 0.0; n -= 1; goto _99; } /* Calculate the final zero or pair of zeros. */ if (n == 2) { quad(p[0],p[1],p[2],&zeror[degree-2],&zeroi[degree-2], &zeror[degree-1],&zeroi[degree-1]); n -= 2; goto _99; } /* Find largest and smallest moduli of coefficients. */ max = 0.0; min = infin; for (i=0;i<=n;i++) { x = fabs(p[i]); if (x > max) max = x; if (x != 0.0 && x < min) min = x; } /* Scale if there are large or very small coefficients. * Computes a scale factor to multiply the coefficients of the * polynomial. The scaling si done to avoid overflow and to * avoid undetected underflow interfering with the convergence * criterion. The factor is a power of the base. */ sc = lo/min; if (sc > 1.0 && infin/sc < max) goto _110; if (sc <= 1.0) { if (max < 10.0) goto _110; if (sc == 0.0) sc = smalno; } l = (int)(log(sc)/log(base) + 0.5); factor = pow(base*1.0,l); if (factor != 1.0) { for (i=0;i<=n;i++) p[i] = factor*p[i]; /* Scale polynomial. */ } _110: /* Compute lower bound on moduli of roots. */ for (i=0;i<=n;i++) { pt[i] = (fabs(p[i])); } pt[n] = - pt[n]; /* Compute upper estimate of bound. */ x = exp((log(-pt[n])-log(pt[0])) / (double)n); /* If Newton step at the origin is better, use it. */ if (pt[n-1] != 0.0) { xm = -pt[n]/pt[n-1]; if (xm < x) x = xm; } /* Chop the interval (0,x) until ff <= 0 */ while (1) { xm = x*0.1; ff = pt[0]; for (i=1;i<=n;i++) ff = ff*xm + pt[i]; if (ff <= 0.0) break; x = xm; } dx = x; /* Do Newton interation until x converges to two * decimal places. */ while (fabs(dx/x) > 0.005) { ff = pt[0]; df = ff; for (i=1;i 0) return; /* Quadratic iteration has failed. Flag that it has * been tried and decrease the convergence criterion. */ vtry = 1; betav *= 0.25; /* Try linear iteration if it has not been tried and * the S sequence is converging. */ if (stry || !spass) goto _50; for (i=0;i 0) return; /* Linear iteration has failed. Flag that it has been * tried and decrease the convergence criterion. */ stry = 1; betas *=0.25; if (iflag == 0) goto _50; /* If linear iteration signals an almost double real * zero attempt quadratic iteration. */ ui = -(s+s); vi = s*s; goto _20; /* Restore variables. */ _50: u = svu; v = svv; for (i=0;i 0.01 * fabs(lzr)) return; /* Evaluate polynomial by quadratic synthetic division. */ quadsd(n,&u,&v,p,qp,&a,&b); mp = fabs(a-szr*b) + fabs(szi*b); /* Compute a rigorous bound on the rounding error in * evaluating p. */ zm = sqrt(fabs(v)); ee = 2.0*fabs(qp[0]); t = -szr*b; for (i=1;i 20) return; if (j < 2) goto _50; if (relstp > 0.01 || mp < omp || tried) goto _50; /* A cluster appears to be stalling the convergence. * Five fixed shift steps are taken with a u,v close * to the cluster. */ if (relstp < eta) relstp = eta; relstp = sqrt(relstp); u = u - u*relstp; v = v + v*relstp; quadsd(n,&u,&v,p,qp,&a,&b); for (i=0;i<5;i++) { calcsc(&type); nextk(&type); } tried = 1; j = 0; _50: omp = mp; /* Calculate next k polynomial and new u and v. */ calcsc(&type); nextk(&type); calcsc(&type); newest(type,&ui,&vi); /* If vi is zero the iteration is not converging. */ if (vi == 0.0) return; relstp = fabs((vi-v)/vi); u = ui; v = vi; goto _10; } /* Variable-shift H polynomial iteration for a real zero. * sss - starting iterate * nz - number of zeros found * iflag - flag to indicate a pair of zeros near real axis. */ void realit(double *sss, int *nz, int *iflag) { double pv,kv,t,s; double ms,mp,omp,ee; int i,j; *nz = 0; s = *sss; *iflag = 0; j = 0; /* Main loop */ while (1) { pv = p[0]; /* Evaluate p at s. */ qp[0] = pv; for (i=1;i<=n;i++) { pv = pv*s + p[i]; qp[i] = pv; } mp = fabs(pv); /* Compute a rigorous bound on the error in evaluating p. */ ms = fabs(s); ee = (mre/(are+mre))*fabs(qp[0]); for (i=1;i<=n;i++) { ee = ee*ms + fabs(qp[i]); } /* Iteration has converged sufficiently if the polynomial * value is less than 20 times this bound. */ if (mp <= 20.0*((are+mre)*ee-mre*mp)) { *nz = 1; szr = s; szi = 0.0; return; } j++; /* Stop iteration after 10 steps. */ if (j > 10) return; if (j < 2) goto _50; if (fabs(t) > 0.001*fabs(s-t) || mp < omp) goto _50; /* A cluster of zeros near the real axis has been * encountered. Return with iflag set to initiate a * quadratic iteration. */ *iflag = 1; *sss = s; return; /* Return if the polynomial value has increased significantly. */ _50: omp = mp; /* Compute t, the next polynomial, and the new iterate. */ kv = k[0]; qk[0] = kv; for (i=1;i fabs(k[n-1]*10.0*eta)) t = -pv/kv; s += t; } } /* This routine calculates scalar quantities used to * compute the next k polynomial and new estimates of * the quadratic coefficients. * type - integer variable set here indicating how the * calculations are normalized to avoid overflow. */ void calcsc(int *type) { /* Synthetic division of k by the quadratic 1,u,v */ quadsd(n-1,&u,&v,k,qk,&c,&d); if (fabs(c) > fabs(k[n-1]*100.0*eta)) goto _10; if (fabs(d) > fabs(k[n-2]*100.0*eta)) goto _10; *type = 3; /* Type=3 indicates the quadratic is almost a factor of k. */ return; _10: if (fabs(d) < fabs(c)) { *type = 1; /* Type=1 indicates that all formulas are divided by c. */ e = a/c; f = d/c; g = u*e; h = v*b; a3 = a*e + (h/c+g)*b; a1 = b - a*(d/c); a7 = a + g*d + h*f; return; } *type = 2; /* Type=2 indicates that all formulas are divided by d. */ e = a/d; f = c/d; g = u*b; h = v*b; a3 = (a+g)*e + h*(b/d); a1 = b*f-a; a7 = (f+u)*a + h; } /* Computes the next k polynomials using scalars * computed in calcsc. */ void nextk(int *type) { double temp; int i; if (*type == 3) { /* Use unscaled form of the recurrence if type is 3. */ k[0] = 0.0; k[1] = 0.0; for (i=2;i= 0.0) /* real zeros. */ d = -d; *lr = (-b+d)/a; *sr = 0.0; if (*lr != 0.0) *sr = (c/ *lr)/a; *si = 0.0; *li = 0.0; } }