Actual source code: nleigs.h

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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:    Private header for NLEIGS
 12: */

 14: #if !defined(SLEPC_NLEIGS_H)
 15: #define SLEPC_NLEIGS_H


 18: #define  LBPOINTS  100   /* default value of the maximum number of Leja-Bagby points */
 19: #define  NDPOINTS  1e4   /* number of discretization points */

 21: typedef struct {
 22:   BV             V;         /* tensor vector basis for the linearization */
 23:   BV             W;         /* tensor vector basis for the linearization */
 24:   PetscInt       nmat;      /* number of interpolation points */
 25:   PetscScalar    *s,*xi;    /* Leja-Bagby points */
 26:   PetscScalar    *beta;     /* scaling factors */
 27:   Mat            *D;        /* divided difference matrices */
 28:   PetscScalar    *coeffD;   /* coefficients for divided differences in split form */
 29:   PetscInt       nshifts;   /* provided number of shifts */
 30:   PetscScalar    *shifts;   /* user-provided shifts for the Rational Krylov variant */
 31:   PetscInt       nshiftsw;  /* actual number of shifts (1 if Krylov-Schur) */
 32:   PetscReal      ddtol;     /* tolerance for divided difference convergence */
 33:   PetscInt       ddmaxit;   /* maximum number of divided difference terms */
 34:   PetscReal      keep;      /* restart parameter */
 35:   PetscBool      lock;      /* locking/non-locking variant */
 36:   PetscInt       idxrk;     /* index of next shift to use */
 37:   KSP            *ksp;      /* ksp array for storing shift factorizations */
 38:   Vec            vrn;       /* random vector with normally distributed value */
 39:   PetscBool      fullbasis; /* use full Krylov basis instead of TOAR basis */
 40:   EPS            eps;       /* eigensolver used in the full basis variant */
 41:   Mat            A;         /* shell matrix used for the eps in full basis */
 42:   Vec            w[6];      /* work vectors */
 43:   void           *singularitiesctx;
 44:   PetscErrorCode (*computesingularities)(NEP,PetscInt*,PetscScalar*,void*);
 45: } NEP_NLEIGS;

 47: typedef struct {
 48:   PetscInt    nmat,maxnmat;
 49:   PetscScalar *coeff;
 50:   Mat         *A;
 51:   Vec         t;
 52: } ShellMatCtx;

 54: PETSC_STATIC_INLINE PetscErrorCode NEPNLEIGSSetShifts(NEP nep,PetscInt *nshiftsw)
 55: {
 56:   NEP_NLEIGS *ctx = (NEP_NLEIGS*)nep->data;

 59:   if (!ctx->nshifts) {
 60:     ctx->shifts = &nep->target;
 61:     *nshiftsw = 1;
 62:   } else *nshiftsw = ctx->nshifts;
 63:   return(0);
 64: }

 66: SLEPC_EXTERN PetscErrorCode NEPSetUp_NLEIGS_FullBasis(NEP);
 67: SLEPC_EXTERN PetscErrorCode NEPNLEIGSSetEPS_NLEIGS(NEP,EPS);
 68: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetEPS_NLEIGS(NEP,EPS*);
 69: SLEPC_EXTERN PetscErrorCode NEPNLEIGSBackTransform(PetscObject,PetscInt,PetscScalar*,PetscScalar *vali);
 70: SLEPC_EXTERN PetscErrorCode NEPNLEIGSEvalNRTFunct(NEP,PetscInt,PetscScalar,PetscScalar*);
 71: SLEPC_EXTERN PetscErrorCode NEPSolve_NLEIGS_FullBasis(NEP);
 72: SLEPC_EXTERN PetscErrorCode NEPSolve_NLEIGS(NEP);

 74: #endif