Actual source code: qepimpl.h
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2011, Universitat 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>
27: extern PetscFList QEPList;
28: extern PetscLogEvent QEP_SetUp, QEP_Solve, QEP_Dense;
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 (*reset)(QEP);
39: PetscErrorCode (*view)(QEP,PetscViewer);
40: };
42: /*
43: Maximum number of monitors you can run with a single QEP
44: */
45: #define MAXQEPMONITORS 5
47: /*
48: Defines the QEP data structure.
49: */
50: struct _p_QEP {
51: PETSCHEADER(struct _QEPOps);
52: /*------------------------- User parameters --------------------------*/
53: PetscInt max_it, /* maximum number of iterations */
54: nev, /* number of eigenvalues to compute */
55: ncv, /* number of basis vectors */
56: mpd, /* maximum dimension of projected problem */
57: nini, ninil, /* number of initial vectors (negative means not copied yet) */
58: allocated_ncv; /* number of basis vectors allocated */
59: PetscReal tol; /* tolerance */
60: PetscReal sfactor; /* scaling factor of the quadratic problem */
61: PetscErrorCode (*conv_func)(QEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
62: void *conv_ctx;
63: QEPWhich which; /* which part of the spectrum to be sought */
64: PetscBool leftvecs; /* if left eigenvectors are requested */
65: PetscErrorCode (*which_func)(QEP,PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
66: void *which_ctx;
67: QEPProblemType problem_type; /* which kind of problem to be solved */
68: PetscBool trackall; /* whether all the residuals must be computed */
70: /*------------------------- Working data --------------------------*/
71: Mat M,C,K; /* problem matrices */
72: Vec *V, /* set of basis vectors and computed eigenvectors */
73: *W, /* set of left basis vectors and computed left eigenvectors */
74: *IS, *ISL; /* placeholder for references to user-provided initial space */
75: PetscScalar *eigr, *eigi, /* real and imaginary parts of eigenvalues */
76: *T; /* matrix for projected eigenproblem */
77: PetscReal *errest; /* error estimates */
78: IP ip; /* innerproduct object */
79: void *data; /* placeholder for misc stuff associated
80: with a particular solver */
81: PetscInt nconv, /* number of converged eigenvalues */
82: its, /* number of iterations so far computed */
83: *perm, /* permutation for eigenvalue ordering */
84: matvecs, linits, /* operation counters */
85: n, nloc; /* problem dimensions (global, local) */
86: PetscRandom rand; /* random number generator */
87: Vec t; /* template vector */
89: /* ---------------- Default work-area and status vars -------------------- */
90: PetscInt nwork;
91: Vec *work;
93: PetscInt setupcalled;
94: QEPConvergedReason reason;
96: PetscErrorCode (*monitor[MAXQEPMONITORS])(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
97: PetscErrorCode (*monitordestroy[MAXQEPMONITORS])(void**);
98: void *monitorcontext[MAXQEPMONITORS];
99: PetscInt numbermonitors;
100: };
102: extern PetscErrorCode QEPMonitor(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt);
104: extern PetscErrorCode QEPDefaultGetWork(QEP,PetscInt);
105: extern PetscErrorCode QEPDefaultFreeWork(QEP);
106: extern PetscErrorCode QEPAllocateSolution(QEP);
107: extern PetscErrorCode QEPFreeSolution(QEP);
108: extern PetscErrorCode QEPComputeVectors_Schur(QEP);
109: extern PetscErrorCode QEPComputeResidualNorm_Private(QEP,PetscScalar,PetscScalar,Vec,Vec,PetscReal*);
110: extern PetscErrorCode QEPComputeRelativeError_Private(QEP,PetscScalar,PetscScalar,Vec,Vec,PetscReal*);
111: extern PetscErrorCode QEPKrylovConvergence(QEP,PetscInt,PetscInt,PetscScalar*,PetscInt,PetscScalar*,PetscInt,PetscReal,PetscInt*,PetscScalar*);
113: #endif