Actual source code: svddense.c

  1: /*                       
  2:      This file contains routines for handling small-size dense problems.
  3:      All routines are simply wrappers to LAPACK routines.

  5:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  6:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  7:    Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain

  9:    This file is part of SLEPc.
 10:       
 11:    SLEPc is free software: you can redistribute it and/or modify it under  the
 12:    terms of version 3 of the GNU Lesser General Public License as published by
 13:    the Free Software Foundation.

 15:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY 
 16:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS 
 17:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for 
 18:    more details.

 20:    You  should have received a copy of the GNU Lesser General  Public  License
 21:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 22:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 23: */

 25: #include <private/svdimpl.h>        /*I "slepcsvd.h" I*/
 26: #include <slepcblaslapack.h>

 30: /*@
 31:    SVDDense - Solves a dense singular value problem.

 33:    Not Collective

 35:    Input Parameters:
 36: +  M  - dimension of the problem (rows)
 37: .  N  - dimension of the problem (colums)
 38: -  A  - pointer to the array containing the matrix values

 40:    Output Parameters:
 41: +  sigma  - pointer to the array to store the computed singular values
 42: .  U  - pointer to the array to store left singular vectors
 43: -  VT  - pointer to the array to store right singular vectors

 45:    Notes:
 46:    Matrix A is overwritten.
 47:    
 48:    This routine uses LAPACK routines xGESDD with JOBZ='O'. Thus, if M>=N
 49:    then U is not referenced and the left singular vectors are returned
 50:    in A, and if M<N then VT is not referenced and the right singular
 51:    vectors are returned in A.

 53:    Level: developer
 54: @*/
 55: PetscErrorCode SVDDense(PetscInt M_,PetscInt N_,PetscScalar* A,PetscReal* sigma,PetscScalar* U,PetscScalar* VT)
 56: {
 57: #if defined(SLEPC_MISSING_LAPACK_GESDD)
 59:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESDD - Lapack routine is unavailable.");
 60: #else
 62:   PetscScalar    qwork,*work;
 63:   PetscBLASInt   n,info,lwork,*iwork,M,N;
 64: #if defined(PETSC_USE_COMPLEX)
 65:   PetscReal      *rwork;
 66: #endif 
 67: 
 69:   /* workspace query & allocation */
 70:   PetscLogEventBegin(SVD_Dense,0,0,0,0);
 71:   M = PetscBLASIntCast(M_);
 72:   N = PetscBLASIntCast(N_);
 73:   n = PetscMin(M,N);
 74:   PetscMalloc(sizeof(PetscInt)*8*n,&iwork);
 75:   lwork = -1;
 76: #if defined(PETSC_USE_COMPLEX)
 77:   PetscMalloc(sizeof(PetscReal)*(5*n*n+7*n),&rwork);
 78:   LAPACKgesdd_("O",&M,&N,A,&M,sigma,U,&M,VT,&N,&qwork,&lwork,rwork,iwork,&info);
 79: #else
 80:   LAPACKgesdd_("O",&M,&N,A,&M,sigma,U,&M,VT,&N,&qwork,&lwork,iwork,&info);
 81: #endif 
 82:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGESDD %d",info);
 83:   lwork = (PetscBLASInt)PetscRealPart(qwork);
 84:   PetscMalloc(sizeof(PetscScalar)*lwork,&work);
 85: 
 86:   /* computation */
 87: #if defined(PETSC_USE_COMPLEX)
 88:   LAPACKgesdd_("O",&M,&N,A,&M,sigma,U,&M,VT,&N,work,&lwork,rwork,iwork,&info);
 89:   PetscFree(rwork);
 90: #else
 91:   LAPACKgesdd_("O",&M,&N,A,&M,sigma,U,&M,VT,&N,work,&lwork,iwork,&info);
 92: #endif
 93:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGESDD %d",info);
 94:   PetscFree(iwork);
 95:   PetscFree(work);
 96:   PetscLogEventEnd(SVD_Dense,0,0,0,0);
 97:   return(0);
 98: #endif 
 99: }