Actual source code: qepimpl.h
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7:
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: #ifndef _QEPIMPL
23: #define _QEPIMPL
25: #include slepcqep.h
30: typedef struct _QEPOps *QEPOps;
32: struct _QEPOps {
33: PetscErrorCode (*solve)(QEP);
34: PetscErrorCode (*setup)(QEP);
35: PetscErrorCode (*setfromoptions)(QEP);
36: PetscErrorCode (*publishoptions)(QEP);
37: PetscErrorCode (*destroy)(QEP);
38: PetscErrorCode (*view)(QEP,PetscViewer);
39: };
41: /*
42: Maximum number of monitors you can run with a single QEP
43: */
44: #define MAXQEPMONITORS 5
46: /*
47: Defines the QEP data structure.
48: */
49: struct _p_QEP {
50: PETSCHEADER(struct _QEPOps);
51: /*------------------------- User parameters --------------------------*/
52: PetscInt max_it, /* maximum number of iterations */
53: nev, /* number of eigenvalues to compute */
54: ncv, /* number of basis vectors */
55: mpd, /* maximum dimension of projected problem */
56: nini, ninil; /* number of initial vectors (negative means not copied yet) */
57: PetscReal tol; /* tolerance */
58: PetscReal sfactor; /* scaling factor of the quadratic problem */
59: PetscErrorCode (*conv_func)(QEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
60: void *conv_ctx;
61: QEPWhich which; /* which part of the spectrum to be sought */
62: PetscTruth leftvecs; /* if left eigenvectors are requested */
63: PetscErrorCode (*which_func)(QEP,PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
64: void *which_ctx;
65: QEPProblemType problem_type; /* which kind of problem to be solved */
66: PetscTruth trackall; /* whether all the residuals must be computed */
68: /*------------------------- Working data --------------------------*/
69: Mat M,C,K; /* problem matrices */
70: Vec *V, /* set of basis vectors and computed eigenvectors */
71: *W, /* set of left basis vectors and computed left eigenvectors */
72: *IS, *ISL; /* placeholder for references to user-provided initial space */
73: PetscScalar *eigr, *eigi, /* real and imaginary parts of eigenvalues */
74: *T; /* matrix for projected eigenproblem */
75: PetscReal *errest; /* error estimates */
76: IP ip; /* innerproduct object */
77: void *data; /* placeholder for misc stuff associated
78: with a particular solver */
79: PetscInt nconv, /* number of converged eigenvalues */
80: its, /* number of iterations so far computed */
81: *perm, /* permutation for eigenvalue ordering */
82: matvecs, linits, /* operation counters */
83: n, nloc; /* problem dimensions (global, local) */
84: PetscRandom rand; /* random number generator */
86: /* ---------------- Default work-area and status vars -------------------- */
87: PetscInt nwork;
88: Vec *work;
90: PetscInt setupcalled;
91: QEPConvergedReason reason;
93: PetscErrorCode (*monitor[MAXQEPMONITORS])(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
94: PetscErrorCode (*monitordestroy[MAXQEPMONITORS])(void*);
95: void *monitorcontext[MAXQEPMONITORS];
96: PetscInt numbermonitors;
97: };
99: #define QEPMonitor(qep,it,nconv,eigr,eigi,errest,nest) \
100: { PetscErrorCode _ierr; PetscInt _i,_im = qep->numbermonitors; \
101: for ( _i=0; _i<_im; _i++ ) {\
102: _ierr=(*qep->monitor[_i])(qep,it,nconv,eigr,eigi,errest,nest,qep->monitorcontext[_i]);\
103: CHKERRQ(_ierr); \
104: } \
105: }
107: /* context for QEPMonitorConverged */
108: typedef struct {
109: PetscViewerASCIIMonitor viewer;
110: PetscInt oldnconv;
111: } QEPMONITOR_CONV;
112: EXTERN PetscErrorCode QEPMonitorDestroy_Converged(QEPMONITOR_CONV*);
114: EXTERN PetscErrorCode QEPRegisterAll(char *);
115: EXTERN PetscErrorCode QEPInitializePackage(char *);
116: EXTERN PetscErrorCode QEPFinalizePackage(void);
118: EXTERN PetscErrorCode QEPDestroy_Default(QEP);
119: EXTERN PetscErrorCode QEPDefaultGetWork(QEP,PetscInt);
120: EXTERN PetscErrorCode QEPDefaultFreeWork(QEP);
121: EXTERN PetscErrorCode QEPComputeVectors_Schur(QEP);
122: EXTERN PetscErrorCode QEPComputeResidualNorm_Private(QEP,PetscScalar,PetscScalar,Vec,Vec,PetscReal*);
123: EXTERN PetscErrorCode QEPComputeRelativeError_Private(QEP,PetscScalar,PetscScalar,Vec,Vec,PetscReal*);
124: EXTERN PetscErrorCode QEPKrylovConvergence(QEP,PetscInt,PetscInt,PetscScalar*,PetscInt,PetscScalar*,PetscInt,PetscReal,PetscInt*,PetscScalar*);
126: #endif