Actual source code: svdscalap.c

slepc-3.18.0 2022-10-01
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:   SVDCheckDefinite(svd);
 28:   MatGetSize(svd->A,&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: {
 44:   SVD_ScaLAPACK  *ctx = (SVD_ScaLAPACK*)svd->data;
 45:   Mat            A = ctx->As,Z,Q,QT,U,V;
 46:   Mat_ScaLAPACK  *a = (Mat_ScaLAPACK*)A->data,*q,*z;
 47:   PetscScalar    *work,minlwork;
 48:   PetscBLASInt   info,lwork=-1,one=1;
 49:   PetscInt       M,N,m,n,mn;
 50: #if defined(PETSC_USE_COMPLEX)
 51:   PetscBLASInt   lrwork;
 52:   PetscReal      *rwork,dummy;
 53: #endif

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

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

 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: {
122:   PetscFree(svd->data);
123:   return 0;
124: }

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

130:   MatDestroy(&ctx->As);
131:   return 0;
132: }

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

138:   PetscNew(&ctx);
139:   svd->data = (void*)ctx;

141:   svd->ops->solve          = SVDSolve_ScaLAPACK;
142:   svd->ops->setup          = SVDSetUp_ScaLAPACK;
143:   svd->ops->destroy        = SVDDestroy_ScaLAPACK;
144:   svd->ops->reset          = SVDReset_ScaLAPACK;
145:   return 0;
146: }