Actual source code: scalapack.c

slepc-3.22.2 2024-12-02
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 eigensolvers in ScaLAPACK.
 12: */

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

 17: typedef struct {
 18:   Mat As,Bs;        /* converted matrices */
 19: } EPS_ScaLAPACK;

 21: static PetscErrorCode EPSSetUp_ScaLAPACK(EPS eps)
 22: {
 23:   EPS_ScaLAPACK  *ctx = (EPS_ScaLAPACK*)eps->data;
 24:   Mat            A,B;
 25:   PetscInt       nmat;
 26:   PetscBool      isshift;
 27:   PetscScalar    shift;

 29:   PetscFunctionBegin;
 30:   EPSCheckHermitianDefinite(eps);
 31:   EPSCheckNotStructured(eps);
 32:   PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift));
 33:   PetscCheck(isshift,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support spectral transformations");
 34:   eps->ncv = eps->n;
 35:   if (eps->mpd!=PETSC_DETERMINE) PetscCall(PetscInfo(eps,"Warning: parameter mpd ignored\n"));
 36:   if (eps->max_it==PETSC_DETERMINE) eps->max_it = 1;
 37:   if (!eps->which) PetscCall(EPSSetWhichEigenpairs_Default(eps));
 38:   PetscCheck(eps->which!=EPS_ALL || eps->inta==eps->intb,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support interval computation");
 39:   EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_STOPPING);
 40:   EPSCheckIgnored(eps,EPS_FEATURE_EXTRACTION | EPS_FEATURE_CONVERGENCE);
 41:   PetscCall(EPSAllocateSolution(eps,0));

 43:   /* convert matrices */
 44:   PetscCall(MatDestroy(&ctx->As));
 45:   PetscCall(MatDestroy(&ctx->Bs));
 46:   PetscCall(STGetNumMatrices(eps->st,&nmat));
 47:   PetscCall(STGetMatrix(eps->st,0,&A));
 48:   PetscCall(MatConvert(A,MATSCALAPACK,MAT_INITIAL_MATRIX,&ctx->As));
 49:   if (nmat>1) {
 50:     PetscCall(STGetMatrix(eps->st,1,&B));
 51:     PetscCall(MatConvert(B,MATSCALAPACK,MAT_INITIAL_MATRIX,&ctx->Bs));
 52:   }
 53:   PetscCall(STGetShift(eps->st,&shift));
 54:   if (shift != 0.0) {
 55:     if (nmat>1) PetscCall(MatAXPY(ctx->As,-shift,ctx->Bs,SAME_NONZERO_PATTERN));
 56:     else PetscCall(MatShift(ctx->As,-shift));
 57:   }
 58:   PetscFunctionReturn(PETSC_SUCCESS);
 59: }

 61: static PetscErrorCode EPSSolve_ScaLAPACK(EPS eps)
 62: {
 63:   EPS_ScaLAPACK  *ctx = (EPS_ScaLAPACK*)eps->data;
 64:   Mat            A = ctx->As,B = ctx->Bs,Q,V;
 65:   Mat_ScaLAPACK  *a = (Mat_ScaLAPACK*)A->data,*b,*q;
 66:   PetscReal      rdummy=0.0,abstol=0.0,*gap=NULL,orfac=-1.0,*w = eps->errest;  /* used to store real eigenvalues */
 67:   PetscScalar    *work,minlwork[3];
 68:   PetscBLASInt   i,m,info,idummy=0,lwork=-1,liwork=-1,minliwork,*iwork,*ifail=NULL,*iclustr=NULL,one=1;
 69: #if defined(PETSC_USE_COMPLEX)
 70:   PetscReal      *rwork,minlrwork[3];
 71:   PetscBLASInt   lrwork=-1;
 72: #endif

 74:   PetscFunctionBegin;
 75:   PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Q));
 76:   PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
 77:   q = (Mat_ScaLAPACK*)Q->data;

 79:   if (B) {

 81:     b = (Mat_ScaLAPACK*)B->data;
 82:     PetscCall(PetscMalloc3(a->grid->nprow*a->grid->npcol,&gap,a->N,&ifail,2*a->grid->nprow*a->grid->npcol,&iclustr));
 83: #if !defined(PETSC_USE_COMPLEX)
 84:     /* allocate workspace */
 85:     PetscCallBLAS("SCALAPACKsygvx",SCALAPACKsygvx_(&one,"V","A","L",&a->N,a->loc,&one,&one,a->desc,b->loc,&one,&one,b->desc,&rdummy,&rdummy,&idummy,&idummy,&abstol,&m,&idummy,w,&orfac,q->loc,&one,&one,q->desc,minlwork,&lwork,&minliwork,&liwork,ifail,iclustr,gap,&info));
 86:     PetscCheckScaLapackInfo("sygvx",info);
 87:     PetscCall(PetscBLASIntCast((PetscInt)minlwork[0],&lwork));
 88:     liwork = minliwork;
 89:     /* call computational routine */
 90:     PetscCall(PetscMalloc2(lwork,&work,liwork,&iwork));
 91:     PetscCallBLAS("SCALAPACKsygvx",SCALAPACKsygvx_(&one,"V","A","L",&a->N,a->loc,&one,&one,a->desc,b->loc,&one,&one,b->desc,&rdummy,&rdummy,&idummy,&idummy,&abstol,&m,&idummy,w,&orfac,q->loc,&one,&one,q->desc,work,&lwork,iwork,&liwork,ifail,iclustr,gap,&info));
 92:     PetscCheckScaLapackInfo("sygvx",info);
 93:     PetscCall(PetscFree2(work,iwork));
 94: #else
 95:     /* allocate workspace */
 96:     PetscCallBLAS("SCALAPACKsygvx",SCALAPACKsygvx_(&one,"V","A","L",&a->N,a->loc,&one,&one,a->desc,b->loc,&one,&one,b->desc,&rdummy,&rdummy,&idummy,&idummy,&abstol,&m,&idummy,w,&orfac,q->loc,&one,&one,q->desc,minlwork,&lwork,minlrwork,&lrwork,&minliwork,&liwork,ifail,iclustr,gap,&info));
 97:     PetscCheckScaLapackInfo("sygvx",info);
 98:     PetscCall(PetscBLASIntCast((PetscInt)PetscRealPart(minlwork[0]),&lwork));
 99:     PetscCall(PetscBLASIntCast((PetscInt)minlrwork[0],&lrwork));
100:     lrwork += a->N*a->N;
101:     liwork = minliwork;
102:     /* call computational routine */
103:     PetscCall(PetscMalloc3(lwork,&work,lrwork,&rwork,liwork,&iwork));
104:     PetscCallBLAS("SCALAPACKsygvx",SCALAPACKsygvx_(&one,"V","A","L",&a->N,a->loc,&one,&one,a->desc,b->loc,&one,&one,b->desc,&rdummy,&rdummy,&idummy,&idummy,&abstol,&m,&idummy,w,&orfac,q->loc,&one,&one,q->desc,work,&lwork,rwork,&lrwork,iwork,&liwork,ifail,iclustr,gap,&info));
105:     PetscCheckScaLapackInfo("sygvx",info);
106:     PetscCall(PetscFree3(work,rwork,iwork));
107: #endif
108:     PetscCall(PetscFree3(gap,ifail,iclustr));

110:   } else {

112: #if !defined(PETSC_USE_COMPLEX)
113:     /* allocate workspace */
114:     PetscCallBLAS("SCALAPACKsyev",SCALAPACKsyev_("V","L",&a->N,a->loc,&one,&one,a->desc,w,q->loc,&one,&one,q->desc,minlwork,&lwork,&info));
115:     PetscCheckScaLapackInfo("syev",info);
116:     PetscCall(PetscBLASIntCast((PetscInt)minlwork[0],&lwork));
117:     PetscCall(PetscMalloc1(lwork,&work));
118:     /* call computational routine */
119:     PetscCallBLAS("SCALAPACKsyev",SCALAPACKsyev_("V","L",&a->N,a->loc,&one,&one,a->desc,w,q->loc,&one,&one,q->desc,work,&lwork,&info));
120:     PetscCheckScaLapackInfo("syev",info);
121:     PetscCall(PetscFree(work));
122: #else
123:     /* allocate workspace */
124:     PetscCallBLAS("SCALAPACKsyev",SCALAPACKsyev_("V","L",&a->N,a->loc,&one,&one,a->desc,w,q->loc,&one,&one,q->desc,minlwork,&lwork,minlrwork,&lrwork,&info));
125:     PetscCheckScaLapackInfo("syev",info);
126:     PetscCall(PetscBLASIntCast((PetscInt)PetscRealPart(minlwork[0]),&lwork));
127:     lrwork = 4*a->N;  /* PetscCall(PetscBLASIntCast((PetscInt)minlrwork[0],&lrwork)); */
128:     PetscCall(PetscMalloc2(lwork,&work,lrwork,&rwork));
129:     /* call computational routine */
130:     PetscCallBLAS("SCALAPACKsyev",SCALAPACKsyev_("V","L",&a->N,a->loc,&one,&one,a->desc,w,q->loc,&one,&one,q->desc,work,&lwork,rwork,&lrwork,&info));
131:     PetscCheckScaLapackInfo("syev",info);
132:     PetscCall(PetscFree2(work,rwork));
133: #endif

135:   }
136:   PetscCall(PetscFPTrapPop());

