Actual source code: elemental.cxx

slepc-3.23.0 2025-03-29
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: static 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:   PetscFunctionBegin;
 30:   EPSCheckHermitianDefinite(eps);
 31:   EPSCheckNotStructured(eps);
 32:   PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift));
 33:   PetscCheck(isshift,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support spectral transformations");
 34:   if (eps->nev==0) eps->nev = 1;
 35:   eps->ncv = eps->n;
 36:   if (eps->mpd!=PETSC_DETERMINE) PetscCall(PetscInfo(eps,"Warning: parameter mpd ignored\n"));
 37:   if (eps->max_it==PETSC_DETERMINE) eps->max_it = 1;
 38:   if (!eps->which) PetscCall(EPSSetWhichEigenpairs_Default(eps));
 39:   PetscCheck(eps->which!=EPS_ALL || eps->inta==eps->intb,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);
 41:   EPSCheckIgnored(eps,EPS_FEATURE_EXTRACTION | EPS_FEATURE_CONVERGENCE | EPS_FEATURE_STOPPING);
 42:   PetscCall(EPSAllocateSolution(eps,0));

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

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

 69:   PetscFunctionBegin;
 70:   El::DistMatrix<PetscReal,El::VR,El::STAR> w(*a->grid);
 71:   PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Q));
 72:   q = (Mat_Elemental*)Q->data;

 74:   if (B) {
 75:     b = (Mat_Elemental*)B->data;
 76:     El::HermitianGenDefEig(El::AXBX,El::LOWER,*a->emat,*b->emat,w,*q->emat);
 77:   } else El::HermitianEig(El::LOWER,*a->emat,w,*q->emat);

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

 89:   eps->nconv  = eps->ncv;
 90:   eps->its    = 1;
 91:   eps->reason = EPS_CONVERGED_TOL;
 92:   PetscFunctionReturn(PETSC_SUCCESS);
 93: }

 95: static PetscErrorCode EPSDestroy_Elemental(EPS eps)
 96: {
 97:   PetscFunctionBegin;
 98:   PetscCall(PetscFree(eps->data));
 99:   PetscFunctionReturn(PETSC_SUCCESS);
100: }

102: static PetscErrorCode EPSReset_Elemental(EPS eps)
103: {
104:   EPS_Elemental  *ctx = (EPS_Elemental*)eps->data;

106:   PetscFunctionBegin;
107:   PetscCall(MatDestroy(&ctx->Ae));
108:   PetscCall(MatDestroy(&ctx->Be));
109:   PetscFunctionReturn(PETSC_SUCCESS);
110: }

112: SLEPC_EXTERN PetscErrorCode EPSCreate_Elemental(EPS eps)
113: {
114:   EPS_Elemental  *ctx;

116:   PetscFunctionBegin;
117:   PetscCall(PetscNew(&ctx));
118:   eps->data = (void*)ctx;

120:   eps->categ = EPS_CATEGORY_OTHER;

122:   eps->ops->solve          = EPSSolve_Elemental;
123:   eps->ops->setup          = EPSSetUp_Elemental;
124:   eps->ops->setupsort      = EPSSetUpSort_Basic;
125:   eps->ops->destroy        = EPSDestroy_Elemental;
126:   eps->ops->reset          = EPSReset_Elemental;
127:   eps->ops->backtransform  = EPSBackTransform_Default;
128:   eps->ops->setdefaultst   = EPSSetDefaultST_NoFactor;
129:   PetscFunctionReturn(PETSC_SUCCESS);
130: }