Actual source code: svdscalap.c

slepc-3.17.1 2022-04-11
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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: {
 23:   SVD_ScaLAPACK  *ctx = (SVD_ScaLAPACK*)svd->data;
 24:   PetscInt       M,N;

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

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

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

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

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

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

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

119: PetscErrorCode SVDDestroy_ScaLAPACK(SVD svd)
120: {
121:   PetscFree(svd->data);
122:   PetscFunctionReturn(0);
123: }

125: PetscErrorCode SVDReset_ScaLAPACK(SVD svd)
126: {
127:   SVD_ScaLAPACK  *ctx = (SVD_ScaLAPACK*)svd->data;

129:   MatDestroy(&ctx->As);
130:   PetscFunctionReturn(0);
131: }

133: SLEPC_EXTERN PetscErrorCode SVDCreate_ScaLAPACK(SVD svd)
134: {
135:   SVD_ScaLAPACK  *ctx;

137:   PetscNewLog(svd,&ctx);
138:   svd->data = (void*)ctx;

140:   svd->ops->solve          = SVDSolve_ScaLAPACK;
141:   svd->ops->setup          = SVDSetUp_ScaLAPACK;
142:   svd->ops->destroy        = SVDDestroy_ScaLAPACK;
143:   svd->ops->reset          = SVDReset_ScaLAPACK;
144:   PetscFunctionReturn(0);
145: }