Actual source code: svdlapack.c

  1: /*                       
  2:        This file implements a wrapper to the LAPACK SVD subroutines.

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

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

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

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

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

 29: PetscErrorCode SVDSetUp_LAPACK(SVD svd)
 30: {
 32:   PetscInt       M,N;

 35:   SVDMatGetSize(svd,&M,&N);
 36:   svd->ncv = N;
 37:   if (svd->mpd) { PetscInfo(svd,"Warning: parameter mpd ignored\n"); }
 38:   svd->max_it = 1;
 39:   if (svd->ncv!=svd->n) {
 40:     VecDuplicateVecs(svd->tl,svd->ncv,&svd->U);
 41:   }
 42:   DSSetType(svd->ds,DSSVD);
 43:   DSAllocate(svd->ds,PetscMax(M,N));
 44:   return(0);
 45: }

 49: PetscErrorCode SVDSolve_LAPACK(SVD svd)
 50: {
 52:   PetscInt       M,N,n,i,j,k,ld;
 53:   Mat            mat;
 54:   PetscScalar    *pU,*pVT,*pmat,*pu,*pv,*A,*w;
 55: 
 57:   DSGetLeadingDimension(svd->ds,&ld);
 58:   MatConvert(svd->OP,MATSEQDENSE,MAT_INITIAL_MATRIX,&mat);
 59:   MatGetSize(mat,&M,&N);
 60:   DSSetDimensions(svd->ds,M,N,0,0);
 61:   MatGetArray(mat,&pmat);
 62:   DSGetArray(svd->ds,DS_MAT_A,&A);
 63:   for (i=0;i<M;i++)
 64:     for (j=0;j<N;j++)
 65:       A[i+j*ld] = pmat[i+j*M];
 66:   DSRestoreArray(svd->ds,DS_MAT_A,&A);
 67:   MatRestoreArray(mat,&pmat);
 68:   DSSetState(svd->ds,DS_STATE_RAW);
 69: 
 70:   n = PetscMin(M,N);
 71:   PetscMalloc(sizeof(PetscScalar)*n,&w);
 72:   DSSolve(svd->ds,w,PETSC_NULL);
 73:   DSSort(svd->ds,w,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);
 74: 
 75:   /* copy singular vectors */
 76:   DSGetArray(svd->ds,DS_MAT_U,&pU);
 77:   DSGetArray(svd->ds,DS_MAT_VT,&pVT);
 78:   for (i=0;i<n;i++) {
 79:     if (svd->which == SVD_SMALLEST) k = n - i - 1;
 80:     else k = i;
 81:     svd->sigma[k] = PetscRealPart(w[i]);
 82:     VecGetArray(svd->U[k],&pu);
 83:     VecGetArray(svd->V[k],&pv);
 84:     if (M>=N) {
 85:       for (j=0;j<M;j++) pu[j] = pU[i*ld+j];
 86:       for (j=0;j<N;j++) pv[j] = PetscConj(pVT[j*ld+i]);
 87:     } else {
 88:       for (j=0;j<N;j++) pu[j] = PetscConj(pVT[j*ld+i]);
 89:       for (j=0;j<M;j++) pv[j] = pU[i*ld+j];
 90:     }
 91:     VecRestoreArray(svd->U[k],&pu);
 92:     VecRestoreArray(svd->V[k],&pv);
 93:   }
 94:   DSRestoreArray(svd->ds,DS_MAT_U,&pU);
 95:   DSRestoreArray(svd->ds,DS_MAT_VT,&pVT);

 97:   svd->nconv = n;
 98:   svd->reason = SVD_CONVERGED_TOL;

100:   MatDestroy(&mat);
101:   PetscFree(w);
102:   return(0);
103: }

107: PetscErrorCode SVDReset_LAPACK(SVD svd)
108: {

112:   VecDestroyVecs(svd->n,&svd->U);
113:   return(0);
114: }

118: PetscErrorCode SVDDestroy_LAPACK(SVD svd)
119: {

123:   PetscFree(svd->data);
124:   return(0);
125: }

127: EXTERN_C_BEGIN
130: PetscErrorCode SVDCreate_LAPACK(SVD svd)
131: {
133:   svd->ops->setup   = SVDSetUp_LAPACK;
134:   svd->ops->solve   = SVDSolve_LAPACK;
135:   svd->ops->destroy = SVDDestroy_LAPACK;
136:   svd->ops->reset   = SVDReset_LAPACK;
137:   if (svd->transmode == PETSC_DECIDE)
138:     svd->transmode = SVD_TRANSPOSE_IMPLICIT; /* don't build the transpose */
139:   return(0);
140: }
141: EXTERN_C_END