Actual source code: elemental.cxx
slepc-3.16.1 2021-11-17
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 Elemental.
12: */
14: #include <slepc/private/epsimpl.h>
15: #include <petsc/private/petscelemental.h>
17: typedef struct {
18: Mat Ae,Be; /* converted matrices */
19: } EPS_Elemental;
21: PetscErrorCode EPSSetUp_Elemental(EPS eps)
22: {
24: EPS_Elemental *ctx = (EPS_Elemental*)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->Ae);
45: MatDestroy(&ctx->Be);
46: STGetNumMatrices(eps->st,&nmat);
47: STGetMatrix(eps->st,0,&A);
48: MatConvert(A,MATELEMENTAL,MAT_INITIAL_MATRIX,&ctx->Ae);
49: if (nmat>1) {
50: STGetMatrix(eps->st,1,&B);
51: MatConvert(B,MATELEMENTAL,MAT_INITIAL_MATRIX,&ctx->Be);
52: }
53: STGetShift(eps->st,&shift);
54: if (shift != 0.0) {
55: if (nmat>1) {
56: MatAXPY(ctx->Ae,-shift,ctx->Be,SAME_NONZERO_PATTERN);
57: } else {
58: MatShift(ctx->Ae,-shift);
59: }
60: }
61: return(0);
62: }
64: PetscErrorCode EPSSolve_Elemental(EPS eps)
65: {
67: EPS_Elemental *ctx = (EPS_Elemental*)eps->data;
68: Mat A = ctx->Ae,B = ctx->Be,Q,V;
69: Mat_Elemental *a = (Mat_Elemental*)A->data,*b,*q;
70: PetscInt i,rrank,ridx,erow;
73: El::DistMatrix<PetscReal,El::VR,El::STAR> w(*a->grid);
74: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Q);
75: q = (Mat_Elemental*)Q->data;
77: if (B) {
78: b = (Mat_Elemental*)B->data;
79: El::HermitianGenDefEig(El::AXBX,El::LOWER,*a->emat,*b->emat,w,*q->emat);
80: } else El::HermitianEig(El::LOWER,*a->emat,w,*q->emat);
82: for (i=0;i<eps->ncv;i++) {
83: P2RO(A,0,i,&rrank,&ridx);
84: RO2E(A,0,rrank,ridx,&erow);
85: eps->eigr[i] = w.Get(erow,0);
86: }
87: BVGetMat(eps->V,&V);
88: MatConvert(Q,MATDENSE,MAT_REUSE_MATRIX,&V);
89: BVRestoreMat(eps->V,&V);
90: MatDestroy(&Q);
92: eps->nconv = eps->ncv;
93: eps->its = 1;
94: eps->reason = EPS_CONVERGED_TOL;
95: return(0);
96: }
98: PetscErrorCode EPSDestroy_Elemental(EPS eps)
99: {
103: PetscFree(eps->data);
104: return(0);
105: }
107: PetscErrorCode EPSReset_Elemental(EPS eps)
108: {
110: EPS_Elemental *ctx = (EPS_Elemental*)eps->data;
113: MatDestroy(&ctx->Ae);
114: MatDestroy(&ctx->Be);
115: return(0);
116: }
118: SLEPC_EXTERN PetscErrorCode EPSCreate_Elemental(EPS eps)
119: {
120: EPS_Elemental *ctx;
124: PetscNewLog(eps,&ctx);
125: eps->data = (void*)ctx;
127: eps->categ = EPS_CATEGORY_OTHER;
129: eps->ops->solve = EPSSolve_Elemental;
130: eps->ops->setup = EPSSetUp_Elemental;
131: eps->ops->setupsort = EPSSetUpSort_Basic;
132: eps->ops->destroy = EPSDestroy_Elemental;
133: eps->ops->reset = EPSReset_Elemental;
134: eps->ops->backtransform = EPSBackTransform_Default;
135: eps->ops->setdefaultst = EPSSetDefaultST_NoFactor;
136: return(0);
137: }