Actual source code: elpa.c
slepc-3.15.1 2021-05-28
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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: typedef struct {
19: Mat As,Bs; /* converted matrices */
20: } EPS_ELPA;
22: PetscErrorCode EPSSetUp_ELPA(EPS eps)
23: {
25: EPS_ELPA *ctx = (EPS_ELPA*)eps->data;
26: Mat A,B;
27: PetscInt nmat;
28: PetscBool isshift;
29: PetscScalar shift;
32: EPSCheckHermitianDefinite(eps);
33: PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift);
34: if (!isshift) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support spectral transformations");
35: eps->ncv = eps->n;
36: if (eps->mpd!=PETSC_DEFAULT) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
37: if (eps->max_it==PETSC_DEFAULT) eps->max_it = 1;
38: if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
39: if (eps->which==EPS_ALL && eps->inta!=eps->intb) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support interval computation");
40: EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_STOPPING);
41: EPSCheckIgnored(eps,EPS_FEATURE_EXTRACTION | EPS_FEATURE_CONVERGENCE);
42: EPSAllocateSolution(eps,0);
44: /* convert matrices */
45: MatDestroy(&ctx->As);
46: MatDestroy(&ctx->Bs);
47: STGetNumMatrices(eps->st,&nmat);
48: STGetMatrix(eps->st,0,&A);
49: MatConvert(A,MATSCALAPACK,MAT_INITIAL_MATRIX,&ctx->As);
50: if (nmat>1) {
51: STGetMatrix(eps->st,1,&B);
52: MatConvert(B,MATSCALAPACK,MAT_INITIAL_MATRIX,&ctx->Bs);
53: }
54: STGetShift(eps->st,&shift);
55: if (shift != 0.0) {
56: if (nmat>1) {
57: MatAXPY(ctx->As,-shift,ctx->Bs,SAME_NONZERO_PATTERN);
58: } else {
59: MatShift(ctx->As,-shift);
60: }
61: }
62: return(0);
63: }
65: PetscErrorCode EPSSolve_ELPA(EPS eps)
66: {
68: EPS_ELPA *ctx = (EPS_ELPA*)eps->data;
69: Mat A = ctx->As,B = ctx->Bs,Q,V;
70: Mat_ScaLAPACK *a = (Mat_ScaLAPACK*)A->data,*b,*q;
71: PetscReal *w = eps->errest; /* used to store real eigenvalues */
72: PetscInt i;
73: elpa_t handle;
76: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Q);
77: q = (Mat_ScaLAPACK*)Q->data;
79: elpa_init(20200417); /* 20171201 */
80: handle = elpa_allocate(&ierr);
82: /* set parameters of the matrix and its MPI distribution */
83: elpa_set(handle,"na",a->N,&ierr); /* matrix size */
84: elpa_set(handle,"nev",a->N,&ierr); /* number of eigenvectors to be computed (1<=nev<=na) */
85: elpa_set(handle,"local_nrows",a->locr,&ierr); /* number of local rows of the distributed matrix on this MPI task */
86: elpa_set(handle,"local_ncols",a->locc,&ierr); /* number of local columns of the distributed matrix on this MPI task */
87: elpa_set(handle,"nblk",a->mb,&ierr); /* size of the BLACS block cyclic distribution */
88: elpa_set(handle,"mpi_comm_parent",MPI_Comm_c2f(PetscObjectComm((PetscObject)eps)),&ierr);
89: elpa_set(handle,"process_row",a->grid->myrow,&ierr); /* row coordinate of MPI process */
90: elpa_set(handle,"process_col",a->grid->mycol,&ierr); /* column coordinate of MPI process */
91: if (B) { elpa_set(handle,"blacs_context",a->grid->ictxt,&ierr); }
93: /* setup and set tunable run-time options */
94: elpa_setup(handle);
95: elpa_set(handle,"solver",ELPA_SOLVER_2STAGE,&ierr);
96: /* elpa_print_settings(handle,&ierr); */
98: /* solve the eigenvalue problem */
99: if (B) {
100: b = (Mat_ScaLAPACK*)B->data;
101: elpa_generalized_eigenvectors(handle,a->loc,b->loc,w,q->loc,0,&ierr);
102: } else {
103: elpa_eigenvectors(handle,a->loc,w,q->loc,&ierr);
104: }
106: /* cleanup */
107: elpa_deallocate(handle,&ierr);
108: elpa_uninit(&ierr);
110: for (i=0;i<eps->ncv;i++) {
111: eps->eigr[i] = eps->errest[i];
112: eps->errest[i] = PETSC_MACHINE_EPSILON;
113: }
115: BVGetMat(eps->V,&V);
116: MatConvert(Q,MATDENSE,MAT_REUSE_MATRIX,&V);
117: BVRestoreMat(eps->V,&V);
118: MatDestroy(&Q);
120: eps->nconv = eps->ncv;
121: eps->its = 1;
122: eps->reason = EPS_CONVERGED_TOL;
123: return(0);
124: }
126: PetscErrorCode EPSDestroy_ELPA(EPS eps)
127: {
131: PetscFree(eps->data);
132: return(0);
133: }
135: PetscErrorCode EPSReset_ELPA(EPS eps)
136: {
138: EPS_ELPA *ctx = (EPS_ELPA*)eps->data;
141: MatDestroy(&ctx->As);
142: MatDestroy(&ctx->Bs);
143: return(0);
144: }
146: SLEPC_EXTERN PetscErrorCode EPSCreate_ELPA(EPS eps)
147: {
148: EPS_ELPA *ctx;
152: PetscNewLog(eps,&ctx);
153: eps->data = (void*)ctx;
155: eps->categ = EPS_CATEGORY_OTHER;
157: eps->ops->solve = EPSSolve_ELPA;
158: eps->ops->setup = EPSSetUp_ELPA;
159: eps->ops->setupsort = EPSSetUpSort_Basic;
160: eps->ops->destroy = EPSDestroy_ELPA;
161: eps->ops->reset = EPSReset_ELPA;
162: eps->ops->backtransform = EPSBackTransform_Default;
163: eps->ops->setdefaultst = EPSSetDefaultST_NoFactor;
164: return(0);
165: }