00001 #ifndef JAMA_EIG_H
00002 #define JAMA_EIG_H
00003
00004
00005 #include "tnt_array1d.h"
00006 #include "tnt_array2d.h"
00007 #include "tnt_math_utils.h"
00008
00009
00010 using namespace TNT;
00011
00012
00013 namespace JAMA
00014 {
00015
00065
00066 template <class Real>
00067 class Eigenvalue
00068 {
00069
00070
00072 int n;
00073
00074 int issymmetric;
00075
00077
00078 TNT::Array1D<Real> d;
00079 TNT::Array1D<Real> e;
00080
00082 TNT::Array2D<Real> V;
00083
00087 TNT::Array2D<Real> H;
00088
00089
00093 TNT::Array1D<Real> ort;
00094
00095
00096
00097
00098 void tred2() {
00099
00100
00101
00102
00103
00104
00105 for (int j = 0; j < n; j++) {
00106 d[j] = V[n-1][j];
00107 }
00108
00109
00110
00111 for (int i = n-1; i > 0; i--) {
00112
00113
00114
00115 Real scale = 0.0;
00116 Real h = 0.0;
00117 for (int k = 0; k < i; k++) {
00118 scale = scale + abs(d[k]);
00119 }
00120 if (scale == 0.0) {
00121 e[i] = d[i-1];
00122 for (int j = 0; j < i; j++) {
00123 d[j] = V[i-1][j];
00124 V[i][j] = 0.0;
00125 V[j][i] = 0.0;
00126 }
00127 } else {
00128
00129
00130
00131 for (int k = 0; k < i; k++) {
00132 d[k] /= scale;
00133 h += d[k] * d[k];
00134 }
00135 Real f = d[i-1];
00136 Real g = sqrt(h);
00137 if (f > 0) {
00138 g = -g;
00139 }
00140 e[i] = scale * g;
00141 h = h - f * g;
00142 d[i-1] = f - g;
00143 for (int j = 0; j < i; j++) {
00144 e[j] = 0.0;
00145 }
00146
00147
00148
00149 for (int j = 0; j < i; j++) {
00150 f = d[j];
00151 V[j][i] = f;
00152 g = e[j] + V[j][j] * f;
00153 for (int k = j+1; k <= i-1; k++) {
00154 g += V[k][j] * d[k];
00155 e[k] += V[k][j] * f;
00156 }
00157 e[j] = g;
00158 }
00159 f = 0.0;
00160 for (int j = 0; j < i; j++) {
00161 e[j] /= h;
00162 f += e[j] * d[j];
00163 }
00164 Real hh = f / (h + h);
00165 for (int j = 0; j < i; j++) {
00166 e[j] -= hh * d[j];
00167 }
00168 for (int j = 0; j < i; j++) {
00169 f = d[j];
00170 g = e[j];
00171 for (int k = j; k <= i-1; k++) {
00172 V[k][j] -= (f * e[k] + g * d[k]);
00173 }
00174 d[j] = V[i-1][j];
00175 V[i][j] = 0.0;
00176 }
00177 }
00178 d[i] = h;
00179 }
00180
00181
00182
00183 for (int i = 0; i < n-1; i++) {
00184 V[n-1][i] = V[i][i];
00185 V[i][i] = 1.0;
00186 Real h = d[i+1];
00187 if (h != 0.0) {
00188 for (int k = 0; k <= i; k++) {
00189 d[k] = V[k][i+1] / h;
00190 }
00191 for (int j = 0; j <= i; j++) {
00192 Real g = 0.0;
00193 for (int k = 0; k <= i; k++) {
00194 g += V[k][i+1] * V[k][j];
00195 }
00196 for (int k = 0; k <= i; k++) {
00197 V[k][j] -= g * d[k];
00198 }
00199 }
00200 }
00201 for (int k = 0; k <= i; k++) {
00202 V[k][i+1] = 0.0;
00203 }
00204 }
00205 for (int j = 0; j < n; j++) {
00206 d[j] = V[n-1][j];
00207 V[n-1][j] = 0.0;
00208 }
00209 V[n-1][n-1] = 1.0;
00210 e[0] = 0.0;
00211 }
00212
00213
00214
00215 void tql2 () {
00216
00217
00218
00219
00220
00221
00222 for (int i = 1; i < n; i++) {
00223 e[i-1] = e[i];
00224 }
00225 e[n-1] = 0.0;
00226
00227 Real f = 0.0;
00228 Real tst1 = 0.0;
00229 Real eps = pow(2.0,-52.0);
00230 for (int l = 0; l < n; l++) {
00231
00232
00233
00234 tst1 = max(tst1,abs(d[l]) + abs(e[l]));
00235 int m = l;
00236
00237
00238 while (m < n) {
00239 if (abs(e[m]) <= eps*tst1) {
00240 break;
00241 }
00242 m++;
00243 }
00244
00245
00246
00247
00248
00249 if (m > l) {
00250 int iter = 0;
00251 do {
00252 iter = iter + 1;
00253
00254
00255
00256 Real g = d[l];
00257 Real p = (d[l+1] - g) / (2.0 * e[l]);
00258 Real r = hypot(p,1.0);
00259 if (p < 0) {
00260 r = -r;
00261 }
00262 d[l] = e[l] / (p + r);
00263 d[l+1] = e[l] * (p + r);
00264 Real dl1 = d[l+1];
00265 Real h = g - d[l];
00266 for (int i = l+2; i < n; i++) {
00267 d[i] -= h;
00268 }
00269 f = f + h;
00270
00271
00272
00273 p = d[m];
00274 Real c = 1.0;
00275 Real c2 = c;
00276 Real c3 = c;
00277 Real el1 = e[l+1];
00278 Real s = 0.0;
00279 Real s2 = 0.0;
00280 for (int i = m-1; i >= l; i--) {
00281 c3 = c2;
00282 c2 = c;
00283 s2 = s;
00284 g = c * e[i];
00285 h = c * p;
00286 r = hypot(p,e[i]);
00287 e[i+1] = s * r;
00288 s = e[i] / r;
00289 c = p / r;
00290 p = c * d[i] - s * g;
00291 d[i+1] = h + s * (c * g + s * d[i]);
00292
00293
00294
00295 for (int k = 0; k < n; k++) {
00296 h = V[k][i+1];
00297 V[k][i+1] = s * V[k][i] + c * h;
00298 V[k][i] = c * V[k][i] - s * h;
00299 }
00300 }
00301 p = -s * s2 * c3 * el1 * e[l] / dl1;
00302 e[l] = s * p;
00303 d[l] = c * p;
00304
00305
00306
00307 } while (abs(e[l]) > eps*tst1);
00308 }
00309 d[l] = d[l] + f;
00310 e[l] = 0.0;
00311 }
00312
00313
00314
00315 for (int i = 0; i < n-1; i++) {
00316 int k = i;
00317 Real p = d[i];
00318 for (int j = i+1; j < n; j++) {
00319 if (d[j] < p) {
00320 k = j;
00321 p = d[j];
00322 }
00323 }
00324 if (k != i) {
00325 d[k] = d[i];
00326 d[i] = p;
00327 for (int j = 0; j < n; j++) {
00328 p = V[j][i];
00329 V[j][i] = V[j][k];
00330 V[j][k] = p;
00331 }
00332 }
00333 }
00334 }
00335
00336
00337
00338 void orthes () {
00339
00340
00341
00342
00343
00344
00345 int low = 0;
00346 int high = n-1;
00347
00348 for (int m = low+1; m <= high-1; m++) {
00349
00350
00351
00352 Real scale = 0.0;
00353 for (int i = m; i <= high; i++) {
00354 scale = scale + abs(H[i][m-1]);
00355 }
00356 if (scale != 0.0) {
00357
00358
00359
00360 Real h = 0.0;
00361 for (int i = high; i >= m; i--) {
00362 ort[i] = H[i][m-1]/scale;
00363 h += ort[i] * ort[i];
00364 }
00365 Real g = sqrt(h);
00366 if (ort[m] > 0) {
00367 g = -g;
00368 }
00369 h = h - ort[m] * g;
00370 ort[m] = ort[m] - g;
00371
00372
00373
00374
00375 for (int j = m; j < n; j++) {
00376 Real f = 0.0;
00377 for (int i = high; i >= m; i--) {
00378 f += ort[i]*H[i][j];
00379 }
00380 f = f/h;
00381 for (int i = m; i <= high; i++) {
00382 H[i][j] -= f*ort[i];
00383 }
00384 }
00385
00386 for (int i = 0; i <= high; i++) {
00387 Real f = 0.0;
00388 for (int j = high; j >= m; j--) {
00389 f += ort[j]*H[i][j];
00390 }
00391 f = f/h;
00392 for (int j = m; j <= high; j++) {
00393 H[i][j] -= f*ort[j];
00394 }
00395 }
00396 ort[m] = scale*ort[m];
00397 H[m][m-1] = scale*g;
00398 }
00399 }
00400
00401
00402
00403 for (int i = 0; i < n; i++) {
00404 for (int j = 0; j < n; j++) {
00405 V[i][j] = (i == j ? 1.0 : 0.0);
00406 }
00407 }
00408
00409 for (int m = high-1; m >= low+1; m--) {
00410 if (H[m][m-1] != 0.0) {
00411 for (int i = m+1; i <= high; i++) {
00412 ort[i] = H[i][m-1];
00413 }
00414 for (int j = m; j <= high; j++) {
00415 Real g = 0.0;
00416 for (int i = m; i <= high; i++) {
00417 g += ort[i] * V[i][j];
00418 }
00419
00420 g = (g / ort[m]) / H[m][m-1];
00421 for (int i = m; i <= high; i++) {
00422 V[i][j] += g * ort[i];
00423 }
00424 }
00425 }
00426 }
00427 }
00428
00429
00430
00431
00432 Real cdivr, cdivi;
00433 void cdiv(Real xr, Real xi, Real yr, Real yi) {
00434 Real r,d;
00435 if (abs(yr) > abs(yi)) {
00436 r = yi/yr;
00437 d = yr + r*yi;
00438 cdivr = (xr + r*xi)/d;
00439 cdivi = (xi - r*xr)/d;
00440 } else {
00441 r = yr/yi;
00442 d = yi + r*yr;
00443 cdivr = (r*xr + xi)/d;
00444 cdivi = (r*xi - xr)/d;
00445 }
00446 }
00447
00448
00449
00450
00451 void hqr2 () {
00452
00453
00454
00455
00456
00457
00458
00459
00460 int nn = this->n;
00461 int n = nn-1;
00462 int low = 0;
00463 int high = nn-1;
00464 Real eps = pow(2.0,-52.0);
00465 Real exshift = 0.0;
00466 Real p=0,q=0,r=0,s=0,z=0,t,w,x,y;
00467
00468
00469
00470 Real norm = 0.0;
00471 for (int i = 0; i < nn; i++) {
00472 if ((i < low) || (i > high)) {
00473 d[i] = H[i][i];
00474 e[i] = 0.0;
00475 }
00476 for (int j = max(i-1,0); j < nn; j++) {
00477 norm = norm + abs(H[i][j]);
00478 }
00479 }
00480
00481
00482
00483 int iter = 0;
00484 while (n >= low) {
00485
00486
00487
00488 int l = n;
00489 while (l > low) {
00490 s = abs(H[l-1][l-1]) + abs(H[l][l]);
00491 if (s == 0.0) {
00492 s = norm;
00493 }
00494 if (abs(H[l][l-1]) < eps * s) {
00495 break;
00496 }
00497 l--;
00498 }
00499
00500
00501
00502
00503 if (l == n) {
00504 H[n][n] = H[n][n] + exshift;
00505 d[n] = H[n][n];
00506 e[n] = 0.0;
00507 n--;
00508 iter = 0;
00509
00510
00511
00512 } else if (l == n-1) {
00513 w = H[n][n-1] * H[n-1][n];
00514 p = (H[n-1][n-1] - H[n][n]) / 2.0;
00515 q = p * p + w;
00516 z = sqrt(abs(q));
00517 H[n][n] = H[n][n] + exshift;
00518 H[n-1][n-1] = H[n-1][n-1] + exshift;
00519 x = H[n][n];
00520
00521
00522
00523 if (q >= 0) {
00524 if (p >= 0) {
00525 z = p + z;
00526 } else {
00527 z = p - z;
00528 }
00529 d[n-1] = x + z;
00530 d[n] = d[n-1];
00531 if (z != 0.0) {
00532 d[n] = x - w / z;
00533 }
00534 e[n-1] = 0.0;
00535 e[n] = 0.0;
00536 x = H[n][n-1];
00537 s = abs(x) + abs(z);
00538 p = x / s;
00539 q = z / s;
00540 r = sqrt(p * p+q * q);
00541 p = p / r;
00542 q = q / r;
00543
00544
00545
00546 for (int j = n-1; j < nn; j++) {
00547 z = H[n-1][j];
00548 H[n-1][j] = q * z + p * H[n][j];
00549 H[n][j] = q * H[n][j] - p * z;
00550 }
00551
00552
00553
00554 for (int i = 0; i <= n; i++) {
00555 z = H[i][n-1];
00556 H[i][n-1] = q * z + p * H[i][n];
00557 H[i][n] = q * H[i][n] - p * z;
00558 }
00559
00560
00561
00562 for (int i = low; i <= high; i++) {
00563 z = V[i][n-1];
00564 V[i][n-1] = q * z + p * V[i][n];
00565 V[i][n] = q * V[i][n] - p * z;
00566 }
00567
00568
00569
00570 } else {
00571 d[n-1] = x + p;
00572 d[n] = x + p;
00573 e[n-1] = z;
00574 e[n] = -z;
00575 }
00576 n = n - 2;
00577 iter = 0;
00578
00579
00580
00581 } else {
00582
00583
00584
00585 x = H[n][n];
00586 y = 0.0;
00587 w = 0.0;
00588 if (l < n) {
00589 y = H[n-1][n-1];
00590 w = H[n][n-1] * H[n-1][n];
00591 }
00592
00593
00594
00595 if (iter == 10) {
00596 exshift += x;
00597 for (int i = low; i <= n; i++) {
00598 H[i][i] -= x;
00599 }
00600 s = abs(H[n][n-1]) + abs(H[n-1][n-2]);
00601 x = y = 0.75 * s;
00602 w = -0.4375 * s * s;
00603 }
00604
00605
00606
00607 if (iter == 30) {
00608 s = (y - x) / 2.0;
00609 s = s * s + w;
00610 if (s > 0) {
00611 s = sqrt(s);
00612 if (y < x) {
00613 s = -s;
00614 }
00615 s = x - w / ((y - x) / 2.0 + s);
00616 for (int i = low; i <= n; i++) {
00617 H[i][i] -= s;
00618 }
00619 exshift += s;
00620 x = y = w = 0.964;
00621 }
00622 }
00623
00624 iter = iter + 1;
00625
00626
00627
00628 int m = n-2;
00629 while (m >= l) {
00630 z = H[m][m];
00631 r = x - z;
00632 s = y - z;
00633 p = (r * s - w) / H[m+1][m] + H[m][m+1];
00634 q = H[m+1][m+1] - z - r - s;
00635 r = H[m+2][m+1];
00636 s = abs(p) + abs(q) + abs(r);
00637 p = p / s;
00638 q = q / s;
00639 r = r / s;
00640 if (m == l) {
00641 break;
00642 }
00643 if (abs(H[m][m-1]) * (abs(q) + abs(r)) <
00644 eps * (abs(p) * (abs(H[m-1][m-1]) + abs(z) +
00645 abs(H[m+1][m+1])))) {
00646 break;
00647 }
00648 m--;
00649 }
00650
00651 for (int i = m+2; i <= n; i++) {
00652 H[i][i-2] = 0.0;
00653 if (i > m+2) {
00654 H[i][i-3] = 0.0;
00655 }
00656 }
00657
00658
00659
00660 for (int k = m; k <= n-1; k++) {
00661 int notlast = (k != n-1);
00662 if (k != m) {
00663 p = H[k][k-1];
00664 q = H[k+1][k-1];
00665 r = (notlast ? H[k+2][k-1] : 0.0);
00666 x = abs(p) + abs(q) + abs(r);
00667 if (x != 0.0) {
00668 p = p / x;
00669 q = q / x;
00670 r = r / x;
00671 }
00672 }
00673 if (x == 0.0) {
00674 break;
00675 }
00676 s = sqrt(p * p + q * q + r * r);
00677 if (p < 0) {
00678 s = -s;
00679 }
00680 if (s != 0) {
00681 if (k != m) {
00682 H[k][k-1] = -s * x;
00683 } else if (l != m) {
00684 H[k][k-1] = -H[k][k-1];
00685 }
00686 p = p + s;
00687 x = p / s;
00688 y = q / s;
00689 z = r / s;
00690 q = q / p;
00691 r = r / p;
00692
00693
00694
00695 for (int j = k; j < nn; j++) {
00696 p = H[k][j] + q * H[k+1][j];
00697 if (notlast) {
00698 p = p + r * H[k+2][j];
00699 H[k+2][j] = H[k+2][j] - p * z;
00700 }
00701 H[k][j] = H[k][j] - p * x;
00702 H[k+1][j] = H[k+1][j] - p * y;
00703 }
00704
00705
00706
00707 for (int i = 0; i <= min(n,k+3); i++) {
00708 p = x * H[i][k] + y * H[i][k+1];
00709 if (notlast) {
00710 p = p + z * H[i][k+2];
00711 H[i][k+2] = H[i][k+2] - p * r;
00712 }
00713 H[i][k] = H[i][k] - p;
00714 H[i][k+1] = H[i][k+1] - p * q;
00715 }
00716
00717
00718
00719 for (int i = low; i <= high; i++) {
00720 p = x * V[i][k] + y * V[i][k+1];
00721 if (notlast) {
00722 p = p + z * V[i][k+2];
00723 V[i][k+2] = V[i][k+2] - p * r;
00724 }
00725 V[i][k] = V[i][k] - p;
00726 V[i][k+1] = V[i][k+1] - p * q;
00727 }
00728 }
00729 }
00730 }
00731 }
00732
00733
00734
00735 if (norm == 0.0) {
00736 return;
00737 }
00738
00739 for (n = nn-1; n >= 0; n--) {
00740 p = d[n];
00741 q = e[n];
00742
00743
00744
00745 if (q == 0) {
00746 int l = n;
00747 H[n][n] = 1.0;
00748 for (int i = n-1; i >= 0; i--) {
00749 w = H[i][i] - p;
00750 r = 0.0;
00751 for (int j = l; j <= n; j++) {
00752 r = r + H[i][j] * H[j][n];
00753 }
00754 if (e[i] < 0.0) {
00755 z = w;
00756 s = r;
00757 } else {
00758 l = i;
00759 if (e[i] == 0.0) {
00760 if (w != 0.0) {
00761 H[i][n] = -r / w;
00762 } else {
00763 H[i][n] = -r / (eps * norm);
00764 }
00765
00766
00767
00768 } else {
00769 x = H[i][i+1];
00770 y = H[i+1][i];
00771 q = (d[i] - p) * (d[i] - p) + e[i] * e[i];
00772 t = (x * s - z * r) / q;
00773 H[i][n] = t;
00774 if (abs(x) > abs(z)) {
00775 H[i+1][n] = (-r - w * t) / x;
00776 } else {
00777 H[i+1][n] = (-s - y * t) / z;
00778 }
00779 }
00780
00781
00782
00783 t = abs(H[i][n]);
00784 if ((eps * t) * t > 1) {
00785 for (int j = i; j <= n; j++) {
00786 H[j][n] = H[j][n] / t;
00787 }
00788 }
00789 }
00790 }
00791
00792
00793
00794 } else if (q < 0) {
00795 int l = n-1;
00796
00797
00798
00799 if (abs(H[n][n-1]) > abs(H[n-1][n])) {
00800 H[n-1][n-1] = q / H[n][n-1];
00801 H[n-1][n] = -(H[n][n] - p) / H[n][n-1];
00802 } else {
00803 cdiv(0.0,-H[n-1][n],H[n-1][n-1]-p,q);
00804 H[n-1][n-1] = cdivr;
00805 H[n-1][n] = cdivi;
00806 }
00807 H[n][n-1] = 0.0;
00808 H[n][n] = 1.0;
00809 for (int i = n-2; i >= 0; i--) {
00810 Real ra,sa,vr,vi;
00811 ra = 0.0;
00812 sa = 0.0;
00813 for (int j = l; j <= n; j++) {
00814 ra = ra + H[i][j] * H[j][n-1];
00815 sa = sa + H[i][j] * H[j][n];
00816 }
00817 w = H[i][i] - p;
00818
00819 if (e[i] < 0.0) {
00820 z = w;
00821 r = ra;
00822 s = sa;
00823 } else {
00824 l = i;
00825 if (e[i] == 0) {
00826 cdiv(-ra,-sa,w,q);
00827 H[i][n-1] = cdivr;
00828 H[i][n] = cdivi;
00829 } else {
00830
00831
00832
00833 x = H[i][i+1];
00834 y = H[i+1][i];
00835 vr = (d[i] - p) * (d[i] - p) + e[i] * e[i] - q * q;
00836 vi = (d[i] - p) * 2.0 * q;
00837 if ((vr == 0.0) && (vi == 0.0)) {
00838 vr = eps * norm * (abs(w) + abs(q) +
00839 abs(x) + abs(y) + abs(z));
00840 }
00841 cdiv(x*r-z*ra+q*sa,x*s-z*sa-q*ra,vr,vi);
00842 H[i][n-1] = cdivr;
00843 H[i][n] = cdivi;
00844 if (abs(x) > (abs(z) + abs(q))) {
00845 H[i+1][n-1] = (-ra - w * H[i][n-1] + q * H[i][n]) / x;
00846 H[i+1][n] = (-sa - w * H[i][n] - q * H[i][n-1]) / x;
00847 } else {
00848 cdiv(-r-y*H[i][n-1],-s-y*H[i][n],z,q);
00849 H[i+1][n-1] = cdivr;
00850 H[i+1][n] = cdivi;
00851 }
00852 }
00853
00854
00855
00856 t = max(abs(H[i][n-1]),abs(H[i][n]));
00857 if ((eps * t) * t > 1) {
00858 for (int j = i; j <= n; j++) {
00859 H[j][n-1] = H[j][n-1] / t;
00860 H[j][n] = H[j][n] / t;
00861 }
00862 }
00863 }
00864 }
00865 }
00866 }
00867
00868
00869
00870 for (int i = 0; i < nn; i++) {
00871 if (i < low || i > high) {
00872 for (int j = i; j < nn; j++) {
00873 V[i][j] = H[i][j];
00874 }
00875 }
00876 }
00877
00878
00879
00880 for (int j = nn-1; j >= low; j--) {
00881 for (int i = low; i <= high; i++) {
00882 z = 0.0;
00883 for (int k = low; k <= min(j,high); k++) {
00884 z = z + V[i][k] * H[k][j];
00885 }
00886 V[i][j] = z;
00887 }
00888 }
00889 }
00890
00891 public:
00892
00893
00897
00898 Eigenvalue(const TNT::Array2D<Real> &A) {
00899 n = A.dim2();
00900 V = Array2D<Real>(n,n);
00901 d = Array1D<Real>(n);
00902 e = Array1D<Real>(n);
00903
00904 issymmetric = 1;
00905 for (int j = 0; (j < n) && issymmetric; j++) {
00906 for (int i = 0; (i < n) && issymmetric; i++) {
00907 issymmetric = (A[i][j] == A[j][i]);
00908 }
00909 }
00910
00911 if (issymmetric) {
00912 for (int i = 0; i < n; i++) {
00913 for (int j = 0; j < n; j++) {
00914 V[i][j] = A[i][j];
00915 }
00916 }
00917
00918
00919 tred2();
00920
00921
00922 tql2();
00923
00924 } else {
00925 H = TNT::Array2D<Real>(n,n);
00926 ort = TNT::Array1D<Real>(n);
00927
00928 for (int j = 0; j < n; j++) {
00929 for (int i = 0; i < n; i++) {
00930 H[i][j] = A[i][j];
00931 }
00932 }
00933
00934
00935 orthes();
00936
00937
00938 hqr2();
00939 }
00940 }
00941
00942
00946
00947 void getV (TNT::Array2D<Real> &V_) {
00948 V_ = V;
00949 return;
00950 }
00951
00955
00956 void getRealEigenvalues (TNT::Array1D<Real> &d_) {
00957 d_ = d;
00958 return ;
00959 }
00960
00966 void getImagEigenvalues (TNT::Array1D<Real> &e_) {
00967 e_ = e;
00968 return;
00969 }
00970
00971
01005 void getD (TNT::Array2D<Real> &D) {
01006 D = Array2D<Real>(n,n);
01007 for (int i = 0; i < n; i++) {
01008 for (int j = 0; j < n; j++) {
01009 D[i][j] = 0.0;
01010 }
01011 D[i][i] = d[i];
01012 if (e[i] > 0) {
01013 D[i][i+1] = e[i];
01014 } else if (e[i] < 0) {
01015 D[i][i-1] = e[i];
01016 }
01017 }
01018 }
01019 };
01020
01021 }
01022
01023
01024 #endif
01025