My Project
Loading...
Searching...
No Matches
tnt_linalg.h
1#ifndef TNT_LINALG_H_
2#define TNT_LINALG_H_
3
4#include <algorithm>
5// for min(), max() below
6
7#include <cmath>
8// for abs(), sqrt() below
9
10
11namespace TNT {
12namespace Linear_Algebra {
13
14 using namespace std;
15
49template <class Real>
51{
52 Matrix<Real> L_; // lower triangular factor
53 int isspd; // 1 if matrix to be factored was SPD
54
55public:
56
57 Cholesky();
58 Cholesky(const Matrix<Real> &A);
59 Matrix<Real> getL() const;
62 int is_spd() const;
63
64};
65
66template <class Real>
67Cholesky<Real>::Cholesky() : L_(Real(0.0)), isspd(0) {}
68
73template <class Real>
75{
76 return isspd;
77}
78
82template <class Real>
84{
85 return L_;
86}
87
94template <class Real>
96{
97
98
99 int m = A.num_rows();
100 int n = A.num_cols();
101
102 isspd = (m == n);
103
104 if (m != n)
105 {
106 L_ = Matrix<Real>(Real(0.0));
107 return;
108 }
109
110 L_ = Matrix<Real>(n,n);
111
112
113 // Main loop.
114 for (int j = 0; j < n; j++)
115 {
116 Real d = Real(0.0);
117 for (int k = 0; k < j; k++)
118 {
119 Real s = Real(0.0);
120 for (int i = 0; i < k; i++)
121 {
122 s += L_[k][i]*L_[j][i];
123 }
124 L_[j][k] = s = (A[j][k] - s)/L_[k][k];
125 d = d + s*s;
126 isspd = isspd && (A[k][j] == A[j][k]);
127 }
128 d = A[j][j] - d;
129 isspd = isspd && (d > Real(0.0));
130 L_[j][j] = sqrt(d > Real(0.0) ? d : Real(0.0));
131 for (int k = j+1; k < n; k++)
132 {
133 L_[j][k] = Real(0.0);
134 }
135 }
136}
137
148template <class Real>
150{
151 int n = L_.num_rows();
152 if (b.dim() != n)
153 return Vector<Real>();
154
155
156 Vector<Real> x = b;
157
158
159 // Solve L*y = b;
160 for (int k = 0; k < n; k++)
161 {
162 for (int i = 0; i < k; i++)
163 x[k] -= x[i]*L_[k][i];
164 x[k] /= L_[k][k];
165
166 }
167
168 // Solve L'*X = Y;
169 for (int k = n-1; k >= 0; k--)
170 {
171 for (int i = k+1; i < n; i++)
172 x[k] -= x[i]*L_[i][k];
173 x[k] /= L_[k][k];
174 }
175
176 return x;
177}
178
179
190template <class Real>
192{
193 int n = L_.num_rows();
194 if (B.num_rows() != n)
195 return Matrix<Real>();
196
197
198 Matrix<Real> X = B;
199 int nx = B.num_cols();
200
201 // Solve L*y = b;
202 for (int j=0; j< nx; j++)
203 {
204 for (int k = 0; k < n; k++)
205 {
206 for (int i = 0; i < k; i++)
207 X[k][j] -= X[i][j]*L_[k][i];
208 X[k][j] /= L_[k][k];
209 }
210 }
211
212 // Solve L'*X = Y;
213 for (int j=0; j<nx; j++)
214 {
215 for (int k = n-1; k >= 0; k--)
216 {
217 for (int i = k+1; i < n; i++)
218 X[k][j] -= X[i][j]*L_[i][k];
219 X[k][j] /= L_[k][k];
220 }
221 }
222
223
224
225 return X;
226}
227
228
229
230
282template <class Real>
284{
285
286
288 int n;
289
290 int issymmetric; /* boolean*/
291
294 Vector<Real> d; /* real part */
295 Vector<Real> e; /* img part */
296
298 Matrix<Real> V;
299
300 /* Array for internal storage of nonsymmetric Hessenberg form.
301 @serial internal storage of nonsymmetric Hessenberg form.
302 */
303 Matrix<Real> H;
304
305
306 /* Working storage for nonsymmetric algorithm.
307 @serial working storage for nonsymmetric algorithm.
308 */
309 Vector<Real> ort;
310
311
312 // Symmetric Householder reduction to tridiagonal form.
313
314 void tred2() {
315
316 // This is derived from the Algol procedures tred2 by
317 // Bowdler, Martin, Reinsch, and Wilkinson, Handbook for
318 // Auto. Comp., Vol.ii-Linear Algebra, and the corresponding
319 // Fortran subroutine in EISPACK.
320
321 for (int j = 0; j < n; j++) {
322 d[j] = V[n-1][j];
323 }
324
325 // Householder reduction to tridiagonal form.
326
327 for (int i = n-1; i > 0; i--) {
328
329 // Scale to avoid under/overflow.
330
331 Real scale = Real(0.0);
332 Real h = Real(0.0);
333 for (int k = 0; k < i; k++) {
334 scale = scale + abs(d[k]);
335 }
336 if (scale == Real(0.0)) {
337 e[i] = d[i-1];
338 for (int j = 0; j < i; j++) {
339 d[j] = V[i-1][j];
340 V[i][j] = Real(0.0);
341 V[j][i] = Real(0.0);
342 }
343 } else {
344
345 // Generate Householder vector.
346
347 for (int k = 0; k < i; k++) {
348 d[k] /= scale;
349 h += d[k] * d[k];
350 }
351 Real f = d[i-1];
352 Real g = sqrt(h);
353 if (f > 0) {
354 g = -g;
355 }
356 e[i] = scale * g;
357 h = h - f * g;
358 d[i-1] = f - g;
359 for (int j = 0; j < i; j++) {
360 e[j] = Real(0.0);
361 }
362
363 // Apply similarity transformation to remaining columns.
364
365 for (int j = 0; j < i; j++) {
366 f = d[j];
367 V[j][i] = f;
368 g = e[j] + V[j][j] * f;
369 for (int k = j+1; k <= i-1; k++) {
370 g += V[k][j] * d[k];
371 e[k] += V[k][j] * f;
372 }
373 e[j] = g;
374 }
375 f = Real(0.0);
376 for (int j = 0; j < i; j++) {
377 e[j] /= h;
378 f += e[j] * d[j];
379 }
380 Real hh = f / (h + h);
381 for (int j = 0; j < i; j++) {
382 e[j] -= hh * d[j];
383 }
384 for (int j = 0; j < i; j++) {
385 f = d[j];
386 g = e[j];
387 for (int k = j; k <= i-1; k++) {
388 V[k][j] -= (f * e[k] + g * d[k]);
389 }
390 d[j] = V[i-1][j];
391 V[i][j] = Real(0.0);
392 }
393 }
394 d[i] = h;
395 }
396
397 // Accumulate transformations.
398
399 for (int i = 0; i < n-1; i++) {
400 V[n-1][i] = V[i][i];
401 V[i][i] = Real(1.0);
402 Real h = d[i+1];
403 if (h != Real(0.0)) {
404 for (int k = 0; k <= i; k++) {
405 d[k] = V[k][i+1] / h;
406 }
407 for (int j = 0; j <= i; j++) {
408 Real g = Real(0.0);
409 for (int k = 0; k <= i; k++) {
410 g += V[k][i+1] * V[k][j];
411 }
412 for (int k = 0; k <= i; k++) {
413 V[k][j] -= g * d[k];
414 }
415 }
416 }
417 for (int k = 0; k <= i; k++) {
418 V[k][i+1] = Real(0.0);
419 }
420 }
421 for (int j = 0; j < n; j++) {
422 d[j] = V[n-1][j];
423 V[n-1][j] = Real(0.0);
424 }
425 V[n-1][n-1] = Real(1.0);
426 e[0] = Real(0.0);
427 }
428
429 // Symmetric tridiagonal QL algorithm.
430
431 void tql2 () {
432
433 // This is derived from the Algol procedures tql2, by
434 // Bowdler, Martin, Reinsch, and Wilkinson, Handbook for
435 // Auto. Comp., Vol.ii-Linear Algebra, and the corresponding
436 // Fortran subroutine in EISPACK.
437
438 for (int i = 1; i < n; i++) {
439 e[i-1] = e[i];
440 }
441 e[n-1] = Real(0.0);
442
443 Real f = Real(0.0);
444 Real tst1 = Real(0.0);
445 Real eps = pow(2.0,-52.0);
446 for (int l = 0; l < n; l++) {
447
448 // Find small subdiagonal element
449
450 tst1 = max(tst1,abs(d[l]) + abs(e[l]));
451 int m = l;
452
453 // Original while-loop from Java code
454 while (m < n) {
455 if (abs(e[m]) <= eps*tst1) {
456 break;
457 }
458 m++;
459 }
460
461
462 // If m == l, d[l] is an eigenvalue,
463 // otherwise, iterate.
464
465 if (m > l) {
466 int iter = 0;
467 do {
468 iter = iter + 1; // (Could check iteration count here.)
469
470 // Compute implicit shift
471
472 Real g = d[l];
473 Real p = (d[l+1] - g) / (2.0 * e[l]);
474 Real r = hypot(p, static_cast<Real>(Real(1.0)));
475 if (p < 0) {
476 r = -r;
477 }
478 d[l] = e[l] / (p + r);
479 d[l+1] = e[l] * (p + r);
480 Real dl1 = d[l+1];
481 Real h = g - d[l];
482 for (int i = l+2; i < n; i++) {
483 d[i] -= h;
484 }
485 f = f + h;
486
487 // Implicit QL transformation.
488
489 p = d[m];
490 Real c = Real(1.0);
491 Real c2 = c;
492 Real c3 = c;
493 Real el1 = e[l+1];
494 Real s = Real(0.0);
495 Real s2 = Real(0.0);
496 for (int i = m-1; i >= l; i--) {
497 c3 = c2;
498 c2 = c;
499 s2 = s;
500 g = c * e[i];
501 h = c * p;
502 r = hypot(p,e[i]);
503 e[i+1] = s * r;
504 s = e[i] / r;
505 c = p / r;
506 p = c * d[i] - s * g;
507 d[i+1] = h + s * (c * g + s * d[i]);
508
509 // Accumulate transformation.
510
511 for (int k = 0; k < n; k++) {
512 h = V[k][i+1];
513 V[k][i+1] = s * V[k][i] + c * h;
514 V[k][i] = c * V[k][i] - s * h;
515 }
516 }
517 p = -s * s2 * c3 * el1 * e[l] / dl1;
518 e[l] = s * p;
519 d[l] = c * p;
520
521 // Check for convergence.
522
523 } while (abs(e[l]) > eps*tst1);
524 }
525 d[l] = d[l] + f;
526 e[l] = Real(0.0);
527 }
528
529 // Sort eigenvalues and corresponding vectors.
530
531 for (int i = 0; i < n-1; i++) {
532 int k = i;
533 Real p = d[i];
534 for (int j = i+1; j < n; j++) {
535 if (d[j] < p) {
536 k = j;
537 p = d[j];
538 }
539 }
540 if (k != i) {
541 d[k] = d[i];
542 d[i] = p;
543 for (int j = 0; j < n; j++) {
544 p = V[j][i];
545 V[j][i] = V[j][k];
546 V[j][k] = p;
547 }
548 }
549 }
550 }
551
552 // Nonsymmetric reduction to Hessenberg form.
553
554 void orthes () {
555
556 // This is derived from the Algol procedures orthes and ortran,
557 // by Martin and Wilkinson, Handbook for Auto. Comp.,
558 // Vol.ii-Linear Algebra, and the corresponding
559 // Fortran subroutines in EISPACK.
560
561 int low = 0;
562 int high = n-1;
563
564 for (int m = low+1; m <= high-1; m++) {
565
566 // Scale column.
567
568 Real scale = Real(0.0);
569 for (int i = m; i <= high; i++) {
570 scale = scale + abs(H[i][m-1]);
571 }
572 if (scale != Real(0.0)) {
573
574 // Compute Householder transformation.
575
576 Real h = Real(0.0);
577 for (int i = high; i >= m; i--) {
578 ort[i] = H[i][m-1]/scale;
579 h += ort[i] * ort[i];
580 }
581 Real g = sqrt(h);
582 if (ort[m] > 0) {
583 g = -g;
584 }
585 h = h - ort[m] * g;
586 ort[m] = ort[m] - g;
587
588 // Apply Householder similarity transformation
589 // H = (I-u*u'/h)*H*(I-u*u')/h)
590
591 for (int j = m; j < n; j++) {
592 Real f = Real(0.0);
593 for (int i = high; i >= m; i--) {
594 f += ort[i]*H[i][j];
595 }
596 f = f/h;
597 for (int i = m; i <= high; i++) {
598 H[i][j] -= f*ort[i];
599 }
600 }
601
602 for (int i = 0; i <= high; i++) {
603 Real f = Real(0.0);
604 for (int j = high; j >= m; j--) {
605 f += ort[j]*H[i][j];
606 }
607 f = f/h;
608 for (int j = m; j <= high; j++) {
609 H[i][j] -= f*ort[j];
610 }
611 }
612 ort[m] = scale*ort[m];
613 H[m][m-1] = scale*g;
614 }
615 }
616
617 // Accumulate transformations (Algol's ortran).
618
619 for (int i = 0; i < n; i++) {
620 for (int j = 0; j < n; j++) {
621 V[i][j] = (i == j ? Real(1.0) : Real(0.0));
622 }
623 }
624
625 for (int m = high-1; m >= low+1; m--) {
626 if (H[m][m-1] != Real(0.0)) {
627 for (int i = m+1; i <= high; i++) {
628 ort[i] = H[i][m-1];
629 }
630 for (int j = m; j <= high; j++) {
631 Real g = Real(0.0);
632 for (int i = m; i <= high; i++) {
633 g += ort[i] * V[i][j];
634 }
635 // Double division avoids possible underflow
636 g = (g / ort[m]) / H[m][m-1];
637 for (int i = m; i <= high; i++) {
638 V[i][j] += g * ort[i];
639 }
640 }
641 }
642 }
643 }
644
645
646 // Complex scalar division.
647
648 Real cdivr, cdivi;
649 void cdiv(Real xr, Real xi, Real yr, Real yi) {
650 Real r,d;
651 if (abs(yr) > abs(yi)) {
652 r = yi/yr;
653 d = yr + r*yi;
654 cdivr = (xr + r*xi)/d;
655 cdivi = (xi - r*xr)/d;
656 } else {
657 r = yr/yi;
658 d = yi + r*yr;
659 cdivr = (r*xr + xi)/d;
660 cdivi = (r*xi - xr)/d;
661 }
662 }
663
664
665 // Nonsymmetric reduction from Hessenberg to real Schur form.
666
667 void hqr2 () {
668
669 // This is derived from the Algol procedure hqr2,
670 // by Martin and Wilkinson, Handbook for Auto. Comp.,
671 // Vol.ii-Linear Algebra, and the corresponding
672 // Fortran subroutine in EISPACK.
673
674 // Initialize
675
676 int nn = this->n;
677 int n = nn-1;
678 int low = 0;
679 int high = nn-1;
680 Real eps = pow(2.0,-52.0);
681 Real exshift = Real(0.0);
682 Real p=0,q=0,r=0,s=0,z=0,t,w,x,y;
683
684 // Store roots isolated by balanc and compute matrix norm
685
686 Real norm = Real(0.0);
687 for (int i = 0; i < nn; i++) {
688 if ((i < low) || (i > high)) {
689 d[i] = H[i][i];
690 e[i] = Real(0.0);
691 }
692 for (int j = max(i-1,0); j < nn; j++) {
693 norm = norm + abs(H[i][j]);
694 }
695 }
696
697 // Outer loop over eigenvalue index
698
699 int iter = 0;
700 while (n >= low) {
701
702 // Look for single small sub-diagonal element
703
704 int l = n;
705 while (l > low) {
706 s = abs(H[l-1][l-1]) + abs(H[l][l]);
707 if (s == Real(0.0)) {
708 s = norm;
709 }
710 if (abs(H[l][l-1]) < eps * s) {
711 break;
712 }
713 l--;
714 }
715
716 // Check for convergence
717 // One root found
718
719 if (l == n) {
720 H[n][n] = H[n][n] + exshift;
721 d[n] = H[n][n];
722 e[n] = Real(0.0);
723 n--;
724 iter = 0;
725
726 // Two roots found
727
728 } else if (l == n-1) {
729 w = H[n][n-1] * H[n-1][n];
730 p = (H[n-1][n-1] - H[n][n]) / 2.0;
731 q = p * p + w;
732 z = sqrt(abs(q));
733 H[n][n] = H[n][n] + exshift;
734 H[n-1][n-1] = H[n-1][n-1] + exshift;
735 x = H[n][n];
736
737 // Real pair
738
739 if (q >= 0) {
740 if (p >= 0) {
741 z = p + z;
742 } else {
743 z = p - z;
744 }
745 d[n-1] = x + z;
746 d[n] = d[n-1];
747 if (z != Real(0.0)) {
748 d[n] = x - w / z;
749 }
750 e[n-1] = Real(0.0);
751 e[n] = Real(0.0);
752 x = H[n][n-1];
753 s = abs(x) + abs(z);
754 p = x / s;
755 q = z / s;
756 r = sqrt(p * p+q * q);
757 p = p / r;
758 q = q / r;
759
760 // Row modification
761
762 for (int j = n-1; j < nn; j++) {
763 z = H[n-1][j];
764 H[n-1][j] = q * z + p * H[n][j];
765 H[n][j] = q * H[n][j] - p * z;
766 }
767
768 // Column modification
769
770 for (int i = 0; i <= n; i++) {
771 z = H[i][n-1];
772 H[i][n-1] = q * z + p * H[i][n];
773 H[i][n] = q * H[i][n] - p * z;
774 }
775
776 // Accumulate transformations
777
778 for (int i = low; i <= high; i++) {
779 z = V[i][n-1];
780 V[i][n-1] = q * z + p * V[i][n];
781 V[i][n] = q * V[i][n] - p * z;
782 }
783
784 // Complex pair
785
786 } else {
787 d[n-1] = x + p;
788 d[n] = x + p;
789 e[n-1] = z;
790 e[n] = -z;
791 }
792 n = n - 2;
793 iter = 0;
794
795 // No convergence yet
796
797 } else {
798
799 // Form shift
800
801 x = H[n][n];
802 y = Real(0.0);
803 w = Real(0.0);
804 if (l < n) {
805 y = H[n-1][n-1];
806 w = H[n][n-1] * H[n-1][n];
807 }
808
809 // Wilkinson's original ad hoc shift
810
811 if (iter == 10) {
812 exshift += x;
813 for (int i = low; i <= n; i++) {
814 H[i][i] -= x;
815 }
816 s = abs(H[n][n-1]) + abs(H[n-1][n-2]);
817 x = y = 0.75 * s;
818 w = -0.4375 * s * s;
819 }
820
821 // MATLAB's new ad hoc shift
822
823 if (iter == 30) {
824 s = (y - x) / 2.0;
825 s = s * s + w;
826 if (s > 0) {
827 s = sqrt(s);
828 if (y < x) {
829 s = -s;
830 }
831 s = x - w / ((y - x) / 2.0 + s);
832 for (int i = low; i <= n; i++) {
833 H[i][i] -= s;
834 }
835 exshift += s;
836 x = y = w = 0.964;
837 }
838 }
839
840 iter = iter + 1; // (Could check iteration count here.)
841
842 // Look for two consecutive small sub-diagonal elements
843
844 int m = n-2;
845 while (m >= l) {
846 z = H[m][m];
847 r = x - z;
848 s = y - z;
849 p = (r * s - w) / H[m+1][m] + H[m][m+1];
850 q = H[m+1][m+1] - z - r - s;
851 r = H[m+2][m+1];
852 s = abs(p) + abs(q) + abs(r);
853 p = p / s;
854 q = q / s;
855 r = r / s;
856 if (m == l) {
857 break;
858 }
859 if (abs(H[m][m-1]) * (abs(q) + abs(r)) <
860 eps * (abs(p) * (abs(H[m-1][m-1]) + abs(z) +
861 abs(H[m+1][m+1])))) {
862 break;
863 }
864 m--;
865 }
866
867 for (int i = m+2; i <= n; i++) {
868 H[i][i-2] = Real(0.0);
869 if (i > m+2) {
870 H[i][i-3] = Real(0.0);
871 }
872 }
873
874 // Double QR step involving rows l:n and columns m:n
875
876 for (int k = m; k <= n-1; k++) {
877 int notlast = (k != n-1);
878 if (k != m) {
879 p = H[k][k-1];
880 q = H[k+1][k-1];
881 r = (notlast ? H[k+2][k-1] : Real(0.0));
882 x = abs(p) + abs(q) + abs(r);
883 if (x != Real(0.0)) {
884 p = p / x;
885 q = q / x;
886 r = r / x;
887 }
888 }
889 if (x == Real(0.0)) {
890 break;
891 }
892 s = sqrt(p * p + q * q + r * r);
893 if (p < 0) {
894 s = -s;
895 }
896 if (s != 0) {
897 if (k != m) {
898 H[k][k-1] = -s * x;
899 } else if (l != m) {
900 H[k][k-1] = -H[k][k-1];
901 }
902 p = p + s;
903 x = p / s;
904 y = q / s;
905 z = r / s;
906 q = q / p;
907 r = r / p;
908
909 // Row modification
910
911 for (int j = k; j < nn; j++) {
912 p = H[k][j] + q * H[k+1][j];
913 if (notlast) {
914 p = p + r * H[k+2][j];
915 H[k+2][j] = H[k+2][j] - p * z;
916 }
917 H[k][j] = H[k][j] - p * x;
918 H[k+1][j] = H[k+1][j] - p * y;
919 }
920
921 // Column modification
922
923 for (int i = 0; i <= min(n,k+3); i++) {
924 p = x * H[i][k] + y * H[i][k+1];
925 if (notlast) {
926 p = p + z * H[i][k+2];
927 H[i][k+2] = H[i][k+2] - p * r;
928 }
929 H[i][k] = H[i][k] - p;
930 H[i][k+1] = H[i][k+1] - p * q;
931 }
932
933 // Accumulate transformations
934
935 for (int i = low; i <= high; i++) {
936 p = x * V[i][k] + y * V[i][k+1];
937 if (notlast) {
938 p = p + z * V[i][k+2];
939 V[i][k+2] = V[i][k+2] - p * r;
940 }
941 V[i][k] = V[i][k] - p;
942 V[i][k+1] = V[i][k+1] - p * q;
943 }
944 } // (s != 0)
945 } // k loop
946 } // check convergence
947 } // while (n >= low)
948
949 // Backsubstitute to find vectors of upper triangular form
950
951 if (norm == Real(0.0)) {
952 return;
953 }
954
955 for (n = nn-1; n >= 0; n--) {
956 p = d[n];
957 q = e[n];
958
959 // Real vector
960
961 if (q == 0) {
962 int l = n;
963 H[n][n] = Real(1.0);
964 for (int i = n-1; i >= 0; i--) {
965 w = H[i][i] - p;
966 r = Real(0.0);
967 for (int j = l; j <= n; j++) {
968 r = r + H[i][j] * H[j][n];
969 }
970 if (e[i] < Real(0.0)) {
971 z = w;
972 s = r;
973 } else {
974 l = i;
975 if (e[i] == Real(0.0)) {
976 if (w != Real(0.0)) {
977 H[i][n] = -r / w;
978 } else {
979 H[i][n] = -r / (eps * norm);
980 }
981
982 // Solve real equations
983
984 } else {
985 x = H[i][i+1];
986 y = H[i+1][i];
987 q = (d[i] - p) * (d[i] - p) + e[i] * e[i];
988 t = (x * s - z * r) / q;
989 H[i][n] = t;
990 if (abs(x) > abs(z)) {
991 H[i+1][n] = (-r - w * t) / x;
992 } else {
993 H[i+1][n] = (-s - y * t) / z;
994 }
995 }
996
997 // Overflow control
998
999 t = abs(H[i][n]);
1000 if ((eps * t) * t > 1) {
1001 for (int j = i; j <= n; j++) {
1002 H[j][n] = H[j][n] / t;
1003 }
1004 }
1005 }
1006 }
1007
1008 // Complex vector
1009
1010 } else if (q < 0) {
1011 int l = n-1;
1012
1013 // Last vector component imaginary so matrix is triangular
1014
1015 if (abs(H[n][n-1]) > abs(H[n-1][n])) {
1016 H[n-1][n-1] = q / H[n][n-1];
1017 H[n-1][n] = -(H[n][n] - p) / H[n][n-1];
1018 } else {
1019 cdiv(Real(0.0),-H[n-1][n],H[n-1][n-1]-p,q);
1020 H[n-1][n-1] = cdivr;
1021 H[n-1][n] = cdivi;
1022 }
1023 H[n][n-1] = Real(0.0);
1024 H[n][n] = Real(1.0);
1025 for (int i = n-2; i >= 0; i--) {
1026 Real ra,sa,vr,vi;
1027 ra = Real(0.0);
1028 sa = Real(0.0);
1029 for (int j = l; j <= n; j++) {
1030 ra = ra + H[i][j] * H[j][n-1];
1031 sa = sa + H[i][j] * H[j][n];
1032 }
1033 w = H[i][i] - p;
1034
1035 if (e[i] < Real(0.0)) {
1036 z = w;
1037 r = ra;
1038 s = sa;
1039 } else {
1040 l = i;
1041 if (e[i] == 0) {
1042 cdiv(-ra,-sa,w,q);
1043 H[i][n-1] = cdivr;
1044 H[i][n] = cdivi;
1045 } else {
1046
1047 // Solve complex equations
1048
1049 x = H[i][i+1];
1050 y = H[i+1][i];
1051 vr = (d[i] - p) * (d[i] - p) + e[i] * e[i] - q * q;
1052 vi = (d[i] - p) * 2.0 * q;
1053 if ((vr == Real(0.0)) && (vi == Real(0.0))) {
1054 vr = eps * norm * (abs(w) + abs(q) +
1055 abs(x) + abs(y) + abs(z));
1056 }
1057 cdiv(x*r-z*ra+q*sa,x*s-z*sa-q*ra,vr,vi);
1058 H[i][n-1] = cdivr;
1059 H[i][n] = cdivi;
1060 if (abs(x) > (abs(z) + abs(q))) {
1061 H[i+1][n-1] = (-ra - w * H[i][n-1] + q * H[i][n]) / x;
1062 H[i+1][n] = (-sa - w * H[i][n] - q * H[i][n-1]) / x;
1063 } else {
1064 cdiv(-r-y*H[i][n-1],-s-y*H[i][n],z,q);
1065 H[i+1][n-1] = cdivr;
1066 H[i+1][n] = cdivi;
1067 }
1068 }
1069
1070 // Overflow control
1071
1072 t = max(abs(H[i][n-1]),abs(H[i][n]));
1073 if ((eps * t) * t > 1) {
1074 for (int j = i; j <= n; j++) {
1075 H[j][n-1] = H[j][n-1] / t;
1076 H[j][n] = H[j][n] / t;
1077 }
1078 }
1079 }
1080 }
1081 }
1082 }
1083
1084 // Vectors of isolated roots
1085
1086 for (int i = 0; i < nn; i++) {
1087 if (i < low || i > high) {
1088 for (int j = i; j < nn; j++) {
1089 V[i][j] = H[i][j];
1090 }
1091 }
1092 }
1093
1094 // Back transformation to get eigenvectors of original matrix
1095
1096 for (int j = nn-1; j >= low; j--) {
1097 for (int i = low; i <= high; i++) {
1098 z = Real(0.0);
1099 for (int k = low; k <= min(j,high); k++) {
1100 z = z + V[i][k] * H[k][j];
1101 }
1102 V[i][j] = z;
1103 }
1104 }
1105 }
1106
1107public:
1108
1109
1115 n = A.num_cols();
1116 V = Matrix<Real>(n,n);
1117 d = Vector<Real>(n);
1118 e = Vector<Real>(n);
1119
1120 issymmetric = 1;
1121 for (int j = 0; (j < n) && issymmetric; j++) {
1122 for (int i = 0; (i < n) && issymmetric; i++) {
1123 issymmetric = (A[i][j] == A[j][i]);
1124 }
1125 }
1126
1127 if (issymmetric) {
1128 for (int i = 0; i < n; i++) {
1129 for (int j = 0; j < n; j++) {
1130 V[i][j] = A[i][j];
1131 }
1132 }
1133
1134 // Tridiagonalize.
1135 tred2();
1136
1137 // Diagonalize.
1138 tql2();
1139
1140 } else {
1141 H = Matrix<Real>(n,n);
1142 ort = Vector<Real>(n);
1143
1144 for (int j = 0; j < n; j++) {
1145 for (int i = 0; i < n; i++) {
1146 H[i][j] = A[i][j];
1147 }
1148 }
1149
1150 // Reduce to Hessenberg form.
1151 orthes();
1152
1153 // Reduce Hessenberg to real Schur form.
1154 hqr2();
1155 }
1156 }
1157
1158
1163 void getV (Matrix<Real> &V_) {
1164 V_ = V;
1165 return;
1166 }
1167
1173 d_ = d;
1174 return ;
1175 }
1176
1183 e_ = e;
1184 return;
1185 }
1186
1187
1221 void getD (Matrix<Real> &D) {
1222 D = Matrix<Real>(n,n);
1223 for (int i = 0; i < n; i++) {
1224 for (int j = 0; j < n; j++) {
1225 D[i][j] = Real(0.0);
1226 }
1227 D[i][i] = d[i];
1228 if (e[i] > 0) {
1229 D[i][i+1] = e[i];
1230 } else if (e[i] < 0) {
1231 D[i][i-1] = e[i];
1232 }
1233 }
1234 }
1235};
1236
1237
1250template <class Real>
1251class LU
1252{
1253
1254 private:
1255
1256 /* Array for internal storage of decomposition. */
1257 Matrix<Real> LU_;
1258 int m, n, pivsign;
1259 Vector<int> piv;
1260
1261
1262 Matrix<Real> permute_copy( const Matrix<Real> &A,
1263 const Vector<int> &piv, int j0, int j1)
1264 {
1265 int piv_length = piv.dim();
1266
1267 Matrix<Real> X(piv_length, j1-j0+1);
1268
1269
1270 for (int i = 0; i < piv_length; i++)
1271 for (int j = j0; j <= j1; j++)
1272 X[i][j-j0] = A[piv[i]][j];
1273
1274 return X;
1275 }
1276
1277 Vector<Real> permute_copy( const Vector<Real> &A,
1278 const Vector<int> &piv)
1279 {
1280 int piv_length = piv.dim();
1281 if (piv_length != A.dim())
1282 return Vector<Real>();
1283
1284 Vector<Real> x(piv_length);
1285
1286
1287 for (int i = 0; i < piv_length; i++)
1288 x[i] = A[piv[i]];
1289
1290 return x;
1291 }
1292
1293
1294 public :
1295
1301 LU (const Matrix<Real> &A) : LU_(A), m(A.num_rows()), n(A.num_cols()),
1302 piv(A.num_rows())
1303
1304 {
1305
1306 // Use a "left-looking", dot-product, Crout/Doolittle algorithm.
1307
1308
1309 for (int i = 0; i < m; i++) {
1310 piv[i] = i;
1311 }
1312 pivsign = 1;
1313 Real *LUrowi = 0;;
1314 Vector<Real> LUcolj(m);
1315
1316 // Outer loop.
1317
1318 for (int j = 0; j < n; j++) {
1319
1320 // Make a copy of the j-th column to localize references.
1321
1322 for (int i = 0; i < m; i++) {
1323 LUcolj[i] = LU_[i][j];
1324 }
1325
1326 // Apply previous transformations.
1327
1328 for (int i = 0; i < m; i++) {
1329 LUrowi = LU_[i];
1330
1331 // Most of the time is spent in the following dot product.
1332
1333 int kmax = min(i,j);
1334 double s = Real(0.0);
1335 for (int k = 0; k < kmax; k++) {
1336 s += LUrowi[k]*LUcolj[k];
1337 }
1338
1339 LUrowi[j] = LUcolj[i] -= s;
1340 }
1341
1342 // Find pivot and exchange if necessary.
1343
1344 int p = j;
1345 for (int i = j+1; i < m; i++) {
1346 if (abs(LUcolj[i]) > abs(LUcolj[p])) {
1347 p = i;
1348 }
1349 }
1350 if (p != j) {
1351 int k=0;
1352 for (k = 0; k < n; k++) {
1353 double t = LU_[p][k];
1354 LU_[p][k] = LU_[j][k];
1355 LU_[j][k] = t;
1356 }
1357 k = piv[p];
1358 piv[p] = piv[j];
1359 piv[j] = k;
1360 pivsign = -pivsign;
1361 }
1362
1363 // Compute multipliers.
1364
1365 if ((j < m) && (LU_[j][j] != Real(0.0))) {
1366 for (int i = j+1; i < m; i++) {
1367 LU_[i][j] /= LU_[j][j];
1368 }
1369 }
1370 }
1371 }
1372
1373
1380 for (int j = 0; j < n; j++) {
1381 if (LU_[j][j] == 0)
1382 return 0;
1383 }
1384 return 1;
1385 }
1386
1392 Matrix<Real> L_(m,n);
1393 for (int i = 0; i < m; i++) {
1394 for (int j = 0; j < n; j++) {
1395 if (i > j) {
1396 L_[i][j] = LU_[i][j];
1397 } else if (i == j) {
1398 L_[i][j] = Real(1.0);
1399 } else {
1400 L_[i][j] = Real(0.0);
1401 }
1402 }
1403 }
1404 return L_;
1405 }
1406
1412 Matrix<Real> U_(n,n);
1413 for (int i = 0; i < n; i++) {
1414 for (int j = 0; j < n; j++) {
1415 if (i <= j) {
1416 U_[i][j] = LU_[i][j];
1417 } else {
1418 U_[i][j] = Real(0.0);
1419 }
1420 }
1421 }
1422 return U_;
1423 }
1424
1430 return piv;
1431 }
1432
1433
1438 Real det () {
1439 if (m != n) {
1440 return Real(0);
1441 }
1442 Real d = Real(pivsign);
1443 for (int j = 0; j < n; j++) {
1444 d *= LU_[j][j];
1445 }
1446 return d;
1447 }
1448
1456 {
1457
1458 /* Dimensions: A is mxn, X is nxk, B is mxk */
1459
1460 if (B.num_rows() != m) {
1461 return Matrix<Real>(Real(0.0));
1462 }
1463 if (!isNonsingular()) {
1464 return Matrix<Real>(Real(0.0));
1465 }
1466
1467 // Copy right hand side with pivoting
1468 int nx = B.num_cols();
1469
1470
1471 Matrix<Real> X = permute_copy(B, piv, 0, nx-1);
1472
1473 // Solve L*Y = B(piv,:)
1474 for (int k = 0; k < n; k++) {
1475 for (int i = k+1; i < n; i++) {
1476 for (int j = 0; j < nx; j++) {
1477 X[i][j] -= X[k][j]*LU_[i][k];
1478 }
1479 }
1480 }
1481 // Solve U*X = Y;
1482 for (int k = n-1; k >= 0; k--) {
1483 for (int j = 0; j < nx; j++) {
1484 X[k][j] /= LU_[k][k];
1485 }
1486 for (int i = 0; i < k; i++) {
1487 for (int j = 0; j < nx; j++) {
1488 X[i][j] -= X[k][j]*LU_[i][k];
1489 }
1490 }
1491 }
1492 return X;
1493 }
1494
1495
1506 {
1507
1508 /* Dimensions: A is mxn, X is nxk, B is mxk */
1509
1510 if (b.dim() != m) {
1511 return Vector<Real>();
1512 }
1513 if (!isNonsingular()) {
1514 return Vector<Real>();
1515 }
1516
1517
1518 Vector<Real> x = permute_copy(b, piv);
1519
1520 // Solve L*Y = B(piv)
1521 for (int k = 0; k < n; k++) {
1522 for (int i = k+1; i < n; i++) {
1523 x[i] -= x[k]*LU_[i][k];
1524 }
1525 }
1526
1527 // Solve U*X = Y;
1528 for (int k = n-1; k >= 0; k--) {
1529 x[k] /= LU_[k][k];
1530 for (int i = 0; i < k; i++)
1531 x[i] -= x[k]*LU_[i][k];
1532 }
1533
1534
1535 return x;
1536 }
1537
1538}; /* class LU */
1539
1540
1564template <class Real>
1565class QR {
1566
1567
1568 /* Array for internal storage of decomposition.
1569 @serial internal array storage.
1570 */
1571
1572 Matrix<Real> QR_;
1573
1574 /* Row and column dimensions.
1575 @serial column dimension.
1576 @serial row dimension.
1577 */
1578 int m, n;
1579
1580 /* Array for internal storage of diagonal of R.
1581 @serial diagonal of R.
1582 */
1583 Vector<Real> Rdiag;
1584
1585
1586public:
1587
1593 QR(const Matrix<Real> &A) /* constructor */
1594 {
1595 QR_ = A;
1596 m = A.num_rows();
1597 n = A.num_cols();
1598 Rdiag = Vector<Real>(n);
1599 int i=0, j=0, k=0;
1600
1601 // Main loop.
1602 for (k = 0; k < n; k++) {
1603 // Compute 2-norm of k-th column without under/overflow.
1604 Real nrm = 0;
1605 for (i = k; i < m; i++) {
1606 nrm = hypot(nrm,QR_[i][k]);
1607 }
1608
1609 if (nrm != Real(0.0)) {
1610 // Form k-th Householder vector.
1611 if (QR_[k][k] < 0) {
1612 nrm = -nrm;
1613 }
1614 for (i = k; i < m; i++) {
1615 QR_[i][k] /= nrm;
1616 }
1617 QR_[k][k] += Real(1.0);
1618
1619 // Apply transformation to remaining columns.
1620 for (j = k+1; j < n; j++) {
1621 Real s = Real(0.0);
1622 for (i = k; i < m; i++) {
1623 s += QR_[i][k]*QR_[i][j];
1624 }
1625 s = -s/QR_[k][k];
1626 for (i = k; i < m; i++) {
1627 QR_[i][j] += s*QR_[i][k];
1628 }
1629 }
1630 }
1631 Rdiag[k] = -nrm;
1632 }
1633 }
1634
1635
1641 int isFullRank() const
1642 {
1643 for (int j = 0; j < n; j++)
1644 {
1645 if (Rdiag[j] == 0)
1646 return 0;
1647 }
1648 return 1;
1649 }
1650
1651
1652
1653
1661 {
1662 Matrix<Real> H(m,n);
1663
1664 /* note: H is completely filled in by algorithm, so
1665 initializaiton of H is not necessary.
1666 */
1667 for (int i = 0; i < m; i++)
1668 {
1669 for (int j = 0; j < n; j++)
1670 {
1671 if (i >= j) {
1672 H[i][j] = QR_[i][j];
1673 } else {
1674 H[i][j] = Real(0.0);
1675 }
1676 }
1677 }
1678 return H;
1679 }
1680
1681
1682
1688 {
1689 Matrix<Real> R(n,n);
1690 for (int i = 0; i < n; i++) {
1691 for (int j = 0; j < n; j++) {
1692 if (i < j) {
1693 R[i][j] = QR_[i][j];
1694 } else if (i == j) {
1695 R[i][j] = Rdiag[i];
1696 } else {
1697 R[i][j] = Real(0.0);
1698 }
1699 }
1700 }
1701 return R;
1702 }
1703
1704
1705
1706
1707
1713 {
1714 int i=0, j=0, k=0;
1715
1716 Matrix<Real> Q(m,n);
1717 for (k = n-1; k >= 0; k--) {
1718 for (i = 0; i < m; i++) {
1719 Q[i][k] = Real(0.0);
1720 }
1721 Q[k][k] = Real(1.0);
1722 for (j = k; j < n; j++) {
1723 if (QR_[k][k] != 0) {
1724 Real s = Real(0.0);
1725 for (i = k; i < m; i++) {
1726 s += QR_[i][k]*Q[i][j];
1727 }
1728 s = -s/QR_[k][k];
1729 for (i = k; i < m; i++) {
1730 Q[i][j] += s*QR_[i][k];
1731 }
1732 }
1733 }
1734 }
1735 return Q;
1736 }
1737
1738
1747 {
1748 if (b.dim() != m) /* arrays must be conformant */
1749 return Vector<Real>();
1750
1751 if ( !isFullRank() ) /* matrix is rank deficient */
1752 {
1753 return Vector<Real>();
1754 }
1755
1756 Vector<Real> x = b;
1757
1758 // Compute Y = transpose(Q)*b
1759 for (int k = 0; k < n; k++)
1760 {
1761 Real s = Real(0.0);
1762 for (int i = k; i < m; i++)
1763 {
1764 s += QR_[i][k]*x[i];
1765 }
1766 s = -s/QR_[k][k];
1767 for (int i = k; i < m; i++)
1768 {
1769 x[i] += s*QR_[i][k];
1770 }
1771 }
1772 // Solve R*X = Y;
1773 for (int k = n-1; k >= 0; k--)
1774 {
1775 x[k] /= Rdiag[k];
1776 for (int i = 0; i < k; i++) {
1777 x[i] -= x[k]*QR_[i][k];
1778 }
1779 }
1780
1781
1782 /* return n x nx portion of X */
1783 Vector<Real> x_(n);
1784 for (int i=0; i<n; i++)
1785 x_[i] = x[i];
1786
1787 return x_;
1788 }
1789
1798 {
1799 if (B.num_rows() != m) /* arrays must be conformant */
1800 return Matrix<Real>(Real(0.0));
1801
1802 if ( !isFullRank() ) /* matrix is rank deficient */
1803 {
1804 return Matrix<Real>(Real(0.0));
1805 }
1806
1807 int nx = B.num_cols();
1808 Matrix<Real> X = B;
1809 int i=0, j=0, k=0;
1810
1811 // Compute Y = transpose(Q)*B
1812 for (k = 0; k < n; k++) {
1813 for (j = 0; j < nx; j++) {
1814 Real s = Real(0.0);
1815 for (i = k; i < m; i++) {
1816 s += QR_[i][k]*X[i][j];
1817 }
1818 s = -s/QR_[k][k];
1819 for (i = k; i < m; i++) {
1820 X[i][j] += s*QR_[i][k];
1821 }
1822 }
1823 }
1824 // Solve R*X = Y;
1825 for (k = n-1; k >= 0; k--) {
1826 for (j = 0; j < nx; j++) {
1827 X[k][j] /= Rdiag[k];
1828 }
1829 for (i = 0; i < k; i++) {
1830 for (j = 0; j < nx; j++) {
1831 X[i][j] -= X[k][j]*QR_[i][k];
1832 }
1833 }
1834 }
1835
1836
1837 /* return n x nx portion of X */
1838 Matrix<Real> X_(n,nx);
1839 for (i=0; i<n; i++)
1840 for (j=0; j<nx; j++)
1841 X_[i][j] = X[i][j];
1842
1843 return X_;
1844 }
1845
1846
1847};
1848
1849
1867template <class Real>
1868class SVD
1869{
1870
1871
1872 Matrix<Real> U, V;
1873 Vector<Real> s;
1874 int m, n;
1875
1876 public:
1877
1878
1879 SVD (const Matrix<Real> &Arg) {
1880
1881
1882 m = Arg.num_rows();
1883 n = Arg.num_cols();
1884 int nu = min(m,n);
1885 s = Vector<Real>(min(m+1,n));
1886 U = Matrix<Real>(m, nu, Real(0));
1887 V = Matrix<Real>(n,n);
1888 Vector<Real> e(n);
1889 Vector<Real> work(m);
1890 Matrix<Real> A(Arg);
1891 int wantu = 1; /* boolean */
1892 int wantv = 1; /* boolean */
1893 int i=0, j=0, k=0;
1894
1895 // Reduce A to bidiagonal form, storing the diagonal elements
1896 // in s and the super-diagonal elements in e.
1897
1898 int nct = min(m-1,n);
1899 int nrt = max(0,min(n-2,m));
1900 for (k = 0; k < max(nct,nrt); k++) {
1901 if (k < nct) {
1902
1903 // Compute the transformation for the k-th column and
1904 // place the k-th diagonal in s[k].
1905 // Compute 2-norm of k-th column without under/overflow.
1906 s[k] = 0;
1907 for (i = k; i < m; i++) {
1908 s[k] = hypot(s[k],A[i][k]);
1909 }
1910 if (s[k] != Real(0.0)) {
1911 if (A[k][k] < Real(0.0)) {
1912 s[k] = -s[k];
1913 }
1914 for (i = k; i < m; i++) {
1915 A[i][k] /= s[k];
1916 }
1917 A[k][k] += Real(1.0);
1918 }
1919 s[k] = -s[k];
1920 }
1921 for (j = k+1; j < n; j++) {
1922 if ((k < nct) && (s[k] != Real(0.0))) {
1923
1924 // Apply the transformation.
1925
1926 double t = 0;
1927 for (i = k; i < m; i++) {
1928 t += A[i][k]*A[i][j];
1929 }
1930 t = -t/A[k][k];
1931 for (i = k; i < m; i++) {
1932 A[i][j] += t*A[i][k];
1933 }
1934 }
1935
1936 // Place the k-th row of A into e for the
1937 // subsequent calculation of the row transformation.
1938
1939 e[j] = A[k][j];
1940 }
1941 if (wantu & (k < nct)) {
1942
1943 // Place the transformation in U for subsequent back
1944 // multiplication.
1945
1946 for (i = k; i < m; i++) {
1947 U[i][k] = A[i][k];
1948 }
1949 }
1950 if (k < nrt) {
1951
1952 // Compute the k-th row transformation and place the
1953 // k-th super-diagonal in e[k].
1954 // Compute 2-norm without under/overflow.
1955 e[k] = 0;
1956 for (i = k+1; i < n; i++) {
1957 e[k] = hypot(e[k],e[i]);
1958 }
1959 if (e[k] != Real(0.0)) {
1960 if (e[k+1] < Real(0.0)) {
1961 e[k] = -e[k];
1962 }
1963 for (i = k+1; i < n; i++) {
1964 e[i] /= e[k];
1965 }
1966 e[k+1] += Real(1.0);
1967 }
1968 e[k] = -e[k];
1969 if ((k+1 < m) & (e[k] != Real(0.0))) {
1970
1971 // Apply the transformation.
1972
1973 for (i = k+1; i < m; i++) {
1974 work[i] = Real(0.0);
1975 }
1976 for (j = k+1; j < n; j++) {
1977 for (i = k+1; i < m; i++) {
1978 work[i] += e[j]*A[i][j];
1979 }
1980 }
1981 for (j = k+1; j < n; j++) {
1982 double t = -e[j]/e[k+1];
1983 for (i = k+1; i < m; i++) {
1984 A[i][j] += t*work[i];
1985 }
1986 }
1987 }
1988 if (wantv) {
1989
1990 // Place the transformation in V for subsequent
1991 // back multiplication.
1992
1993 for (i = k+1; i < n; i++) {
1994 V[i][k] = e[i];
1995 }
1996 }
1997 }
1998 }
1999
2000 // Set up the final bidiagonal matrix or order p.
2001
2002 int p = min(n,m+1);
2003 if (nct < n) {
2004 s[nct] = A[nct][nct];
2005 }
2006 if (m < p) {
2007 s[p-1] = Real(0.0);
2008 }
2009 if (nrt+1 < p) {
2010 e[nrt] = A[nrt][p-1];
2011 }
2012 e[p-1] = Real(0.0);
2013
2014 // If required, generate U.
2015
2016 if (wantu) {
2017 for (j = nct; j < nu; j++) {
2018 for (i = 0; i < m; i++) {
2019 U[i][j] = Real(0.0);
2020 }
2021 U[j][j] = Real(1.0);
2022 }
2023 for (k = nct-1; k >= 0; k--) {
2024 if (s[k] != Real(0.0)) {
2025 for (j = k+1; j < nu; j++) {
2026 double t = 0;
2027 for (i = k; i < m; i++) {
2028 t += U[i][k]*U[i][j];
2029 }
2030 t = -t/U[k][k];
2031 for (i = k; i < m; i++) {
2032 U[i][j] += t*U[i][k];
2033 }
2034 }
2035 for (i = k; i < m; i++ ) {
2036 U[i][k] = -U[i][k];
2037 }
2038 U[k][k] = Real(1.0) + U[k][k];
2039 for (i = 0; i < k-1; i++) {
2040 U[i][k] = Real(0.0);
2041 }
2042 } else {
2043 for (i = 0; i < m; i++) {
2044 U[i][k] = Real(0.0);
2045 }
2046 U[k][k] = Real(1.0);
2047 }
2048 }
2049 }
2050
2051 // If required, generate V.
2052
2053 if (wantv) {
2054 for (k = n-1; k >= 0; k--) {
2055 if ((k < nrt) & (e[k] != Real(0.0))) {
2056 for (j = k+1; j < nu; j++) {
2057 double t = 0;
2058 for (i = k+1; i < n; i++) {
2059 t += V[i][k]*V[i][j];
2060 }
2061 t = -t/V[k+1][k];
2062 for (i = k+1; i < n; i++) {
2063 V[i][j] += t*V[i][k];
2064 }
2065 }
2066 }
2067 for (i = 0; i < n; i++) {
2068 V[i][k] = Real(0.0);
2069 }
2070 V[k][k] = Real(1.0);
2071 }
2072 }
2073
2074 // Main iteration loop for the singular values.
2075
2076 int pp = p-1;
2077 int iter = 0;
2078 double eps = pow(2.0,-52.0);
2079 while (p > 0) {
2080 int k=0;
2081 int kase=0;
2082
2083 // Here is where a test for too many iterations would go.
2084
2085 // This section of the program inspects for
2086 // negligible elements in the s and e arrays. On
2087 // completion the variables kase and k are set as follows.
2088
2089 // kase = 1 if s(p) and e[k-1] are negligible and k<p
2090 // kase = 2 if s(k) is negligible and k<p
2091 // kase = 3 if e[k-1] is negligible, k<p, and
2092 // s(k), ..., s(p) are not negligible (qr step).
2093 // kase = 4 if e(p-1) is negligible (convergence).
2094
2095 for (k = p-2; k >= -1; k--) {
2096 if (k == -1) {
2097 break;
2098 }
2099 if (abs(e[k]) <= eps*(abs(s[k]) + abs(s[k+1]))) {
2100 e[k] = Real(0.0);
2101 break;
2102 }
2103 }
2104 if (k == p-2) {
2105 kase = 4;
2106 } else {
2107 int ks;
2108 for (ks = p-1; ks >= k; ks--) {
2109 if (ks == k) {
2110 break;
2111 }
2112 double t = (ks != p ? abs(e[ks]) : 0.) +
2113 (ks != k+1 ? abs(e[ks-1]) : 0.);
2114 if (abs(s[ks]) <= eps*t) {
2115 s[ks] = Real(0.0);
2116 break;
2117 }
2118 }
2119 if (ks == k) {
2120 kase = 3;
2121 } else if (ks == p-1) {
2122 kase = 1;
2123 } else {
2124 kase = 2;
2125 k = ks;
2126 }
2127 }
2128 k++;
2129
2130 // Perform the task indicated by kase.
2131
2132 switch (kase) {
2133
2134 // Deflate negligible s(p).
2135
2136 case 1: {
2137 double f = e[p-2];
2138 e[p-2] = Real(0.0);
2139 for (j = p-2; j >= k; j--) {
2140 double t = hypot(s[j],f);
2141 double cs = s[j]/t;
2142 double sn = f/t;
2143 s[j] = t;
2144 if (j != k) {
2145 f = -sn*e[j-1];
2146 e[j-1] = cs*e[j-1];
2147 }
2148 if (wantv) {
2149 for (i = 0; i < n; i++) {
2150 t = cs*V[i][j] + sn*V[i][p-1];
2151 V[i][p-1] = -sn*V[i][j] + cs*V[i][p-1];
2152 V[i][j] = t;
2153 }
2154 }
2155 }
2156 }
2157 break;
2158
2159 // Split at negligible s(k).
2160
2161 case 2: {
2162 double f = e[k-1];
2163 e[k-1] = Real(0.0);
2164 for (j = k; j < p; j++) {
2165 double t = hypot(s[j],f);
2166 double cs = s[j]/t;
2167 double sn = f/t;
2168 s[j] = t;
2169 f = -sn*e[j];
2170 e[j] = cs*e[j];
2171 if (wantu) {
2172 for (i = 0; i < m; i++) {
2173 t = cs*U[i][j] + sn*U[i][k-1];
2174 U[i][k-1] = -sn*U[i][j] + cs*U[i][k-1];
2175 U[i][j] = t;
2176 }
2177 }
2178 }
2179 }
2180 break;
2181
2182 // Perform one qr step.
2183
2184 case 3: {
2185
2186 // Calculate the shift.
2187
2188 double scale = max(max(max(max(
2189 abs(s[p-1]),abs(s[p-2])),abs(e[p-2])),
2190 abs(s[k])),abs(e[k]));
2191 double sp = s[p-1]/scale;
2192 double spm1 = s[p-2]/scale;
2193 double epm1 = e[p-2]/scale;
2194 double sk = s[k]/scale;
2195 double ek = e[k]/scale;
2196 double b = ((spm1 + sp)*(spm1 - sp) + epm1*epm1)/2.0;
2197 double c = (sp*epm1)*(sp*epm1);
2198 double shift = Real(0.0);
2199 if ((b != Real(0.0)) || (c != Real(0.0))) {
2200 shift = sqrt(b*b + c);
2201 if (b < Real(0.0)) {
2202 shift = -shift;
2203 }
2204 shift = c/(b + shift);
2205 }
2206 double f = (sk + sp)*(sk - sp) + shift;
2207 double g = sk*ek;
2208
2209 // Chase zeros.
2210
2211 for (j = k; j < p-1; j++) {
2212 double t = hypot(f,g);
2213 double cs = f/t;
2214 double sn = g/t;
2215 if (j != k) {
2216 e[j-1] = t;
2217 }
2218 f = cs*s[j] + sn*e[j];
2219 e[j] = cs*e[j] - sn*s[j];
2220 g = sn*s[j+1];
2221 s[j+1] = cs*s[j+1];
2222 if (wantv) {
2223 for (i = 0; i < n; i++) {
2224 t = cs*V[i][j] + sn*V[i][j+1];
2225 V[i][j+1] = -sn*V[i][j] + cs*V[i][j+1];
2226 V[i][j] = t;
2227 }
2228 }
2229 t = hypot(f,g);
2230 cs = f/t;
2231 sn = g/t;
2232 s[j] = t;
2233 f = cs*e[j] + sn*s[j+1];
2234 s[j+1] = -sn*e[j] + cs*s[j+1];
2235 g = sn*e[j+1];
2236 e[j+1] = cs*e[j+1];
2237 if (wantu && (j < m-1)) {
2238 for (i = 0; i < m; i++) {
2239 t = cs*U[i][j] + sn*U[i][j+1];
2240 U[i][j+1] = -sn*U[i][j] + cs*U[i][j+1];
2241 U[i][j] = t;
2242 }
2243 }
2244 }
2245 e[p-2] = f;
2246 iter = iter + 1;
2247 }
2248 break;
2249
2250 // Convergence.
2251
2252 case 4: {
2253
2254 // Make the singular values positive.
2255
2256 if (s[k] <= Real(0.0)) {
2257 s[k] = (s[k] < Real(0.0) ? -s[k] : Real(0.0));
2258 if (wantv) {
2259 for (i = 0; i <= pp; i++) {
2260 V[i][k] = -V[i][k];
2261 }
2262 }
2263 }
2264
2265 // Order the singular values.
2266
2267 while (k < pp) {
2268 if (s[k] >= s[k+1]) {
2269 break;
2270 }
2271 double t = s[k];
2272 s[k] = s[k+1];
2273 s[k+1] = t;
2274 if (wantv && (k < n-1)) {
2275 for (i = 0; i < n; i++) {
2276 t = V[i][k+1]; V[i][k+1] = V[i][k]; V[i][k] = t;
2277 }
2278 }
2279 if (wantu && (k < m-1)) {
2280 for (i = 0; i < m; i++) {
2281 t = U[i][k+1]; U[i][k+1] = U[i][k]; U[i][k] = t;
2282 }
2283 }
2284 k++;
2285 }
2286 iter = 0;
2287 p--;
2288 }
2289 break;
2290 }
2291 }
2292 }
2293
2294
2295 void getU (Matrix<Real> &A)
2296 {
2297 int minm = min(m+1,n);
2298
2299 A = Matrix<Real>(m, minm);
2300
2301 for (int i=0; i<m; i++)
2302 for (int j=0; j<minm; j++)
2303 A[i][j] = U[i][j];
2304
2305 }
2306
2307 /* Return the right singular vectors */
2308
2309 void getV (Matrix<Real> &A)
2310 {
2311 A = V;
2312 }
2313
2317 {
2318 x = s;
2319 }
2320
2325 void getS (Matrix<Real> &A) {
2326 A = Matrix<Real>(n,n);
2327 for (int i = 0; i < n; i++) {
2328 for (int j = 0; j < n; j++) {
2329 A[i][j] = Real(0.0);
2330 }
2331 A[i][i] = s[i];
2332 }
2333 }
2334
2337 double norm2 () {
2338 return s[0];
2339 }
2340
2343 double cond () {
2344 return s[0]/s[min(m,n)-1];
2345 }
2346
2351 int rank ()
2352 {
2353 double eps = pow(2.0,-52.0);
2354 double tol = max(m,n)*s[0]*eps;
2355 int r = 0;
2356 for (int i = 0; i < s.dim(); i++) {
2357 if (s[i] > tol) {
2358 r++;
2359 }
2360 }
2361 return r;
2362 }
2363};
2364
2365} /* namespace Linear_Algebra */
2366} /* namespace TNT */
2367
2368#endif
2369// TNT_LINALG_H_
Definition tnt_linalg.h:51
Matrix< Real > getL() const
Definition tnt_linalg.h:83
int is_spd() const
Definition tnt_linalg.h:74
Vector< Real > solve(const Vector< Real > &B)
Definition tnt_linalg.h:149
Definition tnt_linalg.h:284
void getD(Matrix< Real > &D)
Definition tnt_linalg.h:1221
void getImagEigenvalues(Vector< Real > &e_)
Definition tnt_linalg.h:1182
void getRealEigenvalues(Vector< Real > &d_)
Definition tnt_linalg.h:1172
void getV(Matrix< Real > &V_)
Definition tnt_linalg.h:1163
Eigenvalue(const Matrix< Real > &A)
Definition tnt_linalg.h:1114
Definition tnt_linalg.h:1252
int isNonsingular()
Definition tnt_linalg.h:1379
Vector< Real > solve(const Vector< Real > &b)
Definition tnt_linalg.h:1505
Matrix< Real > getU()
Definition tnt_linalg.h:1411
Vector< int > getPivot()
Definition tnt_linalg.h:1429
Matrix< Real > getL()
Definition tnt_linalg.h:1391
Matrix< Real > solve(const Matrix< Real > &B)
Definition tnt_linalg.h:1455
LU(const Matrix< Real > &A)
Definition tnt_linalg.h:1301
Real det()
Definition tnt_linalg.h:1438
Definition tnt_linalg.h:1565
Matrix< Real > getQ() const
Definition tnt_linalg.h:1712
Matrix< Real > getR() const
Definition tnt_linalg.h:1687
int isFullRank() const
Definition tnt_linalg.h:1641
Matrix< Real > getHouseholder(void) const
Definition tnt_linalg.h:1660
Matrix< Real > solve(const Matrix< Real > &B) const
Definition tnt_linalg.h:1797
Vector< Real > solve(const Vector< Real > &b) const
Definition tnt_linalg.h:1746
Definition tnt_linalg.h:1869
int rank()
Definition tnt_linalg.h:2351
double cond()
Definition tnt_linalg.h:2343
void getSingularValues(Vector< Real > &x)
Definition tnt_linalg.h:2316
double norm2()
Definition tnt_linalg.h:2337
void getS(Matrix< Real > &A)
Definition tnt_linalg.h:2325
Definition tnt_matrix.h:57
Definition tnt_vector.h:49
Definition tnt_array1d.h:36
double norm(const Matrix< T > &A)
Definition tnt_matrix.h:996
Real hypot(const Real &a, const Real &b)
Definition tnt_math_utils.h:21