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