1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2015, 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 PetscBool PEPRegisterAllCalled;
29: PETSC_EXTERN PetscErrorCode PEPRegisterAll(void);
30: PETSC_EXTERN PetscLogEvent PEP_SetUp,PEP_Solve,PEP_Refine;
32: typedef struct _PEPOps *PEPOps;
34: struct _PEPOps {
35: PetscErrorCode (*solve)(PEP);
36: PetscErrorCode (*setup)(PEP);
37: PetscErrorCode (*setfromoptions)(PetscOptions*,PEP);
38: PetscErrorCode (*publishoptions)(PEP);
39: PetscErrorCode (*destroy)(PEP);
40: PetscErrorCode (*reset)(PEP);
41: PetscErrorCode (*view)(PEP,PetscViewer);
42: PetscErrorCode (*backtransform)(PEP);
43: PetscErrorCode (*computevectors)(PEP);
44: PetscErrorCode (*extractvectors)(PEP);
45: };
47: /*
48: Maximum number of monitors you can run with a single PEP 49: */
50: #define MAXPEPMONITORS 5 52: typedef enum { PEP_STATE_INITIAL,
53: PEP_STATE_SETUP,
54: PEP_STATE_SOLVED,
55: PEP_STATE_EIGENVECTORS } PEPStateType;
57: /*
58: Defines the PEP data structure.
59: */
60: struct _p_PEP {
61: PETSCHEADER(struct _PEPOps);
62: /*------------------------- User parameters ---------------------------*/
63: PetscInt max_it; /* maximum number of iterations */
64: PetscInt nev; /* number of eigenvalues to compute */
65: PetscInt ncv; /* number of basis vectors */
66: PetscInt mpd; /* maximum dimension of projected problem */
67: PetscInt nini; /* number of initial vectors (negative means not copied yet) */
68: PetscScalar target; /* target value */
69: PetscReal tol; /* tolerance */
70: PEPConv conv; /* convergence test */
71: PEPWhich which; /* which part of the spectrum to be sought */
72: PEPBasis basis; /* polynomial basis used to represent the problem */
73: PEPProblemType problem_type; /* which kind of problem to be solved */
74: PEPScale scale; /* scaling strategy to be used */
75: PetscReal sfactor,dsfactor; /* scaling factors */
76: PetscInt sits; /* number of iterations of the scaling method */
77: PetscReal slambda; /* norm eigenvalue approximation for scaling */
78: PEPRefine refine; /* type of refinement to be applied after solve */
79: PetscInt npart; /* number of partitions of the communicator */
80: PetscReal rtol; /* tolerance for refinement */
81: PetscInt rits; /* number of iterations of the refinement method */
82: PetscBool schur; /* use Schur complement in refinement method */
83: PEPExtract extract; /* type of extraction used */
84: PetscBool trackall; /* whether all the residuals must be computed */
86: /*-------------- User-provided functions and contexts -----------------*/
87: PetscErrorCode (*converged)(PEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
88: PetscErrorCode (*convergeddestroy)(void*);
89: void *convergedctx;
90: PetscErrorCode (*monitor[MAXPEPMONITORS])(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
91: PetscErrorCode (*monitordestroy[MAXPEPMONITORS])(void**);
92: void *monitorcontext[MAXPEPMONITORS];
93: PetscInt numbermonitors;
95: /*----------------- Child objects and working data -------------------*/
96: ST st; /* spectral transformation object */
97: DS ds; /* direct solver object */
98: BV V; /* set of basis vectors and computed eigenvectors */
99: RG rg; /* optional region for filtering */
100: PetscRandom rand; /* random number generator */
101: SlepcSC sc; /* sorting criterion data */
102: Mat *A; /* coefficient matrices of the polynomial */
103: PetscInt nmat; /* number of matrices */
104: Vec Dl,Dr; /* diagonal matrices for balancing */
105: Vec *IS; /* references to user-provided initial space */
106: PetscScalar *eigr,*eigi; /* real and imaginary parts of eigenvalues */
107: PetscReal *errest; /* error estimates */
108: PetscInt *perm; /* permutation for eigenvalue ordering */
109: PetscReal *pbc; /* coefficients defining the polynomial basis */
110: PetscScalar *solvematcoeffs; /* coefficients to compute the matrix to be inverted */
111: PetscInt nwork; /* number of work vectors */
112: Vec *work; /* work vectors */
113: KSP refineksp; /* ksp used in refinement */
114: PetscSubcomm refinesubc; /* context for sub-communicators */
115: void *data; /* placeholder for solver-specific stuff */
117: /* ----------------------- Status variables --------------------------*/
118: PEPStateType state; /* initial -> setup -> solved -> eigenvectors */
119: PetscInt nconv; /* number of converged eigenvalues */
120: PetscInt its; /* number of iterations so far computed */
121: PetscInt n,nloc; /* problem dimensions (global, local) */
122: PetscReal *nrma; /* computed matrix norms */
123: PetscReal nrml[2]; /* computed matrix norms for the linearization */
124: PetscBool sfactor_set; /* flag to indicate the user gave sfactor */
125: PetscBool lineariz; /* current solver is based on linearization */
126: PEPConvergedReason reason;
127: };
129: /*
130: Macros to test valid PEP arguments
131: */
132: #if !defined(PETSC_USE_DEBUG)
134: #define PEPCheckSolved(h,arg) do {} while (0)136: #else
138: #define PEPCheckSolved(h,arg) \139: do { \140: if (h->state<PEP_STATE_SOLVED) SETERRQ1(PetscObjectComm((PetscObject)h),PETSC_ERR_ARG_WRONGSTATE,"Must call PEPSolve() first: Parameter #%d",arg); \141: } while (0)143: #endif
145: PETSC_INTERN PetscErrorCode PEPSetDimensions_Default(PEP,PetscInt,PetscInt*,PetscInt*);
146: PETSC_INTERN PetscErrorCode PEPExtractVectors(PEP);
147: PETSC_INTERN PetscErrorCode PEPBackTransform_Default(PEP);
148: PETSC_INTERN PetscErrorCode PEPComputeVectors(PEP);
149: PETSC_INTERN PetscErrorCode PEPComputeVectors_Default(PEP);
150: PETSC_INTERN PetscErrorCode PEPComputeVectors_Indefinite(PEP);
151: PETSC_INTERN PetscErrorCode PEPComputeResidualNorm_Private(PEP,PetscScalar,PetscScalar,Vec,Vec,Vec*,PetscReal*);
152: PETSC_INTERN PetscErrorCode PEPKrylovConvergence(PEP,PetscBool,PetscInt,PetscInt,PetscReal,PetscInt*);
153: PETSC_INTERN PetscErrorCode PEPComputeScaleFactor(PEP);
154: PETSC_INTERN PetscErrorCode PEPBuildDiagonalScaling(PEP);
155: PETSC_INTERN PetscErrorCode PEPBasisCoefficients(PEP,PetscReal*);
156: PETSC_INTERN PetscErrorCode PEPEvaluateBasis(PEP,PetscScalar,PetscScalar,PetscScalar*,PetscScalar*);
157: PETSC_INTERN PetscErrorCode PEPNewtonRefinement_TOAR(PEP,PetscScalar,PetscInt*,PetscReal*,PetscInt,PetscScalar*,PetscInt,PetscInt*);
158: PETSC_INTERN PetscErrorCode PEPNewtonRefinementSimple(PEP,PetscInt*,PetscReal*,PetscInt);
159: PETSC_INTERN PetscErrorCode PEPComputeLinearNorms(PEP);
161: #endif