Actual source code: svdscalap.c

slepc-3.14.0 2020-09-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    This file implements a wrapper to the ScaLAPACK SVD solver
 12: */

 14: #include <slepc/private/svdimpl.h>
 15: #include <slepc/private/slepcscalapack.h>

 17: typedef struct {
 18:   Mat As;        /* converted matrix */
 19: } SVD_ScaLAPACK;

 21: PetscErrorCode SVDSetUp_ScaLAPACK(SVD svd)
 22: {
 24:   SVD_ScaLAPACK  *ctx = (SVD_ScaLAPACK*)svd->data;
 25:   PetscInt       M,N;

 28:   SVDMatGetSize(svd,&M,&N);
 29:   svd->ncv = N;
 30:   if (svd->mpd!=PETSC_DEFAULT) { PetscInfo(svd,"Warning: parameter mpd ignored\n"); }
 31:   if (svd->max_it==PETSC_DEFAULT) svd->max_it = 1;
 32:   svd->leftbasis = PETSC_TRUE;
 33:   SVDCheckUnsupported(svd,SVD_FEATURE_STOPPING);
 34:   SVDAllocateSolution(svd,0);

 36:   /* convert matrix */
 37:   MatDestroy(&ctx->As);
 38:   MatConvert(svd->OP,MATSCALAPACK,MAT_INITIAL_MATRIX,&ctx->As);
 39:   return(0);
 40: }

 42: PetscErrorCode SVDSolve_ScaLAPACK(SVD svd)
 43: {
 45:   SVD_ScaLAPACK  *ctx = (SVD_ScaLAPACK*)svd->data;
 46:   Mat            A = ctx->As,Z,Q,QT,U,V;
 47:   Mat_ScaLAPACK  *a = (Mat_ScaLAPACK*)A->data,*q,*z;
 48:   PetscScalar    *work,minlwork;
 49:   PetscBLASInt   info,lwork=-1,one=1;
 50: #if defined(PETSC_USE_COMPLEX)
 51:   PetscBLASInt   lrwork;
 52:   PetscReal      *rwork,dummy;
 53: #endif

 56:   MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Z);
 57:   z = (Mat_ScaLAPACK*)Z->data;
 58:   MatCreate(PetscObjectComm((PetscObject)A),&QT);
 59:   MatSetSizes(QT,A->cmap->n,A->cmap->n,PETSC_DECIDE,PETSC_DECIDE);
 60:   MatSetType(QT,MATSCALAPACK);
 61:   MatSetUp(QT);
 62:   MatAssemblyBegin(QT,MAT_FINAL_ASSEMBLY);
 63:   MatAssemblyEnd(QT,MAT_FINAL_ASSEMBLY);
 64:   q = (Mat_ScaLAPACK*)QT->data;

 66: #if !defined(PETSC_USE_COMPLEX)
 67:   /* allocate workspace */
 68:   PetscStackCallBLAS("SCALAPACKgesvd",SCALAPACKgesvd_("V","V",&a->M,&a->N,a->loc,&one,&one,a->desc,svd->sigma,z->loc,&one,&one,z->desc,q->loc,&one,&one,q->desc,&minlwork,&lwork,&info));
 70:   PetscBLASIntCast(minlwork,&lwork);
 71:   PetscMalloc1(lwork,&work);
 72:   /* call computational routine */
 73:   PetscStackCallBLAS("SCALAPACKgesvd",SCALAPACKgesvd_("V","V",&a->M,&a->N,a->loc,&one,&one,a->desc,svd->sigma,z->loc,&one,&one,z->desc,q->loc,&one,&one,q->desc,work,&lwork,&info));
 75:   PetscFree(work);
 76: #else
 77:   /* allocate workspace */
 78:   PetscStackCallBLAS("SCALAPACKgesvd",SCALAPACKgesvd_("V","V",&a->M,&a->N,a->loc,&one,&one,a->desc,svd->sigma,z->loc,&one,&one,z->desc,q->loc,&one,&one,q->desc,&minlwork,&lwork,&dummy,&info));
 80:   PetscBLASIntCast(minlwork,&lwork);
 81:   lrwork = 1+4*PetscMax(a->M,a->N);
 82:   PetscMalloc2(lwork,&work,lrwork,&rwork);
 83:   /* call computational routine */
 84:   PetscStackCallBLAS("SCALAPACKgesvd",SCALAPACKgesvd_("V","V",&a->M,&a->N,a->loc,&one,&one,a->desc,svd->sigma,z->loc,&one,&one,z->desc,q->loc,&one,&one,q->desc,work,&lwork,rwork,&info));
 86:   PetscFree2(work,rwork);
 87: #endif

 89:   BVGetMat(svd->U,&U);
 90:   MatConvert(Z,MATDENSE,MAT_REUSE_MATRIX,&U);
 91:   BVRestoreMat(svd->U,&U);
 92:   MatDestroy(&Z);
 93:   MatHermitianTranspose(QT,MAT_INITIAL_MATRIX,&Q);
 94:   MatDestroy(&QT);
 95:   BVGetMat(svd->V,&V);
 96:   MatConvert(Q,MATDENSE,MAT_REUSE_MATRIX,&V);
 97:   BVRestoreMat(svd->V,&V);
 98:   MatDestroy(&Q);

100:   svd->nconv  = svd->ncv;
101:   svd->its    = 1;
102:   svd->reason = SVD_CONVERGED_TOL;
103:   return(0);
104: }

106: PetscErrorCode SVDDestroy_ScaLAPACK(SVD svd)
107: {

111:   PetscFree(svd->data);
112:   return(0);
113: }

115: PetscErrorCode SVDReset_ScaLAPACK(SVD svd)
116: {
118:   SVD_ScaLAPACK  *ctx = (SVD_ScaLAPACK*)svd->data;

121:   MatDestroy(&ctx->As);
122:   return(0);
123: }

125: SLEPC_EXTERN PetscErrorCode SVDCreate_ScaLAPACK(SVD svd)
126: {
128:   SVD_ScaLAPACK  *ctx;

131:   PetscNewLog(svd,&ctx);
132:   svd->data = (void*)ctx;

134:   svd->ops->solve          = SVDSolve_ScaLAPACK;
135:   svd->ops->setup          = SVDSetUp_ScaLAPACK;
136:   svd->ops->destroy        = SVDDestroy_ScaLAPACK;
137:   svd->ops->reset          = SVDReset_ScaLAPACK;
138:   return(0);
139: }