/**************************************************************************** * VCGLib o o * * Visual and Computer Graphics Library o o * * _ O _ * * Copyright(C) 2006 \/)\/ * * Visual Computing Lab /\/| * * ISTI - Italian National Research Council | * * \ * * All rights reserved. * * * * This program is free software; you can redistribute it and/or modify * * it under the terms of the GNU General Public License as published by * * the Free Software Foundation; either version 2 of the License, or * * (at your option) any later version. * * * * This program is distributed in the hope that it will be useful, * * but WITHOUT ANY WARRANTY; without even the implied warranty of * * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the * * GNU General Public License (http://www.gnu.org/licenses/gpl.txt) * * for more details. * * * ****************************************************************************/ /**************************************************************************** History $Log: not supported by cvs2svn $ Revision 1.14 2006/09/28 22:49:49 fiorin Removed some warnings Revision 1.13 2006/07/28 12:39:05 zifnab1974 added some typename directives Revision 1.12 2006/07/24 07:26:47 fiorin Changed the template argument in JacobiRotate and added method for sorting eigenvalues and eigenvectors (SortEigenvaluesAndEigenvectors) Revision 1.11 2006/05/25 09:35:55 cignoni added missing internal prototype to Sort function Revision 1.10 2006/05/17 09:26:35 cignoni Added initial disclaimer ****************************************************************************/ #ifndef __VCGLIB_LINALGEBRA_H #define __VCGLIB_LINALGEBRA_H #include namespace vcg { /** \addtogroup math */ /* @{ */ /*! * */ template< typename MATRIX_TYPE > static void JacobiRotate(MATRIX_TYPE &A, typename MATRIX_TYPE::ScalarType s, typename MATRIX_TYPE::ScalarType tau, int i,int j,int k,int l) { typename MATRIX_TYPE::ScalarType g=A[i][j]; typename MATRIX_TYPE::ScalarType h=A[k][l]; A[i][j]=g-s*(h+g*tau); A[k][l]=h+s*(g-h*tau); }; /*! * Computes all eigenvalues and eigenvectors of a real symmetric matrix . * On output, elements of the input matrix above the diagonal are destroyed. * \param d returns the eigenvalues of a. * \param v is a matrix whose columns contain, the normalized eigenvectors * \param nrot returns the number of Jacobi rotations that were required. */ template static void Jacobi(MATRIX_TYPE &w, POINT_TYPE &d, MATRIX_TYPE &v, int &nrot) { typedef typename MATRIX_TYPE::ScalarType ScalarType; assert(w.RowsNumber()==w.ColumnsNumber()); int dimension = w.RowsNumber(); int j,iq,ip,i; //assert(w.IsSymmetric()); typename MATRIX_TYPE::ScalarType tresh, theta, tau, t, sm, s, h, g, c; POINT_TYPE b, z; v.SetIdentity(); for (ip=0;ip4 && (float)(fabs(d[ip])+g) == (float)fabs(d[ip]) && (float)(fabs(d[iq])+g) == (float)fabs(d[iq])) w[ip][iq]=ScalarType(0.0); else if (fabs(w[ip][iq]) > tresh) { h=d[iq]-d[ip]; if ((float)(fabs(h)+g) == (float)fabs(h)) t=(w[ip][iq])/h; //t =1/(2#) else { theta=ScalarType(0.5)*h/(w[ip][iq]); //Equation (11.1.10). t=ScalarType(1.0)/(fabs(theta)+sqrt(ScalarType(1.0)+theta*theta)); if (theta < ScalarType(0.0)) t = -t; } c=ScalarType(1.0)/sqrt(ScalarType(1.0)+t*t); s=t*c; tau=s/(ScalarType(1.0)+c); h=t*w[ip][iq]; z[ip] -= h; z[iq] += h; d[ip] -= h; d[iq] += h; w[ip][iq]=ScalarType(0.0); for (j=0;j<=ip-1;j++) { //Case of rotations 1 <= j < p. JacobiRotate(w,s,tau,j,ip,j,iq) ; } for (j=ip+1;j<=iq-1;j++) { //Case of rotations p < j < q. JacobiRotate(w,s,tau,ip,j,j,iq); } for (j=iq+1;j(w,s,tau,ip,j,iq,j); } for (j=0;j(v,s,tau,j,ip,j,iq); } ++nrot; } } } for (ip=0;ip void SortEigenvaluesAndEigenvectors(POINT_TYPE &eigenvalues, MATRIX_TYPE &eigenvectors) { assert(eigenvectors.ColumnsNumber()==eigenvectors.RowsNumber()); int dimension = eigenvectors.ColumnsNumber(); int i, j, k; float p; for (i=0; i= p) p = eigenvalues[ k=j ]; if (k != i) { eigenvalues[k] = eigenvalues[i]; // i.e. eigenvalues[i] = p; // swaps the value of the elements i-th and k-th for (j=0; j inline static TYPE pythagora(TYPE a, TYPE b) { TYPE abs_a = fabs(a); TYPE abs_b = fabs(b); if (abs_a > abs_b) return abs_a*sqrt(1.0+sqr(abs_b/abs_a)); else return (abs_b == 0.0 ? 0.0 : abs_b*sqrt(1.0+sqr(abs_a/abs_b))); }; template inline static TYPE sign(TYPE a, TYPE b) { return (b >= 0.0 ? fabs(a) : -fabs(a)); }; template inline static TYPE sqr(TYPE a) { TYPE sqr_arg = a; return (sqr_arg == 0 ? 0 : sqr_arg*sqr_arg); } /*! * */ enum SortingStrategy {LeaveUnsorted=0, SortAscending=1, SortDescending=2}; template< typename MATRIX_TYPE > void Sort(MATRIX_TYPE &U, typename MATRIX_TYPE::ScalarType W[], MATRIX_TYPE &V, const SortingStrategy sorting) ; /*! * Given a matrix Am×n, this routine computes its singular value decomposition, * i.e. A=U·W·VT. The matrix A will be destroyed! * (This is the implementation described in Numerical Recipies). * \param A the matrix to be decomposed * \param W the diagonal matrix of singular values W, stored as a vector W[1...N] * \param V the matrix V (not the transpose VT) * \param max_iters max iteration number (default = 30). * \return */ template static bool SingularValueDecomposition(MATRIX_TYPE &A, typename MATRIX_TYPE::ScalarType *W, MATRIX_TYPE &V, const SortingStrategy sorting=LeaveUnsorted, const int max_iters=30) { typedef typename MATRIX_TYPE::ScalarType ScalarType; int m = (int) A.RowsNumber(); int n = (int) A.ColumnsNumber(); int flag,i,its,j,jj,k,l,nm; double anorm, c, f, g, h, s, scale, x, y, z, *rv1; bool convergence = true; rv1 = new double[n]; g = scale = anorm = 0; // Householder reduction to bidiagonal form. for (i=0; i( sqrt(s), f ); h = f*g - s; A[i][i]=f-g; for (j=l; j(sqrt(s),f); h = f*g - s; A[i][l] = f-g; for (k=l; k=0; i--) { //Accumulation of right-hand transformations. if (i < (n-1)) { if (g) { for (j=l; j=0; i--) { l = i+1; g = W[i]; for (j=l; j=0; k--) { for (its=1; its<=max_iters; its++) { flag=1; for (l=k; l>=0; l--) { // Test for splitting. nm=l-1; // Note that rv1[1] is always zero. if ((double)(fabs(rv1[l])+anorm) == anorm) { flag=0; break; } if ((double)(fabs(W[nm])+anorm) == anorm) break; } if (flag) { c=0.0; //Cancellation of rv1[l], if l > 1. s=1.0; for (i=l ;i<=k; i++) { f = s*rv1[i]; rv1[i] = c*rv1[i]; if ((double)(fabs(f)+anorm) == anorm) break; g = W[i]; h = pythagora(f,g); W[i] = h; h = 1.0/h; c = g*h; s = -f*h; for (j=0; j(f,1.0); f=((x-z)*(x+z) + h*((y/(f+sign(g,f)))-h))/x; c=s=1.0; //Next QR transformation: for (j=l; j<= nm;j++) { i = j+1; g = rv1[i]; y = W[i]; h = s*g; g = c*g; z = pythagora(f,h); rv1[j] = z; c = f/z; s = h/z; f = x*c + g*s; g = g*c - x*s; h = y*s; y *= c; for (jj=0; jj(f,h); W[j] = z; // Rotation can be arbitrary if z = 0. if (z) { z = 1.0/z; c = f*z; s = h*z; } f = c*g + s*y; x = c*y - s*g; for (jj=0; jj(A, W, V, sorting); return convergence; }; /*! * Sort the singular values computed by the SingularValueDecomposition procedure and * modify the matrices U and V accordingly. */ // TODO modify the last parameter type template< typename MATRIX_TYPE > void Sort(MATRIX_TYPE &U, typename MATRIX_TYPE::ScalarType W[], MATRIX_TYPE &V, const SortingStrategy sorting) { typedef typename MATRIX_TYPE::ScalarType ScalarType; assert(U.ColumnsNumber()==V.ColumnsNumber()); int mu = U.RowsNumber(); int mv = V.RowsNumber(); int n = U.ColumnsNumber(); //ScalarType* u = &U[0][0]; //ScalarType* v = &V[0][0]; for (int i=0; i p) { k = j; p = W[j]; } } break; } } if (k != i) { W[k] = W[i]; // i.e. W[i] = p; // swaps the i-th and the k-th elements int j = mu; //ScalarType* uji = u + i; // uji = &U[0][i] //ScalarType* ujk = u + k; // ujk = &U[0][k] //ScalarType* vji = v + i; // vji = &V[0][i] //ScalarType* vjk = v + k; // vjk = &V[0][k] //if (j) //{ // for(;;) for( ; j!=0; --j, uji+=n, ujk+=n) // { { // p = *uji; p = *uji; // i.e. // *uji = *ujk; *uji = *ujk; // swap( U[s][i], U[s][k] ) // *ujk = p; *ujk = p; // // if (!(--j)) } // break; // uji += n; // ujk += n; // } //} for(int s=0; j!=0; ++s, --j) std::swap(U[s][i], U[s][k]); j = mv; //if (j!=0) //{ // for(;;) for ( ; j!=0; --j, vji+=n, ujk+=n) // { { // p = *vji; p = *vji; // i.e. // *vji = *vjk; *vji = *vjk; // swap( V[s][i], V[s][k] ) // *vjk = p; *vjk = p; // // if (!(--j)) } // break; // vji += n; // vjk += n; // } //} for (int s=0; j!=0; ++s, --j) std::swap(V[s][i], V[s][k]); } } } /*! * Solves A·X = B for a vector X, where A is specified by the matrices Um×n, * Wn×1 and Vn×n as returned by SingularValueDecomposition. * No input quantities are destroyed, so the routine may be called sequentially with different b’s. * \param x is the output solution vector (xn×1) * \param b is the input right-hand side (bn×1) */ template static void SingularValueBacksubstitution(const MATRIX_TYPE &U, const typename MATRIX_TYPE::ScalarType *W, const MATRIX_TYPE &V, typename MATRIX_TYPE::ScalarType *x, const typename MATRIX_TYPE::ScalarType *b) { typedef typename MATRIX_TYPE::ScalarType ScalarType; unsigned int jj, j, i; unsigned int columns_number = U.ColumnsNumber(); unsigned int rows_number = U.RowsNumber(); ScalarType s; ScalarType *tmp = new ScalarType[columns_number]; for (j=0; j