PLaSK library
Loading...
Searching...
No Matches
jama_svd.h
Go to the documentation of this file.
1#ifndef JAMA_SVD_H
2#define JAMA_SVD_H
3
4
5#include "../tnt/tnt_array1d.h"
6#include "../tnt/tnt_array1d_utils.h"
7#include "../tnt/tnt_array2d.h"
8#include "../tnt/tnt_array2d_utils.h"
9#include "../tnt/tnt_math_utils.h"
10
11#include <algorithm>
12// for min(), max() below
13#include <cmath>
14// for abs() below
15#if(defined(_MSC_VER) && _MSC_VER < 1600)
16#define abs fabs
17#endif
18
19
20using namespace TNT;
21//using namespace std;
22
23namespace JAMA
24{
42template <class Real>
43class SVD
44{
45
46
47 Array2D<Real> U, V;
49 int m, n;
50
51 public:
52
53
54 SVD (const Array2D<Real> &Arg) {
55
56
57 m = Arg.dim1();
58 n = Arg.dim2();
59 int nu = std::min(m,n);
60 s = Array1D<Real>(std::min(m+1,n));
61 U = Array2D<Real>(m, nu, Real(0));
62 V = Array2D<Real>(n,n);
64 Array1D<Real> work(m);
65 Array2D<Real> A(Arg.copy());
66 int wantu = 1; /* boolean */
67 int wantv = 1; /* boolean */
68 int i=0, j=0, k=0;
69
70 // Reduce A to bidiagonal form, storing the diagonal elements
71 // in s and the super-diagonal elements in e.
72
73 int nct = std::min(m-1,n);
74 int nrt = std::max(0,std::min(n-2,m));
75 for (k = 0; k < std::max(nct,nrt); k++) {
76 if (k < nct) {
77
78 // Compute the transformation for the k-th column and
79 // place the k-th diagonal in s[k].
80 // Compute 2-norm of k-th column without under/overflow.
81 s[k] = 0;
82 for (i = k; i < m; i++) {
83 s[k] = hypot(s[k],A[i][k]);
84 }
85 if (s[k] != 0.0) {
86 if (A[k][k] < 0.0) {
87 s[k] = -s[k];
88 }
89 for (i = k; i < m; i++) {
90 A[i][k] /= s[k];
91 }
92 A[k][k] += 1.0;
93 }
94 s[k] = -s[k];
95 }
96 for (j = k+1; j < n; j++) {
97 if ((k < nct) && (s[k] != 0.0)) {
98
99 // Apply the transformation.
100
101 double t = 0;
102 for (i = k; i < m; i++) {
103 t += A[i][k]*A[i][j];
104 }
105 t = -t/A[k][k];
106 for (i = k; i < m; i++) {
107 A[i][j] += t*A[i][k];
108 }
109 }
110
111 // Place the k-th row of A into e for the
112 // subsequent calculation of the row transformation.
113
114 e[j] = A[k][j];
115 }
116 if (wantu & (k < nct)) {
117
118 // Place the transformation in U for subsequent back
119 // multiplication.
120
121 for (i = k; i < m; i++) {
122 U[i][k] = A[i][k];
123 }
124 }
125 if (k < nrt) {
126
127 // Compute the k-th row transformation and place the
128 // k-th super-diagonal in e[k].
129 // Compute 2-norm without under/overflow.
130 e[k] = 0;
131 for (i = k+1; i < n; i++) {
132 e[k] = hypot(e[k],e[i]);
133 }
134 if (e[k] != 0.0) {
135 if (e[k+1] < 0.0) {
136 e[k] = -e[k];
137 }
138 for (i = k+1; i < n; i++) {
139 e[i] /= e[k];
140 }
141 e[k+1] += 1.0;
142 }
143 e[k] = -e[k];
144 if ((k+1 < m) & (e[k] != 0.0)) {
145
146 // Apply the transformation.
147
148 for (i = k+1; i < m; i++) {
149 work[i] = 0.0;
150 }
151 for (j = k+1; j < n; j++) {
152 for (i = k+1; i < m; i++) {
153 work[i] += e[j]*A[i][j];
154 }
155 }
156 for (j = k+1; j < n; j++) {
157 double t = -e[j]/e[k+1];
158 for (i = k+1; i < m; i++) {
159 A[i][j] += t*work[i];
160 }
161 }
162 }
163 if (wantv) {
164
165 // Place the transformation in V for subsequent
166 // back multiplication.
167
168 for (i = k+1; i < n; i++) {
169 V[i][k] = e[i];
170 }
171 }
172 }
173 }
174
175 // Set up the final bidiagonal matrix or order p.
176
177 int p = std::min(n,m+1);
178 if (nct < n) {
179 s[nct] = A[nct][nct];
180 }
181 if (m < p) {
182 s[p-1] = 0.0;
183 }
184 if (nrt+1 < p) {
185 e[nrt] = A[nrt][p-1];
186 }
187 e[p-1] = 0.0;
188
189 // If required, generate U.
190
191 if (wantu) {
192 for (j = nct; j < nu; j++) {
193 for (i = 0; i < m; i++) {
194 U[i][j] = 0.0;
195 }
196 U[j][j] = 1.0;
197 }
198 for (k = nct-1; k >= 0; k--) {
199 if (s[k] != 0.0) {
200 for (j = k+1; j < nu; j++) {
201 double t = 0;
202 for (i = k; i < m; i++) {
203 t += U[i][k]*U[i][j];
204 }
205 t = -t/U[k][k];
206 for (i = k; i < m; i++) {
207 U[i][j] += t*U[i][k];
208 }
209 }
210 for (i = k; i < m; i++ ) {
211 U[i][k] = -U[i][k];
212 }
213 U[k][k] = 1.0 + U[k][k];
214 for (i = 0; i < k-1; i++) {
215 U[i][k] = 0.0;
216 }
217 } else {
218 for (i = 0; i < m; i++) {
219 U[i][k] = 0.0;
220 }
221 U[k][k] = 1.0;
222 }
223 }
224 }
225
226 // If required, generate V.
227
228 if (wantv) {
229 for (k = n-1; k >= 0; k--) {
230 if ((k < nrt) & (e[k] != 0.0)) {
231 for (j = k+1; j < nu; j++) {
232 double t = 0;
233 for (i = k+1; i < n; i++) {
234 t += V[i][k]*V[i][j];
235 }
236 t = -t/V[k+1][k];
237 for (i = k+1; i < n; i++) {
238 V[i][j] += t*V[i][k];
239 }
240 }
241 }
242 for (i = 0; i < n; i++) {
243 V[i][k] = 0.0;
244 }
245 V[k][k] = 1.0;
246 }
247 }
248
249 // Main iteration loop for the singular values.
250
251 int pp = p-1;
252 int iter = 0;
253 double eps = pow(2.0,-52.0);
254 while (p > 0) {
255 int k=0;
256 int kase=0;
257
258 // Here is where a test for too many iterations would go.
259
260 // This section of the program inspects for
261 // negligible elements in the s and e arrays. On
262 // completion the variables kase and k are set as follows.
263
264 // kase = 1 if s(p) and e[k-1] are negligible and k<p
265 // kase = 2 if s(k) is negligible and k<p
266 // kase = 3 if e[k-1] is negligible, k<p, and
267 // s(k), ..., s(p) are not negligible (qr step).
268 // kase = 4 if e(p-1) is negligible (convergence).
269
270 for (k = p-2; k >= -1; k--) {
271 if (k == -1) {
272 break;
273 }
274 if (abs(e[k]) <= eps*(abs(s[k]) + abs(s[k+1]))) {
275 e[k] = 0.0;
276 break;
277 }
278 }
279 if (k == p-2) {
280 kase = 4;
281 } else {
282 int ks;
283 for (ks = p-1; ks >= k; ks--) {
284 if (ks == k) {
285 break;
286 }
287 double t = (ks != p ? abs(e[ks]) : 0.) +
288 (ks != k+1 ? abs(e[ks-1]) : 0.);
289 if (abs(s[ks]) <= eps*t) {
290 s[ks] = 0.0;
291 break;
292 }
293 }
294 if (ks == k) {
295 kase = 3;
296 } else if (ks == p-1) {
297 kase = 1;
298 } else {
299 kase = 2;
300 k = ks;
301 }
302 }
303 k++;
304
305 // Perform the task indicated by kase.
306
307 switch (kase) {
308
309 // Deflate negligible s(p).
310
311 case 1: {
312 double f = e[p-2];
313 e[p-2] = 0.0;
314 for (j = p-2; j >= k; j--) {
315 double t = hypot(s[j],f);
316 double cs = s[j]/t;
317 double sn = f/t;
318 s[j] = t;
319 if (j != k) {
320 f = -sn*e[j-1];
321 e[j-1] = cs*e[j-1];
322 }
323 if (wantv) {
324 for (i = 0; i < n; i++) {
325 t = cs*V[i][j] + sn*V[i][p-1];
326 V[i][p-1] = -sn*V[i][j] + cs*V[i][p-1];
327 V[i][j] = t;
328 }
329 }
330 }
331 }
332 break;
333
334 // Split at negligible s(k).
335
336 case 2: {
337 double f = e[k-1];
338 e[k-1] = 0.0;
339 for (j = k; j < p; j++) {
340 double t = hypot(s[j],f);
341 double cs = s[j]/t;
342 double sn = f/t;
343 s[j] = t;
344 f = -sn*e[j];
345 e[j] = cs*e[j];
346 if (wantu) {
347 for (i = 0; i < m; i++) {
348 t = cs*U[i][j] + sn*U[i][k-1];
349 U[i][k-1] = -sn*U[i][j] + cs*U[i][k-1];
350 U[i][j] = t;
351 }
352 }
353 }
354 }
355 break;
356
357 // Perform one qr step.
358
359 case 3: {
360
361 // Calculate the shift.
362
363 double scale = std::max(std::max(std::max(std::max(
364 abs(s[p-1]),abs(s[p-2])),abs(e[p-2])),
365 abs(s[k])),abs(e[k]));
366 double sp = s[p-1]/scale;
367 double spm1 = s[p-2]/scale;
368 double epm1 = e[p-2]/scale;
369 double sk = s[k]/scale;
370 double ek = e[k]/scale;
371 double b = ((spm1 + sp)*(spm1 - sp) + epm1*epm1)/2.0;
372 double c = (sp*epm1)*(sp*epm1);
373 double shift = 0.0;
374 if ((b != 0.0) || (c != 0.0)) {
375 shift = sqrt(b*b + c);
376 if (b < 0.0) {
377 shift = -shift;
378 }
379 shift = c/(b + shift);
380 }
381 double f = (sk + sp)*(sk - sp) + shift;
382 double g = sk*ek;
383
384 // Chase zeros.
385
386 for (j = k; j < p-1; j++) {
387 double t = hypot(f,g);
388 double cs = f/t;
389 double sn = g/t;
390 if (j != k) {
391 e[j-1] = t;
392 }
393 f = cs*s[j] + sn*e[j];
394 e[j] = cs*e[j] - sn*s[j];
395 g = sn*s[j+1];
396 s[j+1] = cs*s[j+1];
397 if (wantv) {
398 for (i = 0; i < n; i++) {
399 t = cs*V[i][j] + sn*V[i][j+1];
400 V[i][j+1] = -sn*V[i][j] + cs*V[i][j+1];
401 V[i][j] = t;
402 }
403 }
404 t = hypot(f,g);
405 cs = f/t;
406 sn = g/t;
407 s[j] = t;
408 f = cs*e[j] + sn*s[j+1];
409 s[j+1] = -sn*e[j] + cs*s[j+1];
410 g = sn*e[j+1];
411 e[j+1] = cs*e[j+1];
412 if (wantu && (j < m-1)) {
413 for (i = 0; i < m; i++) {
414 t = cs*U[i][j] + sn*U[i][j+1];
415 U[i][j+1] = -sn*U[i][j] + cs*U[i][j+1];
416 U[i][j] = t;
417 }
418 }
419 }
420 e[p-2] = f;
421 iter = iter + 1;
422 }
423 break;
424
425 // Convergence.
426
427 case 4: {
428
429 // Make the singular values positive.
430
431 if (s[k] <= 0.0) {
432 s[k] = (s[k] < 0.0 ? -s[k] : 0.0);
433 if (wantv) {
434 for (i = 0; i <= pp; i++) {
435 V[i][k] = -V[i][k];
436 }
437 }
438 }
439
440 // Order the singular values.
441
442 while (k < pp) {
443 if (s[k] >= s[k+1]) {
444 break;
445 }
446 double t = s[k];
447 s[k] = s[k+1];
448 s[k+1] = t;
449 if (wantv && (k < n-1)) {
450 for (i = 0; i < n; i++) {
451 t = V[i][k+1]; V[i][k+1] = V[i][k]; V[i][k] = t;
452 }
453 }
454 if (wantu && (k < m-1)) {
455 for (i = 0; i < m; i++) {
456 t = U[i][k+1]; U[i][k+1] = U[i][k]; U[i][k] = t;
457 }
458 }
459 k++;
460 }
461 iter = 0;
462 p--;
463 }
464 break;
465 }
466 }
467 }
468
469
471 {
472 int minm = std::min(m+1,n);
473
474 A = Array2D<Real>(m, minm);
475
476 for (int i=0; i<m; i++)
477 for (int j=0; j<minm; j++)
478 A[i][j] = U[i][j];
479
480 }
481
482 /* Return the right singular vectors */
483
485 {
486 A = V;
487 }
488
492 {
493 x = s;
494 }
495
500 void getS (Array2D<Real> &A) {
501 A = Array2D<Real>(n,n);
502 for (int i = 0; i < n; i++) {
503 for (int j = 0; j < n; j++) {
504 A[i][j] = 0.0;
505 }
506 A[i][i] = s[i];
507 }
508 }
509
512 double norm2 () {
513 return s[0];
514 }
515
518 double cond () {
519 return s[0]/s[std::min(m,n)-1];
520 }
521
526 int rank ()
527 {
528 double eps = pow(2.0,-52.0);
529 double tol = std::max(m,n)*s[0]*eps;
530 int r = 0;
531 for (int i = 0; i < s.dim(); i++) {
532 if (s[i] > tol) {
533 r++;
534 }
535 }
536 return r;
537 }
538};
539
540}
541#endif
542// JAMA_SVD_H