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