OPAL (Object Oriented Parallel Accelerator Library)  2.2.0
OPAL
DoubleEigen.cpp
Go to the documentation of this file.
1 // ------------------------------------------------------------------------
2 // $RCSfile: DoubleEigen.cpp,v $
3 // ------------------------------------------------------------------------
4 // $Revision: 1.1.1.1 $
5 // ------------------------------------------------------------------------
6 // Copyright: see Copyright.readme
7 // ------------------------------------------------------------------------
8 //
9 // Class ComplexEigen:
10 // Eigenvector methods for class: Matrix<double>.
11 //
12 // ------------------------------------------------------------------------
13 // Class category: Algebra
14 // ------------------------------------------------------------------------
15 //
16 // $Date: 2000/03/27 09:32:32 $
17 // $Author: fci $
18 //
19 // ------------------------------------------------------------------------
20 
21 #include "Algebra/DoubleEigen.h"
22 #include "Algebra/Vector.h"
24 #include "Utilities/SizeError.h"
25 #include <cmath>
26 
27 using std::abs;
28 using std::imag;
29 using std::max;
30 using std::real;
31 using std::swap;
32 
33 
34 // Helper function.
35 // ------------------------------------------------------------------------
36 
37 namespace {
38  void cdiv(double ar, double ai, double br, double bi,
39  double &cr, double &ci)
40  // Complex division: (cr,ci) = (ar,ai) / (br,bi).
41  // Adapted March 1997 by F. Christoph Iselin, CERN, SL/AP.
42  {
43  double s = std::abs(br) + std::abs(bi);
44  double ars = ar / s;
45  double ais = ai / s;
46  double brs = br / s;
47  double bis = bi / s;
48 
49  s = brs * brs + bis * bis;
50 
51  cr = (ars * brs + ais * bis) / s;
52  ci = (ais * brs - ars * bis) / s;
53  }
54 }
55 
56 
57 // Class DoubleEigen; public methods.
58 // --------------------------------------------------------------------------
59 
61  lambda(), vectors()
62 {}
63 
64 
66  lambda(M.lambda), vectors(M.vectors)
67 {}
68 
69 
71  lambda(M.nrows()), vectors(M.nrows(), M.nrows(), 1.0)
72  // Before execution:
73  // M is the matrix whose eigenvalues are requested.
74  // After execution:
75  // lambda contains the eigenvalues.
76  // If vect is true, then:
77  // vectors contains the eigenvectors, packed in a real matrix.
78 {
79  int n = M.nrows();
80  int m = M.ncols();
81  if(m != n) {
82  throw SizeError("DoubleEigen::DoubleEigen()", "Matrix is not square.");
83  }
84  int low, upp;
85 
86  // Copy original matrix.
87  Matrix<double> h(M);
88 
89  // Balance the copy.
90  Vector<double> scale(n);
91  balance(h, low, upp, scale);
92 
93  // Reduce the copy to upper Hessenberg form.
94  Array1D<int> index(n);
95  elmhes(h, low, upp, index);
96 
97  if(vec) {
98  elmtran(h, low, upp, index);
99 
100  if(hqr2(h, low, upp) != 0) {
101  throw EigenvalueError("DoubleEigen::DoubleEigen()",
102  "Unable to find all eigenvalues.");
103  } else {
104  balbak(low, upp, scale);
105  }
106  } else {
107  if(hqr(h, low, upp) != 0) {
108  throw EigenvalueError("DoubleEigen::DoubleEigen()",
109  "Unable to find all eigenvalues.");
110  }
111  }
112 }
113 
114 
116 {}
117 
118 
120 // Return Eigenvalues as complex vector.
121 {
122  return lambda;
123 }
124 
125 
127 // Return eigenvectors as a complex matrix.
128 {
129  int n = vectors.nrows();
130  Matrix<complex<double> > R(n, n);
131 
132  for(int i = 0; i < n; i++) {
133  if(imag(lambda[i]) == 0.0) {
134  // One real eigenvector.
135  for(int j = 0; j < n; j++) {
136  R[j][i] = complex<double>(vectors[j][i]);
137  }
138  } else {
139  // Two complex eigenvectors.
140  for(int j = 0; j < n; j++) {
141  R[j][i] = complex<double>(vectors[j][i], + vectors[j][i+1]);
142  R[j][i+1] = complex<double>(vectors[j][i], - vectors[j][i+1]);
143  }
144  i++;
145  }
146  }
147 
148  return R;
149 }
150 
151 
153 // Return eigenvectors packed in a real matrix.
154 {
155  return vectors;
156 }
157 
158 
159 // Class DoubleEigen; private methods.
160 // --------------------------------------------------------------------------
161 
163 (Matrix<double> &copy, int &low, int &upp, Vector<double> &scale)
164 
165 // This method is a translation of the Algol procedure balance,
166 // Num. Math. 13, 293-304(1969) by Parlett and Reinsch.
167 // Handbook for Auto. Comp., Vol.ii-Linear Algebra, 315-326(1971).
168 //
169 // It balances a real matrix and isolates eigenvalues whenever possible.
170 //
171 // On input:
172 // copy Contains a copy of the original matrix A.
173 //
174 // On output:
175 // copy Contains the balanced matrix A.
176 //
177 // low, upp Are two indices such that A[i][j] is equal to zero if
178 // (1) i is greater than j and
179 // (2) j = 0, ..., low - 1 or i = upp + 1, ..., n - 1.
180 //
181 // scale Contains information determining the
182 // permutations and scaling factors used.
183 //
184 // Suppose that the principal submatrix in rows "low" through "upp"
185 // has been balanced, that p[j] denotes the index interchanged with j
186 // during the permutation step, and that the elements of the diagonal
187 // matrix used are denoted by d[j][j]. Then
188 // scale[j] = p[j], for j = 0, ..., low - 1,
189 // = d[j,j], j = low, ..., upp,
190 // = p[j] j = upp + 1, ..., n - 1.
191 // The order in which the interchanges are made is n to upp + 1,
192 // then 0 to low - 1.
193 //
194 // Translated to FORTRAN by Burton S. Garbow, ANL
195 // Adapted March 1997 by F. Christoph Iselin, CERN, SL/AP
196 {
197  // Tentatively set the limits "low" and "upp".
198  int n = copy.nrows();
199  low = 0;
200  upp = n - 1;
201 
202  // Search for rows isolating an eigenvalue and push them down.
203  for(int j = upp; ; j--) {
204  for(int i = 0; i <= upp; i++) {
205  if(i != j && copy[j][i] != 0.0) goto next_column;
206  }
207 
208  // Column "j" qualifies for exchange.
209  exchange(copy, j, upp, low, upp);
210  scale[upp] = double(j);
211  j = --upp;
212 
213  // Restart search for previous column.
214 next_column:
215  if(j == 0) break;
216  }
217 
218  // Search for columns isolating an eigenvalue and push them left.
219  for(int j = low; ; j++) {
220  for(int i = low; i <= upp; i++) {
221  if(i != j && copy[i][j] != 0.0) goto next_row;
222  }
223 
224  // Row "j" qualifies for exchange.
225  exchange(copy, j, low, low, upp);
226  scale[low] = double(j);
227  j = ++low;
228 
229  // Restart search for next column.
230 next_row:
231  if(j == upp) break;
232  }
233 
234  // Now balance the submatrix in rows low, ..., upp.
235  for(int i = low; i <= upp; i++) scale[i] = 1.0;
236 
237  // Iterative loop for norm reduction.
238  const double radix = 16.0;
239  const double b2 = radix * radix;
240  bool noconv = false;
241 
242  do {
243  noconv = false;
244 
245  for(int i = low; i <= upp; i++) {
246  double r = 0.0;
247  double c = 0.0;
248 
249  for(int j = low; j <= upp; j++) {
250  if(j != i) {
251  c += std::abs(copy[j][i]);
252  r += std::abs(copy[i][j]);
253  }
254  }
255 
256  // Guard against zero c or r due to underflow.
257  if(c != 0.0 && r != 0.0) {
258  double g = r / radix;
259  double f = 1.0;
260  double s = c + r;
261 
262  while(c < g) {
263  f *= radix;
264  c *= b2;
265  }
266 
267  g = r * radix;
268 
269  while(c >= g) {
270  f /= radix;
271  c /= b2;
272  }
273 
274  // Now balance.
275  if((c + r) / f < s * .95) {
276  g = 1.0 / f;
277  scale[i] *= f;
278  noconv = true;
279  for(int j = low; j < n; j++) copy[i][j] *= g;
280  for(int j = 0; j <= upp; j++) copy[j][i] *= f;
281  }
282  }
283  }
284  } while(noconv);
285 }
286 
287 
288 void DoubleEigen::balbak(int low, int upp, Vector<double> &scale)
289 // This method is a translation of the Algol procedure balbak,
290 // Num. Math. 13, 293-304(1969) by Parlett and Reinsch.
291 // Handbook for Auto. Comp., Vol.ii-Linear Algebra, 315-326(1971).
292 //
293 // It forms the eigenvectors of a real general matrix by back transforming
294 // those of the corresponding balanced matrix determined by "balance".
295 //
296 // On input:
297 // vectors Contains the real and imaginary parts of the eigenvectors
298 // to be back transformed.
299 //
300 // low, upp Are indices determined by "balance".
301 //
302 // scale Contains information determining the permutations and
303 // scaling factors used by "balance".
304 //
305 // On output:
306 // vectors Contains the real and imaginary parts of the transformed
307 // eigenvectors.
308 //
309 // Translated to FORTRAN by Burton S. Garbow, ANL.
310 // Adapted March 1997 by F. Christoph Iselin, CERN, SL/AP.
311 {
312  int n = vectors.nrows();
313 
314  // Apply scale factors found by "balance" to rows low, ..., upp.
315  if(upp != low) {
316  for(int i = low; i <= upp; i++) {
317  double s = scale[i];
318  for(int j = 0; j < n; j++) vectors[i][j] *= s;
319  }
320  }
321 
322  // Exchange rows which were interchanged by "balance".
323  for(int i = low; i-- > 0;) {
324  int k = int(scale[i]);
325  if(k != i) vectors.swapRows(i, k);
326  }
327 
328  for(int i = upp + 1; i < n; i++) {
329  int k = int(scale[i]);
330  if(k != i) vectors.swapRows(i, k);
331  }
332 }
333 
334 
335 void DoubleEigen::elmhes(Matrix<double> &copy, int low, int upp,
336  Array1D<int> &index)
337 // This method is a translation of the Algol procedure elmhes,
338 // Num. Math. 12, 349-368(1968) by Martin and Wilkinson.
339 // Handbook for Auto. Comp., Vol.ii-Linear Algebra, 339-358(1971).
340 //
341 // Given a real general matrix, it reduces a submatrix situated in
342 // rows and columns low, ..., upp to upper Hessenberg form by
343 // stabilized elementary similarity transformations.
344 //
345 // On input:
346 // copy Contains the balanced matrix.
347 //
348 // low, upp Are indices determined by "balance". If "balance"has not
349 // been used, set low = 1, upp = n.
350 //
351 // On output:
352 // copy Contains the Hessenberg matrix. The multipliers
353 // which were used in the reduction are stored in the
354 // remaining triangle under the Hessenberg matrix.
355 //
356 // index Contains information on the rows and columns
357 // interchanged in the reduction.
358 // Only elements low, ..., upp are used.
359 //
360 // Translated to FORTRAN by Burton S. Garbow, ANL.
361 // Adapted March 1997 by F. Christoph Iselin, CERN, SL/AP.
362 {
363  int n = copy.nrows();
364 
365  for(int m = low + 1; m < upp; m++) {
366  double x = 0.0;
367  int i = m;
368 
369  for(int j = m; j <= upp; j++) {
370  if(std::abs(copy[j][m-1]) > std::abs(x)) {
371  x = copy[j][m-1];
372  i = j;
373  }
374  }
375 
376  // Interchange rows and columns of A.
377  index[m] = i;
378 
379  if(i != m) {
380  for(int j = m - 1; j < n; j++) swap(copy[i][j], copy[m][j]);
381  for(int j = 0; j <= upp; j++) swap(copy[j][i], copy[j][m]);
382  }
383 
384  if(x != 0.0) {
385  for(i = m + 1; i <= upp; i++) {
386  double y = copy[i][m-1];
387  if(y != 0.0) {
388  y /= x;
389  copy[i][m-1] = y;
390  for(int j = m; j < n; j++) copy[i][j] -= y * copy[m][j];
391  for(int j = 0; j <= upp; j++) copy[j][m] += y * copy[j][i];
392  }
393  }
394  }
395  }
396 }
397 
398 
399 void DoubleEigen::elmtran(Matrix<double> &copy, int low, int upp,
400  Array1D<int> &index)
401 // This method is a translation of the Algol procedure elmtrans,
402 // Num. Math. 16, 181-204(1970) by Peters and Wilkinson.
403 // Handbook for Auto. Comp., Vol.ii-Linear Algebra, 372-395(1971).
404 //
405 // It accumulates the stabilized elementary similarity transformations
406 // used in the reduction of a real general matrix to upper Hessenberg
407 // form by "elmhes".
408 //
409 // On input:
410 // copy Contains the matrix A to be transformed.
411 //
412 // low, upp Are indices determined by the balancing subroutine
413 // "balance". If "balance" has not been used,
414 // set low=1, upp=n.
415 //
416 // copy Contains the multipliers which were used in the reduction
417 // by "elmhes" in its lower triangle below the subdiagonal.
418 //
419 // index Contains information on the rows and columns interchanged
420 // in the reduction by "elmhes". Only elements low, ..., upp
421 // are used.
422 //
423 // On output:
424 // vectors Contains the transformation matrix produced in the
425 // reduction by elmhes.
426 //
427 // Translated to FORTRAN by Burton S. Garbow, ANL.
428 // Adapted March 1997 by F. Christoph Iselin, CERN, SL/AP.
429 {
430  for(int mp = upp; --mp > low;) {
431  int i;
432  for(i = mp + 1; i <= upp; i++) vectors[i][mp] = copy[i][mp-1];
433 
434  i = index[mp];
435  if(i != mp) {
436  for(int j = mp; j <= upp; j++) {
437  vectors[mp][j] = vectors[i][j];
438  vectors[i][j] = 0.0;
439  }
440 
441  vectors[i][mp] = 1.0;
442  }
443  }
444 }
445 
446 
448  int j, int m, int low, int upp) {
449  int n = copy.nrows();
450 
451  if(j != m) {
452  for(int i = 0; i <= upp; i++) swap(copy[i][j], copy[i][m]);
453  for(int i = low; i < n; i++) swap(copy[j][i], copy[m][i]);
454  }
455 }
456 
457 
458 int DoubleEigen::hqr(Matrix<double> &h, int low, int upp)
459 // This method is a translation of the Algol procedure hqr,
460 // Num. Math. 14, 219-231(1970) by Martin, Peters, and Wilkinson.
461 // Handbook for Auto. Comp., Vol.ii-Linear Algebra, 359-371(1971).
462 //
463 // It finds the eigenvalues of a real
464 // upper Hessenberg matrix by the QR method.
465 //
466 // On input:
467 // h Contains the upper Hessenberg matrix. Information about
468 // the transformations used in the reduction to Hessenberg
469 // form by "elmhes" or "orthes", if performed, is stored
470 // in the remaining triangle under the Hessenberg matrix.
471 //
472 // low, upp Are indices determined by the balancing subroutine
473 // "balance". If "balance" has not been used,
474 // set low=1, upp=n.
475 //
476 // On output:
477 // h Is destroyed. Therefore, it must be saved before calling
478 // "hqr" if subsequent calculation and back transformation
479 // of eigenvectors is to be performed.
480 //
481 // lambda Contains the eigenvalues. The eigenvalues are unordered
482 // except that complex conjugate pairs of values appear
483 // consecutively with the eigenvalue having the positive
484 // imaginary part first. If an error exit is made, the
485 // eigenvalues should be correct for indices ierr, ..., n - 1.
486 //
487 // The result value is set to:
488 // zero For normal return,
489 // j + 1 If the limit of 30*n iterations is exhausted
490 // while the j-th eigenvalue is being sought.
491 //
492 // Translated to FORTRAN by Burton S. Garbow, ANL.
493 // Adapted March 1997 by F. Christoph Iselin, CERN, SL/AP.
494 {
495  int n = h.nrows();
496 
497  // Store roots isolated by "balance".
498  for(int i = 0; i < low; i++) lambda[i] = complex<double>(h[i][i]);
499  for(int i = upp + 1; i < n; i++) lambda[i] = complex<double>(h[i][i]);
500 
501  // Compute matrix norm.
502  double norm = 0.0;
503  int k = 0;
504  for(int i = 0; i < n; i++) {
505  for(int j = k; j < n; j++) norm += std::abs(h[i][j]);
506  k = i;
507  }
508 
509  // Search for next eigenvalues.
510  double t = 0.0;
511  double p = 0.0, q = 0.0, r = 0.0, s = 0.0;
512  double w = 0.0, x = 0.0, y = 0.0, z = 0.0;
513  int itn = n * 30;
514  for(int en = upp + 1; en-- > low;) {
515  int its = 0;
516 
517  // Look for single small sub-diagonal element.
518  while(true) {
519  int l;
520 
521  for(l = en; l > low; l--) {
522  s = std::abs(h[l-1][l-1]) + std::abs(h[l][l]);
523  if(s == 0.0) s = norm;
524  if((s + std::abs(h[l][l-1])) == s) break;
525  }
526 
527  // Form shift.
528  x = h[en][en];
529  if(l == en) goto one_root;
530 
531  y = h[en-1][en-1];
532  w = h[en][en-1] * h[en-1][en];
533  if(l == en - 1) goto two_roots;
534 
535  // Set error, if all eigenvalues have not converged.
536  if(itn == 0) return (en + 1);
537 
538  // Form exceptional shift.
539  if(its == 10 || its == 20) {
540  t += x;
541  for(int i = low; i <= en; i++) h[i][i] -= x;
542  s = std::abs(h[en][en-1]) + std::abs(h[en-1][en-2]);
543  x = y = 0.75 * s;
544  w = -.4375 * s * s;
545  }
546 
547  its++;
548  itn--;
549 
550  // Look for two consecutive small sub-diagonal elements.
551  int m;
552  for(m = en - 1; m-- > l;) {
553  z = h[m][m];
554  r = x - z;
555  s = y - z;
556  p = (r * s - w) / h[m+1][m] + h[m][m+1];
557  q = h[m+1][m+1] - z - r - s;
558  r = h[m+2][m+1];
559  s = std::abs(p) + std::abs(q) + std::abs(r);
560  p /= s;
561  q /= s;
562  r /= s;
563  if(m == l) break;
564  double tst1 = std::abs(p) * (std::abs(h[m-1][m-1]) + std::abs(z) + std::abs(h[m+1][m+1]));
565  double tst2 = tst1 + std::abs(h[m][m-1]) * (std::abs(q) + std::abs(r));
566  if(tst2 == tst1) break;
567  }
568 
569  h[m+2][m] = 0.0;
570  for(int i = m + 3; i <= en; i++) h[i][i-2] = h[i][i-3] = 0.0;
571 
572  // Double QR step involving rows l, ..., en and columns m, ..., en.
573  for(k = m; k < en; k++) {
574  bool notLast = (k != en - 1);
575 
576  if(k != m) {
577  p = h[k][k-1];
578  q = h[k+1][k-1];
579  r = notLast ? h[k+2][k-1] : 0.0;
580  x = std::abs(p) + std::abs(q) + std::abs(r);
581  if(x == 0.0) continue;
582  p /= x;
583  q /= x;
584  r /= x;
585  }
586 
587  s = sqrt(p * p + q * q + r * r);
588  if(p < 0.0) s = - s;
589 
590  if(k != m)
591  h[k][k-1] = - s * x;
592  else if(l != m)
593  h[k][k-1] = - h[k][k-1];
594 
595  p += s;
596  x = p / s;
597  y = q / s;
598  z = r / s;
599  q /= p;
600  r /= p;
601 
602  // Row modification.
603  for(int j = k; j <= en; j++) {
604  p = h[k][j] + q * h[k+1][j];
605  if(notLast) {
606  p += r * h[k+2][j];
607  h[k+2][j] -= p * z;
608  }
609  h[k][j] -= p * x;
610  h[k+1][j] -= p * y;
611  }
612 
613  // Column modification.
614  int j = (en < k + 3) ? en : (k + 3);
615  for(int i = l; i <= j; i++) {
616  p = x * h[i][k] + y * h[i][k+1];
617  if(notLast) {
618  p += z * h[i][k+2];
619  h[i][k+2] -= p * r;
620  }
621  h[i][k] -= p;
622  h[i][k+1] -= p * q;
623  }
624  }
625  }
626 
627 one_root:
628  lambda[en] = complex<double>(x + t);
629  continue;
630 
631 two_roots:
632  p = (y - x) / 2.0;
633  q = p * p + w;
634  z = sqrt(std::abs(q));
635  x += t;
636 
637  if(q >= 0.0) {
638  // Real pair.
639  z = (p > 0.0) ? (p + z) : (p - z);
640  lambda[en-1] = complex<double>(x + z);
641  lambda[en] = complex<double>((z != 0.0) ? (x - w / z) : x + z);
642  } else {
643  // Complex pair.
644  lambda[en-1] = complex<double>(x + p, + z);
645  lambda[en] = complex<double>(x + p, - z);
646  }
647 
648  en--;
649  }
650 
651  return 0;
652 }
653 
654 
655 int DoubleEigen::hqr2(Matrix<double> &h, int low, int upp)
656 // This method is a translation of the Algol procedure hqr2,
657 // Num. Math. 16, 181-204(1970) by Peters and Wilkinson.
658 // Handbook for Auto. Comp., Vol.ii-Linear Algebra, 372-395(1971).
659 //
660 // It finds the eigenvalues and eigenvectors of a real upper Hessenberg
661 // matrix by the QR method. The eigenvectors of a real general matrix can
662 // also be found if elmhes and elmtran or orthes and ortran have been used
663 // to reduce this general matrix to Hessenberg form and to accumulate the
664 // similarity transformations.
665 //
666 // On input:
667 // h Contains the upper Hessenberg matrix. Information about
668 // the transformations used in the reduction to Hessenberg
669 // form by elmhes or orthes, if performed, is stored in the
670 // remaining triangle under the Hessenberg matrix.
671 //
672 // low, upp Are indices determined by "balance". If "balance" has not
673 // been used, set low = 1, upp = n.
674 //
675 // vectors Contains the transformation matrix produced by elmtran
676 // after the reduction by elmhes, or by ortran after the
677 // reduction by orthes, if performed. If the eigenvectors
678 // of the Hessenberg matrix are desired, vectors must contain
679 // the identity matrix.
680 //
681 // On output:
682 // h Has been destroyed.
683 //
684 // lambda Contains the eigenvalues. The eigenvalues are unordered
685 // except that complex conjugate pairs of values appear
686 // consecutively with the eigenvalue having the positive
687 // imaginary part first. If an error exit is made, the
688 // eigenvalues should be correct for indices ierr, ..., n - 1.
689 //
690 // vectors Contains the real and imaginary parts of the eigenvectors.
691 // If the i-th eigenvalue is real, the i-th column of z
692 // contains its eigenvector. If the i-th eigenvalue is complex
693 // with positive imaginary part, the i-th and (i+1)-th
694 // columns of z contain the real and imaginary parts of its
695 // eigenvector. The eigenvectors are unnormalized. If an
696 // error exit is made, none of the eigenvectors has been found.
697 //
698 // The return value is set to:
699 // zero For normal return,
700 // j+1 If the limit of 30*n iterations is exhausted
701 // while the j-th eigenvalue is being sought.
702 //
703 // Translated to FORTRAN by Burton S. Garbow, ANL.
704 // Adapted March 1997 by F. Christoph Iselin, CERN, SL/AP.
705 {
706  int n = h.nrows();
707  double tst1, tst2;
708 
709  // Store roots isolated by "balance".
710  for(int i = 0; i < low; i++) lambda[i] = complex<double>(h[i][i]);
711  for(int i = upp + 1; i < n; i++) lambda[i] = complex<double>(h[i][i]);
712 
713  // Compute matrix norm.
714  double norm = 0.0;
715  int k = 0;
716  for(int i = 0; i < n; i++) {
717  for(int j = k; j < n; j++) norm += std::abs(h[i][j]);
718  k = i;
719  }
720 
721  double t = 0.0;
722  double p = 0.0, q = 0.0, r = 0.0, s = 0.0;
723  double w = 0.0, x = 0.0, y = 0.0, z = 0.0;
724  int itn = n * 30;
725 
726  // Search for eigenvalues.
727  for(int en = upp + 1; en-- > low;) {
728  int its = 0;
729 
730  // Look for single small sub-diagonal element.
731  while(true) {
732  int l;
733  for(l = en; l > low; l--) {
734  s = std::abs(h[l-1][l-1]) + std::abs(h[l][l]);
735  if(s == 0.0) s = norm;
736  if((s + std::abs(h[l][l-1])) == s) break;
737  }
738 
739  // Form shift.
740  x = h[en][en];
741  if(l == en) goto one_root;
742 
743  y = h[en-1][en-1];
744  w = h[en][en-1] * h[en-1][en];
745  if(l == en - 1) goto two_roots;
746 
747  // Set error, if all eigenvalues have not converged.
748  if(itn == 0) return (en + 1);
749 
750  // Form exceptional shift.
751  if(its == 10 || its == 20) {
752  t += x;
753  for(int i = low; i <= en; i++) h[i][i] -= x;
754  s = std::abs(h[en][en-1]) + std::abs(h[en-1][en-2]);
755  x = y = 0.75 * s;
756  w = -0.4375 * s * s;
757  }
758 
759  its++;
760  itn--;
761 
762  // Look for two consecutive small sub-diagonal elements.
763  int m;
764  for(m = en - 1; m-- > l;) {
765  z = h[m][m];
766  r = x - z;
767  s = y - z;
768  p = (r * s - w) / h[m+1][m] + h[m][m+1];
769  q = h[m+1][m+1] - z - r - s;
770  r = h[m+2][m+1];
771  s = std::abs(p) + std::abs(q) + std::abs(r);
772  p /= s;
773  q /= s;
774  r /= s;
775  if(m == l) break;
776  tst1 = std::abs(p) * (std::abs(h[m-1][m-1]) + std::abs(z) + std::abs(h[m+1][m+1]));
777  tst2 = tst1 + std::abs(h[m][m-1]) * (std::abs(q) + std::abs(r));
778  if(tst2 == tst1) break;
779  }
780 
781  h[m+2][m] = 0.0;
782  for(int i = m + 3; i <= en; i++) h[i][i-2] = h[i][i-3] = 0.0;
783 
784  // Double QR step involving rows l to en and columns m to en.
785  for(k = m; k < en; k++) {
786  bool notLast = (k != en - 1);
787 
788  if(k != m) {
789  p = h[k][k-1];
790  q = h[k+1][k-1];
791  r = notLast ? h[k+2][k-1] : 0.0;
792  x = std::abs(p) + std::abs(q) + std::abs(r);
793  if(x == 0.0) continue;
794  p /= x;
795  q /= x;
796  r /= x;
797  }
798 
799  s = sqrt(p * p + q * q + r * r);
800  if(p < 0.0) s = - s;
801 
802  if(k != m)
803  h[k][k-1] = - s * x;
804  else if(l != m)
805  h[k][k-1] = - h[k][k-1];
806 
807  p += s;
808  x = p / s;
809  y = q / s;
810  z = r / s;
811  q /= p;
812  r /= p;
813 
814  // Row modification.
815  for(int j = k; j < n; j++) {
816  p = h[k][j] + q * h[k+1][j];
817  if(notLast) {
818  p += r * h[k+2][j];
819  h[k+2][j] -= p * z;
820  }
821  h[k][j] -= p * x;
822  h[k+1][j] -= p * y;
823  }
824 
825  // Column modification.
826  int j = (en < k + 3) ? en : (k + 3);
827  for(int i = 0; i <= j; i++) {
828  p = x * h[i][k] + y * h[i][k+1];
829  if(notLast) {
830  p += z * h[i][k+2];
831  h[i][k+2] -= p * r;
832  }
833  h[i][k] -= p;
834  h[i][k+1] -= p * q;
835  }
836 
837  // Accumulate transformations.
838  for(int i = low; i <= upp; i++) {
839  p = x * vectors[i][k] + y * vectors[i][k+1];
840  if(notLast) {
841  p += z * vectors[i][k+2];
842  vectors[i][k+2] -= p * r;
843  }
844  vectors[i][k] -= p;
845  vectors[i][k+1] -= p * q;
846  }
847  }
848  }
849 
850 one_root:
851  lambda[en] = complex<double>(h[en][en] = x + t);
852  continue;
853 
854 two_roots:
855  p = (y - x) / 2.0;
856  q = p * p + w;
857  z = sqrt(std::abs(q));
858  h[en][en] = x = x + t;
859  h[en-1][en-1] = y + t;
860 
861  if(q >= 0.0) {
862  // Real pair.
863  z = (p > 0.0) ? (p + z) : (p - z);
864  lambda[en-1] = complex<double>(x + z);
865  lambda[en] = complex<double>((z != 0.0) ? (x - w / z) : (x + z));
866  x = h[en][en-1];
867  s = std::abs(x) + std::abs(z);
868  p = x / s;
869  q = z / s;
870  r = sqrt(p * p + q * q);
871  p /= r;
872  q /= r;
873 
874  // Row modification.
875  for(int j = en - 1; j < n; j++) {
876  z = h[en-1][j];
877  h[en-1][j] = q * z + p * h[en][j];
878  h[en][j] = q * h[en][j] - p * z;
879  }
880 
881  // Column modification.
882  for(int i = 0; i <= en; i++) {
883  z = h[i][en-1];
884  h[i][en-1] = q * z + p * h[i][en];
885  h[i][en] = q * h[i][en] - p * z;
886  }
887 
888  // Accumulate transformations.
889  for(int i = low; i <= upp; i++) {
890  z = vectors[i][en-1];
891  vectors[i][en-1] = q * z + p * vectors[i][en];
892  vectors[i][en] = q * vectors[i][en] - p * z;
893  }
894  } else {
895  // Complex pair.
896  lambda[en-1] = complex<double>(x + p, + z);
897  lambda[en] = complex<double>(x + p, - z);
898  }
899 
900  en--;
901  }
902 
903  // All roots found;
904  // Backsubstitute to find vectors of upper triangular form.
905  if(norm == 0.0) return 0;
906 
907  for(int en = n; en-- > 0;) {
908  p = real(lambda[en]);
909  q = imag(lambda[en]);
910 
911  if(q < 0.0) {
912  // Complex vector.
913  int m = en - 1;
914 
915  // Last vector component chosen imaginary so that eigenvector matrix
916  // is triangular.
917  if(std::abs(h[en][en-1]) > std::abs(h[en-1][en])) {
918  h[en-1][en-1] = q / h[en][en-1];
919  h[en-1][en] = - (h[en][en] - p) / h[en][en-1];
920  } else {
921  cdiv(0.0, - h[en-1][en], h[en-1][en-1] - p, q,
922  h[en-1][en-1], h[en-1][en]);
923  }
924 
925  h[en][en-1] = 0.0;
926  h[en][en] = 1.0;
927 
928  if(en > 1) {
929  for(int i = en - 1; i-- > 0;) {
930  w = h[i][i] - p;
931  double ra = 0.0;
932  double sa = 0.0;
933 
934  for(int j = m; j <= en; j++) {
935  ra += h[i][j] * h[j][en-1];
936  sa += h[i][j] * h[j][en];
937  }
938 
939  if(imag(lambda[i]) < 0.0) {
940  z = w;
941  r = ra;
942  s = sa;
943  } else {
944  m = i;
945  if(imag(lambda[i]) == 0.0) {
946  cdiv(- ra, -sa, w, q, h[i][en-1], h[i][en]);
947  } else {
948  // Solve complex equations.
949  x = h[i][i+1];
950  y = h[i+1][i];
951  double vr = (real(lambda[i]) - p) * (real(lambda[i]) - p) +
952  imag(lambda[i]) * imag(lambda[i]) - q * q;
953  double vi = real(lambda[i] - p) * 2.0 * q;
954 
955  if(vr == 0.0 && vi == 0.0) {
956  tst1 =
957  norm * (std::abs(w) + std::abs(q) + std::abs(x) + std::abs(y) + std::abs(z));
958  vr = tst1;
959 
960  do {
961  vr *= .01;
962  tst2 = tst1 + vr;
963  } while(tst2 > tst1);
964  }
965 
966  cdiv(x * r - z * ra + q * sa, x * s - z * sa - q * ra, vr, vi,
967  h[i][en-1], h[i][en]);
968 
969  if(std::abs(x) > std::abs(z) + std::abs(q)) {
970  h[i+1][en-1] = (- ra - w * h[i][en-1] + q * h[i][en]) / x;
971  h[i+1][en] = (- sa - w * h[i][en] - q * h[i][en-1]) / x;
972  } else {
973  cdiv(- r - y * h[i][en-1], - s - y * h[i][en], z, q,
974  h[i+1][en-1], h[i+1][en]);
975  }
976  }
977 
978  // Overflow control.
979  t = max(std::abs(h[i][en-1]), std::abs(h[i][en]));
980 
981  if(t != 0.0 && (t + 1.0 / t) <= t) {
982  for(int j = i; j <= en; j++) {
983  h[j][en-1] /= t;
984  h[j][en] /= t;
985  }
986  }
987  }
988  }
989  }
990  // End complex vector.
991  } else if(q == 0) {
992  // Real vector.
993  int m = en;
994  h[en][en] = 1.0;
995 
996  if(en > 0) {
997  for(int i = en; i-- > 0;) {
998  w = h[i][i] - p;
999  r = 0.0;
1000  for(int j = m; j <= en; j++) r += h[i][j] * h[j][en];
1001 
1002  if(imag(lambda[i]) < 0.0) {
1003  z = w;
1004  s = r;
1005  } else {
1006  m = i;
1007  if(imag(lambda[i]) == 0.0) {
1008  t = w;
1009  if(t == 0.0) {
1010  t = tst1 = norm;
1011 
1012  do {
1013  t *= .01;
1014  tst2 = norm + t;
1015  } while(tst2 > tst1);
1016  }
1017 
1018  h[i][en] = - r / t;
1019  } else {
1020  // Solve real equations.
1021  x = h[i][i+1];
1022  y = h[i+1][i];
1023  q = (real(lambda[i]) - p) * (real(lambda[i]) - p) +
1024  imag(lambda[i]) * imag(lambda[i]);
1025  t = (x * s - z * r) / q;
1026  h[i][en] = t;
1027  h[i+1][en] = (std::abs(x) > std::abs(z)) ?
1028  (- r - w * t) / x : (- s - y * t) / z;
1029  }
1030 
1031  // Overflow control.
1032  t = std::abs(h[i][en]);
1033  if(t != 0.0 && (t + 1.0 / t) <= t) {
1034  for(int j = i; j <= en; j++) h[j][en] /= t;
1035  }
1036  }
1037  }
1038  }
1039  // End real vector.
1040  }
1041  }
1042 
1043  // End back substitution; vectors of isolated roots.
1044  for(int i = 0; i < low; i++) {
1045  for(int j = i; j < n; j++) vectors[i][j] = h[i][j];
1046  }
1047 
1048  for(int i = upp + 1; i < n; i++) {
1049  for(int j = i; j < n; j++) vectors[i][j] = h[i][j];
1050  }
1051 
1052  // Multiply by transformation matrix to give vectors of original
1053  // full matrix.
1054  for(int j = n; j-- > low;) {
1055  int m = (j < upp) ? j : upp;
1056 
1057  for(int i = low; i <= upp; i++) {
1058  z = 0.0;
1059  for(k = low; k <= m; k++) z += vectors[i][k] * h[k][j];
1060  vectors[i][j] = z;
1061  }
1062  }
1063 
1064  return 0;
1065 }
static void balance(Matrix< double > &, int &low, int &high, Vector< double > &)
PETE_TUTree< FnAbs, typename T::PETE_Expr_t > abs(const PETE_Expr< T > &l)
Double eigenvector routines.
Definition: DoubleEigen.h:36
void cdiv(double ar, double ai, double br, double bi, double &cr, double &ci)
Definition: FDoubleEigen.h:103
Size error exception.
Definition: SizeError.h:33
static void exchange(Matrix< double > &, int j, int m, int low, int high)
T::PETE_Expr_t::PETE_Return_t max(const PETE_Expr< T > &expr, NDIndex< D > &loc)
Definition: ReductionLoc.h:123
Matrix< double > vectors
Definition: DoubleEigen.h:92
Matrix< double > packedEigenVectors() const
Get eigenvectors.
int ncols() const
Get number of columns.
Definition: Array2D.h:307
void elmtran(Matrix< double > &, int low, int high, Array1D< int > &index)
Vector< complex< double > > lambda
Definition: DoubleEigen.h:91
FLieGenerator< T, N > real(const FLieGenerator< std::complex< T >, N > &)
Take real part of a complex generator.
Matrix< complex< double > > eigenVectors() const
Get eigenvectors.
constexpr double c
The velocity of light in m/s.
Definition: Physics.h:52
int hqr(Matrix< double > &, int low, int high)
int hqr2(Matrix< double > &, int low, int high)
Eigenvalue error exception.
Tps< T > sqrt(const Tps< T > &x)
Square root.
Definition: TpsMath.h:91
Definition: Vec.h:21
static void elmhes(Matrix< double > &, int low, int high, Array1D< int > &index)
FLieGenerator< T, N > imag(const FLieGenerator< std::complex< T >, N > &)
Take imaginary part of a complex generator.
int nrows() const
Get number of rows.
Definition: Array2D.h:301
Vector< complex< double > > eigenValues() const
Get eigenvalues.
void balbak(int low, int high, Vector< double > &scale)
void swapRows(int row1, int row2)
Exchange rows.
Definition: Array2D.h:465