Actual source code: elemental.cxx

slepc-3.18.1 2022-11-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 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: {
 23:   EPS_Elemental  *ctx = (EPS_Elemental*)eps->data;
 24:   Mat            A,B;
 25:   PetscInt       nmat;
 26:   PetscBool      isshift;
 27:   PetscScalar    shift;

 29:   EPSCheckHermitianDefinite(eps);
 30:   PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift);
 32:   eps->ncv = eps->n;
 33:   if (eps->mpd!=PETSC_DEFAULT) PetscInfo(eps,"Warning: parameter mpd ignored\n");
 34:   if (eps->max_it==PETSC_DEFAULT) eps->max_it = 1;
 35:   if (!eps->which) EPSSetWhichEigenpairs_Default(eps);
 37:   EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_STOPPING);
 38:   EPSCheckIgnored(eps,EPS_FEATURE_EXTRACTION | EPS_FEATURE_CONVERGENCE);
 39:   EPSAllocateSolution(eps,0);

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

 59: PetscErrorCode EPSSolve_Elemental(EPS eps)
 60: {
 61:   EPS_Elemental  *ctx = (EPS_Elemental*)eps->data;
 62:   Mat            A = ctx->Ae,B = ctx->Be,Q,V;
 63:   Mat_Elemental  *a = (Mat_Elemental*)A->data,*b,*q;
 64:   PetscInt       i,rrank,ridx,erow;

 66:   El::DistMatrix<PetscReal,El::VR,El::STAR> w(*a->grid);
 67:   MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Q);
 68:   q = (Mat_Elemental*)Q->data;

 70:   if (B) {
 71:     b = (Mat_Elemental*)B->data;
 72:     El::HermitianGenDefEig(El::AXBX,El::LOWER,*a->emat,*b->emat,w,*q->emat);
 73:   } else El::HermitianEig(El::LOWER,*a->emat,w,*q->emat);

 75:   for (i=0;i<eps->ncv;i++) {
 76:     P2RO(A,0,i,&rrank,&ridx);
 77:     RO2E(A,0,rrank,ridx,&erow);
 78:     eps->eigr[i] = w.Get(erow,0);
 79:   }
 80:   BVGetMat(eps->V,&V);
 81:   MatConvert(Q,MATDENSE,MAT_REUSE_MATRIX,&V);
 82:   BVRestoreMat(eps->V,&V);
 83:   MatDestroy(&Q);

 85:   eps->nconv  = eps->ncv;
 86:   eps->its    = 1;
 87:   eps->reason = EPS_CONVERGED_TOL;
 88:   return 0;
 89: }

 91: PetscErrorCode EPSDestroy_Elemental(EPS eps)
 92: {
 93:   PetscFree(eps->data);
 94:   return 0;
 95: }

 97: PetscErrorCode EPSReset_Elemental(EPS eps)
 98: {
 99:   EPS_Elemental  *ctx = (EPS_Elemental*)eps->data;

101:   MatDestroy(&ctx->Ae);
102:   MatDestroy(&ctx->Be);
103:   return 0;
104: }

106: SLEPC_EXTERN PetscErrorCode EPSCreate_Elemental(EPS eps)
107: {
108:   EPS_Elemental  *ctx;

110:   PetscNew(&ctx);
111:   eps->data = (void*)ctx;

113:   eps->categ = EPS_CATEGORY_OTHER;

115:   eps->ops->solve          = EPSSolve_Elemental;
116:   eps->ops->setup          = EPSSetUp_Elemental;
117:   eps->ops->setupsort      = EPSSetUpSort_Basic;
118:   eps->ops->destroy        = EPSDestroy_Elemental;
119:   eps->ops->reset          = EPSReset_Elemental;
120:   eps->ops->backtransform  = EPSBackTransform_Default;
121:   eps->ops->setdefaultst   = EPSSetDefaultST_NoFactor;
122:   return 0;
123: }