|
bpp-core
2.1.0
|
00001 // 00002 // File: EigenValue.h 00003 // Created by: Julien Dutheil 00004 // Created on: Tue Apr 7 16:24 2004 00005 // 00006 00007 /* 00008 Copyright or © or Copr. Bio++ Development Team, (November 17, 2004) 00009 00010 This software is a computer program whose purpose is to provide classes 00011 for numerical calculus. This file is part of the Bio++ project. 00012 00013 This software is governed by the CeCILL license under French law and 00014 abiding by the rules of distribution of free software. You can use, 00015 modify and/ or redistribute the software under the terms of the CeCILL 00016 license as circulated by CEA, CNRS and INRIA at the following URL 00017 "http://www.cecill.info". 00018 00019 As a counterpart to the access to the source code and rights to copy, 00020 modify and redistribute granted by the license, users are provided only 00021 with a limited warranty and the software's author, the holder of the 00022 economic rights, and the successive licensors have only limited 00023 liability. 00024 00025 In this respect, the user's attention is drawn to the risks associated 00026 with loading, using, modifying and/or developing or reproducing the 00027 software by the user in light of its specific status of free software, 00028 that may mean that it is complicated to manipulate, and that also 00029 therefore means that it is reserved for developers and experienced 00030 professionals having in-depth computer knowledge. Users are therefore 00031 encouraged to load and test the software's suitability as regards their 00032 requirements in conditions enabling the security of their systems and/or 00033 data to be ensured and, more generally, to use and operate it in the 00034 same conditions as regards security. 00035 00036 The fact that you are presently reading this means that you have had 00037 knowledge of the CeCILL license and that you accept its terms. 00038 */ 00039 00040 00041 00042 #ifndef _EIGENVALUE_H_ 00043 #define _EIGENVALUE_H_ 00044 00045 #include <algorithm> 00046 // for min(), max() below 00047 00048 #include <cmath> 00049 // for abs() below 00050 #include <climits> 00051 00052 #include "Matrix.h" 00053 #include "../NumTools.h" 00054 00055 namespace bpp 00056 { 00057 00103 template <class Real> 00104 class EigenValue 00105 { 00106 00107 private: 00111 size_t n_; 00112 00116 bool issymmetric_; 00117 00123 std::vector<Real> d_; /* real part */ 00124 std::vector<Real> e_; /* img part */ 00130 RowMatrix<Real> V_; 00131 00137 RowMatrix<Real> H_; 00138 00139 00145 mutable RowMatrix<Real> D_; 00146 00152 std::vector<Real> ort_; 00153 00162 void tred2() 00163 { 00164 for (size_t j = 0; j < n_; j++) 00165 { 00166 d_[j] = V_(n_-1,j); 00167 } 00168 00169 // Householder reduction to tridiagonal form. 00170 00171 for (int i = static_cast<int>(n_)-1; i > 0; i--) 00172 { 00173 // Scale to avoid under/overflow. 00174 00175 Real scale = 0.0; 00176 Real h = 0.0; 00177 for (int k = 0; k < i; k++) 00178 { 00179 scale = scale + NumTools::abs<Real>(d_[k]); 00180 } 00181 if (scale == 0.0) 00182 { 00183 e_[i] = d_[i-1]; 00184 for (int j = 0; j < i; j++) 00185 { 00186 d_[j] = V_(i-1,j); 00187 V_(i,j) = 0.0; 00188 V_(j,i) = 0.0; 00189 } 00190 } 00191 else 00192 { 00193 // Generate Householder vector. 00194 00195 for (int k = 0; k < i; k++) 00196 { 00197 d_[k] /= scale; 00198 h += d_[k] * d_[k]; 00199 } 00200 Real f = d_[i-1]; 00201 Real g = sqrt(h); 00202 if (f > 0) 00203 { 00204 g = -g; 00205 } 00206 e_[i] = scale * g; 00207 h = h - f * g; 00208 d_[i-1] = f - g; 00209 for (int j = 0; j < i; j++) 00210 { 00211 e_[j] = 0.0; 00212 } 00213 00214 // Apply similarity transformation to remaining columns. 00215 00216 for (int j = 0; j < i; j++) 00217 { 00218 f = d_[j]; 00219 V_(j,i) = f; 00220 g = e_[j] + V_(j,j) * f; 00221 for (int k = j+1; k <= i-1; k++) 00222 { 00223 g += V_(k,j) * d_[k]; 00224 e_[k] += V_(k,j) * f; 00225 } 00226 e_[j] = g; 00227 } 00228 f = 0.0; 00229 for (int j = 0; j < i; j++) 00230 { 00231 e_[j] /= h; 00232 f += e_[j] * d_[j]; 00233 } 00234 Real hh = f / (h + h); 00235 for (int j = 0; j < i; j++) 00236 { 00237 e_[j] -= hh * d_[j]; 00238 } 00239 for (int j = 0; j < i; j++) 00240 { 00241 f = d_[j]; 00242 g = e_[j]; 00243 for (int k = j; k <= i-1; k++) 00244 { 00245 V_(k,j) -= (f * e_[k] + g * d_[k]); 00246 } 00247 d_[j] = V_(i-1,j); 00248 V_(i,j) = 0.0; 00249 } 00250 } 00251 d_[i] = h; 00252 } 00253 00254 // Accumulate transformations. 00255 00256 for (size_t i = 0; i < n_-1; i++) 00257 { 00258 V_(n_-1,i) = V_(i,i); 00259 V_(i,i) = 1.0; 00260 Real h = d_[i+1]; 00261 if (h != 0.0) 00262 { 00263 for (size_t k = 0; k <= i; k++) 00264 { 00265 d_[k] = V_(k,i+1) / h; 00266 } 00267 for (size_t j = 0; j <= i; j++) 00268 { 00269 Real g = 0.0; 00270 for (size_t k = 0; k <= i; k++) 00271 { 00272 g += V_(k,i+1) * V_(k,j); 00273 } 00274 for (size_t k = 0; k <= i; k++) 00275 { 00276 V_(k,j) -= g * d_[k]; 00277 } 00278 } 00279 } 00280 for (size_t k = 0; k <= i; k++) 00281 { 00282 V_(k,i+1) = 0.0; 00283 } 00284 } 00285 for (size_t j = 0; j < n_; j++) 00286 { 00287 d_[j] = V_(n_-1,j); 00288 V_(n_-1,j) = 0.0; 00289 } 00290 V_(n_-1,n_-1) = 1.0; 00291 e_[0] = 0.0; 00292 } 00293 00302 void tql2 () 00303 { 00304 for (size_t i = 1; i < n_; i++) 00305 { 00306 e_[i-1] = e_[i]; 00307 } 00308 e_[n_-1] = 0.0; 00309 00310 Real f = 0.0; 00311 Real tst1 = 0.0; 00312 Real eps = pow(2.0,-52.0); 00313 for (int l = 0; l < static_cast<int>(n_); l++) 00314 { 00315 // Find small subdiagonal element 00316 00317 tst1 = std::max(tst1, NumTools::abs<Real>(d_[l]) + NumTools::abs<Real>(e_[l])); 00318 int m = l; 00319 00320 // Original while-loop from Java code 00321 while (m < static_cast<int>(n_)) 00322 { 00323 if (NumTools::abs<Real>(e_[m]) <= eps*tst1) 00324 { 00325 break; 00326 } 00327 m++; 00328 } 00329 00330 // If m == l, d_[l] is an eigenvalue, 00331 // otherwise, iterate. 00332 00333 if (m > l) 00334 { 00335 int iter = 0; 00336 do 00337 { 00338 iter = iter + 1; // (Could check iteration count here.) 00339 00340 // Compute implicit shift 00341 00342 Real g = d_[l]; 00343 Real p = (d_[l+1] - g) / (2.0 * e_[l]); 00344 Real r = hypot(p,1.0); 00345 if (p < 0) 00346 { 00347 r = -r; 00348 } 00349 d_[l] = e_[l] / (p + r); 00350 d_[l+1] = e_[l] * (p + r); 00351 Real dl1 = d_[l+1]; 00352 Real h = g - d_[l]; 00353 for (size_t i = l+2; i < n_; i++) 00354 { 00355 d_[i] -= h; 00356 } 00357 f = f + h; 00358 00359 // Implicit QL transformation. 00360 00361 p = d_[m]; 00362 Real c = 1.0; 00363 Real c2 = c; 00364 Real c3 = c; 00365 Real el1 = e_[l+1]; 00366 Real s = 0.0; 00367 Real s2 = 0.0; 00368 for (int i = m-1; i >= l; i--) 00369 { 00370 c3 = c2; 00371 c2 = c; 00372 s2 = s; 00373 g = c * e_[i]; 00374 h = c * p; 00375 r = hypot(p,e_[i]); 00376 e_[i+1] = s * r; 00377 s = e_[i] / r; 00378 c = p / r; 00379 p = c * d_[i] - s * g; 00380 d_[i+1] = h + s * (c * g + s * d_[i]); 00381 00382 // Accumulate transformation. 00383 00384 for (size_t k = 0; k < n_; k++) 00385 { 00386 h = V_(k,i+1); 00387 V_(k,i+1) = s * V_(k,i) + c * h; 00388 V_(k,i) = c * V_(k,i) - s * h; 00389 } 00390 } 00391 p = -s * s2 * c3 * el1 * e_[l] / dl1; 00392 e_[l] = s * p; 00393 d_[l] = c * p; 00394 00395 // Check for convergence. 00396 00397 } while (NumTools::abs<Real>(e_[l]) > eps*tst1); 00398 } 00399 d_[l] = d_[l] + f; 00400 e_[l] = 0.0; 00401 } 00402 00403 // Sort eigenvalues and corresponding vectors. 00404 00405 for (size_t i = 0; n_ > 0 && i < n_-1; i++) 00406 { 00407 size_t k = i; 00408 Real p = d_[i]; 00409 for (size_t j = i+1; j < n_; j++) 00410 { 00411 if (d_[j] < p) 00412 { 00413 k = j; 00414 p = d_[j]; 00415 } 00416 } 00417 if (k != i) 00418 { 00419 d_[k] = d_[i]; 00420 d_[i] = p; 00421 for (size_t j = 0; j < n_; j++) 00422 { 00423 p = V_(j,i); 00424 V_(j,i) = V_(j,k); 00425 V_(j,k) = p; 00426 } 00427 } 00428 } 00429 } 00430 00439 void orthes () 00440 { 00441 int low = 0; 00442 int high = static_cast<int>(n_)-1; 00443 00444 for (int m = low+1; m <= high-1; m++) 00445 { 00446 // Scale column. 00447 00448 Real scale = 0.0; 00449 for (int i = m; i <= high; i++) 00450 { 00451 scale = scale + NumTools::abs<Real>(H_(i,m-1)); 00452 } 00453 if (scale != 0.0) 00454 { 00455 // Compute Householder transformation. 00456 00457 Real h = 0.0; 00458 for (int i = high; i >= m; i--) 00459 { 00460 ort_[i] = H_(i,m-1)/scale; 00461 h += ort_[i] * ort_[i]; 00462 } 00463 Real g = sqrt(h); 00464 if (ort_[m] > 0) 00465 { 00466 g = -g; 00467 } 00468 h = h - ort_[m] * g; 00469 ort_[m] = ort_[m] - g; 00470 00471 // Apply Householder similarity transformation 00472 // H = (I-u*u'/h)*H*(I-u*u')/h) 00473 00474 for (int j = m; j < static_cast<int>(n_); j++) 00475 { 00476 Real f = 0.0; 00477 for (int i = high; i >= m; i--) 00478 { 00479 f += ort_[i]*H_(i,j); 00480 } 00481 f = f/h; 00482 for (int i = m; i <= high; i++) 00483 { 00484 H_(i,j) -= f*ort_[i]; 00485 } 00486 } 00487 00488 for (int i = 0; i <= high; i++) 00489 { 00490 Real f = 0.0; 00491 for (int j = high; j >= m; j--) 00492 { 00493 f += ort_[j]*H_(i,j); 00494 } 00495 f = f/h; 00496 for (int j = m; j <= high; j++) 00497 { 00498 H_(i,j) -= f*ort_[j]; 00499 } 00500 } 00501 ort_[m] = scale*ort_[m]; 00502 H_(m,m-1) = scale*g; 00503 } 00504 } 00505 00506 // Accumulate transformations (Algol's ortran). 00507 00508 for (size_t i = 0; i < n_; i++) 00509 { 00510 for (size_t j = 0; j < n_; j++) 00511 { 00512 V_(i,j) = (i == j ? 1.0 : 0.0); 00513 } 00514 } 00515 00516 for (int m = high-1; m >= low+1; m--) 00517 { 00518 if (H_(m,m-1) != 0.0) 00519 { 00520 for (int i = m+1; i <= high; i++) 00521 { 00522 ort_[i] = H_(i,m-1); 00523 } 00524 for (int j = m; j <= high; j++) 00525 { 00526 Real g = 0.0; 00527 for (int i = m; i <= high; i++) 00528 { 00529 g += ort_[i] * V_(i,j); 00530 } 00531 // Double division avoids possible underflow 00532 g = (g / ort_[m]) / H_(m,m-1); 00533 for (int i = m; i <= high; i++) 00534 { 00535 V_(i,j) += g * ort_[i]; 00536 } 00537 } 00538 } 00539 } 00540 } 00541 00542 00543 // Complex scalar division. 00544 00545 Real cdivr, cdivi; 00546 void cdiv(Real xr, Real xi, Real yr, Real yi) 00547 { 00548 Real r,d; 00549 if (NumTools::abs<Real>(yr) > NumTools::abs<Real>(yi)) 00550 { 00551 r = yi/yr; 00552 d = yr + r*yi; 00553 cdivr = (xr + r*xi)/d; 00554 cdivi = (xi - r*xr)/d; 00555 } 00556 else 00557 { 00558 r = yr/yi; 00559 d = yi + r*yr; 00560 cdivr = (r*xr + xi)/d; 00561 cdivi = (r*xi - xr)/d; 00562 } 00563 } 00564 00565 00566 // Nonsymmetric reduction from Hessenberg to real Schur form. 00567 00568 void hqr2 () 00569 { 00570 // This is derived from the Algol procedure hqr2, 00571 // by Martin and Wilkinson, Handbook for Auto. Comp., 00572 // Vol.ii-Linear Algebra, and the corresponding 00573 // Fortran subroutine in EISPACK. 00574 00575 // Initialize 00576 00577 int nn = static_cast<int>(this->n_); 00578 int n = nn-1; 00579 int low = 0; 00580 int high = nn-1; 00581 Real eps = pow(2.0,-52.0); 00582 Real exshift = 0.0; 00583 Real p=0,q=0,r=0,s=0,z=0,t,w,x,y; 00584 00585 // Store roots isolated by balanc and compute matrix norm 00586 00587 Real norm = 0.0; 00588 for (int i = 0; i < nn; i++) 00589 { 00590 if ((i < low) || (i > high)) 00591 { 00592 d_[i] = H_(i,i); 00593 e_[i] = 0.0; 00594 } 00595 for (int j = std::max(i-1,0); j < nn; j++) 00596 { 00597 norm = norm + NumTools::abs<Real>(H_(i,j)); 00598 } 00599 } 00600 00601 // Outer loop over eigenvalue index 00602 00603 int iter = 0; 00604 while (n >= low) 00605 { 00606 // Look for single small sub-diagonal element 00607 00608 int l = n; 00609 while (l > low) 00610 { 00611 s = NumTools::abs<Real>(H_(l-1,l-1)) + NumTools::abs<Real>(H_(l,l)); 00612 if (s == 0.0) 00613 { 00614 s = norm; 00615 } 00616 if (NumTools::abs<Real>(H_(l,l-1)) < eps * s) 00617 { 00618 break; 00619 } 00620 l--; 00621 } 00622 00623 // Check for convergence 00624 // One root found 00625 00626 if (l == n) 00627 { 00628 H_(n,n) = H_(n,n) + exshift; 00629 d_[n] = H_(n,n); 00630 e_[n] = 0.0; 00631 n--; 00632 iter = 0; 00633 00634 // Two roots found 00635 00636 } 00637 else if (l == n-1) 00638 { 00639 w = H_(n,n-1) * H_(n-1,n); 00640 p = (H_(n-1,n-1) - H_(n,n)) / 2.0; 00641 q = p * p + w; 00642 z = sqrt(NumTools::abs<Real>(q)); 00643 H_(n,n) = H_(n,n) + exshift; 00644 H_(n-1,n-1) = H_(n-1,n-1) + exshift; 00645 x = H_(n,n); 00646 00647 // Real pair 00648 00649 if (q >= 0) 00650 { 00651 if (p >= 0) 00652 { 00653 z = p + z; 00654 } 00655 else 00656 { 00657 z = p - z; 00658 } 00659 d_[n-1] = x + z; 00660 d_[n] = d_[n-1]; 00661 if (z != 0.0) 00662 { 00663 d_[n] = x - w / z; 00664 } 00665 e_[n-1] = 0.0; 00666 e_[n] = 0.0; 00667 x = H_(n,n-1); 00668 s = NumTools::abs<Real>(x) + NumTools::abs<Real>(z); 00669 p = x / s; 00670 q = z / s; 00671 r = sqrt(p * p+q * q); 00672 p = p / r; 00673 q = q / r; 00674 00675 // Row modification 00676 00677 for (int j = n-1; j < nn; j++) 00678 { 00679 z = H_(n-1,j); 00680 H_(n-1,j) = q * z + p * H_(n,j); 00681 H_(n,j) = q * H_(n,j) - p * z; 00682 } 00683 00684 // Column modification 00685 00686 for (int i = 0; i <= n; i++) 00687 { 00688 z = H_(i,n-1); 00689 H_(i,n-1) = q * z + p * H_(i,n); 00690 H_(i,n) = q * H_(i,n) - p * z; 00691 } 00692 00693 // Accumulate transformations 00694 00695 for (int i = low; i <= high; i++) 00696 { 00697 z = V_(i,n-1); 00698 V_(i,n-1) = q * z + p * V_(i,n); 00699 V_(i,n) = q * V_(i,n) - p * z; 00700 } 00701 00702 // Complex pair 00703 00704 } 00705 else 00706 { 00707 d_[n-1] = x + p; 00708 d_[n] = x + p; 00709 e_[n-1] = z; 00710 e_[n] = -z; 00711 } 00712 n = n - 2; 00713 iter = 0; 00714 00715 // No convergence yet 00716 00717 } 00718 else 00719 { 00720 // Form shift 00721 00722 x = H_(n,n); 00723 y = 0.0; 00724 w = 0.0; 00725 if (l < n) 00726 { 00727 y = H_(n-1,n-1); 00728 w = H_(n,n-1) * H_(n-1,n); 00729 } 00730 00731 // Wilkinson's original ad hoc shift 00732 00733 if (iter == 10) 00734 { 00735 exshift += x; 00736 for (int i = low; i <= n; i++) 00737 { 00738 H_(i,i) -= x; 00739 } 00740 s = NumTools::abs<Real>(H_(n,n-1)) + NumTools::abs<Real>(H_(n-1,n-2)); 00741 x = y = 0.75 * s; 00742 w = -0.4375 * s * s; 00743 } 00744 00745 // MATLAB's new ad hoc shift 00746 if (iter == 30) 00747 { 00748 s = (y - x) / 2.0; 00749 s = s * s + w; 00750 if (s > 0) 00751 { 00752 s = sqrt(s); 00753 if (y < x) 00754 { 00755 s = -s; 00756 } 00757 s = x - w / ((y - x) / 2.0 + s); 00758 for (int i = low; i <= n; i++) 00759 { 00760 H_(i,i) -= s; 00761 } 00762 exshift += s; 00763 x = y = w = 0.964; 00764 } 00765 } 00766 00767 iter = iter + 1; // (Could check iteration count here.) 00768 00769 // Look for two consecutive small sub-diagonal elements 00770 00771 int m = n-2; 00772 while (m >= l) 00773 { 00774 z = H_(m,m); 00775 r = x - z; 00776 s = y - z; 00777 p = (r * s - w) / H_(m+1,m) + H_(m,m+1); 00778 q = H_(m+1,m+1) - z - r - s; 00779 r = H_(m+2,m+1); 00780 s = NumTools::abs<Real>(p) + NumTools::abs<Real>(q) + NumTools::abs<Real>(r); 00781 p = p / s; 00782 q = q / s; 00783 r = r / s; 00784 if (m == l) 00785 { 00786 break; 00787 } 00788 if (NumTools::abs<Real>(H_(m,m-1)) * (NumTools::abs<Real>(q) + NumTools::abs<Real>(r)) < 00789 eps * (NumTools::abs<Real>(p) * (NumTools::abs<Real>(H_(m-1,m-1)) + NumTools::abs<Real>(z) + 00790 NumTools::abs<Real>(H_(m+1,m+1))))) 00791 { 00792 break; 00793 } 00794 m--; 00795 } 00796 00797 for (int i = m+2; i <= n; i++) 00798 { 00799 H_(i,i-2) = 0.0; 00800 if (i > m+2) 00801 { 00802 H_(i,i-3) = 0.0; 00803 } 00804 } 00805 00806 // Double QR step involving rows l:n and columns m:n 00807 00808 for (int k = m; k <= n-1; k++) 00809 { 00810 int notlast = (k != n-1); 00811 if (k != m) 00812 { 00813 p = H_(k,k-1); 00814 q = H_(k+1,k-1); 00815 r = (notlast ? H_(k+2,k-1) : 0.0); 00816 x = NumTools::abs<Real>(p) + NumTools::abs<Real>(q) + NumTools::abs<Real>(r); 00817 if (x != 0.0) 00818 { 00819 p = p / x; 00820 q = q / x; 00821 r = r / x; 00822 } 00823 } 00824 if (x == 0.0) 00825 { 00826 break; 00827 } 00828 s = sqrt(p * p + q * q + r * r); 00829 if (p < 0) 00830 { 00831 s = -s; 00832 } 00833 if (s != 0) 00834 { 00835 if (k != m) 00836 { 00837 H_(k,k-1) = -s * x; 00838 } 00839 else if (l != m) 00840 { 00841 H_(k,k-1) = -H_(k,k-1); 00842 } 00843 p = p + s; 00844 x = p / s; 00845 y = q / s; 00846 z = r / s; 00847 q = q / p; 00848 r = r / p; 00849 00850 // Row modification 00851 00852 for (int j = k; j < nn; j++) 00853 { 00854 p = H_(k,j) + q * H_(k+1,j); 00855 if (notlast) 00856 { 00857 p = p + r * H_(k+2,j); 00858 H_(k+2,j) = H_(k+2,j) - p * z; 00859 } 00860 H_(k,j) = H_(k,j) - p * x; 00861 H_(k+1,j) = H_(k+1,j) - p * y; 00862 } 00863 00864 // Column modification 00865 00866 for (int i = 0; i <= std::min(n,k+3); i++) 00867 { 00868 p = x * H_(i,k) + y * H_(i,k+1); 00869 if (notlast) 00870 { 00871 p = p + z * H_(i,k+2); 00872 H_(i,k+2) = H_(i,k+2) - p * r; 00873 } 00874 H_(i,k) = H_(i,k) - p; 00875 H_(i,k+1) = H_(i,k+1) - p * q; 00876 } 00877 00878 // Accumulate transformations 00879 00880 for (int i = low; i <= high; i++) 00881 { 00882 p = x * V_(i,k) + y * V_(i,k+1); 00883 if (notlast) 00884 { 00885 p = p + z * V_(i,k+2); 00886 V_(i,k+2) = V_(i,k+2) - p * r; 00887 } 00888 V_(i,k) = V_(i,k) - p; 00889 V_(i,k+1) = V_(i,k+1) - p * q; 00890 } 00891 } // (s != 0) 00892 } // k loop 00893 } // check convergence 00894 } // while (n >= low) 00895 00896 // Backsubstitute to find vectors of upper triangular form 00897 00898 if (norm == 0.0) 00899 { 00900 return; 00901 } 00902 00903 for (n = nn-1; n >= 0; n--) 00904 { 00905 p = d_[n]; 00906 q = e_[n]; 00907 00908 // Real vector 00909 00910 if (q == 0) 00911 { 00912 int l = n; 00913 H_(n,n) = 1.0; 00914 for (int i = n-1; i >= 0; i--) 00915 { 00916 w = H_(i,i) - p; 00917 r = 0.0; 00918 for (int j = l; j <= n; j++) 00919 { 00920 r = r + H_(i,j) * H_(j,n); 00921 } 00922 if (e_[i] < 0.0) 00923 { 00924 z = w; 00925 s = r; 00926 } 00927 else 00928 { 00929 l = i; 00930 if (e_[i] == 0.0) 00931 { 00932 if (w != 0.0) 00933 { 00934 H_(i,n) = -r / w; 00935 } 00936 else 00937 { 00938 H_(i,n) = -r / (eps * norm); 00939 } 00940 00941 // Solve real equations 00942 00943 } 00944 else 00945 { 00946 x = H_(i,i+1); 00947 y = H_(i+1,i); 00948 q = (d_[i] - p) * (d_[i] - p) + e_[i] * e_[i]; 00949 t = (x * s - z * r) / q; 00950 H_(i,n) = t; 00951 if (NumTools::abs<Real>(x) > NumTools::abs<Real>(z)) 00952 { 00953 H_(i+1,n) = (-r - w * t) / x; 00954 } 00955 else 00956 { 00957 H_(i+1,n) = (-s - y * t) / z; 00958 } 00959 } 00960 00961 // Overflow control 00962 00963 t = NumTools::abs<Real>(H_(i,n)); 00964 if ((eps * t) * t > 1) 00965 { 00966 for (int j = i; j <= n; j++) 00967 { 00968 H_(j,n) = H_(j,n) / t; 00969 } 00970 } 00971 } 00972 } 00973 00974 // Complex vector 00975 00976 } 00977 else if (q < 0) 00978 { 00979 int l = n-1; 00980 00981 // Last vector component imaginary so matrix is triangular 00982 00983 if (NumTools::abs<Real>(H_(n,n-1)) > NumTools::abs<Real>(H_(n-1,n))) 00984 { 00985 H_(n-1,n-1) = q / H_(n,n-1); 00986 H_(n-1,n) = -(H_(n,n) - p) / H_(n,n-1); 00987 } 00988 else 00989 { 00990 cdiv(0.0,-H_(n-1,n),H_(n-1,n-1)-p,q); 00991 H_(n-1,n-1) = cdivr; 00992 H_(n-1,n) = cdivi; 00993 } 00994 H_(n,n-1) = 0.0; 00995 H_(n,n) = 1.0; 00996 for (int i = n-2; i >= 0; i--) 00997 { 00998 Real ra,sa,vr,vi; 00999 ra = 0.0; 01000 sa = 0.0; 01001 for (int j = l; j <= n; j++) 01002 { 01003 ra = ra + H_(i,j) * H_(j,n-1); 01004 sa = sa + H_(i,j) * H_(j,n); 01005 } 01006 w = H_(i,i) - p; 01007 01008 if (e_[i] < 0.0) 01009 { 01010 z = w; 01011 r = ra; 01012 s = sa; 01013 } 01014 else 01015 { 01016 l = i; 01017 if (e_[i] == 0) 01018 { 01019 cdiv(-ra,-sa,w,q); 01020 H_(i,n-1) = cdivr; 01021 H_(i,n) = cdivi; 01022 } 01023 else 01024 { 01025 // Solve complex equations 01026 01027 x = H_(i,i+1); 01028 y = H_(i+1,i); 01029 vr = (d_[i] - p) * (d_[i] - p) + e_[i] * e_[i] - q * q; 01030 vi = (d_[i] - p) * 2.0 * q; 01031 if ((vr == 0.0) && (vi == 0.0)) 01032 { 01033 vr = eps * norm * (NumTools::abs<Real>(w) + NumTools::abs<Real>(q) + 01034 NumTools::abs<Real>(x) + NumTools::abs<Real>(y) + NumTools::abs<Real>(z)); 01035 } 01036 cdiv(x*r-z*ra+q*sa,x*s-z*sa-q*ra,vr,vi); 01037 H_(i,n-1) = cdivr; 01038 H_(i,n) = cdivi; 01039 if (NumTools::abs<Real>(x) > (NumTools::abs<Real>(z) + NumTools::abs<Real>(q))) 01040 { 01041 H_(i+1,n-1) = (-ra - w * H_(i,n-1) + q * H_(i,n)) / x; 01042 H_(i+1,n) = (-sa - w * H_(i,n) - q * H_(i,n-1)) / x; 01043 } 01044 else 01045 { 01046 cdiv(-r-y*H_(i,n-1),-s-y*H_(i,n),z,q); 01047 H_(i+1,n-1) = cdivr; 01048 H_(i+1,n) = cdivi; 01049 } 01050 } 01051 01052 // Overflow control 01053 t = std::max(NumTools::abs<Real>(H_(i,n-1)),NumTools::abs<Real>(H_(i,n))); 01054 if ((eps * t) * t > 1) 01055 { 01056 for (int j = i; j <= n; j++) 01057 { 01058 H_(j,n-1) = H_(j,n-1) / t; 01059 H_(j,n) = H_(j,n) / t; 01060 } 01061 } 01062 } 01063 } 01064 } 01065 } 01066 01067 // Vectors of isolated roots 01068 01069 for (int i = 0; i < nn; i++) 01070 { 01071 if (i < low || i > high) 01072 { 01073 for (int j = i; j < nn; j++) 01074 { 01075 V_(i,j) = H_(i,j); 01076 } 01077 } 01078 } 01079 01080 // Back transformation to get eigenvectors of original matrix 01081 01082 for (int j = nn-1; j >= low; j--) 01083 { 01084 for (int i = low; i <= high; i++) 01085 { 01086 z = 0.0; 01087 for (int k = low; k <= std::min(j,high); k++) 01088 { 01089 z = z + V_(i,k) * H_(k,j); 01090 } 01091 V_(i,j) = z; 01092 } 01093 } 01094 } 01095 01096 public: 01097 01098 bool isSymmetric() const { return issymmetric_; } 01099 01100 01106 EigenValue(const Matrix<Real>& A) : 01107 n_(A.getNumberOfColumns()), 01108 issymmetric_(true), 01109 d_(n_), 01110 e_(n_), 01111 V_(n_,n_), 01112 H_(), 01113 D_(n_, n_), 01114 ort_(), 01115 cdivr(), cdivi() 01116 { 01117 if (n_ > INT_MAX) 01118 throw Exception("EigenValue: can only be computed for matrices <= " + TextTools::toString(INT_MAX)); 01119 for (size_t j = 0; (j < n_) && issymmetric_; j++) 01120 { 01121 for (size_t i = 0; (i < n_) && issymmetric_; i++) 01122 { 01123 issymmetric_ = (A(i,j) == A(j,i)); 01124 } 01125 } 01126 01127 if (issymmetric_) 01128 { 01129 for (size_t i = 0; i < n_; i++) 01130 { 01131 for (size_t j = 0; j < n_; j++) 01132 { 01133 V_(i,j) = A(i,j); 01134 } 01135 } 01136 01137 // Tridiagonalize. 01138 tred2(); 01139 01140 // Diagonalize. 01141 tql2(); 01142 } 01143 else 01144 { 01145 H_.resize(n_,n_); 01146 ort_.resize(n_); 01147 01148 for (size_t j = 0; j < n_; j++) 01149 { 01150 for (size_t i = 0; i < n_; i++) 01151 { 01152 H_(i,j) = A(i,j); 01153 } 01154 } 01155 01156 // Reduce to Hessenberg form. 01157 orthes(); 01158 01159 // Reduce Hessenberg to real Schur form. 01160 hqr2(); 01161 } 01162 } 01163 01164 01170 const RowMatrix<Real>& getV() const { return V_; } 01171 01177 const std::vector<Real>& getRealEigenValues() const { return d_; } 01178 01184 const std::vector<Real>& getImagEigenValues() const { return e_; } 01185 01219 const RowMatrix<Real>& getD() const 01220 { 01221 for (size_t i = 0; i < n_; i++) 01222 { 01223 for (size_t j = 0; j < n_; j++) 01224 { 01225 D_(i,j) = 0.0; 01226 } 01227 D_(i,i) = d_[i]; 01228 if (e_[i] > 0) 01229 { 01230 D_(i,i+1) = e_[i]; 01231 } 01232 else if (e_[i] < 0) 01233 { 01234 D_(i,i-1) = e_[i]; 01235 } 01236 } 01237 return D_; 01238 } 01239 }; 01240 01241 } //end of namespace bpp. 01242 01243 #endif //_EIGENVALUE_H_ 01244