00001 #ifndef JAMA_SVD_H
00002 #define JAMA_SVD_H
00003
00004
00005 #include "tnt_array1d.h"
00006 #include "tnt_array1d_utils.h"
00007 #include "tnt_array2d.h"
00008 #include "tnt_array2d_utils.h"
00009 #include "tnt_math_utils.h"
00010
00011 #include <algorithm>
00012
00013 #include <cmath>
00014
00015
00016 using namespace TNT;
00017 using namespace std;
00018
00019 namespace JAMA
00020 {
00038 template <class Real>
00039 class SVD
00040 {
00041
00042
00043 Array2D<Real> U, V;
00044 Array1D<Real> s;
00045 int m, n;
00046
00047 public:
00048
00049
00050 SVD (const Array2D<Real> &Arg) {
00051
00052
00053 m = Arg.dim1();
00054 n = Arg.dim2();
00055 int nu = min(m,n);
00056 s = Array1D<Real>(min(m+1,n));
00057 U = Array2D<Real>(m, nu, Real(0));
00058 V = Array2D<Real>(n,n);
00059 Array1D<Real> e(n);
00060 Array1D<Real> work(m);
00061 Array2D<Real> A(Arg.copy());
00062 int wantu = 1;
00063 int wantv = 1;
00064 int i=0, j=0, k=0;
00065
00066
00067
00068
00069 int nct = min(m-1,n);
00070 int nrt = max(0,min(n-2,m));
00071 for (k = 0; k < max(nct,nrt); k++) {
00072 if (k < nct) {
00073
00074
00075
00076
00077 s[k] = 0;
00078 for (i = k; i < m; i++) {
00079 s[k] = hypot(s[k],A[i][k]);
00080 }
00081 if (s[k] != 0.0) {
00082 if (A[k][k] < 0.0) {
00083 s[k] = -s[k];
00084 }
00085 for (i = k; i < m; i++) {
00086 A[i][k] /= s[k];
00087 }
00088 A[k][k] += 1.0;
00089 }
00090 s[k] = -s[k];
00091 }
00092 for (j = k+1; j < n; j++) {
00093 if ((k < nct) && (s[k] != 0.0)) {
00094
00095
00096
00097 Real t(0.0);
00098 for (i = k; i < m; i++) {
00099 t += A[i][k]*A[i][j];
00100 }
00101 t = -t/A[k][k];
00102 for (i = k; i < m; i++) {
00103 A[i][j] += t*A[i][k];
00104 }
00105 }
00106
00107
00108
00109
00110 e[j] = A[k][j];
00111 }
00112 if (wantu & (k < nct)) {
00113
00114
00115
00116
00117 for (i = k; i < m; i++) {
00118 U[i][k] = A[i][k];
00119 }
00120 }
00121 if (k < nrt) {
00122
00123
00124
00125
00126 e[k] = 0;
00127 for (i = k+1; i < n; i++) {
00128 e[k] = hypot(e[k],e[i]);
00129 }
00130 if (e[k] != 0.0) {
00131 if (e[k+1] < 0.0) {
00132 e[k] = -e[k];
00133 }
00134 for (i = k+1; i < n; i++) {
00135 e[i] /= e[k];
00136 }
00137 e[k+1] += 1.0;
00138 }
00139 e[k] = -e[k];
00140 if ((k+1 < m) & (e[k] != 0.0)) {
00141
00142
00143
00144 for (i = k+1; i < m; i++) {
00145 work[i] = 0.0;
00146 }
00147 for (j = k+1; j < n; j++) {
00148 for (i = k+1; i < m; i++) {
00149 work[i] += e[j]*A[i][j];
00150 }
00151 }
00152 for (j = k+1; j < n; j++) {
00153 Real t(-e[j]/e[k+1]);
00154 for (i = k+1; i < m; i++) {
00155 A[i][j] += t*work[i];
00156 }
00157 }
00158 }
00159 if (wantv) {
00160
00161
00162
00163
00164 for (i = k+1; i < n; i++) {
00165 V[i][k] = e[i];
00166 }
00167 }
00168 }
00169 }
00170
00171
00172
00173 int p = min(n,m+1);
00174 if (nct < n) {
00175 s[nct] = A[nct][nct];
00176 }
00177 if (m < p) {
00178 s[p-1] = 0.0;
00179 }
00180 if (nrt+1 < p) {
00181 e[nrt] = A[nrt][p-1];
00182 }
00183 e[p-1] = 0.0;
00184
00185
00186
00187 if (wantu) {
00188 for (j = nct; j < nu; j++) {
00189 for (i = 0; i < m; i++) {
00190 U[i][j] = 0.0;
00191 }
00192 U[j][j] = 1.0;
00193 }
00194 for (k = nct-1; k >= 0; k--) {
00195 if (s[k] != 0.0) {
00196 for (j = k+1; j < nu; j++) {
00197 Real t(0.0);
00198 for (i = k; i < m; i++) {
00199 t += U[i][k]*U[i][j];
00200 }
00201 t = -t/U[k][k];
00202 for (i = k; i < m; i++) {
00203 U[i][j] += t*U[i][k];
00204 }
00205 }
00206 for (i = k; i < m; i++ ) {
00207 U[i][k] = -U[i][k];
00208 }
00209 U[k][k] = 1.0 + U[k][k];
00210 for (i = 0; i < k-1; i++) {
00211 U[i][k] = 0.0;
00212 }
00213 } else {
00214 for (i = 0; i < m; i++) {
00215 U[i][k] = 0.0;
00216 }
00217 U[k][k] = 1.0;
00218 }
00219 }
00220 }
00221
00222
00223
00224 if (wantv) {
00225 for (k = n-1; k >= 0; k--) {
00226 if ((k < nrt) & (e[k] != 0.0)) {
00227 for (j = k+1; j < nu; j++) {
00228 Real t(0.0);
00229 for (i = k+1; i < n; i++) {
00230 t += V[i][k]*V[i][j];
00231 }
00232 t = -t/V[k+1][k];
00233 for (i = k+1; i < n; i++) {
00234 V[i][j] += t*V[i][k];
00235 }
00236 }
00237 }
00238 for (i = 0; i < n; i++) {
00239 V[i][k] = 0.0;
00240 }
00241 V[k][k] = 1.0;
00242 }
00243 }
00244
00245
00246
00247 int pp = p-1;
00248 int iter = 0;
00249 Real eps(pow(2.0,-52.0));
00250 while (p > 0) {
00251 int k=0;
00252 int kase=0;
00253
00254
00255
00256
00257
00258
00259
00260
00261
00262
00263
00264
00265
00266 for (k = p-2; k >= -1; k--) {
00267 if (k == -1) {
00268 break;
00269 }
00270 if (abs(e[k]) <= eps*(abs(s[k]) + abs(s[k+1]))) {
00271 e[k] = 0.0;
00272 break;
00273 }
00274 }
00275 if (k == p-2) {
00276 kase = 4;
00277 } else {
00278 int ks;
00279 for (ks = p-1; ks >= k; ks--) {
00280 if (ks == k) {
00281 break;
00282 }
00283 Real t( (ks != p ? abs(e[ks]) : 0.) +
00284 (ks != k+1 ? abs(e[ks-1]) : 0.));
00285 if (abs(s[ks]) <= eps*t) {
00286 s[ks] = 0.0;
00287 break;
00288 }
00289 }
00290 if (ks == k) {
00291 kase = 3;
00292 } else if (ks == p-1) {
00293 kase = 1;
00294 } else {
00295 kase = 2;
00296 k = ks;
00297 }
00298 }
00299 k++;
00300
00301
00302
00303 switch (kase) {
00304
00305
00306
00307 case 1: {
00308 Real f(e[p-2]);
00309 e[p-2] = 0.0;
00310 for (j = p-2; j >= k; j--) {
00311 Real t( hypot(s[j],f));
00312 Real cs(s[j]/t);
00313 Real sn(f/t);
00314 s[j] = t;
00315 if (j != k) {
00316 f = -sn*e[j-1];
00317 e[j-1] = cs*e[j-1];
00318 }
00319 if (wantv) {
00320 for (i = 0; i < n; i++) {
00321 t = cs*V[i][j] + sn*V[i][p-1];
00322 V[i][p-1] = -sn*V[i][j] + cs*V[i][p-1];
00323 V[i][j] = t;
00324 }
00325 }
00326 }
00327 }
00328 break;
00329
00330
00331
00332 case 2: {
00333 Real f(e[k-1]);
00334 e[k-1] = 0.0;
00335 for (j = k; j < p; j++) {
00336 Real t(hypot(s[j],f));
00337 Real cs( s[j]/t);
00338 Real sn(f/t);
00339 s[j] = t;
00340 f = -sn*e[j];
00341 e[j] = cs*e[j];
00342 if (wantu) {
00343 for (i = 0; i < m; i++) {
00344 t = cs*U[i][j] + sn*U[i][k-1];
00345 U[i][k-1] = -sn*U[i][j] + cs*U[i][k-1];
00346 U[i][j] = t;
00347 }
00348 }
00349 }
00350 }
00351 break;
00352
00353
00354
00355 case 3: {
00356
00357
00358
00359 Real scale = max(max(max(max(
00360 abs(s[p-1]),abs(s[p-2])),abs(e[p-2])),
00361 abs(s[k])),abs(e[k]));
00362 Real sp = s[p-1]/scale;
00363 Real spm1 = s[p-2]/scale;
00364 Real epm1 = e[p-2]/scale;
00365 Real sk = s[k]/scale;
00366 Real ek = e[k]/scale;
00367 Real b = ((spm1 + sp)*(spm1 - sp) + epm1*epm1)/2.0;
00368 Real c = (sp*epm1)*(sp*epm1);
00369 Real shift = 0.0;
00370 if ((b != 0.0) || (c != 0.0)) {
00371 shift = sqrt(b*b + c);
00372 if (b < 0.0) {
00373 shift = -shift;
00374 }
00375 shift = c/(b + shift);
00376 }
00377 Real f = (sk + sp)*(sk - sp) + shift;
00378 Real g = sk*ek;
00379
00380
00381
00382 for (j = k; j < p-1; j++) {
00383 Real t = hypot(f,g);
00384 Real cs = f/t;
00385 Real sn = g/t;
00386 if (j != k) {
00387 e[j-1] = t;
00388 }
00389 f = cs*s[j] + sn*e[j];
00390 e[j] = cs*e[j] - sn*s[j];
00391 g = sn*s[j+1];
00392 s[j+1] = cs*s[j+1];
00393 if (wantv) {
00394 for (i = 0; i < n; i++) {
00395 t = cs*V[i][j] + sn*V[i][j+1];
00396 V[i][j+1] = -sn*V[i][j] + cs*V[i][j+1];
00397 V[i][j] = t;
00398 }
00399 }
00400 t = hypot(f,g);
00401 cs = f/t;
00402 sn = g/t;
00403 s[j] = t;
00404 f = cs*e[j] + sn*s[j+1];
00405 s[j+1] = -sn*e[j] + cs*s[j+1];
00406 g = sn*e[j+1];
00407 e[j+1] = cs*e[j+1];
00408 if (wantu && (j < m-1)) {
00409 for (i = 0; i < m; i++) {
00410 t = cs*U[i][j] + sn*U[i][j+1];
00411 U[i][j+1] = -sn*U[i][j] + cs*U[i][j+1];
00412 U[i][j] = t;
00413 }
00414 }
00415 }
00416 e[p-2] = f;
00417 iter = iter + 1;
00418 }
00419 break;
00420
00421
00422
00423 case 4: {
00424
00425
00426
00427 if (s[k] <= 0.0) {
00428 s[k] = (s[k] < 0.0 ? -s[k] : 0.0);
00429 if (wantv) {
00430 for (i = 0; i <= pp; i++) {
00431 V[i][k] = -V[i][k];
00432 }
00433 }
00434 }
00435
00436
00437
00438 while (k < pp) {
00439 if (s[k] >= s[k+1]) {
00440 break;
00441 }
00442 Real t = s[k];
00443 s[k] = s[k+1];
00444 s[k+1] = t;
00445 if (wantv && (k < n-1)) {
00446 for (i = 0; i < n; i++) {
00447 t = V[i][k+1]; V[i][k+1] = V[i][k]; V[i][k] = t;
00448 }
00449 }
00450 if (wantu && (k < m-1)) {
00451 for (i = 0; i < m; i++) {
00452 t = U[i][k+1]; U[i][k+1] = U[i][k]; U[i][k] = t;
00453 }
00454 }
00455 k++;
00456 }
00457 iter = 0;
00458 p--;
00459 }
00460 break;
00461 }
00462 }
00463 }
00464
00465
00466 void getU (Array2D<Real> &A)
00467 {
00468 int minm = min(m+1,n);
00469
00470 A = Array2D<Real>(m, minm);
00471
00472 for (int i=0; i<m; i++)
00473 for (int j=0; j<minm; j++)
00474 A[i][j] = U[i][j];
00475
00476 }
00477
00478
00479
00480 void getV (Array2D<Real> &A)
00481 {
00482 A = V;
00483 }
00484
00487 void getSingularValues (Array1D<Real> &x)
00488 {
00489 x = s;
00490 }
00491
00496 void getS (Array2D<Real> &A) {
00497 A = Array2D<Real>(n,n);
00498 for (int i = 0; i < n; i++) {
00499 for (int j = 0; j < n; j++) {
00500 A[i][j] = 0.0;
00501 }
00502 A[i][i] = s[i];
00503 }
00504 }
00505
00508 Real norm2 () {
00509 return s[0];
00510 }
00511
00514 Real cond () {
00515 return s[0]/s[min(m,n)-1];
00516 }
00517
00522 int rank ()
00523 {
00524 Real eps = pow(2.0,-52.0);
00525 Real tol = max(m,n)*s[0]*eps;
00526 int r = 0;
00527 for (int i = 0; i < s.dim(); i++) {
00528 if (s[i] > tol) {
00529 r++;
00530 }
00531 }
00532 return r;
00533 }
00534 };
00535
00536 }
00537 #endif
00538