138:   for (i=0;i<eps->ncv;i++) {
139:     eps->eigr[i]   = eps->errest[i];
140:     eps->errest[i] = PETSC_MACHINE_EPSILON;
141:   }

143:   PetscCall(BVGetMat(eps->V,&V));
144:   PetscCall(MatConvert(Q,MATDENSE,MAT_REUSE_MATRIX,&V));
145:   PetscCall(BVRestoreMat(eps->V,&V));
146:   PetscCall(MatDestroy(&Q));

148:   eps->nconv  = eps->ncv;
149:   eps->its    = 1;
150:   eps->reason = EPS_CONVERGED_TOL;
151:   PetscFunctionReturn(PETSC_SUCCESS);
152: }

154: static PetscErrorCode EPSDestroy_ScaLAPACK(EPS eps)
155: {
156:   PetscFunctionBegin;
157:   PetscCall(PetscFree(eps->data));
158:   PetscFunctionReturn(PETSC_SUCCESS);
159: }

161: static PetscErrorCode EPSReset_ScaLAPACK(EPS eps)
162: {
163:   EPS_ScaLAPACK  *ctx = (EPS_ScaLAPACK*)eps->data;

165:   PetscFunctionBegin;
166:   PetscCall(MatDestroy(&ctx->As));
167:   PetscCall(MatDestroy(&ctx->Bs));
168:   PetscFunctionReturn(PETSC_SUCCESS);
169: }

171: SLEPC_EXTERN PetscErrorCode EPSCreate_ScaLAPACK(EPS eps)
172: {
173:   EPS_ScaLAPACK  *ctx;

175:   PetscFunctionBegin;
176:   PetscCall(PetscNew(&ctx));
177:   eps->data = (void*)ctx;

179:   eps->categ = EPS_CATEGORY_OTHER;

181:   eps->ops->solve          = EPSSolve_ScaLAPACK;
182:   eps->ops->setup          = EPSSetUp_ScaLAPACK;
183:   eps->ops->setupsort      = EPSSetUpSort_Basic;
184:   eps->ops->destroy        = EPSDestroy_ScaLAPACK;
185:   eps->ops->reset          = EPSReset_ScaLAPACK;
186:   eps->ops->backtransform  = EPSBackTransform_Default;
187:   eps->ops->setdefaultst   = EPSSetDefaultST_NoFactor;
188:   PetscFunctionReturn(PETSC_SUCCESS);
189: }