/** Singular Value Decomposition. U,S,V=svd(A);get...(): A=U*S*V^_T_
For an m-by-n matrix A with m >= n, the singular value decomposition is an m-by-n orthogonal matrix U, an n-by-n diagonal matrix S, and an n-by-n orthogonal matrix V so that A = U*S*V'.
The singular values, sigma[k] = S[k][k], are ordered so that sigma[0] >= sigma[1] >= ... >= sigma[n-1].
The singular value decompostion always exists, so the constructor will never fail. The matrix condition number and the effective numerical rank can be computed from this decomposition.
(Adapted from JAMA, a Java Matrix Library, developed by jointly
by the Mathworks and NIST; see http://math.nist.gov/javanumerics/jama).
*/
#ifdef JAMA_DYNAMIC
int m, n;
int nu,minm;
#else
#define nu MIN(m,n)
#define minm MIN(m+1,n)
#endif
JamaArray2(Real,m,nu) U;
JamaArray2(Real,n,n) V;
JamaArray1(Real,minm) s;
PUBLIC
member(SVD) (const JamaArray2(Real,m,n) &Arg) {
SVD0(Arg,1,1);
}
member(SVD) (const JamaArray2(Real,m,n) &Arg,int wantu,int wantv) {
SVD0(Arg,wantu,wantv);
}
void member(SVD0) (const JamaArray2(Real,m,n) &Arg,int wantu,int wantv) {
//allocation
#ifdef JAMA_DYNAMIC
m = Arg.NROWS();
n = Arg.NCOLS();
nu = MIN(m,n);//length of diagonal
minm=MIN(m+1,n);//width of bidiagonal matrix.
//s = JamaArray1(Real,minm)(minm);
//U = JamaArray2(Real,m,nu)(m, nu);//, (0)
//V = JamaArray2(Real,n,n)(n,n);
RESIZE1(Real, s, minm);
RESIZE2(Real, U, m ,nu);
RESIZE2(Real, V, n, n);
JamaArray1(Real,n) e(n); //PZ: e serves as storage for upper diagonal elements and also temporarily storage of last portion of rows.
//PZ: thus, minm-1 elements would be not enough.
//e[0]=D[0][1]
JamaArray1(Real,m) work(m);
JamaArray2(Real,m,n) A(COPY(Arg));
//end of allocation
for (int i_=0;i_ = -1; k--) {
if (k == -1) {
break;
}
if (ABS($ e[k]) <= eps*(ABS($ s[k]) + ABS($ s[k+1]))) {
$ e[k] = 0.0;
break;
}
}
if (k == p-2) {
kase = 4;
} else {
int ks;
for (ks = p-1; ks >= k; ks--) {
if (ks == k) {
break;
}
Real t=(Real)( (ks != p ? ABS($ e[ks]) : 0.) +
(ks != k+1 ? ABS($ e[ks-1]) : 0.));
if (ABS($ s[ks]) <= eps*t) {
$ s[ks] = 0.0;
break;
}
}
if (ks == k) {
kase = 3;
} else if (ks == p-1) {
kase = 1;
} else {
kase = 2;
k = ks;
}
}
k++;
// Perform the task indicated by kase.
switch (kase) {
// Deflate negligible s(p).
case 1: {
Real f($ e[p-2]);
$ e[p-2] = 0.0;
for (j = p-2; j >= k; j--) {
Real t( hypot_($ s[j],f));
Real cs($ s[j]/t);
Real sn(f/t);
$ s[j] = t;
if (j != k) {
f = -sn*$ e[j-1];
$ e[j-1] = cs*$ e[j-1];
}
if (wantv) {
for (i = 0; i < n; i++) {
t = cs*$ V[i][j] + sn*$ V[i][p-1];
$ V[i][p-1] = -sn*$ V[i][j] + cs*$ V[i][p-1];
$ V[i][j] = t;
}
}
}
}
break;
// Split at negligible s(k).
case 2: {
Real f($ e[k-1]);
$ e[k-1] = 0.0;
for (j = k; j < p; j++) {
Real t(hypot_($ s[j],f));
Real cs( $ s[j]/t);
Real sn(f/t);
$ s[j] = t;
f = -sn*$ e[j];
$ e[j] = cs*$ e[j];
if (wantu) {
for (i = 0; i < m; i++) {
t = cs*$ U[i][j] + sn*$ U[i][k-1];
$ U[i][k-1] = -sn*$ U[i][j] + cs*$ U[i][k-1];
$ U[i][j] = t;
}
}
}
}
break;
// Perform one qr step.
case 3: {
// Calculate the shift.
Real scale = MAX(MAX(MAX(MAX(
ABS($ s[p-1]),ABS($ s[p-2])),ABS($ e[p-2])),
ABS($ s[k])),ABS($ e[k]));
Real sp = $ s[p-1]/scale;//PZ: calculating 1/scale would not include numeric overflow....?
Real spm1 = $ s[p-2]/scale;
Real epm1 = $ e[p-2]/scale;
Real sk = $ s[k]/scale;
Real ek = $ e[k]/scale;
Real b = (Real)(((spm1 + sp)*(spm1 - sp) + epm1*epm1)/2.0);
Real c = (sp*epm1)*(sp*epm1);
Real shift = 0.0;
if ((b != 0.0) || (c != 0.0)) {
shift = sqrt(b*b + c);
if (b < 0.0) {
shift = -shift;
}
shift = c/(b + shift);
}
Real f = (sk + sp)*(sk - sp) + shift;
Real g = sk*ek;
// Chase zeros.
for (j = k; j < p-1; j++) {
Real t = hypot_(f,g);
Real cs = f/t;
Real sn = g/t;
if (j != k) {
$ e[j-1] = t;
}
f = cs*$ s[j] + sn*$ e[j];
$ e[j] = cs*$ e[j] - sn*$ s[j];
g = sn*$ s[j+1];
$ s[j+1] = cs*$ s[j+1];
if (wantv) {
for (i = 0; i < n; i++) {
t = cs*$ V[i][j] + sn*$ V[i][j+1];
$ V[i][j+1] = -sn*$ V[i][j] + cs*$ V[i][j+1];
$ V[i][j] = t;
}
}
t = hypot_(f,g);
cs = f/t;
sn = g/t;
$ s[j] = t;
f = cs*$ e[j] + sn*$ s[j+1];
$ s[j+1] = -sn*$ e[j] + cs*$ s[j+1];
g = sn*$ e[j+1];
$ e[j+1] = cs*$ e[j+1];
if (wantu && (j < m-1)) {
for (i = 0; i < m; i++) {
t = cs*$ U[i][j] + sn*$ U[i][j+1];
$ U[i][j+1] = -sn*$ U[i][j] + cs*$ U[i][j+1];
$ U[i][j] = t;
}
}
}
$ e[p-2] = f;
iter = iter + 1;
}
break;
// Convergence.
case 4: {
// Make the singular values positive.
if ($ s[k] <= 0.0) {
$ s[k] = (Real)($ s[k] < 0.0 ? -$ s[k] : 0.0);
if (wantv) {
for (i = 0; i <= pp; i++) {
$ V[i][k] = -$ V[i][k];
}
}
}
// Order the singular values by bubble sort
while (k < pp) {
if ($ s[k] >= $ s[k+1]) {
break;//nothing to do
}
//PZ exchange s[k <> k+1], U[:,k<>k+1], V[:,k<>k+1],
Real t = $ s[k];
$ s[k] = $ s[k+1];
$ s[k+1] = t;
if (wantv && (k < n-1)) {
for (i = 0; i < n; i++) {
t = $ V[i][k+1]; $ V[i][k+1] = $ V[i][k]; $ V[i][k] = t;
}
}
if (wantu && (k < m-1)) {
for (i = 0; i < m; i++) {
t = $ U[i][k+1]; $ U[i][k+1] = $ U[i][k]; $ U[i][k] = t;
}
}
k++;
}
iter = 0;
p--;
}
break;
}
}
}
// void member(getU) (JamaArray2(Real,m,minm) &A)
void member(getU) (JamaArray2(Real,m,nu) &A) //PZ
{
// int minm = MIN(m+1,n);
// RESIZE2(Real,A,m,minm);
RESIZE2(Real,A,m,nu);//PZ
//#ifdef __cplusplus
// A = JamaArray2(Real,m,minm)(m, minm);
//#endif
for (int i=0; i