Actual source code: elpa.c

slepc-3.17.0 2022-03-31
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 ELPA.
 12: */

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

 18: #define PetscCallELPA(func, ...) do {                                                   \
 19:     PetscErrorCode elpa_ierr_; \
 20:     PetscStackPush(PetscStringize(func));                                   \
 21:     func(__VA_ARGS__,&elpa_ierr_);                                              \
 22:     PetscStackPop;                                                                             \
 24:   } while (0)

 26: #define PetscCallELPARET(func, ...) do {                                                   \
 27:     PetscStackPush(PetscStringize(func));                                                      \
 28:     PetscErrorCode elpa_ierr_ = func(__VA_ARGS__);                                              \
 29:     PetscStackPop;                                                                             \
 31:   } while (0)

 33: #define PetscCallELPANOARG(func) do {                                                   \
 34:     PetscErrorCode elpa_ierr_; \
 35:     PetscStackPush(PetscStringize(func));                                   \
 36:     func(&elpa_ierr_);                                              \
 37:     PetscStackPop;                                                                             \
 39:   } while (0)

 41: typedef struct {
 42:   Mat As,Bs;        /* converted matrices */
 43: } EPS_ELPA;

 45: PetscErrorCode EPSSetUp_ELPA(EPS eps)
 46: {
 47:   EPS_ELPA       *ctx = (EPS_ELPA*)eps->data;
 48:   Mat            A,B;
 49:   PetscInt       nmat;
 50:   PetscBool      isshift;
 51:   PetscScalar    shift;

 53:   EPSCheckHermitianDefinite(eps);
 54:   PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift);
 56:   eps->ncv = eps->n;
 57:   if (eps->mpd!=PETSC_DEFAULT) PetscInfo(eps,"Warning: parameter mpd ignored\n");
 58:   if (eps->max_it==PETSC_DEFAULT) eps->max_it = 1;
 59:   if (!eps->which) EPSSetWhichEigenpairs_Default(eps);
 61:   EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_STOPPING);
 62:   EPSCheckIgnored(eps,EPS_FEATURE_EXTRACTION | EPS_FEATURE_CONVERGENCE);
 63:   EPSAllocateSolution(eps,0);

 65:   /* convert matrices */
 66:   MatDestroy(&ctx->As);
 67:   MatDestroy(&ctx->Bs);
 68:   STGetNumMatrices(eps->st,&nmat);
 69:   STGetMatrix(eps->st,0,&A);
 70:   MatConvert(A,MATSCALAPACK,MAT_INITIAL_MATRIX,&ctx->As);
 71:   if (nmat>1) {
 72:     STGetMatrix(eps->st,1,&B);
 73:     MatConvert(B,MATSCALAPACK,MAT_INITIAL_MATRIX,&ctx->Bs);
 74:   }
 75:   STGetShift(eps->st,&shift);
 76:   if (shift != 0.0) {
 77:     if (nmat>1) MatAXPY(ctx->As,-shift,ctx->Bs,SAME_NONZERO_PATTERN);
 78:     else MatShift(ctx->As,-shift);
 79:   }
 80:   PetscFunctionReturn(0);
 81: }

 83: PetscErrorCode EPSSolve_ELPA(EPS eps)
 84: {
 85:   EPS_ELPA       *ctx = (EPS_ELPA*)eps->data;
 86:   Mat            A = ctx->As,B = ctx->Bs,Q,V;
 87:   Mat_ScaLAPACK  *a = (Mat_ScaLAPACK*)A->data,*b,*q;
 88:   PetscReal      *w = eps->errest;  /* used to store real eigenvalues */
 89:   PetscInt       i;
 90:   elpa_t         handle;

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

 95:   elpa_init,20200417;    /* 20171201 */
 96:   handle = elpa_allocate;

 98:   /* set parameters of the matrix and its MPI distribution */
 99:   elpa_set,handle,"na",a->N;                         /* matrix size */
100:   elpa_set,handle,"nev",a->N;                        /* number of eigenvectors to be computed (1<=nev<=na) */
101:   elpa_set,handle,"local_nrows",a->locr;             /* number of local rows of the distributed matrix on this MPI task  */
102:   elpa_set,handle,"local_ncols",a->locc;             /* number of local columns of the distributed matrix on this MPI task */
103:   elpa_set,handle,"nblk",a->mb;                      /* size of the BLACS block cyclic distribution */
104:   elpa_set,handle,"mpi_comm_parent",MPI_Comm_c2f(PetscObjectComm((PetscObject)eps));
105:   elpa_set,handle,"process_row",a->grid->myrow;      /* row coordinate of MPI process */
106:   elpa_set,handle,"process_col",a->grid->mycol;      /* column coordinate of MPI process */
107:   if (B) elpa_set,handle,"blacs_context",a->grid->ictxt;

109:   /* setup and set tunable run-time options */
110:   elpa_setup,handle;
111:   elpa_set,handle,"solver",ELPA_SOLVER_2STAGE;
112:   /* elpa_print_settings,handle; */

114:   /* solve the eigenvalue problem */
115:   if (B) {
116:     b = (Mat_ScaLAPACK*)B->data;
117:     elpa_generalized_eigenvectors,handle,a->loc,b->loc,w,q->loc,0;
118:   } else elpa_eigenvectors,handle,a->loc,w,q->loc;

120:   /* cleanup */
121:   elpa_deallocate,handle;
122:   elpa_uninit;

124:   for (i=0;i<eps->ncv;i++) {
125:     eps->eigr[i]   = eps->errest[i];
126:     eps->errest[i] = PETSC_MACHINE_EPSILON;
127:   }

129:   BVGetMat(eps->V,&V);
130:   MatConvert(Q,MATDENSE,MAT_REUSE_MATRIX,&V);
131:   BVRestoreMat(eps->V,&V);
132:   MatDestroy(&Q);

134:   eps->nconv  = eps->ncv;
135:   eps->its    = 1;
136:   eps->reason = EPS_CONVERGED_TOL;
137:   PetscFunctionReturn(0);
138: }

140: PetscErrorCode EPSDestroy_ELPA(EPS eps)
141: {
142:   PetscFree(eps->data);
143:   PetscFunctionReturn(0);
144: }

146: PetscErrorCode EPSReset_ELPA(EPS eps)
147: {
148:   EPS_ELPA       *ctx = (EPS_ELPA*)eps->data;

150:   MatDestroy(&ctx->As);
151:   MatDestroy(&ctx->Bs);
152:   PetscFunctionReturn(0);
153: }

155: SLEPC_EXTERN PetscErrorCode EPSCreate_ELPA(EPS eps)
156: {
157:   EPS_ELPA       *ctx;

159:   PetscNewLog(eps,&ctx);
160:   eps->data = (void*)ctx;

162:   eps->categ = EPS_CATEGORY_OTHER;

164:   eps->ops->solve          = EPSSolve_ELPA;
165:   eps->ops->setup          = EPSSetUp_ELPA;
166:   eps->ops->setupsort      = EPSSetUpSort_Basic;
167:   eps->ops->destroy        = EPSDestroy_ELPA;
168:   eps->ops->reset          = EPSReset_ELPA;
169:   eps->ops->backtransform  = EPSBackTransform_Default;
170:   eps->ops->setdefaultst   = EPSSetDefaultST_NoFactor;
171:   PetscFunctionReturn(0);
172: }