| [833] | 1 | //
|
|---|
| 2 | // ********************************************************************
|
|---|
| 3 | // * License and Disclaimer *
|
|---|
| 4 | // * *
|
|---|
| 5 | // * The Geant4 software is copyright of the Copyright Holders of *
|
|---|
| 6 | // * the Geant4 Collaboration. It is provided under the terms and *
|
|---|
| 7 | // * conditions of the Geant4 Software License, included in the file *
|
|---|
| 8 | // * LICENSE and available at http://cern.ch/geant4/license . These *
|
|---|
| 9 | // * include a list of copyright holders. *
|
|---|
| 10 | // * *
|
|---|
| 11 | // * Neither the authors of this software system, nor their employing *
|
|---|
| 12 | // * institutes,nor the agencies providing financial support for this *
|
|---|
| 13 | // * work make any representation or warranty, express or implied, *
|
|---|
| 14 | // * regarding this software system or assume any liability for its *
|
|---|
| 15 | // * use. Please see the license in the file LICENSE and URL above *
|
|---|
| 16 | // * for the full disclaimer and the limitation of liability. *
|
|---|
| 17 | // * *
|
|---|
| 18 | // * This code implementation is the result of the scientific and *
|
|---|
| 19 | // * technical work of the GEANT4 collaboration. *
|
|---|
| 20 | // * By using, copying, modifying or distributing the software (or *
|
|---|
| 21 | // * any work based on the software) you agree to acknowledge its *
|
|---|
| 22 | // * use in resulting scientific publications, and indicate your *
|
|---|
| 23 | // * acceptance of all terms of the Geant4 Software license. *
|
|---|
| 24 | // ********************************************************************
|
|---|
| 25 | //
|
|---|
| [850] | 26 | // $Id: G4JTPolynomialSolver.cc,v 1.7 2008/03/13 09:35:57 gcosmo Exp $
|
|---|
| 27 | // GEANT4 tag $Name: HEAD $
|
|---|
| [833] | 28 | //
|
|---|
| 29 | // --------------------------------------------------------------------
|
|---|
| 30 | // GEANT 4 class source file
|
|---|
| 31 | //
|
|---|
| 32 | // G4JTPolynomialSolver
|
|---|
| 33 | //
|
|---|
| 34 | // Implementation based on Jenkins-Traub algorithm.
|
|---|
| 35 | // --------------------------------------------------------------------
|
|---|
| 36 |
|
|---|
| 37 | #include "G4JTPolynomialSolver.hh"
|
|---|
| 38 |
|
|---|
| 39 | const G4double G4JTPolynomialSolver::base = 2;
|
|---|
| 40 | const G4double G4JTPolynomialSolver::eta = DBL_EPSILON;
|
|---|
| 41 | const G4double G4JTPolynomialSolver::infin = DBL_MAX;
|
|---|
| 42 | const G4double G4JTPolynomialSolver::smalno = DBL_MIN;
|
|---|
| 43 | const G4double G4JTPolynomialSolver::are = DBL_EPSILON;
|
|---|
| 44 | const G4double G4JTPolynomialSolver::mre = DBL_EPSILON;
|
|---|
| 45 | const G4double G4JTPolynomialSolver::lo = DBL_MIN/DBL_EPSILON ;
|
|---|
| 46 |
|
|---|
| 47 | G4JTPolynomialSolver::G4JTPolynomialSolver()
|
|---|
| 48 | {
|
|---|
| 49 | }
|
|---|
| 50 |
|
|---|
| 51 | G4JTPolynomialSolver::~G4JTPolynomialSolver()
|
|---|
| 52 | {
|
|---|
| 53 | }
|
|---|
| 54 |
|
|---|
| 55 | G4int G4JTPolynomialSolver::FindRoots(G4double *op, G4int degr,
|
|---|
| 56 | G4double *zeror, G4double *zeroi)
|
|---|
| 57 | {
|
|---|
| 58 | G4double t=0.0, aa=0.0, bb=0.0, cc=0.0, factor=1.0;
|
|---|
| 59 | G4double max=0.0, min=infin, xxx=0.0, x=0.0, sc=0.0, bnd=0.0;
|
|---|
| 60 | G4double xm=0.0, ff=0.0, df=0.0, dx=0.0;
|
|---|
| 61 | G4int cnt=0, nz=0, i=0, j=0, jj=0, l=0, nm1=0, zerok=0;
|
|---|
| 62 |
|
|---|
| 63 | // Initialization of constants for shift rotation.
|
|---|
| 64 | //
|
|---|
| 65 | G4double xx = std::sqrt(0.5);
|
|---|
| 66 | G4double yy = -xx,
|
|---|
| 67 | rot = 94.0*deg;
|
|---|
| 68 | G4double cosr = std::cos(rot),
|
|---|
| 69 | sinr = std::sin(rot);
|
|---|
| 70 | n = degr;
|
|---|
| 71 |
|
|---|
| 72 | // Algorithm fails if the leading coefficient is zero.
|
|---|
| 73 | //
|
|---|
| 74 | if (!(op[0] != 0.0)) { return -1; }
|
|---|
| 75 |
|
|---|
| 76 | // Remove the zeros at the origin, if any.
|
|---|
| 77 | //
|
|---|
| 78 | while (!(op[n] != 0.0))
|
|---|
| 79 | {
|
|---|
| 80 | j = degr - n;
|
|---|
| 81 | zeror[j] = 0.0;
|
|---|
| 82 | zeroi[j] = 0.0;
|
|---|
| 83 | n--;
|
|---|
| 84 | }
|
|---|
| 85 | if (n < 1) { return -1; }
|
|---|
| 86 |
|
|---|
| 87 | // Allocate buffers here
|
|---|
| 88 | //
|
|---|
| 89 | std::vector<G4double> temp(degr+1) ;
|
|---|
| 90 | std::vector<G4double> pt(degr+1) ;
|
|---|
| 91 |
|
|---|
| 92 | p.assign(degr+1,0) ;
|
|---|
| 93 | qp.assign(degr+1,0) ;
|
|---|
| 94 | k.assign(degr+1,0) ;
|
|---|
| 95 | qk.assign(degr+1,0) ;
|
|---|
| 96 | svk.assign(degr+1,0) ;
|
|---|
| 97 |
|
|---|
| 98 | // Make a copy of the coefficients.
|
|---|
| 99 | //
|
|---|
| 100 | for (i=0;i<=n;i++)
|
|---|
| 101 | { p[i] = op[i]; }
|
|---|
| 102 |
|
|---|
| 103 | do
|
|---|
| 104 | {
|
|---|
| 105 | if (n == 1) // Start the algorithm for one zero.
|
|---|
| 106 | {
|
|---|
| 107 | zeror[degr-1] = -p[1]/p[0];
|
|---|
| 108 | zeroi[degr-1] = 0.0;
|
|---|
| 109 | n -= 1;
|
|---|
| 110 | return degr - n ;
|
|---|
| 111 | }
|
|---|
| 112 | if (n == 2) // Calculate the final zero or pair of zeros.
|
|---|
| 113 | {
|
|---|
| 114 | Quadratic(p[0],p[1],p[2],&zeror[degr-2],&zeroi[degr-2],
|
|---|
| 115 | &zeror[degr-1],&zeroi[degr-1]);
|
|---|
| 116 | n -= 2;
|
|---|
| 117 | return degr - n ;
|
|---|
| 118 | }
|
|---|
| 119 |
|
|---|
| 120 | // Find largest and smallest moduli of coefficients.
|
|---|
| 121 | //
|
|---|
| 122 | max = 0.0;
|
|---|
| 123 | min = infin;
|
|---|
| 124 | for (i=0;i<=n;i++)
|
|---|
| 125 | {
|
|---|
| 126 | x = std::fabs(p[i]);
|
|---|
| 127 | if (x > max) { max = x; }
|
|---|
| 128 | if (x != 0.0 && x < min) { min = x; }
|
|---|
| 129 | }
|
|---|
| 130 |
|
|---|
| 131 | // Scale if there are large or very small coefficients.
|
|---|
| 132 | // Computes a scale factor to multiply the coefficients of the
|
|---|
| 133 | // polynomial. The scaling is done to avoid overflow and to
|
|---|
| 134 | // avoid undetected underflow interfering with the convergence
|
|---|
| 135 | // criterion. The factor is a power of the base.
|
|---|
| 136 | //
|
|---|
| 137 | sc = lo/min;
|
|---|
| 138 |
|
|---|
| 139 | if ( ((sc <= 1.0) && (max >= 10.0))
|
|---|
| 140 | || ((sc > 1.0) && (infin/sc >= max))
|
|---|
| 141 | || ((infin/sc >= max) && (max >= 10)) )
|
|---|
| 142 | {
|
|---|
| 143 | if (!( sc != 0.0 ))
|
|---|
| 144 | { sc = smalno ; }
|
|---|
| 145 | l = (G4int)(std::log(sc)/std::log(base) + 0.5);
|
|---|
| 146 | factor = std::pow(base*1.0,l);
|
|---|
| 147 | if (factor != 1.0)
|
|---|
| 148 | {
|
|---|
| 149 | for (i=0;i<=n;i++)
|
|---|
| 150 | { p[i] = factor*p[i]; } // Scale polynomial.
|
|---|
| 151 | }
|
|---|
| 152 | }
|
|---|
| 153 |
|
|---|
| 154 | // Compute lower bound on moduli of roots.
|
|---|
| 155 | //
|
|---|
| 156 | for (i=0;i<=n;i++)
|
|---|
| 157 | {
|
|---|
| 158 | pt[i] = (std::fabs(p[i]));
|
|---|
| 159 | }
|
|---|
| 160 | pt[n] = - pt[n];
|
|---|
| 161 |
|
|---|
| 162 | // Compute upper estimate of bound.
|
|---|
| 163 | //
|
|---|
| 164 | x = std::exp((std::log(-pt[n])-std::log(pt[0])) / (G4double)n);
|
|---|
| 165 |
|
|---|
| 166 | // If Newton step at the origin is better, use it.
|
|---|
| 167 | //
|
|---|
| 168 | if (pt[n-1] != 0.0)
|
|---|
| 169 | {
|
|---|
| 170 | xm = -pt[n]/pt[n-1];
|
|---|
| 171 | if (xm < x) { x = xm; }
|
|---|
| 172 | }
|
|---|
| 173 |
|
|---|
| 174 | // Chop the interval (0,x) until ff <= 0
|
|---|
| 175 | //
|
|---|
| 176 | while (1)
|
|---|
| 177 | {
|
|---|
| 178 | xm = x*0.1;
|
|---|
| 179 | ff = pt[0];
|
|---|
| 180 | for (i=1;i<=n;i++)
|
|---|
| 181 | { ff = ff*xm + pt[i]; }
|
|---|
| 182 | if (ff <= 0.0) { break; }
|
|---|
| 183 | x = xm;
|
|---|
| 184 | }
|
|---|
| 185 | dx = x;
|
|---|
| 186 |
|
|---|
| 187 | // Do Newton interation until x converges to two decimal places.
|
|---|
| 188 | //
|
|---|
| 189 | while (std::fabs(dx/x) > 0.005)
|
|---|
| 190 | {
|
|---|
| 191 | ff = pt[0];
|
|---|
| 192 | df = ff;
|
|---|
| 193 | for (i=1;i<n;i++)
|
|---|
| 194 | {
|
|---|
| 195 | ff = ff*x + pt[i];
|
|---|
| 196 | df = df*x + ff;
|
|---|
| 197 | }
|
|---|
| 198 | ff = ff*x + pt[n];
|
|---|
| 199 | dx = ff/df;
|
|---|
| 200 | x -= dx;
|
|---|
| 201 | }
|
|---|
| 202 | bnd = x;
|
|---|
| 203 |
|
|---|
| 204 | // Compute the derivative as the initial k polynomial
|
|---|
| 205 | // and do 5 steps with no shift.
|
|---|
| 206 | //
|
|---|
| 207 | nm1 = n - 1;
|
|---|
| 208 | for (i=1;i<n;i++)
|
|---|
| 209 | { k[i] = (G4double)(n-i)*p[i]/(G4double)n; }
|
|---|
| 210 | k[0] = p[0];
|
|---|
| 211 | aa = p[n];
|
|---|
| 212 | bb = p[n-1];
|
|---|
| 213 | zerok = (k[n-1] == 0);
|
|---|
| 214 | for(jj=0;jj<5;jj++)
|
|---|
| 215 | {
|
|---|
| 216 | cc = k[n-1];
|
|---|
| 217 | if (!zerok) // Use a scaled form of recurrence if k at 0 is nonzero.
|
|---|
| 218 | {
|
|---|
| 219 | // Use a scaled form of recurrence if value of k at 0 is nonzero.
|
|---|
| 220 | //
|
|---|
| 221 | t = -aa/cc;
|
|---|
| 222 | for (i=0;i<nm1;i++)
|
|---|
| 223 | {
|
|---|
| 224 | j = n-i-1;
|
|---|
| 225 | k[j] = t*k[j-1]+p[j];
|
|---|
| 226 | }
|
|---|
| 227 | k[0] = p[0];
|
|---|
| 228 | zerok = (std::fabs(k[n-1]) <= std::fabs(bb)*eta*10.0);
|
|---|
| 229 | }
|
|---|
| 230 | else // Use unscaled form of recurrence.
|
|---|
| 231 | {
|
|---|
| 232 | for (i=0;i<nm1;i++)
|
|---|
| 233 | {
|
|---|
| 234 | j = n-i-1;
|
|---|
| 235 | k[j] = k[j-1];
|
|---|
| 236 | }
|
|---|
| 237 | k[0] = 0.0;
|
|---|
| 238 | zerok = (!(k[n-1] != 0.0));
|
|---|
| 239 | }
|
|---|
| 240 | }
|
|---|
| 241 |
|
|---|
| 242 | // Save k for restarts with new shifts.
|
|---|
| 243 | //
|
|---|
| 244 | for (i=0;i<n;i++)
|
|---|
| 245 | { temp[i] = k[i]; }
|
|---|
| 246 |
|
|---|
| 247 | // Loop to select the quadratic corresponding to each new shift.
|
|---|
| 248 | //
|
|---|
| 249 | for (cnt = 0;cnt < 20;cnt++)
|
|---|
| 250 | {
|
|---|
| 251 | // Quadratic corresponds to a double shift to a
|
|---|
| 252 | // non-real point and its complex conjugate. The point
|
|---|
| 253 | // has modulus bnd and amplitude rotated by 94 degrees
|
|---|
| 254 | // from the previous shift.
|
|---|
| 255 | //
|
|---|
| 256 | xxx = cosr*xx - sinr*yy;
|
|---|
| 257 | yy = sinr*xx + cosr*yy;
|
|---|
| 258 | xx = xxx;
|
|---|
| 259 | sr = bnd*xx;
|
|---|
| 260 | si = bnd*yy;
|
|---|
| 261 | u = -2.0 * sr;
|
|---|
| 262 | v = bnd;
|
|---|
| 263 | ComputeFixedShiftPolynomial(20*(cnt+1),&nz);
|
|---|
| 264 | if (nz != 0)
|
|---|
| 265 | {
|
|---|
| 266 | // The second stage jumps directly to one of the third
|
|---|
| 267 | // stage iterations and returns here if successful.
|
|---|
| 268 | // Deflate the polynomial, store the zero or zeros and
|
|---|
| 269 | // return to the main algorithm.
|
|---|
| 270 | //
|
|---|
| 271 | j = degr - n;
|
|---|
| 272 | zeror[j] = szr;
|
|---|
| 273 | zeroi[j] = szi;
|
|---|
| 274 | n -= nz;
|
|---|
| 275 | for (i=0;i<=n;i++)
|
|---|
| 276 | { p[i] = qp[i]; }
|
|---|
| 277 | if (nz != 1)
|
|---|
| 278 | {
|
|---|
| 279 | zeror[j+1] = lzr;
|
|---|
| 280 | zeroi[j+1] = lzi;
|
|---|
| 281 | }
|
|---|
| 282 | break;
|
|---|
| 283 | }
|
|---|
| 284 | else
|
|---|
| 285 | {
|
|---|
| 286 | // If the iteration is unsuccessful another quadratic
|
|---|
| 287 | // is chosen after restoring k.
|
|---|
| 288 | //
|
|---|
| 289 | for (i=0;i<n;i++)
|
|---|
| 290 | {
|
|---|
| 291 | k[i] = temp[i];
|
|---|
| 292 | }
|
|---|
| 293 | }
|
|---|
| 294 | }
|
|---|
| 295 | }
|
|---|
| 296 | while (nz != 0); // End of initial DO loop
|
|---|
| 297 |
|
|---|
| 298 | // Return with failure if no convergence with 20 shifts.
|
|---|
| 299 | //
|
|---|
| 300 | return degr - n;
|
|---|
| 301 | }
|
|---|
| 302 |
|
|---|
| 303 | void G4JTPolynomialSolver::ComputeFixedShiftPolynomial(G4int l2, G4int *nz)
|
|---|
| 304 | {
|
|---|
| 305 | // Computes up to L2 fixed shift k-polynomials, testing for convergence
|
|---|
| 306 | // in the linear or quadratic case. Initiates one of the variable shift
|
|---|
| 307 | // iterations and returns with the number of zeros found.
|
|---|
| 308 |
|
|---|
| 309 | G4double svu=0.0, svv=0.0, ui=0.0, vi=0.0, xs=0.0;
|
|---|
| 310 | G4double betas=0.25, betav=0.25, oss=sr, ovv=v,
|
|---|
| 311 | ss=0.0, vv=0.0, ts=1.0, tv=1.0;
|
|---|
| 312 | G4double ots=0.0, otv=0.0;
|
|---|
| 313 | G4double tvv=1.0, tss=1.0;
|
|---|
| 314 | G4int type=0, i=0, j=0, iflag=0, vpass=0, spass=0, vtry=0, stry=0;
|
|---|
| 315 |
|
|---|
| 316 | *nz = 0;
|
|---|
| 317 |
|
|---|
| 318 | // Evaluate polynomial by synthetic division.
|
|---|
| 319 | //
|
|---|
| 320 | QuadraticSyntheticDivision(n,&u,&v,p,qp,&a,&b);
|
|---|
| 321 | ComputeScalarFactors(&type);
|
|---|
| 322 | for (j=0;j<l2;j++)
|
|---|
| 323 | {
|
|---|
| 324 | // Calculate next k polynomial and estimate v.
|
|---|
| 325 | //
|
|---|
| 326 | ComputeNextPolynomial(&type);
|
|---|
| 327 | ComputeScalarFactors(&type);
|
|---|
| 328 | ComputeNewEstimate(type,&ui,&vi);
|
|---|
| 329 | vv = vi;
|
|---|
| 330 |
|
|---|
| 331 | // Estimate xs.
|
|---|
| 332 | //
|
|---|
| 333 | ss = 0.0;
|
|---|
| 334 | if (k[n-1] != 0.0) { ss = -p[n]/k[n-1]; }
|
|---|
| 335 | tv = 1.0;
|
|---|
| 336 | ts = 1.0;
|
|---|
| 337 | if (j == 0 || type == 3)
|
|---|
| 338 | {
|
|---|
| 339 | ovv = vv;
|
|---|
| 340 | oss = ss;
|
|---|
| 341 | otv = tv;
|
|---|
| 342 | ots = ts;
|
|---|
| 343 | continue;
|
|---|
| 344 | }
|
|---|
| 345 |
|
|---|
| 346 | // Compute relative measures of convergence of xs and v sequences.
|
|---|
| 347 | //
|
|---|
| 348 | if (vv != 0.0) { tv = std::fabs((vv-ovv)/vv); }
|
|---|
| 349 | if (ss != 0.0) { ts = std::fabs((ss-oss)/ss); }
|
|---|
| 350 |
|
|---|
| 351 | // If decreasing, multiply two most recent convergence measures.
|
|---|
| 352 | tvv = 1.0;
|
|---|
| 353 | if (tv < otv) { tvv = tv*otv; }
|
|---|
| 354 | tss = 1.0;
|
|---|
| 355 | if (ts < ots) { tss = ts*ots; }
|
|---|
| 356 |
|
|---|
| 357 | // Compare with convergence criteria.
|
|---|
| 358 | vpass = (tvv < betav);
|
|---|
| 359 | spass = (tss < betas);
|
|---|
| 360 | if (!(spass || vpass))
|
|---|
| 361 | {
|
|---|
| 362 | ovv = vv;
|
|---|
| 363 | oss = ss;
|
|---|
| 364 | otv = tv;
|
|---|
| 365 | ots = ts;
|
|---|
| 366 | continue;
|
|---|
| 367 | }
|
|---|
| 368 |
|
|---|
| 369 | // At least one sequence has passed the convergence test.
|
|---|
| 370 | // Store variables before iterating.
|
|---|
| 371 | //
|
|---|
| 372 | svu = u;
|
|---|
| 373 | svv = v;
|
|---|
| 374 | for (i=0;i<n;i++)
|
|---|
| 375 | {
|
|---|
| 376 | svk[i] = k[i];
|
|---|
| 377 | }
|
|---|
| 378 | xs = ss;
|
|---|
| 379 |
|
|---|
| 380 | // Choose iteration according to the fastest converging sequence.
|
|---|
| 381 | //
|
|---|
| 382 | vtry = 0;
|
|---|
| 383 | stry = 0;
|
|---|
| 384 | if ((spass && (!vpass)) || (tss < tvv))
|
|---|
| 385 | {
|
|---|
| 386 | RealPolynomialIteration(&xs,nz,&iflag);
|
|---|
| 387 | if (*nz > 0) { return; }
|
|---|
| 388 |
|
|---|
| 389 | // Linear iteration has failed. Flag that it has been
|
|---|
| 390 | // tried and decrease the convergence criterion.
|
|---|
| 391 | //
|
|---|
| 392 | stry = 1;
|
|---|
| 393 | betas *=0.25;
|
|---|
| 394 | if (iflag == 0) { goto _restore_variables; }
|
|---|
| 395 |
|
|---|
| 396 | // If linear iteration signals an almost double real
|
|---|
| 397 | // zero attempt quadratic iteration.
|
|---|
| 398 | //
|
|---|
| 399 | ui = -(xs+xs);
|
|---|
| 400 | vi = xs*xs;
|
|---|
| 401 | }
|
|---|
| 402 |
|
|---|
| 403 | _quadratic_iteration:
|
|---|
| 404 |
|
|---|
| 405 | do
|
|---|
| 406 | {
|
|---|
| 407 | QuadraticPolynomialIteration(&ui,&vi,nz);
|
|---|
| 408 | if (*nz > 0) { return; }
|
|---|
| 409 |
|
|---|
| 410 | // Quadratic iteration has failed. Flag that it has
|
|---|
| 411 | // been tried and decrease the convergence criterion.
|
|---|
| 412 | //
|
|---|
| 413 | vtry = 1;
|
|---|
| 414 | betav *= 0.25;
|
|---|
| 415 |
|
|---|
| 416 | // Try linear iteration if it has not been tried and
|
|---|
| 417 | // the S sequence is converging.
|
|---|
| 418 | //
|
|---|
| 419 | if (stry || !spass) { break; }
|
|---|
| 420 | for (i=0;i<n;i++)
|
|---|
| 421 | {
|
|---|
| 422 | k[i] = svk[i];
|
|---|
| 423 | }
|
|---|
| 424 | RealPolynomialIteration(&xs,nz,&iflag);
|
|---|
| 425 | if (*nz > 0) { return; }
|
|---|
| 426 |
|
|---|
| 427 | // Linear iteration has failed. Flag that it has been
|
|---|
| 428 | // tried and decrease the convergence criterion.
|
|---|
| 429 | //
|
|---|
| 430 | stry = 1;
|
|---|
| 431 | betas *=0.25;
|
|---|
| 432 | if (iflag == 0) { break; }
|
|---|
| 433 |
|
|---|
| 434 | // If linear iteration signals an almost double real
|
|---|
| 435 | // zero attempt quadratic iteration.
|
|---|
| 436 | //
|
|---|
| 437 | ui = -(xs+xs);
|
|---|
| 438 | vi = xs*xs;
|
|---|
| 439 | }
|
|---|
| 440 | while (iflag != 0);
|
|---|
| 441 |
|
|---|
| 442 | // Restore variables.
|
|---|
| 443 |
|
|---|
| 444 | _restore_variables:
|
|---|
| 445 |
|
|---|
| 446 | u = svu;
|
|---|
| 447 | v = svv;
|
|---|
| 448 | for (i=0;i<n;i++)
|
|---|
| 449 | {
|
|---|
| 450 | k[i] = svk[i];
|
|---|
| 451 | }
|
|---|
| 452 |
|
|---|
| 453 | // Try quadratic iteration if it has not been tried
|
|---|
| 454 | // and the V sequence is converging.
|
|---|
| 455 | //
|
|---|
| 456 | if (vpass && !vtry) { goto _quadratic_iteration; }
|
|---|
| 457 |
|
|---|
| 458 | // Recompute QP and scalar values to continue the
|
|---|
| 459 | // second stage.
|
|---|
| 460 | //
|
|---|
| 461 | QuadraticSyntheticDivision(n,&u,&v,p,qp,&a,&b);
|
|---|
| 462 | ComputeScalarFactors(&type);
|
|---|
| 463 |
|
|---|
| 464 | ovv = vv;
|
|---|
| 465 | oss = ss;
|
|---|
| 466 | otv = tv;
|
|---|
| 467 | ots = ts;
|
|---|
| 468 | }
|
|---|
| 469 | }
|
|---|
| 470 |
|
|---|
| 471 | void G4JTPolynomialSolver::
|
|---|
| 472 | QuadraticPolynomialIteration(G4double *uu, G4double *vv, G4int *nz)
|
|---|
| 473 | {
|
|---|
| 474 | // Variable-shift k-polynomial iteration for a
|
|---|
| 475 | // quadratic factor converges only if the zeros are
|
|---|
| 476 | // equimodular or nearly so.
|
|---|
| 477 | // uu, vv - coefficients of starting quadratic.
|
|---|
| 478 | // nz - number of zeros found.
|
|---|
| 479 | //
|
|---|
| 480 | G4double ui=0.0, vi=0.0;
|
|---|
| 481 | G4double omp=0.0;
|
|---|
| 482 | G4double relstp=0.0;
|
|---|
| 483 | G4double mp=0.0, ee=0.0, t=0.0, zm=0.0;
|
|---|
| 484 | G4int type=0, i=1, j=0, tried=0;
|
|---|
| 485 |
|
|---|
| 486 | *nz = 0;
|
|---|
| 487 | tried = 0;
|
|---|
| 488 | u = *uu;
|
|---|
| 489 | v = *vv;
|
|---|
| 490 |
|
|---|
| 491 | // Main loop.
|
|---|
| 492 |
|
|---|
| 493 | while (1)
|
|---|
| 494 | {
|
|---|
| 495 | Quadratic(1.0,u,v,&szr,&szi,&lzr,&lzi);
|
|---|
| 496 |
|
|---|
| 497 | // Return if roots of the quadratic are real and not
|
|---|
| 498 | // close to multiple or nearly equal and of opposite
|
|---|
| 499 | // sign.
|
|---|
| 500 | //
|
|---|
| 501 | if (std::fabs(std::fabs(szr)-std::fabs(lzr)) > 0.01 * std::fabs(lzr))
|
|---|
| 502 | { return; }
|
|---|
| 503 |
|
|---|
| 504 | // Evaluate polynomial by quadratic synthetic division.
|
|---|
| 505 | //
|
|---|
| 506 | QuadraticSyntheticDivision(n,&u,&v,p,qp,&a,&b);
|
|---|
| 507 | mp = std::fabs(a-szr*b) + std::fabs(szi*b);
|
|---|
| 508 |
|
|---|
| 509 | // Compute a rigorous bound on the rounding error in evaluating p.
|
|---|
| 510 | //
|
|---|
| 511 | zm = std::sqrt(std::fabs(v));
|
|---|
| 512 | ee = 2.0*std::fabs(qp[0]);
|
|---|
| 513 | t = -szr*b;
|
|---|
| 514 | for (i=1;i<n;i++)
|
|---|
| 515 | {
|
|---|
| 516 | ee = ee*zm + std::fabs(qp[i]);
|
|---|
| 517 | }
|
|---|
| 518 | ee = ee*zm + std::fabs(a+t);
|
|---|
| 519 | ee *= (5.0 *mre + 4.0*are);
|
|---|
| 520 | ee = ee - (5.0*mre+2.0*are)*(std::fabs(a+t)+std::fabs(b)*zm)
|
|---|
| 521 | + 2.0*are*std::fabs(t);
|
|---|
| 522 |
|
|---|
| 523 | // Iteration has converged sufficiently if the
|
|---|
| 524 | // polynomial value is less than 20 times this bound.
|
|---|
| 525 | //
|
|---|
| 526 | if (mp <= 20.0*ee)
|
|---|
| 527 | {
|
|---|
| 528 | *nz = 2;
|
|---|
| 529 | return;
|
|---|
| 530 | }
|
|---|
| 531 | j++;
|
|---|
| 532 |
|
|---|
| 533 | // Stop iteration after 20 steps.
|
|---|
| 534 | //
|
|---|
| 535 | if (j > 20) { return; }
|
|---|
| 536 | if (j >= 2)
|
|---|
| 537 | {
|
|---|
| 538 | if (!(relstp > 0.01 || mp < omp || tried))
|
|---|
| 539 | {
|
|---|
| 540 | // A cluster appears to be stalling the convergence.
|
|---|
| 541 | // Five fixed shift steps are taken with a u,v close to the cluster.
|
|---|
| 542 | //
|
|---|
| 543 | if (relstp < eta) { relstp = eta; }
|
|---|
| 544 | relstp = std::sqrt(relstp);
|
|---|
| 545 | u = u - u*relstp;
|
|---|
| 546 | v = v + v*relstp;
|
|---|
| 547 | QuadraticSyntheticDivision(n,&u,&v,p,qp,&a,&b);
|
|---|
| 548 | for (i=0;i<5;i++)
|
|---|
| 549 | {
|
|---|
| 550 | ComputeScalarFactors(&type);
|
|---|
| 551 | ComputeNextPolynomial(&type);
|
|---|
| 552 | }
|
|---|
| 553 | tried = 1;
|
|---|
| 554 | j = 0;
|
|---|
| 555 | }
|
|---|
| 556 | }
|
|---|
| 557 | omp = mp;
|
|---|
| 558 |
|
|---|
| 559 | // Calculate next k polynomial and new u and v.
|
|---|
| 560 | //
|
|---|
| 561 | ComputeScalarFactors(&type);
|
|---|
| 562 | ComputeNextPolynomial(&type);
|
|---|
| 563 | ComputeScalarFactors(&type);
|
|---|
| 564 | ComputeNewEstimate(type,&ui,&vi);
|
|---|
| 565 |
|
|---|
| 566 | // If vi is zero the iteration is not converging.
|
|---|
| 567 | //
|
|---|
| 568 | if (!(vi != 0.0)) { return; }
|
|---|
| 569 | relstp = std::fabs((vi-v)/vi);
|
|---|
| 570 | u = ui;
|
|---|
| 571 | v = vi;
|
|---|
| 572 | }
|
|---|
| 573 | }
|
|---|
| 574 |
|
|---|
| 575 | void G4JTPolynomialSolver::
|
|---|
| 576 | RealPolynomialIteration(G4double *sss, G4int *nz, G4int *iflag)
|
|---|
| 577 | {
|
|---|
| 578 | // Variable-shift H polynomial iteration for a real zero.
|
|---|
| 579 | // sss - starting iterate
|
|---|
| 580 | // nz - number of zeros found
|
|---|
| 581 | // iflag - flag to indicate a pair of zeros near real axis.
|
|---|
| 582 |
|
|---|
| 583 | G4double t=0.;
|
|---|
| 584 | G4double omp=0.;
|
|---|
| 585 | G4double pv=0.0, kv=0.0, xs= *sss;
|
|---|
| 586 | G4double mx=0.0, mp=0.0, ee=0.0;
|
|---|
| 587 | G4int i=1, j=0;
|
|---|
| 588 |
|
|---|
| 589 | *nz = 0;
|
|---|
| 590 | *iflag = 0;
|
|---|
| 591 |
|
|---|
| 592 | // Main loop
|
|---|
| 593 | //
|
|---|
| 594 | while (1)
|
|---|
| 595 | {
|
|---|
| 596 | pv = p[0];
|
|---|
| 597 |
|
|---|
| 598 | // Evaluate p at xs.
|
|---|
| 599 | //
|
|---|
| 600 | qp[0] = pv;
|
|---|
| 601 | for (i=1;i<=n;i++)
|
|---|
| 602 | {
|
|---|
| 603 | pv = pv*xs + p[i];
|
|---|
| 604 | qp[i] = pv;
|
|---|
| 605 | }
|
|---|
| 606 | mp = std::fabs(pv);
|
|---|
| 607 |
|
|---|
| 608 | // Compute a rigorous bound on the error in evaluating p.
|
|---|
| 609 | //
|
|---|
| 610 | mx = std::fabs(xs);
|
|---|
| 611 | ee = (mre/(are+mre))*std::fabs(qp[0]);
|
|---|
| 612 | for (i=1;i<=n;i++)
|
|---|
| 613 | {
|
|---|
| 614 | ee = ee*mx + std::fabs(qp[i]);
|
|---|
| 615 | }
|
|---|
| 616 |
|
|---|
| 617 | // Iteration has converged sufficiently if the polynomial
|
|---|
| 618 | // value is less than 20 times this bound.
|
|---|
| 619 | //
|
|---|
| 620 | if (mp <= 20.0*((are+mre)*ee-mre*mp))
|
|---|
| 621 | {
|
|---|
| 622 | *nz = 1;
|
|---|
| 623 | szr = xs;
|
|---|
| 624 | szi = 0.0;
|
|---|
| 625 | return;
|
|---|
| 626 | }
|
|---|
| 627 | j++;
|
|---|
| 628 |
|
|---|
| 629 | // Stop iteration after 10 steps.
|
|---|
| 630 | //
|
|---|
| 631 | if (j > 10) { return; }
|
|---|
| 632 | if (j >= 2)
|
|---|
| 633 | {
|
|---|
| 634 | if (!(std::fabs(t) > 0.001*std::fabs(xs-t) || mp < omp))
|
|---|
| 635 | {
|
|---|
| 636 | // A cluster of zeros near the real axis has been encountered.
|
|---|
| 637 | // Return with iflag set to initiate a quadratic iteration.
|
|---|
| 638 | //
|
|---|
| 639 | *iflag = 1;
|
|---|
| 640 | *sss = xs;
|
|---|
| 641 | return;
|
|---|
| 642 | } // Return if the polynomial value has increased significantly.
|
|---|
| 643 | }
|
|---|
| 644 |
|
|---|
| 645 | omp = mp;
|
|---|
| 646 |
|
|---|
| 647 | // Compute t, the next polynomial, and the new iterate.
|
|---|
| 648 | //
|
|---|
| 649 | kv = k[0];
|
|---|
| 650 | qk[0] = kv;
|
|---|
| 651 | for (i=1;i<n;i++)
|
|---|
| 652 | {
|
|---|
| 653 | kv = kv*xs + k[i];
|
|---|
| 654 | qk[i] = kv;
|
|---|
| 655 | }
|
|---|
| 656 | if (std::fabs(kv) <= std::fabs(k[n-1])*10.0*eta) // Use unscaled form.
|
|---|
| 657 | {
|
|---|
| 658 | k[0] = 0.0;
|
|---|
| 659 | for (i=1;i<n;i++)
|
|---|
| 660 | {
|
|---|
| 661 | k[i] = qk[i-1];
|
|---|
| 662 | }
|
|---|
| 663 | }
|
|---|
| 664 | else // Use the scaled form of the recurrence if k at xs is nonzero.
|
|---|
| 665 | {
|
|---|
| 666 | t = -pv/kv;
|
|---|
| 667 | k[0] = qp[0];
|
|---|
| 668 | for (i=1;i<n;i++)
|
|---|
| 669 | {
|
|---|
| 670 | k[i] = t*qk[i-1] + qp[i];
|
|---|
| 671 | }
|
|---|
| 672 | }
|
|---|
| 673 | kv = k[0];
|
|---|
| 674 | for (i=1;i<n;i++)
|
|---|
| 675 | {
|
|---|
| 676 | kv = kv*xs + k[i];
|
|---|
| 677 | }
|
|---|
| 678 | t = 0.0;
|
|---|
| 679 | if (std::fabs(kv) > std::fabs(k[n-1]*10.0*eta)) { t = -pv/kv; }
|
|---|
| 680 | xs += t;
|
|---|
| 681 | }
|
|---|
| 682 | }
|
|---|
| 683 |
|
|---|
| 684 | void G4JTPolynomialSolver::ComputeScalarFactors(G4int *type)
|
|---|
| 685 | {
|
|---|
| 686 | // This function calculates scalar quantities used to
|
|---|
| 687 | // compute the next k polynomial and new estimates of
|
|---|
| 688 | // the quadratic coefficients.
|
|---|
| 689 | // type - integer variable set here indicating how the
|
|---|
| 690 | // calculations are normalized to avoid overflow.
|
|---|
| 691 |
|
|---|
| 692 | // Synthetic division of k by the quadratic 1,u,v
|
|---|
| 693 | //
|
|---|
| 694 | QuadraticSyntheticDivision(n-1,&u,&v,k,qk,&c,&d);
|
|---|
| 695 | if (std::fabs(c) <= std::fabs(k[n-1]*100.0*eta))
|
|---|
| 696 | {
|
|---|
| 697 | if (std::fabs(d) <= std::fabs(k[n-2]*100.0*eta))
|
|---|
| 698 | {
|
|---|
| 699 | *type = 3; // Type=3 indicates the quadratic is almost a factor of k.
|
|---|
| 700 | return;
|
|---|
| 701 | }
|
|---|
| 702 | }
|
|---|
| 703 |
|
|---|
| 704 | if (std::fabs(d) < std::fabs(c))
|
|---|
| 705 | {
|
|---|
| 706 | *type = 1; // Type=1 indicates that all formulas are divided by c.
|
|---|
| 707 | e = a/c;
|
|---|
| 708 | f = d/c;
|
|---|
| 709 | g = u*e;
|
|---|
| 710 | h = v*b;
|
|---|
| 711 | a3 = a*e + (h/c+g)*b;
|
|---|
| 712 | a1 = b - a*(d/c);
|
|---|
| 713 | a7 = a + g*d + h*f;
|
|---|
| 714 | return;
|
|---|
| 715 | }
|
|---|
| 716 | *type = 2; // Type=2 indicates that all formulas are divided by d.
|
|---|
| 717 | e = a/d;
|
|---|
| 718 | f = c/d;
|
|---|
| 719 | g = u*b;
|
|---|
| 720 | h = v*b;
|
|---|
| 721 | a3 = (a+g)*e + h*(b/d);
|
|---|
| 722 | a1 = b*f-a;
|
|---|
| 723 | a7 = (f+u)*a + h;
|
|---|
| 724 | }
|
|---|
| 725 |
|
|---|
| 726 | void G4JTPolynomialSolver::ComputeNextPolynomial(G4int *type)
|
|---|
| 727 | {
|
|---|
| 728 | // Computes the next k polynomials using scalars
|
|---|
| 729 | // computed in ComputeScalarFactors.
|
|---|
| 730 |
|
|---|
| 731 | G4int i=2;
|
|---|
| 732 |
|
|---|
| 733 | if (*type == 3) // Use unscaled form of the recurrence if type is 3.
|
|---|
| 734 | {
|
|---|
| 735 | k[0] = 0.0;
|
|---|
| 736 | k[1] = 0.0;
|
|---|
| 737 | for (i=2;i<n;i++)
|
|---|
| 738 | {
|
|---|
| 739 | k[i] = qk[i-2];
|
|---|
| 740 | }
|
|---|
| 741 | return;
|
|---|
| 742 | }
|
|---|
| 743 | G4double temp = a;
|
|---|
| 744 | if (*type == 1) { temp = b; }
|
|---|
| 745 | if (std::fabs(a1) <= std::fabs(temp)*eta*10.0)
|
|---|
| 746 | {
|
|---|
| 747 | // If a1 is nearly zero then use a special form of the recurrence.
|
|---|
| 748 | //
|
|---|
| 749 | k[0] = 0.0;
|
|---|
| 750 | k[1] = -a7*qp[0];
|
|---|
| 751 | for(i=2;i<n;i++)
|
|---|
| 752 | {
|
|---|
| 753 | k[i] = a3*qk[i-2] - a7*qp[i-1];
|
|---|
| 754 | }
|
|---|
| 755 | return;
|
|---|
| 756 | }
|
|---|
| 757 |
|
|---|
| 758 | // Use scaled form of the recurrence.
|
|---|
| 759 | //
|
|---|
| 760 | a7 /= a1;
|
|---|
| 761 | a3 /= a1;
|
|---|
| 762 | k[0] = qp[0];
|
|---|
| 763 | k[1] = qp[1] - a7*qp[0];
|
|---|
| 764 | for (i=2;i<n;i++)
|
|---|
| 765 | {
|
|---|
| 766 | k[i] = a3*qk[i-2] - a7*qp[i-1] + qp[i];
|
|---|
| 767 | }
|
|---|
| 768 | }
|
|---|
| 769 |
|
|---|
| 770 | void G4JTPolynomialSolver::
|
|---|
| 771 | ComputeNewEstimate(G4int type, G4double *uu, G4double *vv)
|
|---|
| 772 | {
|
|---|
| 773 | // Compute new estimates of the quadratic coefficients
|
|---|
| 774 | // using the scalars computed in calcsc.
|
|---|
| 775 |
|
|---|
| 776 | G4double a4=0.0, a5=0.0, b1=0.0, b2=0.0,
|
|---|
| 777 | c1=0.0, c2=0.0, c3=0.0, c4=0.0, temp=0.0;
|
|---|
| 778 |
|
|---|
| 779 | // Use formulas appropriate to setting of type.
|
|---|
| 780 | //
|
|---|
| 781 | if (type == 3) // If type=3 the quadratic is zeroed.
|
|---|
| 782 | {
|
|---|
| 783 | *uu = 0.0;
|
|---|
| 784 | *vv = 0.0;
|
|---|
| 785 | return;
|
|---|
| 786 | }
|
|---|
| 787 | if (type == 2)
|
|---|
| 788 | {
|
|---|
| 789 | a4 = (a+g)*f + h;
|
|---|
| 790 | a5 = (f+u)*c + v*d;
|
|---|
| 791 | }
|
|---|
| 792 | else
|
|---|
| 793 | {
|
|---|
| 794 | a4 = a + u*b +h*f;
|
|---|
| 795 | a5 = c + (u+v*f)*d;
|
|---|
| 796 | }
|
|---|
| 797 |
|
|---|
| 798 | // Evaluate new quadratic coefficients.
|
|---|
| 799 | //
|
|---|
| 800 | b1 = -k[n-1]/p[n];
|
|---|
| 801 | b2 = -(k[n-2]+b1*p[n-1])/p[n];
|
|---|
| 802 | c1 = v*b2*a1;
|
|---|
| 803 | c2 = b1*a7;
|
|---|
| 804 | c3 = b1*b1*a3;
|
|---|
| 805 | c4 = c1 - c2 - c3;
|
|---|
| 806 | temp = a5 + b1*a4 - c4;
|
|---|
| 807 | if (!(temp != 0.0))
|
|---|
| 808 | {
|
|---|
| 809 | *uu = 0.0;
|
|---|
| 810 | *vv = 0.0;
|
|---|
| 811 | return;
|
|---|
| 812 | }
|
|---|
| 813 | *uu = u - (u*(c3+c2)+v*(b1*a1+b2*a7))/temp;
|
|---|
| 814 | *vv = v*(1.0+c4/temp);
|
|---|
| 815 | return;
|
|---|
| 816 | }
|
|---|
| 817 |
|
|---|
| 818 | void G4JTPolynomialSolver::
|
|---|
| 819 | QuadraticSyntheticDivision(G4int nn, G4double *uu, G4double *vv,
|
|---|
| 820 | std::vector<G4double> &pp, std::vector<G4double> &qq,
|
|---|
| 821 | G4double *aa, G4double *bb)
|
|---|
| 822 | {
|
|---|
| 823 | // Divides pp by the quadratic 1,uu,vv placing the quotient
|
|---|
| 824 | // in qq and the remainder in aa,bb.
|
|---|
| 825 |
|
|---|
| 826 | G4double cc=0.0;
|
|---|
| 827 | *bb = pp[0];
|
|---|
| 828 | qq[0] = *bb;
|
|---|
| 829 | *aa = pp[1] - (*bb)*(*uu);
|
|---|
| 830 | qq[1] = *aa;
|
|---|
| 831 | for (G4int i=2;i<=nn;i++)
|
|---|
| 832 | {
|
|---|
| 833 | cc = pp[i] - (*aa)*(*uu) - (*bb)*(*vv);
|
|---|
| 834 | qq[i] = cc;
|
|---|
| 835 | *bb = *aa;
|
|---|
| 836 | *aa = cc;
|
|---|
| 837 | }
|
|---|
| 838 | }
|
|---|
| 839 |
|
|---|
| 840 | void G4JTPolynomialSolver::Quadratic(G4double aa,G4double b1,
|
|---|
| 841 | G4double cc,G4double *ssr,G4double *ssi,
|
|---|
| 842 | G4double *lr,G4double *li)
|
|---|
| 843 | {
|
|---|
| 844 |
|
|---|
| 845 | // Calculate the zeros of the quadratic aa*z^2 + b1*z + cc.
|
|---|
| 846 | // The quadratic formula, modified to avoid overflow, is used
|
|---|
| 847 | // to find the larger zero if the zeros are real and both
|
|---|
| 848 | // are complex. The smaller real zero is found directly from
|
|---|
| 849 | // the product of the zeros c/a.
|
|---|
| 850 |
|
|---|
| 851 | G4double bb=0.0, dd=0.0, ee=0.0;
|
|---|
| 852 |
|
|---|
| 853 | if (!(aa != 0.0)) // less than two roots
|
|---|
| 854 | {
|
|---|
| 855 | if (b1 != 0.0)
|
|---|
| 856 | { *ssr = -cc/b1; }
|
|---|
| 857 | else
|
|---|
| 858 | { *ssr = 0.0; }
|
|---|
| 859 | *lr = 0.0;
|
|---|
| 860 | *ssi = 0.0;
|
|---|
| 861 | *li = 0.0;
|
|---|
| 862 | return;
|
|---|
| 863 | }
|
|---|
| 864 | if (!(cc != 0.0)) // one real root, one zero root
|
|---|
| 865 | {
|
|---|
| 866 | *ssr = 0.0;
|
|---|
| 867 | *lr = -b1/aa;
|
|---|
| 868 | *ssi = 0.0;
|
|---|
| 869 | *li = 0.0;
|
|---|
| 870 | return;
|
|---|
| 871 | }
|
|---|
| 872 |
|
|---|
| 873 | // Compute discriminant avoiding overflow.
|
|---|
| 874 | //
|
|---|
| 875 | bb = b1/2.0;
|
|---|
| 876 | if (std::fabs(bb) < std::fabs(cc))
|
|---|
| 877 | {
|
|---|
| 878 | if (cc < 0.0)
|
|---|
| 879 | { ee = -aa; }
|
|---|
| 880 | else
|
|---|
| 881 | { ee = aa; }
|
|---|
| 882 | ee = bb*(bb/std::fabs(cc)) - ee;
|
|---|
| 883 | dd = std::sqrt(std::fabs(ee))*std::sqrt(std::fabs(cc));
|
|---|
| 884 | }
|
|---|
| 885 | else
|
|---|
| 886 | {
|
|---|
| 887 | ee = 1.0 - (aa/bb)*(cc/bb);
|
|---|
| 888 | dd = std::sqrt(std::fabs(ee))*std::fabs(bb);
|
|---|
| 889 | }
|
|---|
| 890 | if (ee < 0.0) // complex conjugate zeros
|
|---|
| 891 | {
|
|---|
| 892 | *ssr = -bb/aa;
|
|---|
| 893 | *lr = *ssr;
|
|---|
| 894 | *ssi = std::fabs(dd/aa);
|
|---|
| 895 | *li = -(*ssi);
|
|---|
| 896 | }
|
|---|
| 897 | else
|
|---|
| 898 | {
|
|---|
| 899 | if (bb >= 0.0) // real zeros.
|
|---|
| 900 | { dd = -dd; }
|
|---|
| 901 | *lr = (-bb+dd)/aa;
|
|---|
| 902 | *ssr = 0.0;
|
|---|
| 903 | if (*lr != 0.0)
|
|---|
| 904 | { *ssr = (cc/ *lr)/aa; }
|
|---|
| 905 | *ssi = 0.0;
|
|---|
| 906 | *li = 0.0;
|
|---|
| 907 | }
|
|---|
| 908 | }
|
|---|