Actual source code: svdscalap.c

slepc-3.14.2 2021-02-01
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:   PetscInt       M,N,m,n,mn;
 51: #if defined(PETSC_USE_COMPLEX)
 52:   PetscBLASInt   lrwork;
 53:   PetscReal      *rwork,dummy;
 54: #endif

 57:   MatGetSize(A,&M,&N);
 58:   MatGetLocalSize(A,&m,&n);
 59:   mn = (M>=N)? n: m;
 60:   MatCreate(PetscObjectComm((PetscObject)A),&Z);
 61:   MatSetSizes(Z,m,mn,PETSC_DECIDE,PETSC_DECIDE);
 62:   MatSetType(Z,MATSCALAPACK);
 63:   MatSetUp(Z);
 64:   MatAssemblyBegin(Z,MAT_FINAL_ASSEMBLY);
 65:   MatAssemblyEnd(Z,MAT_FINAL_ASSEMBLY);
 66:   z = (Mat_ScaLAPACK*)Z->data;
 67:   MatCreate(PetscObjectComm((PetscObject)A),&QT);
 68:   MatSetSizes(QT,mn,n,PETSC_DECIDE,PETSC_DECIDE);
 69:   MatSetType(QT,MATSCALAPACK);
 70:   MatSetUp(QT);
 71:   MatAssemblyBegin(QT,MAT_FINAL_ASSEMBLY);
 72:   MatAssemblyEnd(QT,MAT_FINAL_ASSEMBLY);
 73:   q = (Mat_ScaLAPACK*)QT->data;

 75: #if !defined(PETSC_USE_COMPLEX)
 76:   /* allocate workspace */
 77:   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));
 79:   PetscBLASIntCast((PetscInt)minlwork,&lwork);
 80:   PetscMalloc1(lwork,&work);
 81:   /* call computational routine */
 82:   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));
 84:   PetscFree(work);
 85: #else
 86:   /* allocate workspace */
 87:   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));
 89:   PetscBLASIntCast((PetscInt)PetscRealPart(minlwork),&lwork);
 90:   lrwork = 1+4*PetscMax(a->M,a->N);
 91:   PetscMalloc2(lwork,&work,lrwork,&rwork);
 92:   /* call computational routine */
 93:   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));
 95:   PetscFree2(work,rwork);
 96: #endif

 98:   MatHermitianTranspose(QT,MAT_INITIAL_MATRIX,&Q);
 99:   MatDestroy(&QT);
100:   BVGetMat(svd->U,&U);
101:   BVGetMat(svd->V,&V);
102:   if (M>=N) {
103:     MatConvert(Z,MATDENSE,MAT_REUSE_MATRIX,&U);
104:     MatConvert(Q,MATDENSE,MAT_REUSE_MATRIX,&V);
105:   } else {
106:     MatConvert(Q,MATDENSE,MAT_REUSE_MATRIX,&U);
107:     MatConvert(Z,MATDENSE,MAT_REUSE_MATRIX,&V);
108:   }
109:   BVRestoreMat(svd->U,&U);
110:   BVRestoreMat(svd->V,&V);
111:   MatDestroy(&Z);
112:   MatDestroy(&Q);

114:   svd->nconv  = svd->ncv;
115:   svd->its    = 1;
116:   svd->reason = SVD_CONVERGED_TOL;
117:   return(0);
118: }

120: PetscErrorCode SVDDestroy_ScaLAPACK(SVD svd)
121: {

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

129: PetscErrorCode SVDReset_ScaLAPACK(SVD svd)
130: {
132:   SVD_ScaLAPACK  *ctx = (SVD_ScaLAPACK*)svd->data;

135:   MatDestroy(&ctx->As);
136:   return(0);
137: }

139: SLEPC_EXTERN PetscErrorCode SVDCreate_ScaLAPACK(SVD svd)
140: {
142:   SVD_ScaLAPACK  *ctx;

145:   PetscNewLog(svd,&ctx);
146:   svd->data = (void*)ctx;

148:   svd->ops->solve          = SVDSolve_ScaLAPACK;
149:   svd->ops->setup          = SVDSetUp_ScaLAPACK;
150:   svd->ops->destroy        = SVDDestroy_ScaLAPACK;
151:   svd->ops->reset          = SVDReset_ScaLAPACK;
152:   return(0);
153: }