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-2013, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

 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:     PetscLogObjectParents(svd,svd->ncv,svd->U);
 42:   }
 43:   DSSetType(svd->ds,DSSVD);
 44:   DSAllocate(svd->ds,PetscMax(M,N));
 45:   return(0);
 46: }

 50: PetscErrorCode SVDSolve_LAPACK(SVD svd)
 51: {
 53:   PetscInt       M,N,n,i,j,k,ld;
 54:   Mat            mat;
 55:   PetscScalar    *pU,*pVT,*pmat,*pu,*pv,*A,*w;

 58:   DSGetLeadingDimension(svd->ds,&ld);
 59:   MatConvert(svd->OP,MATSEQDENSE,MAT_INITIAL_MATRIX,&mat);
 60:   MatGetSize(mat,&M,&N);
 61:   DSSetDimensions(svd->ds,M,N,0,0);
 62:   MatDenseGetArray(mat,&pmat);
 63:   DSGetArray(svd->ds,DS_MAT_A,&A);
 64:   for (i=0;i<M;i++)
 65:     for (j=0;j<N;j++)
 66:       A[i+j*ld] = pmat[i+j*M];
 67:   DSRestoreArray(svd->ds,DS_MAT_A,&A);
 68:   MatDenseRestoreArray(mat,&pmat);
 69:   DSSetState(svd->ds,DS_STATE_RAW);

 71:   n = PetscMin(M,N);
 72:   PetscMalloc(sizeof(PetscScalar)*n,&w);
 73:   DSSolve(svd->ds,w,NULL);
 74:   DSSort(svd->ds,w,NULL,NULL,NULL,NULL);

 76:   /* copy singular vectors */
 77:   DSGetArray(svd->ds,DS_MAT_U,&pU);
 78:   DSGetArray(svd->ds,DS_MAT_VT,&pVT);
 79:   for (i=0;i<n;i++) {
 80:     if (svd->which == SVD_SMALLEST) k = n - i - 1;
 81:     else k = i;
 82:     svd->sigma[k] = PetscRealPart(w[i]);
 83:     VecGetArray(svd->U[k],&pu);
 84:     VecGetArray(svd->V[k],&pv);
 85:     if (M>=N) {
 86:       for (j=0;j<M;j++) pu[j] = pU[i*ld+j];
 87:       for (j=0;j<N;j++) pv[j] = PetscConj(pVT[j*ld+i]);
 88:     } else {
 89:       for (j=0;j<N;j++) pu[j] = PetscConj(pVT[j*ld+i]);
 90:       for (j=0;j<M;j++) pv[j] = pU[i*ld+j];
 91:     }
 92:     VecRestoreArray(svd->U[k],&pu);
 93:     VecRestoreArray(svd->V[k],&pv);
 94:   }
 95:   DSRestoreArray(svd->ds,DS_MAT_U,&pU);
 96:   DSRestoreArray(svd->ds,DS_MAT_VT,&pVT);

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

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

108: PetscErrorCode SVDReset_LAPACK(SVD svd)
109: {

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

119: PetscErrorCode SVDDestroy_LAPACK(SVD svd)
120: {

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

130: PETSC_EXTERN 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: }