Actual source code: pepkrylov.h
slepc-3.10.0 2018-09-18
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2018, 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 TOAR and STOAR
12: */
17: PETSC_INTERN PetscErrorCode PEPExtractVectors_TOAR(PEP);
18: PETSC_INTERN PetscErrorCode PEPSTOARSetUpInnerMatrix(PEP,Mat*);
19: PETSC_INTERN PetscErrorCode PEPSolve_STOAR(PEP);
20: PETSC_INTERN PetscErrorCode PEPSolve_STOAR_QSlice(PEP);
21: PETSC_INTERN PetscErrorCode PEPSetUp_STOAR_QSlice(PEP);
22: PETSC_INTERN PetscErrorCode PEPReset_STOAR_QSlice(PEP);
24: /* Structure characterizing a shift in spectrum slicing */
25: typedef struct _n_shift *PEP_shift;
26: struct _n_shift {
27: PetscReal value;
28: PetscInt inertia;
29: PetscBool comp[2]; /* Shows completion of subintervals (left and right) */
30: PEP_shift neighb[2]; /* Adjacent shifts */
31: PetscInt index; /* Index in eig where found values are stored */
32: PetscInt neigs; /* Number of values found */
33: PetscReal ext[2]; /* Limits for accepted values */
34: PetscInt nsch[2]; /* Number of missing values for each subinterval */
35: PetscInt nconv[2]; /* Converged on each side (accepted or not) */
36: };
38: /* Identifies the TOAR vectors for each Pseudo-Lanczos vector in the global array */
39: typedef struct {
40: PetscInt nq;
41: PetscInt *q;
42: } PEP_QInfo;
44: /* Structure for storing the state of spectrum slicing */
45: struct _n_SR {
46: PetscReal int0,int1; /* Extremes of the interval */
47: PetscInt dir; /* Determines the order of values in eig (+1 incr, -1 decr) */
48: PetscBool hasEnd; /* Tells whether the interval has an end */
49: PetscBool dirch; /* Tells if dir has been changed */
50: PetscInt inertia0,inertia1;
51: PetscScalar *back;
52: PetscInt numEigs; /* Number of eigenvalues in the interval */
53: PetscInt indexEig;
54: PEP_shift sPres; /* Present shift */
55: PEP_shift *pending; /* Pending shifts array */
56: PetscInt nPend; /* Number of pending shifts */
57: PetscInt maxPend; /* Size of "pending" array */
58: PetscInt *idxDef0,*idxDef1; /* For deflation */
59: PetscInt ndef0,ndef1; /* Index in deflation arrays */
60: PetscInt nMAXCompl;
61: PetscInt iterCompl;
62: PetscInt itsKs; /* Krylovschur restarts */
63: PetscInt nleap;
64: PEP_shift s0; /* Initial shift */
65: PEP_shift sPrev;
66: PetscInt nv; /* position of restart vector */
67: BV V; /* full TOAR basis */
68: PetscScalar *S; /* TOAR coefficients */
69: PetscInt ld; /* Leading dimension for each block of S */
70: BV Vnext; /* temporary working basis during change of shift */
71: PetscScalar *eigr,*eigi; /* eigenvalues */
72: PetscReal *errest; /* error estimates */
73: PetscInt *perm; /* permutation */
74: PEP_QInfo *qinfo; /* TOAR vectors for each pseudo-Lanczos vector */
75: PetscInt intcorr; /* Global inertia correction */
76: PetscInt symmlost; /* Counter for symmetry lost */
77: Vec v[3];
78: EPS eps;
79: };
80: typedef struct _n_SR *PEP_SR;
82: typedef struct {
83: PetscReal keep; /* restart parameter */
84: PetscBool lock; /* locking/non-locking variant */
85: BV V; /* tensor basis vectors object for the linearization */
86: PEP_SR sr; /* spectrum slicing context */
87: PetscReal *shifts; /* array containing global shifts */
88: PetscInt *inertias; /* array containing global inertias */
89: PetscInt nshifts; /* elements in the arrays of shifts and inertias */
90: PetscInt nev; /* number of eigenvalues to compute */
91: PetscInt ncv; /* number of basis vectors */
92: PetscInt mpd; /* maximum dimension of projected problem */
93: PetscBool detect; /* check for zeros during factorizations */
94: PetscBool hyperbolic; /* hyperbolic problem flag */
95: PetscReal alpha,beta; /* coefficients defining the linearization */
96: } PEP_TOAR;
98: typedef struct {
99: PetscReal scal[2];
100: Mat A[2];
101: Vec t;
102: } ShellMatCtx;
103: #endif