Actual source code: nleigs.h
slepc-3.22.2 2024-12-02
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: Private header for NLEIGS
12: */
14: #pragma once
16: #define LBPOINTS 100 /* default value of the maximum number of Leja-Bagby points */
17: #define NDPOINTS 1e4 /* number of discretization points */
19: typedef struct {
20: BV V; /* tensor vector basis for the linearization */
21: BV W; /* tensor vector basis for the linearization */
22: PetscInt nmat; /* number of interpolation points */
23: PetscScalar *s,*xi; /* Leja-Bagby points */
24: PetscScalar *beta; /* scaling factors */
25: Mat *D; /* divided difference matrices */
26: PetscScalar *coeffD; /* coefficients for divided differences in split form */
27: PetscInt nshifts; /* provided number of shifts */
28: PetscScalar *shifts; /* user-provided shifts for the Rational Krylov variant */
29: PetscInt nshiftsw; /* actual number of shifts (1 if Krylov-Schur) */
30: PetscReal ddtol; /* tolerance for divided difference convergence */
31: PetscInt ddmaxit; /* maximum number of divided difference terms */
32: PetscReal keep; /* restart parameter */
33: PetscBool lock; /* locking/non-locking variant */
34: PetscInt idxrk; /* index of next shift to use */
35: KSP *ksp; /* ksp array for storing shift factorizations */
36: Vec vrn; /* random vector with normally distributed value */
37: PetscBool fullbasis; /* use full Krylov basis instead of TOAR basis */
38: EPS eps; /* eigensolver used in the full basis variant */
39: Mat A; /* shell matrix used for the eps in full basis */
40: Vec w[6]; /* work vectors */
41: void *singularitiesctx;
42: PetscErrorCode (*computesingularities)(NEP,PetscInt*,PetscScalar*,void*);
43: } NEP_NLEIGS;
45: typedef struct {
46: PetscInt nmat,maxnmat;
47: PetscScalar *coeff;
48: Mat *A;
49: Vec t;
50: } NEP_NLEIGS_MATSHELL;
52: static inline PetscErrorCode NEPNLEIGSSetShifts(NEP nep,PetscInt *nshiftsw)
53: {
54: NEP_NLEIGS *ctx = (NEP_NLEIGS*)nep->data;
56: PetscFunctionBegin;
57: if (!ctx->nshifts) {
58: ctx->shifts = &nep->target;
59: *nshiftsw = 1;
60: } else *nshiftsw = ctx->nshifts;
61: PetscFunctionReturn(PETSC_SUCCESS);
62: }
64: SLEPC_INTERN PetscErrorCode NEPSetUp_NLEIGS_FullBasis(NEP);
65: SLEPC_INTERN PetscErrorCode NEPNLEIGSSetEPS_NLEIGS(NEP,EPS);
66: SLEPC_INTERN PetscErrorCode NEPNLEIGSGetEPS_NLEIGS(NEP,EPS*);
67: SLEPC_INTERN PetscErrorCode NEPNLEIGSBackTransform(PetscObject,PetscInt,PetscScalar*,PetscScalar *vali);
68: SLEPC_INTERN PetscErrorCode NEPNLEIGSEvalNRTFunct(NEP,PetscInt,PetscScalar,PetscScalar*);
69: SLEPC_INTERN PetscErrorCode NEPSolve_NLEIGS_FullBasis(NEP);
70: SLEPC_INTERN PetscErrorCode NEPSolve_NLEIGS(NEP);