1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
22: #if !defined(_PEPIMPL)
23: #define _PEPIMPL 25: #include <slepcpep.h>
26: #include <slepc-private/slepcimpl.h>
28: PETSC_EXTERN PetscLogEvent PEP_SetUp,PEP_Solve,PEP_Refine;
30: typedef struct _PEPOps *PEPOps;
32: struct _PEPOps {
33: PetscErrorCode (*solve)(PEP);
34: PetscErrorCode (*setup)(PEP);
35: PetscErrorCode (*setfromoptions)(PEP);
36: PetscErrorCode (*publishoptions)(PEP);
37: PetscErrorCode (*destroy)(PEP);
38: PetscErrorCode (*reset)(PEP);
39: PetscErrorCode (*view)(PEP,PetscViewer);
40: PetscErrorCode (*computevectors)(PEP);
41: };
43: /*
44: Maximum number of monitors you can run with a single PEP 45: */
46: #define MAXPEPMONITORS 5 48: typedef enum { PEP_STATE_INITIAL,
49: PEP_STATE_SETUP,
50: PEP_STATE_SOLVED,
51: PEP_STATE_EIGENVECTORS } PEPStateType;
53: /*
54: Defines the PEP data structure.
55: */
56: struct _p_PEP {
57: PETSCHEADER(struct _PEPOps);
58: /*------------------------- User parameters ---------------------------*/
59: PetscInt max_it; /* maximum number of iterations */
60: PetscInt nev; /* number of eigenvalues to compute */
61: PetscInt ncv; /* number of basis vectors */
62: PetscInt mpd; /* maximum dimension of projected problem */
63: PetscInt nini; /* number of initial vectors (negative means not copied yet) */
64: PetscScalar target; /* target value */
65: PetscReal tol; /* tolerance */
66: PEPConv conv; /* convergence test */
67: PEPWhich which; /* which part of the spectrum to be sought */
68: PEPBasis basis; /* polynomial basis used to represent the problem */
69: PEPProblemType problem_type; /* which kind of problem to be solved */
70: PEPScale scale; /* scaling strategy to be used */
71: PetscReal sfactor; /* scaling factor */
72: PetscInt sits; /* number of iterations of the scaling method */
73: PetscReal slambda; /* norm eigenvalue approximation for scaling */
74: PEPRefine refine; /* type of refinement to be applied after solve */
75: PetscInt npart; /* number of partitions of the communicator */
76: PetscReal rtol; /* tolerance for refinement */
77: PetscInt rits; /* number of iterations of the refinement method */
78: PetscBool schur; /* use Schur complement in refinement method */
79: PEPExtract extract; /* type of extraction used */
80: PetscBool trackall; /* whether all the residuals must be computed */
82: /*-------------- User-provided functions and contexts -----------------*/
83: PetscErrorCode (*converged)(PEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
84: PetscErrorCode (*convergeddestroy)(void*);
85: void *convergedctx;
86: PetscErrorCode (*monitor[MAXPEPMONITORS])(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
87: PetscErrorCode (*monitordestroy[MAXPEPMONITORS])(void**);
88: void *monitorcontext[MAXPEPMONITORS];
89: PetscInt numbermonitors;
91: /*----------------- Child objects and working data -------------------*/
92: ST st; /* spectral transformation object */
93: DS ds; /* direct solver object */
94: BV V; /* set of basis vectors and computed eigenvectors */
95: RG rg; /* optional region for filtering */
96: PetscRandom rand; /* random number generator */
97: SlepcSC sc; /* sorting criterion data */
98: Mat *A; /* coefficient matrices of the polynomial */
99: PetscInt nmat; /* number of matrices */
100: Vec Dl,Dr; /* diagonal matrices for balancing */
101: Vec *IS; /* references to user-provided initial space */
102: PetscScalar *eigr,*eigi; /* real and imaginary parts of eigenvalues */
103: PetscReal *errest; /* error estimates */
104: PetscInt *perm; /* permutation for eigenvalue ordering */
105: PetscReal *pbc; /* coefficients defining the polynomial basis */
106: PetscScalar *solvematcoeffs; /* coefficients to compute the matrix to be inverted */
107: PetscInt nwork; /* number of work vectors */
108: Vec *work; /* work vectors */
109: void *data; /* placeholder for solver-specific stuff */
111: /* ----------------------- Status variables --------------------------*/
112: PEPStateType state; /* initial -> setup -> solved -> eigenvectors */
113: PetscInt nconv; /* number of converged eigenvalues */
114: PetscInt its; /* number of iterations so far computed */
115: PetscInt n,nloc; /* problem dimensions (global, local) */
116: PetscReal *nrma; /* computed matrix norms */
117: PetscBool sfactor_set; /* flag to indicate the user gave sfactor */
118: PEPConvergedReason reason;
119: };
121: /*
122: Macros to test valid PEP arguments
123: */
124: #if !defined(PETSC_USE_DEBUG)
126: #define PEPCheckSolved(h,arg) do {} while (0)128: #else
130: #define PEPCheckSolved(h,arg) \131: do { \132: if (h->state<PEP_STATE_SOLVED) SETERRQ1(PetscObjectComm((PetscObject)h),PETSC_ERR_ARG_WRONGSTATE,"Must call PEPSolve() first: Parameter #%d",arg); \133: } while (0)135: #endif
137: PETSC_INTERN PetscErrorCode PEPSetDimensions_Default(PEP,PetscInt,PetscInt*,PetscInt*);
138: PETSC_INTERN PetscErrorCode PEPComputeVectors_Schur(PEP);
139: PETSC_INTERN PetscErrorCode PEPComputeVectors_Indefinite(PEP);
140: PETSC_INTERN PetscErrorCode PEPComputeResidualNorm_Private(PEP,PetscScalar,PetscScalar,Vec,Vec,PetscReal*);
141: PETSC_INTERN PetscErrorCode PEPComputeRelativeError_Private(PEP,PetscScalar,PetscScalar,Vec,Vec,PetscReal*);
142: PETSC_INTERN PetscErrorCode PEPKrylovConvergence(PEP,PetscBool,PetscInt,PetscInt,PetscReal,PetscInt*);
143: PETSC_INTERN PetscErrorCode PEPComputeScaleFactor(PEP);
144: PETSC_INTERN PetscErrorCode PEPBuildDiagonalScaling(PEP);
145: PETSC_INTERN PetscErrorCode PEPBasisCoefficients(PEP,PetscReal*);
146: PETSC_INTERN PetscErrorCode PEPEvaluateBasis(PEP,PetscScalar,PetscScalar,PetscScalar*,PetscScalar*);
147: PETSC_INTERN PetscErrorCode PEPNewtonRefinement_TOAR(PEP,PetscScalar,PetscInt*,PetscReal*,PetscInt,PetscScalar*,PetscInt,PetscInt*);
148: PETSC_INTERN PetscErrorCode PEPNewtonRefinementSimple(PEP,PetscInt*,PetscReal*,PetscInt);
150: #endif