Image Component Library (ICL)
PolynomialSolver.h
Go to the documentation of this file.
1 /********************************************************************
2 ** Image Component Library (ICL) **
3 ** **
4 ** Copyright (C) 2006-2013 CITEC, University of Bielefeld **
5 ** Neuroinformatics Group **
6 ** Website: www.iclcv.org and **
7 ** http://opensource.cit-ec.de/projects/icl **
8 ** **
9 ** File : ICLUtils/src/ICLMath/PolynomialSolver.h **
10 ** Module : ICLMath **
11 ** Authors: Sergius Gaulik **
12 ** **
13 ** **
14 ** GNU LESSER GENERAL PUBLIC LICENSE **
15 ** This file may be used under the terms of the GNU Lesser General **
16 ** Public License version 3.0 as published by the **
17 ** **
18 ** Free Software Foundation and appearing in the file LICENSE.LGPL **
19 ** included in the packaging of this file. Please review the **
20 ** following information to ensure the license requirements will **
21 ** be met: http://www.gnu.org/licenses/lgpl-3.0.txt **
22 ** **
23 ** The development of this software was supported by the **
24 ** Excellence Cluster EXC 277 Cognitive Interaction Technology. **
25 ** The Excellence Cluster EXC 277 is a grant of the Deutsche **
26 ** Forschungsgemeinschaft (DFG) in the context of the German **
27 ** Excellence Initiative. **
28 ** **
29 ********************************************************************/
30 
31 #pragma once
32 
33 #include <ICLUtils/BasicTypes.h>
34 #include <float.h>
35 
36 #if (defined _MSC_VER && _MSC_VER < 1800)
37 #pragma WARNING("This compiler does not support the functions ilogb and scalbln.")
38 #else
39 
40 namespace icl{
41  namespace math{
42 
43  typedef std::complex<double> xcomplex;
44 
45  /****************************************************************************
46  *
47  * cpoly.C René Hartung
48  * Laurent Bartholdi
49  *
50  * @(#)$id: fr_dll.c,v 1.18 2010/10/26 05:19:40 gap exp $
51  *
52  * Copyright (c) 2011, Laurent Bartholdi
53  *
54  ****************************************************************************
55  *
56  * template for root-finding of complex polynomial
57  *
58  * The ICL developers applied only minor changes to this file
59  ****************************************************************************/
60 
61  typedef double xreal;
62  static const struct { double ZERO, INFIN; int MIN_EXP, MAX_EXP; }
63  xdata = { 0.0, DBL_MAX, DBL_MIN_EXP, DBL_MAX_EXP };
64  static xreal xnorm(xcomplex z) { return std::norm(z); }
65  static xreal xabs(xcomplex z) { return std::abs(z); }
66  static xreal xroot(xreal x, int n) { return pow(x, 1.0 / n); }
67  static int xlogb(xcomplex z) { return ilogb(xnorm(z)) / 2; }
68  #define xbits(z) DBL_MANT_DIG
69  #define xeta(z) DBL_EPSILON
70  static void xscalbln(xcomplex *z, int e) {
71  z->real(scalbln(z->real(), e));
72  z->imag(scalbln(z->imag(), e));
73  }
74 
75  // CAUCHY COMPUTES A LOWER BOUND ON THE MODULI OF THE ZEROS OF A
76  // POLYNOMIAL - PT IS THE MODULUS OF THE COEFFICIENTS.
77  //
78  static xcomplex cauchy(const int deg, xcomplex *P)
79  {
80  xreal x, xm, f, dx, df;
81  xreal *tmp = new xreal[deg + 1];
82 
83  for(int i = 0; i<=deg; i++){ tmp[i] = xabs(P[i]); };
84 
85  // Compute upper estimate bound
86  x = xroot(tmp[deg],deg) / xroot(tmp[0],deg);
87  if(tmp[deg - 1] != 0.0) {
88  // Newton step at the origin is better, use it
89  xm = tmp[deg] / tmp[deg-1];
90  if (xm < x) x = xm;
91  }
92 
93  tmp[deg] = -tmp[deg];
94 
95  // Chop the interval (0,x) until f < 0
96  while(1) {
97  xm = x * 0.1;
98  // Evaluate the polynomial <tmp> at <xm>
99  f = tmp[0];
100  for(int i = 1; i <= deg; i++)
101  f = f * xm + tmp[i];
102 
103  if(f <= 0.0) break;
104  x = xm;
105  }
106  dx = x;
107 
108  // Do Newton iteration until x converges to two decimal places
109  while(fabs(dx / x) > 0.005) {
110  f = tmp[0];
111  df = 0.0;
112  for(int i = 1; i <= deg; i++){
113  df = df * x + f;
114  f = f * x + tmp[i];
115  }
116 
117  dx = f / df;
118  x -= dx; // Newton step
119  }
120 
121  delete tmp;
122 
123  return (xcomplex)(x);
124  }
125 
126  // RETURNS A SCALE FACTOR TO MULTIPLY THE COEFFICIENTS OF THE POLYNOMIAL.
127  // THE SCALING IS DONE TO AVOID OVERFLOW AND TO AVOID UNDETECTED UNDERFLOW
128  // INTERFERING WITH THE CONVERGENCE CRITERION. THE FACTOR IS A POWER OF THE
129  //int BASE.
130  // PT - MODULUS OF COEFFICIENTS OF P
131  // ETA, INFIN, SMALNO, BASE - CONSTANTS DESCRIBING THE FLOATING POINT ARITHMETIC.
132  //
133  static void scale(const int deg, xcomplex* P)
134  {
135  int hi, lo, max, min, x, sc;
136 
137  // Find largest and smallest moduli of coefficients
138  hi = (int)(xdata.MAX_EXP / 2.0);
139  lo = (int)(xdata.MIN_EXP - xbits(P[0]));
140  max = xlogb(P[0]); // leading coefficient does not vanish!
141  min = xlogb(P[0]);
142 
143  for(int i = 0; i <= deg; i++) {
144  if (P[i] != xdata.ZERO){
145  x = xlogb(P[i]);
146  if(x > max) max = x;
147  if(x < min) min = x;
148  }
149  }
150 
151  // Scale only if there are very large or very small components
152  if(min >= lo && max <= hi) return;
153 
154  x = lo - min;
155  if(x <= 0)
156  sc = -(max+min) / 2;
157  else {
158  sc = x;
159  if(xdata.MAX_EXP - sc > max) sc = 0;
160  }
161 
162  // Scale the polynomial
163  for(int i = 0; i<= deg; i++){ xscalbln(&P[i],sc); }
164  }
165 
166  // COMPUTES THE DERIVATIVE POLYNOMIAL AS THE INITIAL H
167  // POLYNOMIAL AND COMPUTES L1 NO-SHIFT H POLYNOMIALS.
168  //
169  static void noshft(const int l1, int deg, xcomplex *P, xcomplex *H)
170  {
171  int i, j, jj;
172  //#ifndef ICL_SYSTEM_WINDOWS removed by CE
173  //#warning "why is tmp not iniitialized here?"
174  //#endif
175  xcomplex t, tmp;
176 
177  // compute the first H-polynomial as the (normed) derivative of P
178  for (i = 0; i < deg; i++) {
179  // CE: shouldn't it be this here (I add it since it seems to not make sense
180  // to use the uninitiylized tmp here a multipliler ...
181  //
182  tmp = P[i]; // added by CE
183  tmp *= (deg - i);
184  tmp /= deg;
185  H[i] = tmp;
186  }
187 #if 0
188  // possible fix? found here: https://github.com/gap-packages/float/blob/master/src/cpoly.C
189  for(i = 0; i < deg; i++)
190  H[i] = (P[i] * (deg-i)) / deg;
191 #endif
192 
193  for(jj = 1; jj <= l1; jj++) {
194  if(xnorm(H[deg - 1]) > xeta(P[deg-1])*xeta(P[deg-1])* 10*10 * xnorm(P[deg - 1])) {
195  t = -P[deg] / H[deg-1];
196  for(i = 0; i < deg-1; i++){
197  j = deg - i - 1;
198  H[j] = t * H[j-1] + P[j];
199  }
200  H[0] = P[0];
201  } else {
202  // if the constant term is essentially zero, shift H coefficients
203  for(i = 0; i < deg-1; i++) {
204  j = deg - i - 1;
205  H[j] = H[j - 1];
206  }
207  H[0] = xdata.ZERO;
208  }
209  }
210  }
211 
212  // EVALUATES A POLYNOMIAL P AT S BY THE HORNER RECURRENCE
213  // PLACING THE PARTIAL SUMS IN Q AND THE COMPUTED VALUE IN PV.
214  //
215  static xcomplex polyev(const int deg, const xcomplex s, const xcomplex *Q, xcomplex *q) {
216  q[0] = Q[0];
217  for(int i = 1; i <= deg; i++){ q[i] = q[i-1] * s + Q[i]; };
218  return q[deg];
219  }
220 
221  // COMPUTES T = -P(S)/H(S).
222  // BOOL - LOGICAL, SET TRUE IF H(S) IS ESSENTIALLY ZERO.
223  //
224  static xcomplex calct(bool *bol, int deg, xcomplex Ps, xcomplex *H, xcomplex *h, xcomplex s){
225  xcomplex Hs;
226  Hs = polyev(deg-1, s, H, h);
227  *bol = xnorm(Hs) <= xeta(H[deg-1])*xeta(H[deg-1]) * 10*10 * xnorm(H[deg-1]);
228  if(!*bol)
229  return -Ps / Hs;
230  else
231  return xdata.ZERO;
232  }
233 
234  // CALCULATES THE NEXT SHIFTED H POLYNOMIAL.
235  // BOOL - LOGICAL, IF .TRUE. H(S) IS ESSENTIALLY ZERO
236  //
237  static void nexth(const bool bol, int deg, xcomplex t, xcomplex *H, xcomplex *h, xcomplex *p){
238  if(!bol){
239  for(int j = 1; j < deg; j++)
240  H[j] = t * h[j-1] + p[j];
241  H[0] = p[0];
242  }
243  else {
244  // If h[s] is zero replace H with qh
245  for(int j = 1; j < deg; j++)
246  H[j] = h[j - 1];
247  h[0] = xdata.ZERO;
248  }
249  }
250 
251  // BOUNDS THE ERROR IN EVALUATING THE POLYNOMIAL BY THE HORNER RECURRENCE.
252  // QR,QI - THE PARTIAL SUMS
253  // MS -MODULUS OF THE POINT
254  // MP -MODULUS OF POLYNOMIAL VALUE
255  // ARE, MRE -ERROR BOUNDS ON COMPLEX ADDITION AND MULTIPLICATION
256  //
257  static xreal errev(const int deg, const xcomplex *p, const xreal ms, const xreal mp){
258  xreal MRE = 2.0 * sqrt(2.0) * xeta(p[0]);
259  xreal e = xabs(p[0]) * MRE / (xeta(p[0]) + MRE);
260 
261  for(int i = 0; i <= deg; i++)
262  e = e * ms + xabs(p[i]);
263 
264  return e * (xeta(p[0]) + MRE) - MRE * mp;
265  }
266 
267  // CARRIES OUT THE THIRD STAGE ITERATION.
268  // L3 - LIMIT OF STEPS IN STAGE 3.
269  // ZR,ZI - ON ENTRY CONTAINS THE INITIAL ITERATE, IF THE
270  // ITERATION CONVERGES IT CONTAINS THE FINAL ITERATE ON EXIT.
271  // CONV - .TRUE. IF ITERATION CONVERGES
272  //
273  static bool vrshft(const int l3, int deg, xcomplex *P, xcomplex *p, xcomplex *H, xcomplex *h, xcomplex *zero, xcomplex *s){
274  bool bol, conv, b;
275  int i, j;
276  xcomplex Ps, t;
277  xreal mp, ms, omp = 0.0, relstp = 0.0, tp;
278 
279  conv = b = false;
280  *s = *zero;
281 
282  // Main loop for stage three
283  for(i = 1; i <= l3; i++) {
284  // Evaluate P at S and test for convergence
285  Ps = polyev(deg, *s, P, p);
286  mp = xabs(Ps);
287  ms = xabs(*s);
288  if(mp <= 20 * errev(deg, p, ms, mp)) {
289  // Polynomial value is smaller in value than a bound on the error
290  // in evaluating P, terminate the iteration
291  conv = true;
292  *zero = *s;
293  return conv;
294  }
295 
296  if(i != 1) {
297  if(!(b || mp < omp || relstp >= 0.05)){
298  // if(!(b || xlogb(mp) < omp || real(relstp) >= 0.05)){
299  // Iteration has stalled. Probably a cluster of zeros. Do 5 fixed
300  // shift steps into the cluster to force one zero to dominate
301  tp = relstp;
302  b = true;
303  if(relstp < xeta(P[0])) tp = xeta(P[0]);
304 
305  *s *= 1.0 + xcomplex(1.0, 1.0)*sqrt(tp);
306 
307  Ps = polyev(deg, *s, P, p);
308  for(j = 1; j <= 5; j++){
309  t = calct(&bol, deg, Ps, H, h, *s);
310  nexth(bol, deg, t, H, h, p);
311  }
312  omp = xdata.INFIN;
313  goto _20;
314  }
315 
316  // Exit if polynomial value increase significantly
317  if(mp * 0.1 > omp) return conv;
318  }
319 
320  omp = mp;
321 
322  // Calculate next iterate
323  _20: t = calct(&bol, deg, Ps, H, h, *s);
324  nexth(bol, deg, t, H, h, p);
325  t = calct(&bol, deg, Ps, H, h, *s);
326  if(!bol) {
327  relstp = xabs(t) / xabs(*s);
328  *s += t;
329  }
330  } // end for
331 
332  return conv;
333  }
334 
335  // COMPUTES L2 FIXED-SHIFT H POLYNOMIALS AND TESTS FOR CONVERGENCE.
336  // INITIATES A VARIABLE-SHIFT ITERATION AND RETURNS WITH THE
337  // APPROXIMATE ZERO IF SUCCESSFUL.
338  // L2 - LIMIT OF FIXED SHIFT STEPS
339  // ZR,ZI - APPROXIMATE ZERO IF CONV IS .TRUE.
340  // CONV - LOGICAL INDICATING CONVERGENCE OF STAGE 3 ITERATION
341  //
342  static bool fxshft(const int l2, int deg, xcomplex *P, xcomplex *p, xcomplex *H, xcomplex *h, xcomplex *zero, xcomplex *s){
343  bool bol, conv; // boolean for convergence of stage 2
344  bool test, pasd;
345  xcomplex old_T, old_S, Ps, t;
346  xcomplex *Tmp = new xcomplex[deg+1];
347 
348  Ps = polyev(deg, *s, P, p);
349  test = true;
350  pasd = false;
351 
352  // Calculate first T = -P(S)/H(S)
353  t = calct(&bol, deg, Ps, H, h, *s);
354 
355  // Main loop for second stage
356  for(int j = 1; j <= l2; j++){
357  old_T = t;
358 
359  // Compute the next H Polynomial and new t
360  nexth(bol, deg, t, H, h, p);
361  t = calct(&bol, deg, Ps, H, h, *s);
362  *zero = *s + t;
363 
364  // Test for convergence unless stage 3 has failed once or this
365  // is the last H Polynomial
366  if(!(bol || !test || j == l2)){
367  if(xabs(t - old_T) < 0.5 * xabs(*zero)) {
368  if(pasd) {
369  // The weak convergence test has been passwed twice, start the third stage
370  // Iteration, after saving the current H polynomial and shift
371  for(int i = 0; i < deg; i++)
372  Tmp[i] = H[i];
373  old_S = *s;
374 
375  conv = vrshft(10, deg, P, p, H, h, zero, s);
376  if(conv) return conv;
377 
378  //The iteration failed to converge. Turn off testing and restore h,s,pv and T
379  test = false;
380  for(int i = 0; i < deg; i++)
381  H[i] = Tmp[i];
382  *s = old_S;
383 
384  Ps = polyev(deg, *s, P, p);
385  t = calct(&bol, deg, Ps, H, h, *s);
386  continue;
387  }
388  pasd = true;
389  }
390  else
391  pasd = false;
392  }
393  }
394 
395  // Attempt an iteration with final H polynomial from second stage
396  conv = vrshft(10, deg, P, p, H, h, zero, s);
397 
398  delete Tmp;
399 
400  return conv;
401  }
402 
403  // Main function
404  //
405  int cpoly(int degree, const xcomplex poly[], xcomplex Roots[])
406  {
407  xcomplex PhiDiff = xcomplex(-0.069756473, 0.99756405);
408  xcomplex PhiRand = xcomplex(1.0, -1.0) / sqrt(2.0);
409  xcomplex *P = new xcomplex[degree + 1], *H = new xcomplex[degree + 1];
410  xcomplex *h = new xcomplex[degree + 1], *p = new xcomplex[degree + 1];
411  xcomplex zero, s, bnd;
412  unsigned int conv = 0;
413 
414  while(poly[0] == xdata.ZERO) {
415  poly++;
416  degree--;
417  if (degree < 0)
418  return -1;
419  };
420 
421  int deg = degree;
422 
423  // Remove the zeros at the origin if any
424  while(poly[deg] == xdata.ZERO){
425  Roots[degree - deg] = xdata.ZERO;
426  deg--;
427  }
428 
429  if (deg == 0) return degree;
430 
431  // Make a copy of the coefficients
432  for(int i = 0; i <= deg; i++) { P[i] = poly[i]; }
433 
434  scale(deg, P);
435 
436  search:
437 
438  if(deg <= 1){
439  Roots[degree-1] = - P[1] / P[0];
440  return degree;
441  }
442 
443  // compute a bound of the moduli of the roots (Newton-Raphson)
444  bnd = cauchy(deg, P);
445 
446  // Outer loop to control 2 Major passes with different sequences of shifts
447  for(int cnt1 = 1; cnt1 <= 2; cnt1++) {
448  // First stage calculation , no shift
449  noshft(5, deg, P, H);
450 
451  // Inner loop to select a shift
452  for(int cnt2 = 1; cnt2 <= 9; cnt2++) {
453  // Shift is chosen with modulus bnd and amplitude rotated by 94 degree from the previous shif
454  PhiRand = PhiDiff * PhiRand;
455  s = bnd * PhiRand;
456 
457  // Second stage calculation, fixed shift
458  conv = fxshft(10 * cnt2, deg, P, p, H, h, &zero, &s);
459  if(conv) {
460  // The second stage jumps directly to the third stage iteration
461  // If successful the zero is stored and the polynomial deflated
462  Roots[degree - deg] = zero;
463 
464  // continue with the remaining polynomial
465  deg--;
466  for(int i = 0; i <= deg; i++){ P[i] = p[i]; };
467  goto search;
468  }
469  // if the iteration is unsuccessful another shift is chosen
470  }
471 
472  // if 9 shifts fail, the outer loop is repeated with another sequence of shifts
473  }
474 
475  delete P;
476  delete H;
477  delete h;
478  delete p;
479 
480  // The zerofinder has failed on two major passes
481  // return empty handed with the number of roots found (less than the original degree)
482  return degree - deg;
483  }
484 
485  // Wrapper function:
486  // This function calls the one made by Laurent Bartholdi
487  int solve_poly(int degree, const double *in_real, const double *in_imag, double *out_real, double *out_imag)
488  {
489  int nroots;
490  xcomplex *poly = new xcomplex[degree+1];
491  xcomplex *roots = new xcomplex[degree+1];
492 
493  for (int i = 0; i <= degree; ++i) {
494  poly[i].real(in_real[i]);
495  poly[i].imag(in_imag[i]);
496  }
497 
498  nroots = cpoly(degree, poly, roots);
499 
500  for (int i = 0; i < nroots; ++i) {
501  out_real[i] = roots[i].real();
502  out_imag[i] = roots[i].imag();
503  }
504 
505  delete[] poly;
506  delete[] roots;
507 
508  return nroots;
509  }
510  }
511 }
512 
513 #endif
static xcomplex calct(bool *bol, int deg, xcomplex Ps, xcomplex *H, xcomplex *h, xcomplex s)
Definition: PolynomialSolver.h:224
ICLQt_API ImgQ sqrt(const ImgQ &image)
calls sqrt( each pixel)
static bool vrshft(const int l3, int deg, xcomplex *P, xcomplex *p, xcomplex *H, xcomplex *h, xcomplex *zero, xcomplex *s)
Definition: PolynomialSolver.h:273
undocument this line if you encounter any issues!
Definition: Any.h:37
double ZERO
Definition: PolynomialSolver.h:62
ICLQt_API core::Img< T > norm(const core::Img< T > &image)
normalize an images range to [0,255]
static void xscalbln(xcomplex *z, int e)
Definition: PolynomialSolver.h:70
int MIN_EXP
Definition: PolynomialSolver.h:62
static const struct icl::math::@18 xdata
static xreal errev(const int deg, const xcomplex *p, const xreal ms, const xreal mp)
Definition: PolynomialSolver.h:257
static bool fxshft(const int l2, int deg, xcomplex *P, xcomplex *p, xcomplex *H, xcomplex *h, xcomplex *zero, xcomplex *s)
Definition: PolynomialSolver.h:342
double xreal
Definition: PolynomialSolver.h:61
static void nexth(const bool bol, int deg, xcomplex t, xcomplex *H, xcomplex *h, xcomplex *p)
Definition: PolynomialSolver.h:237
static xreal xnorm(xcomplex z)
Definition: PolynomialSolver.h:64
ICLQt_API ImgQ abs(const ImgQ &image)
calls abs ( each pixel)
std::complex< double > xcomplex
Definition: PolynomialSolver.h:43
double INFIN
Definition: PolynomialSolver.h:62
static xreal xroot(xreal x, int n)
Definition: PolynomialSolver.h:66
static int xlogb(xcomplex z)
Definition: PolynomialSolver.h:67
static xcomplex cauchy(const int deg, xcomplex *P)
Definition: PolynomialSolver.h:78
int solve_poly(int degree, const double *in_real, const double *in_imag, double *out_real, double *out_imag)
Definition: PolynomialSolver.h:487
static void scale(const int deg, xcomplex *P)
Definition: PolynomialSolver.h:133
static xreal xabs(xcomplex z)
Definition: PolynomialSolver.h:65
int MAX_EXP
Definition: PolynomialSolver.h:62
static void noshft(const int l1, int deg, xcomplex *P, xcomplex *H)
Definition: PolynomialSolver.h:169
#define xeta(z)
Definition: PolynomialSolver.h:69
static xcomplex polyev(const int deg, const xcomplex s, const xcomplex *Q, xcomplex *q)
Definition: PolynomialSolver.h:215
int cpoly(int degree, const xcomplex poly[], xcomplex Roots[])
Definition: PolynomialSolver.h:405
#define xbits(z)
Definition: PolynomialSolver.h:68