Actual source code: scalapack.c

slepc-3.14.1 2020-12-08
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 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: PetscErrorCode EPSSetUp_ScaLAPACK(EPS eps)
 22: {
 24:   EPS_ScaLAPACK  *ctx = (EPS_ScaLAPACK*)eps->data;
 25:   Mat            A,B;
 26:   PetscInt       nmat;
 27:   PetscBool      isshift;
 28:   PetscScalar    shift;

 31:   EPSCheckHermitianDefinite(eps);
 32:   PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift);
 33:   if (!isshift) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support spectral transformations");
 34:   eps->ncv = eps->n;
 35:   if (eps->mpd!=PETSC_DEFAULT) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
 36:   if (eps->max_it==PETSC_DEFAULT) eps->max_it = 1;
 37:   if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
 38:   if (eps->which==EPS_ALL && eps->inta!=eps->intb) SETERRQ(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:   EPSAllocateSolution(eps,0);

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

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

 79:   MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Q);
 80:   q = (Mat_ScaLAPACK*)Q->data;

 82:   if (B) {

 84:     b = (Mat_ScaLAPACK*)B->data;
 85:     PetscMalloc3(a->grid->nprow*a->grid->npcol,&gap,a->N,&ifail,2*a->grid->nprow*a->grid->npcol,&iclustr);
 86: #if !defined(PETSC_USE_COMPLEX)
 87:     /* allocate workspace */
 88:     PetscStackCallBLAS("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));
 90:     PetscBLASIntCast(minlwork[0],&lwork);
 91:     liwork = minliwork;
 92:     /* call computational routine */
 93:     PetscMalloc2(lwork,&work,liwork,&iwork);
 94:     PetscStackCallBLAS("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));
 96:     PetscFree2(work,iwork);
 97: #else
 98:     /* allocate workspace */
 99:     PetscStackCallBLAS("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));
101:     PetscBLASIntCast(minlwork[0],&lwork);
102:     PetscBLASIntCast(minlrwork[0],&lrwork);
103:     lrwork += a->N*a->N;
104:     liwork = minliwork;
105:     /* call computational routine */
106:     PetscMalloc3(lwork,&work,lrwork,&rwork,liwork,&iwork);
107:     PetscStackCallBLAS("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));
109:     PetscFree3(work,rwork,iwork);
110: #endif
111:     PetscFree3(gap,ifail,iclustr);

113:   } else {

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

138:   }

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

145:   BVGetMat(eps->V,&V);
146:   MatConvert(Q,MATDENSE,MAT_REUSE_MATRIX,&V);
147:   BVRestoreMat(eps->V,&V);
148:   MatDestroy(&Q);

150:   eps->nconv  = eps->ncv;
151:   eps->its    = 1;
152:   eps->reason = EPS_CONVERGED_TOL;
153:   return(0);
154: }

156: PetscErrorCode EPSDestroy_ScaLAPACK(EPS eps)
157: {

161:   PetscFree(eps->data);
162:   return(0);
163: }

165: PetscErrorCode EPSReset_ScaLAPACK(EPS eps)
166: {
168:   EPS_ScaLAPACK  *ctx = (EPS_ScaLAPACK*)eps->data;

171:   MatDestroy(&ctx->As);
172:   MatDestroy(&ctx->Bs);
173:   return(0);
174: }

176: SLEPC_EXTERN PetscErrorCode EPSCreate_ScaLAPACK(EPS eps)
177: {
178:   EPS_ScaLAPACK  *ctx;

182:   PetscNewLog(eps,&ctx);
183:   eps->data = (void*)ctx;

185:   eps->categ = EPS_CATEGORY_OTHER;

187:   eps->ops->solve          = EPSSolve_ScaLAPACK;
188:   eps->ops->setup          = EPSSetUp_ScaLAPACK;
189:   eps->ops->setupsort      = EPSSetUpSort_Basic;
190:   eps->ops->destroy        = EPSDestroy_ScaLAPACK;
191:   eps->ops->reset          = EPSReset_ScaLAPACK;
192:   eps->ops->backtransform  = EPSBackTransform_Default;
193:   eps->ops->setdefaultst   = EPSSetDefaultST_NoFactor;
194:   return(0);
195: }