File Coverage

lib/PDL/Math/cpoly.c
Criterion Covered Total %
statement 214 247 86.6
branch 114 166 68.6
condition n/a
subroutine n/a
pod n/a
total 328 413 79.4


line stmt bran cond sub pod time code
1             /* Translated from F77 to C, rjrw 10/04/2000 */
2             /* replaced 'bool' by 'boolvar' to get it to compile on my
3             linux machine, DJB Aug 02 2000 */
4              
5             /* algorithm 419 collected algorithms from acm.
6             algorithm appeared in comm. acm, vol. 15, no. 02, p. 097. */
7             /* available in 2024 from https://calgo.acm.org/
8             see https://en.wikipedia.org/wiki/Jenkins%E2%80%93Traub_algorithm */
9              
10             #include
11             #include
12             #include
13             #include
14             /*
15             #if !defined(WIN32) && !defined(_WIN32) && !defined(__APPLE__) && !defined(__CYGWIN__)
16             #include
17             #endif
18             */
19             #include
20             /* #define DEBUGMAIN */ /* Set up debugging main, etc. */
21             #include "cpoly.h"
22              
23             /* Internal routines */
24             static complex double noshft(int l1, int nn, complex double tc, complex double hc[], complex double pc[]);
25             static int fxshft(int l2, int nn, complex double shc[], complex double qpc[], complex double hc[], complex double pc[], complex double qhc[], complex double *tc, complex double *sc, complex double *pvc, complex double *zc);
26             static int vrshft(int l3, int nn, complex double qpc[], complex double pc[], complex double qhc[], complex double hc[], complex double *tc, complex double *sc, complex double *pvc, complex double *zc);
27             static int calct(int nn, complex double sc, complex double pvc, complex double qhc[], complex double hc[], complex double *tc);
28             static void nexth(int boolvar, int nn, complex double tc, complex double qhc[], complex double qpc[], complex double hc[]);
29             static complex double polyev(int nn, complex double sc, complex double pc[], complex double qc[]);
30             static double errev(int nn, double ms, double mp, complex double qc[]);
31             static double cauchy(int nn, complex double pc[]);
32             static double scale(int nn, complex double pc[]);
33             static complex double cdivid(complex double a, complex double b);
34             static double cmod(complex double a);
35             static void mcon(void);
36             static void init(int nncr);
37              
38             /* Internal global variables */
39             static double are,mre,eta,infin,smalno,base;
40              
41             #ifdef DEBUGMAIN
42             /* driver to test cpoly */
43             int main()
44             {
45             char *fail = NULL;
46             double p[50],pi[50],zr[50],zi[50];
47              
48             int i;
49              
50             printf("Example 1. polynomial with zeros 1,2,...,10.\n");
51             p[0]=1L;
52             p[1]=-55L;
53             p[2]=1320L;
54             p[3]=-18150L;
55             p[4]=157773L;
56             p[5]=-902055L;
57             p[6] = 3416930L;
58             p[7]=-8409500L;
59             p[8]=12753576L;
60             p[9]=-10628640L;
61             p[10]=3628800L;
62             for (i=0;i<11;i++)
63             pi[i]=0;
64             prtc(11,p,pi);
65             fail = cpoly(p,pi,10,zr,zi);
66             if(fail)
67             printf("cpoly has failed on this example\n");
68             prtz (10,zr,zi);
69             printf("Example 2. zeros on imaginary axis degree 3.\n");
70             p[0]=1;
71             p[1]=0;
72             p[2]=-10001.0001L;
73             p[3]=0;
74             pi[0]=0;
75             pi[1]=-10001.0001L;
76             pi[2]=0;
77             pi[3]=1;
78             prtc(4,p,pi);
79             fail = cpoly(p,pi,3,zr,zi);
80             if (fail)
81             printf("cpoly has failed on this example\n");
82             prtz(3,zr,zi);
83             printf("Example 3. zeros at 1+i,1/2*(1+i)....1/(2**-9)*(1+i)\n");
84             p[0]=1.0;
85             p[1]=-1.998046875L;
86             p[2]=0.0;
87             p[3]=.7567065954208374L;
88             p[4]=-.2002119533717632L;
89             p[5]=1.271507365163416e-2L;
90             p[6]=0;
91             p[7]=-1.154642632172909e-5L;
92             p[8]=1.584803612786345e-7L;
93             p[9]=-4.652065399568528e-10L;
94             p[10]=0;
95             pi[0]=0;
96             pi[1]=p[1];
97             pi[2]=2.658859252929688L;
98             pi[3]=-7.567065954208374e-1L;
99             pi[4]=0;
100             pi[5]=p[5];
101             pi[6]=-7.820779428584501e-4L;
102             pi[7]=-p[7];
103             pi[8]=0;
104             pi[9]=p[9];
105             pi[10]=9.094947017729282e-13L;
106             prtc(11,p,pi);
107             fail = cpoly(p,pi,10,zr,zi);
108             if (fail)
109             printf("cpoly has failed on this example\n");
110             prtz(10,zr,zi);
111             printf("Example 4. multiple zeros\n");
112             p[0]=1L;
113             p[1]=-10L;
114             p[2]=3L;
115             p[3]=284L;
116             p[4]=-1293L;
117             p[5]=2374L;
118             p[6]=-1587L;
119             p[7]=-920L;
120             p[8]=2204L;
121             p[9]=-1344L;
122             p[10]=288L;
123             pi[0]=0;
124             pi[1]=-10L;
125             pi[2]=100L;
126             pi[3]=-334L;
127             pi[4]=200L;
128             pi[5]=1394L;
129             pi[6] =-3836L;
130             pi[7]=4334L;
131             pi[8]=-2352L;
132             pi[9]=504L;
133             pi[10]=0;
134             prtc(11,p,pi);
135             fail = cpoly(p,pi,10,zr,zi);
136             if (fail)
137             printf("cpoly has failed on this example\n");
138             prtz(10,zr,zi);
139             printf("Example 5. 12 zeros evenly distributed on a circle of radius 1. centered at 0+2i.\n");
140             p[0]=1L;
141             p[1]=0;
142             p[2]=-264L;
143             p[3]=0;
144             p[4]=7920L;
145             p[5]=0;
146             p[6]=-59136L;
147             p[7]=0;
148             p[8]=126720L;
149             p[9]=0;
150             p[10]=-67584L;
151             p[11]=0;
152             p[12]=4095L;
153             pi[0]=0;
154             pi[1]=-24L;
155             pi[2]=0;
156             pi[3]=1760L;
157             pi[4]=0;
158             pi[5]=-25344L;
159             pi[6]=0;
160             pi[7]=101376L;
161             pi[8]=0;
162             pi[9]=-112640L;
163             pi[10]=0;
164             pi[11]=24576L;
165             pi[12]=0;
166             prtc(13,p,pi);
167             fail = cpoly(p,pi,12,zr,zi);
168             if(fail)
169             printf("cpoly has failed on this example\n");
170             prtz(12,zr,zi);
171             return 0;
172             }
173             void prtc(int n, double p[], double q[])
174             {
175             int i;
176             printf("Coefficients\n");
177             for (i=0;i
178             printf("%26.16g %26.16g\n",p[i],q[i]);
179             }
180              
181             void prtz(int n,double zr[], double zi[])
182             {
183             int i;
184             printf("Zeroes\n");
185             for (i=0;i
186             printf("%26.16g %26.16g\n",zr[i],zi[i]);
187             }
188             #endif
189              
190             #define COSR (-0.0697564737441253008L)
191             #define SINR (0.997564050259824248L)
192              
193 11           char *cpoly(double opr[], double opi[], int degree,
194             double zeror[], double zeroi[])
195             {
196             /* Finds the zeros of a complex polynomial.
197              
198             opr, opi - double precision vectors of real and imaginary parts
199             of the coefficients in order of decreasing powers.
200             degree - integer degree of polynomial
201             zeror, zeroi
202             - output double precision vectors of real and imaginary
203             parts of the zeros.
204             fail - output logical parameter, TRUE if leading coefficient
205             is zero, if cpoly has found fewer than degree zeros,
206             or if there is another internal error.
207              
208             The program has been written to reduce the chance of overflow
209             occurring. If it does occur, there is still a possibility that
210             the zerofinder will work provided the overflowed quantity is
211             replaced by a large number. */
212              
213 11 100         if (degree < 1)
214 1           return "algorithm only works for degree >= 1.";
215              
216 10 50         if (opr[0] == 0.0 && opi[0] == 0.0)
    0          
217 0           return "algorithm fails if the leading coefficient is zero.";
218              
219 10           double xx=0.0,yy=0.0,xxx=0.0,bnd=0.0;
220 10           complex double zc=0.0,tc=0.0,sc=0.0,pvc=0.0;
221 10           char *failreason = NULL;
222             int cnt1,cnt2,i,idnn2;
223              
224             /* initialization of constants */
225 10           int nn = degree+1;
226 10           init(nn);
227 10           complex double *pc = malloc(nn*sizeof(complex double));
228 10           complex double *hc = malloc(nn*sizeof(complex double));
229 10           complex double *qpc = malloc(nn*sizeof(complex double));
230 10           complex double *qhc = malloc(nn*sizeof(complex double));
231 10           complex double *shc = malloc(nn*sizeof(complex double));
232              
233 10 50         if (!(pc && hc && qpc && qhc && shc)) {
    50          
    50          
    50          
    50          
234 0           failreason = "Couldn't allocate space for cpoly";
235 0           goto returnlab;
236             }
237              
238 10           xx = 0.70710678118654752438L;
239 10           yy = -xx;
240              
241             /* Remove the zeros at the origin if any */
242 10 50         while (opr[nn-1] == 0.0 && opi[nn-1] == 0.0) {
    0          
243 0           idnn2 = degree+1-nn;
244 0           zeror[idnn2] = 0.0;
245 0           zeroi[idnn2] = 0.0;
246 0           nn--;
247             }
248              
249             /* Make a copy of the coefficients */
250 74 100         for (i=0;i
251 64           pc[i] = opr[i] + I*opi[i];
252 64           shc[i] = cmod(pc[i]);
253             }
254              
255             /* Scale the polynomial */
256 10           bnd = scale(nn,shc);
257 10 50         if (bnd != 1.0) {
258 0 0         for (i=0;i
259 0           pc[i] *= bnd;
260             }
261             }
262              
263 54 50         while (!failreason) {
264              
265             /* Start the algorithm for one zero */
266 54 100         if (nn < 3) {
267             /* Calculate the final zero and return */
268 10           complex double zeroc = cdivid(-pc[1], pc[0]);
269 10           zeror[degree-1] = creal(zeroc); zeroi[degree-1] = cimag(zeroc);
270 10           goto returnlab;
271             }
272              
273             /* Calculate bnd, a lower bound on the modulus of the zeros */
274 327 100         for (i=0;i
275 283           shc[i] = cmod(pc[i]) + I*cimag(shc[i]);
276             }
277 44           bnd = cauchy(nn,shc);
278              
279             /* Outer loop to control 2 major passes with different sequences
280             of shifts */
281 44           failreason = "found fewer than degree zeros";
282 88 100         for(cnt1=1;failreason && (cnt1<=2);cnt1++) {
    50          
283              
284             /* First stage calculation, no shift */
285 44           tc = noshft(5,nn,tc,hc,pc);
286 44 50         if (isnan(creal(tc)) || isnan(cimag(tc)))
    50          
287 0           return "noshft returned NaN";
288              
289             /* Inner loop to select a shift. */
290 88 100         for (cnt2=1;failreason && (cnt2<10);cnt2++) {
    50          
291             /* Shift is chosen with modulus bnd and amplitude rotated by
292             94 degrees from the previous shift */
293 44           xxx = COSR*xx-SINR*yy;
294 44           yy = SINR*xx+COSR*yy;
295 44           xx = xxx;
296 44           sc = bnd*xx + I*bnd*yy;
297              
298             /* Second stage calculation, fixed shift */
299 44 50         if (fxshft(10*cnt2,nn,shc,qpc,hc,pc,qhc,&tc,&sc,&pvc,&zc)) {
300              
301             /* The second stage jumps directly to the third stage iteration
302             If successful the zero is stored and the polynomial deflated */
303 44           idnn2 = degree+1-nn;
304 44           zeror[idnn2] = creal(zc);
305 44           zeroi[idnn2] = cimag(zc);
306 44           nn--;
307 283 100         for(i=0;i
308 239           pc[i] = qpc[i];
309 44           failreason = NULL;
310             }
311             /* If the iteration is unsuccessful another shift is chosen */
312             }
313             /* If 9 shifts fail, the outer loop is repeated with another
314             sequence of shifts */
315             }
316             }
317              
318 0           returnlab:
319 10 50         if (shc) free(shc);
320 10 50         if (qhc) free(qhc);
321 10 50         if (qpc) free(qpc);
322 10 50         if (hc) free(hc);
323 10 50         if (pc) free(pc);
324             /* The zerofinder has failed on two major passes
325             Return empty handed */
326 10           return failreason;
327             }
328              
329 44           static complex double noshft(int l1, int nn, complex double tc, complex double hc[], complex double pc[])
330             {
331             /* Computes the derivative polynomial as the initial h
332             polynomial and computes l1 no-shift h polynomials. */
333 44           int i,jj,n = nn-1,nm1 = n-1;
334 283 100         for (i=0;i
335 239           hc[i] = (n-i)*pc[i] / n;
336 264 100         for (jj=0;jj
337 220 100         if (cmod(hc[nm1]) > eta*10.0*cmod(pc[nm1])) {
338 214           tc = cdivid(-pc[n], hc[nm1]);
339 214 50         if (isnan(creal(tc)) || isnan(cimag(tc)))
    50          
340 0           return tc; /* error, stop now */
341 1174 100         for (i=0;i
342 960           int j = nm1-i;
343 960           hc[j] = pc[j] + tc * hc[j-1];
344             }
345 214           hc[0] = pc[0];
346             } else {
347             /* If the constant term is essentially zero, shift h coefficients */
348 21 100         for (i=0;i
349 15           int j = nm1-i;
350 15           hc[j] = hc[j-1];
351             }
352 6           hc[0] = 0.0;
353             }
354             }
355 44           return tc;
356             }
357              
358 44           static int fxshft(int l2, int nn, complex double shc[], complex double qpc[], complex double hc[], complex double pc[], complex double qhc[], complex double *tc, complex double *sc, complex double *pvc, complex double *zc)
359             /* Computes l2 fixed-shift h polynomials and tests for convergence
360              
361             Initiates a variable-shift iteration and returns with the
362             approximate zero if successful.
363              
364             l2 - Limit of fixed shift steps
365              
366             output:
367             zc - Approximate zero if true return
368             returns true if convergence of stage 3 iteration
369             */
370             {
371             int test,pasd,boolvar;
372 44           int i,j,n = nn-1;
373              
374             /* Evaluate p at s */
375 44           *pvc = polyev(nn,*sc,pc,qpc);
376 44           test = TRUE;
377 44           pasd = FALSE;
378              
379             /* Calculate first t = -p(s)/h(s) */
380 44           boolvar = calct(nn,*sc,*pvc,qhc,hc,tc);
381              
382             /* Main loop for one second stage step */
383 103 100         for (j=0;j
384 102           complex double otc = *tc;
385              
386             /* Compute next h polynomial and new t */
387 102           nexth(boolvar,nn, *tc, qhc, qpc, hc);
388 102           boolvar = calct(nn,*sc,*pvc,qhc,hc,tc);
389 102           *zc = *sc+*tc;
390              
391             /* Test for convergence unless stage 3 has failed once or
392             this is the last h polynomial */
393 102 50         if (boolvar || !test || j == l2-1) continue;
    50          
    100          
394 101 100         if (cmod(*tc-otc) >= .5*cmod(*zc)) { pasd = FALSE; continue; }
395 87 100         if (!pasd) { pasd = TRUE; continue; }
396             /* The weak convergence test has been passed twice, start the
397             third stage iteration, after saving the current h polynomial
398             and shift */
399 278 100         for (i=0;i
400 235           shc[i] = hc[i];
401             }
402 43           complex double svsc = *sc;
403 43 50         if (vrshft(10,nn,qpc,pc,qhc,hc,tc,sc,pvc,zc))
404 43           return TRUE;
405              
406             /* The iteration failed to converge
407             Turn off testing and restore h,s,pv and t */
408 0           test = FALSE;
409 0 0         for (i=0;i
410 0           hc[i] = shc[i];
411             }
412 0           *sc = svsc;
413 0           *pvc = polyev(nn,*sc,pc,qpc);
414 0           boolvar = calct(nn,*sc,*pvc,qhc,hc,tc);
415             }
416              
417             /* Attempt an iteration with final h polynomial from second stage */
418 1           return vrshft(10,nn,qpc,pc,qhc,hc,tc,sc,pvc,zc);
419             }
420              
421 44           static int vrshft(int l3, int nn, complex double qpc[], complex double pc[], complex double qhc[], complex double hc[], complex double *tc, complex double *sc, complex double *pvc, complex double *zc)
422             /* Carries out the third stage iteration
423              
424             l3 - Limit of steps in stage 3
425             zc - On entry contains the initial iterate,
426             On exit, it contains the final iterate (if it converges).
427             returns TRUE if iteration converges
428             */
429             {
430 44           double omp = 0,relstp = 0;
431             int i,j,boolvar;
432 44           int b = FALSE;
433 44           *sc = *zc;
434              
435             /* Main loop for stage three */
436 182 50         for (i=0; i
437              
438             /* Evaluate p at s and test for convergence */
439 182           *pvc = polyev(nn,*sc,pc,qpc);
440 182           double mp = cmod(*pvc), ms = cmod(*sc);
441 182 100         if (mp <= 20.0L*errev(nn,ms,mp,qpc)) {
442             /* Polynomial value is smaller in value than a bound on the error
443             in evaluating p, terminate the iteration */
444 44           *zc = *sc;
445 44           return TRUE;
446             } else {
447 138 100         if (i!=0) {
448 97 100         if (!b && mp>=omp && relstp < .05L) {
    100          
    50          
449             /* Iteration has stalled, probably a cluster of zeros
450             Do 5 fixed shift steps into the cluster to force one zero
451             to dominate */
452 4           b = TRUE;
453 4 50         double tp = (relstp < eta) ? eta : relstp;
454 4           *sc *= 1.0L + sqrt(tp)*(1 + I);
455 4           *pvc = polyev(nn,*sc,pc,qpc);
456 24 100         for (j=0;j<5;j++) {
457 20           boolvar = calct(nn,*sc,*pvc,qhc,hc,tc);
458 20           nexth(boolvar,nn, *tc, qhc, qpc, hc);
459             }
460 4           omp = infin;
461             } else {
462             /* Exit if polynomial value increases significantly */
463 93 50         if (mp*0.1L > omp)
464 0           return FALSE;
465 93           omp = mp;
466             }
467             } else {
468 41           omp = mp;
469             }
470             }
471              
472             /* Calculate next iterate. */
473 138           boolvar = calct(nn,*sc,*pvc,qhc,hc,tc);
474 138           nexth(boolvar,nn, *tc, qhc, qpc, hc);
475 138           boolvar = calct(nn,*sc,*pvc,qhc,hc,tc);
476 138 50         if (!boolvar) {
477 138           relstp = cmod(*tc)/cmod(*sc);
478 138           *sc += *tc;
479             }
480             }
481 0           return FALSE;
482             }
483              
484 442           static int calct(int nn, complex double sc, complex double pvc, complex double qhc[], complex double hc[], complex double *tc)
485             /* Computes t = -p(s)/h(s)
486             Returns TRUE if h(s) is essentially zero
487             */
488             {
489 442           complex double hvc = polyev(nn-1,sc,hc,qhc); /* Evaluate h(s) */
490 442           int boolvar = (cmod(hvc) <= are*10.0*cmod(hc[nn-2]));
491 442 50         *tc = boolvar ? 0.0 : cdivid(-pvc, hvc);
492 442           return boolvar;
493             }
494              
495 260           static void nexth(int boolvar, int nn, complex double tc, complex double qhc[], complex double qpc[], complex double hc[])
496             /* Calculates the next shifted h polynomial
497             boolvar - TRUE if h(s) is essentially zero
498             */
499             {
500 260           int j,n = nn-1;
501              
502 260 50         if (!boolvar) {
503 1354 100         for (j=1;j
504 1094           hc[j] = qpc[j] + tc*qhc[j-1];
505             }
506 260           hc[0] = qpc[0];
507             } else {
508             /* If h(s) is zero, replace h with qh */
509 0 0         for (j=1;j
510 0           hc[j] = qhc[j-1];
511             }
512 0           hc[0] = 0.0;
513             }
514 260           }
515              
516 672           static complex double polyev(int nn, complex double sc, complex double pc[],
517             complex double qc[])
518             /* Evaluates a polynomial p at s by the Horner recurrence,
519             placing the partial sums in q, returns the computed value
520             */
521             {
522             int i;
523 672           complex double vc = pc[0];
524 672           qc[0] = vc;
525 3783 100         for (i=1;i
526 3111           qc[i] = vc = vc*sc + pc[i];
527 672           return vc;
528             }
529              
530 182           static double errev(int nn, double ms, double mp, complex double qc[])
531             /* Bounds the error in evaluating the polynomial by the Horner recurrence
532            
533             qr,qi - The partial sums
534             ms - Modulus of the point
535             mp - Modulus of polynomial value
536             */
537             {
538             double e;
539             int i;
540              
541 182           e = cmod(qc[0])*mre/(are+mre);
542 1334 100         for (i=0;i
543 1152           e = e*ms+cmod(qc[i]);
544 182           return e*(are+mre)-mp*mre;
545             }
546              
547 44           static double cauchy(int nn, complex double pc[])
548             /* Cauchy computes a lower bound on the moduli of the zeros of a
549             polynomial - the real component of pc is the modulus of the coefficients
550             */
551             {
552             double x,xm,f,dx,df;
553 44           int n=nn-1, nm=nn-2, i;
554              
555 44           pc[n] = -creal(pc[n]) + I*cimag(pc[n]);
556              
557             /* Compute upper estimate of bound */
558 44           x = exp( (log(-creal(pc[n])) - log(creal(pc[0])))/n );
559 44 100         if (creal(pc[nm]) != 0.0) {
560             /* If Newton step at the origin is better, use it */
561 42           xm = -creal(pc[n])/creal(pc[nm]);
562 42 100         if (xm < x)
563 41           x = xm;
564             }
565              
566             /* Chop the interval (0,x) until f <= 0 */
567             while (1) {
568 44           xm = x*.1;
569 44           f = creal(pc[0]);
570 283 100         for (i=1;i
571 239           f = f*xm+creal(pc[i]);
572 44 50         if (f <= 0.) break;
573 0           x = xm;
574             }
575 44           dx = x;
576              
577             /* Do Newton iteration until x converges to two decimal places */
578 175 100         while (fabs(dx/x) > .005L) {
579 131           pc[0] = creal(pc[0]) + I*creal(pc[0]);
580 843 100         for(i=1;i
581 712           pc[i] = creal(pc[i]) + I*(cimag(pc[i-1])*x+creal(pc[i]));
582 131           f = cimag(pc[n]);
583 131           df = cimag(pc[0]);
584 712 100         for (i=1;i
585 581           df = df*x+cimag(pc[i]);
586 131           dx = f/df;
587 131           x -= dx;
588             }
589 44           return x;
590             }
591              
592 10           static double scale(int nn, complex double pc[])
593             /* Returns a scale factor to multiply the coefficients of the
594             polynomial. The scaling is done to avoid overflow and to avoid
595             undetected underflow interfering with the convergence
596             criterion. The factor is a power of the base.
597            
598             pc - the real components are modulus of coefficients of p
599             */
600             {
601             double hi,lo,max,min,x,sc;
602             int i,l;
603              
604             /* Find largest and smallest moduli of coefficients */
605 10           hi = sqrt(infin);
606 10           lo = smalno/eta;
607 10           max = 0.0;
608 10           min = infin;
609 74 100         for (i=0;i
610 64           x = creal(pc[i]);
611 64 100         if (x > max)
612 46           max = x;
613 64 100         if (x != 0.0 && x < min)
    100          
614 11           min = x;
615             }
616              
617             /* Scale only if there are very large or very small components */
618 10 50         if (min >= lo && max <= hi)
    50          
619 10           return 1.0;
620 0           x = lo/min;
621 0 0         if (x <= 1.0L) {
622 0           sc = 1.0L/(sqrt(max)*sqrt(min));
623             } else {
624 0           sc = x;
625 0 0         if (infin/sc > max)
626 0           sc = 1.0;
627             }
628 0           l = log(sc)/log(base) + .500;
629 0           return pow(base,l);
630             }
631              
632 666           static complex double cdivid(complex double a, complex double b)
633             /* Complex division c = a/b, avoiding overflow */
634             {
635 666           double ar=creal(a),ai=cimag(a),br=creal(b),bi=cimag(b);
636 666 50         if (br == 0.0 && bi == 0.0) {
    0          
637             /* division by zero, c = infinity. */
638 0           return infin*(1 + I);
639 666 100         } else if (fabs(br) < fabs(bi)) {
640 73           double r = br/bi;
641 73           double d = bi+r*br;
642 73           return (ar*r+ai + I*(ai*r-ar))/d;
643             } else {
644 593           double r = bi/br;
645 593           double d = br+r*bi;
646 593 50         if (isinf(d)) return 0 + 0*I; /* clang 3.4.1 has 1/inf = NaN so we dodge */
647 593           return (ar+ai*r + I*(ai-ar*r))/d;
648             }
649             }
650              
651 3847           static double cmod(complex double a)
652             /* Modulus of a complex number avoiding overflow */
653             {
654 3847           double ar = fabs(creal(a)), ai = fabs(cimag(a)), f;
655 3847 100         if (ar < ai) {
656 368           f = ar/ai;
657 368           return ai*sqrt(1.0+f*f);
658 3479 100         } else if (ar > ai) {
659 3443           f = ai/ar;
660 3443           return ar*sqrt(1.0+f*f);
661             } else {
662 36           return ar*sqrt(2.0);
663             }
664             }
665              
666 3           static void mcon(void)
667             /* mcon provides machine constants used in various parts of the
668             program. The user may either set them directly or use the
669             statements below to compute them. The meaning of the four
670             constants are -
671            
672             eta the maximum relative representation error
673             which can be described as the smallest positive
674             floating-point number such that 1.0d0 + eta is
675             greater than 1.0d0.
676             infin the largest floating-point number
677             smalno the smallest positive floating-point number
678             base the base of the floating-point number system used
679              
680             Let t be the number of base-digits in each floating-point
681             number (double precision). Then eta is either .5*b**(1-t)
682             or b**(1-t) depending on whether rounding or truncation
683             is used.
684              
685             Let m be the largest exponent and n the smallest exponent
686             in the number system. Then infiny is (1-base**(-t))*base**m
687             and smalno is base**n.
688             */
689             {
690              
691             /*
692             #if !defined(WIN32) && !defined(_WIN32) && !defined(__APPLE__) && !defined(__CYGWIN__)
693             base = 2;
694             eta = DBL_EPSILON;
695             smalno = MINDOUBLE;
696             infin = MAXDOUBLE;
697             #else
698             */
699 3           base = 2;
700 3           eta = DBL_EPSILON;
701 3           smalno = DBL_MIN;
702 3           infin = DBL_MAX;
703             /* #endif */
704              
705             #ifdef IBM360
706             /* These values for base,t,m,n correspond to the ibm/360. */
707             int m,n,t;
708             base = 16.0;
709             t = 14;
710             m = 63;
711             n = -65;
712             eta = pow(base,1-t);
713             infin = (base)*(1.0-pow(base,-t))*pow(base,m-1);
714             smalno = pow(base,n+3)/pow(base,3);
715             #endif
716 3           }
717              
718 10           static void init(int nncr)
719             {
720             static int nmax=0;
721              
722 10 100         if (nmax == 0) {
723             /* Set up once-off constants */
724 3           mcon();
725              
726             /* are, mre - Error bounds on complex addition and multiplication,
727             cf e.g. errev() above */
728 3           are = eta;
729 3           mre = 2.0L*sqrt(2.0L)*eta;
730 3           nmax = 1;
731             }
732 10           }