Actual source code: slepcqep.h
1: /*
2: User interface for SLEPc's quadratic eigenvalue solvers.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
26: #include <slepceps.h>
28: PETSC_EXTERN PetscErrorCode QEPInitializePackage(void);
30: /*S
31: QEP - Abstract SLEPc object that manages all the quadratic eigenvalue
32: problem solvers.
34: Level: beginner
36: .seealso: QEPCreate()
37: S*/
38: typedef struct _p_QEP* QEP;
40: /*J
41: QEPType - String with the name of a quadratic eigensolver
43: Level: beginner
45: .seealso: QEPSetType(), QEP
46: J*/
47: typedef const char* QEPType;
48: #define QEPLINEAR "linear"
49: #define QEPQARNOLDI "qarnoldi"
50: #define QEPQLANCZOS "qlanczos"
52: /* Logging support */
53: PETSC_EXTERN PetscClassId QEP_CLASSID;
55: /*E
56: QEPProblemType - Determines the type of the quadratic eigenproblem
58: Level: intermediate
60: .seealso: QEPSetProblemType(), QEPGetProblemType()
61: E*/
62: typedef enum { QEP_GENERAL=1,
63: QEP_HERMITIAN, /* M, C, K Hermitian */
64: QEP_GYROSCOPIC /* M, K Hermitian, M>0, C skew-Hermitian */
65: } QEPProblemType;
67: /*E
68: QEPWhich - Determines which part of the spectrum is requested
70: Level: intermediate
72: .seealso: QEPSetWhichEigenpairs(), QEPGetWhichEigenpairs()
73: E*/
74: typedef enum { QEP_LARGEST_MAGNITUDE=1,
75: QEP_SMALLEST_MAGNITUDE,
76: QEP_LARGEST_REAL,
77: QEP_SMALLEST_REAL,
78: QEP_LARGEST_IMAGINARY,
79: QEP_SMALLEST_IMAGINARY,
80: QEP_TARGET_MAGNITUDE,
81: QEP_TARGET_REAL,
82: QEP_TARGET_IMAGINARY} QEPWhich;
84: PETSC_EXTERN PetscErrorCode QEPCreate(MPI_Comm,QEP*);
85: PETSC_EXTERN PetscErrorCode QEPDestroy(QEP*);
86: PETSC_EXTERN PetscErrorCode QEPReset(QEP);
87: PETSC_EXTERN PetscErrorCode QEPSetType(QEP,QEPType);
88: PETSC_EXTERN PetscErrorCode QEPGetType(QEP,QEPType*);
89: PETSC_EXTERN PetscErrorCode QEPSetProblemType(QEP,QEPProblemType);
90: PETSC_EXTERN PetscErrorCode QEPGetProblemType(QEP,QEPProblemType*);
91: PETSC_EXTERN PetscErrorCode QEPSetOperators(QEP,Mat,Mat,Mat);
92: PETSC_EXTERN PetscErrorCode QEPGetOperators(QEP,Mat*,Mat*,Mat*);
93: PETSC_EXTERN PetscErrorCode QEPSetTarget(QEP,PetscScalar);
94: PETSC_EXTERN PetscErrorCode QEPGetTarget(QEP,PetscScalar*);
95: PETSC_EXTERN PetscErrorCode QEPSetST(QEP,ST);
96: PETSC_EXTERN PetscErrorCode QEPGetST(QEP,ST*);
97: PETSC_EXTERN PetscErrorCode QEPSetFromOptions(QEP);
98: PETSC_EXTERN PetscErrorCode QEPSetUp(QEP);
99: PETSC_EXTERN PetscErrorCode QEPSolve(QEP);
100: PETSC_EXTERN PetscErrorCode QEPView(QEP,PetscViewer);
101: PETSC_EXTERN PetscErrorCode QEPPrintSolution(QEP,PetscViewer);
103: PETSC_EXTERN PetscErrorCode QEPSetIP(QEP,IP);
104: PETSC_EXTERN PetscErrorCode QEPGetIP(QEP,IP*);
105: PETSC_EXTERN PetscErrorCode QEPSetDS(QEP,DS);
106: PETSC_EXTERN PetscErrorCode QEPGetDS(QEP,DS*);
107: PETSC_EXTERN PetscErrorCode QEPSetTolerances(QEP,PetscReal,PetscInt);
108: PETSC_EXTERN PetscErrorCode QEPGetTolerances(QEP,PetscReal*,PetscInt*);
109: PETSC_EXTERN PetscErrorCode QEPSetConvergenceTest(QEP,PetscErrorCode (*)(QEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void*);
110: PETSC_EXTERN PetscErrorCode QEPConvergedDefault(QEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
111: PETSC_EXTERN PetscErrorCode QEPConvergedAbsolute(QEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
112: PETSC_EXTERN PetscErrorCode QEPSetDimensions(QEP,PetscInt,PetscInt,PetscInt);
113: PETSC_EXTERN PetscErrorCode QEPGetDimensions(QEP,PetscInt*,PetscInt*,PetscInt*);
114: PETSC_EXTERN PetscErrorCode QEPSetScaleFactor(QEP,PetscReal);
115: PETSC_EXTERN PetscErrorCode QEPGetScaleFactor(QEP,PetscReal*);
117: PETSC_EXTERN PetscErrorCode QEPGetConverged(QEP,PetscInt*);
118: PETSC_EXTERN PetscErrorCode QEPGetEigenpair(QEP,PetscInt,PetscScalar*,PetscScalar*,Vec,Vec);
119: PETSC_EXTERN PetscErrorCode QEPComputeRelativeError(QEP,PetscInt,PetscReal*);
120: PETSC_EXTERN PetscErrorCode QEPComputeResidualNorm(QEP,PetscInt,PetscReal*);
121: PETSC_EXTERN PetscErrorCode QEPGetErrorEstimate(QEP,PetscInt,PetscReal*);
123: PETSC_EXTERN PetscErrorCode QEPMonitor(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt);
124: PETSC_EXTERN PetscErrorCode QEPMonitorSet(QEP,PetscErrorCode (*)(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*),void*,PetscErrorCode (*)(void**));
125: PETSC_EXTERN PetscErrorCode QEPMonitorCancel(QEP);
126: PETSC_EXTERN PetscErrorCode QEPGetMonitorContext(QEP,void **);
127: PETSC_EXTERN PetscErrorCode QEPGetIterationNumber(QEP,PetscInt*);
128: PETSC_EXTERN PetscErrorCode QEPGetOperationCounters(QEP,PetscInt*,PetscInt*,PetscInt*);
130: PETSC_EXTERN PetscErrorCode QEPSetInitialSpace(QEP,PetscInt,Vec*);
131: PETSC_EXTERN PetscErrorCode QEPSetInitialSpaceLeft(QEP,PetscInt,Vec*);
132: PETSC_EXTERN PetscErrorCode QEPSetWhichEigenpairs(QEP,QEPWhich);
133: PETSC_EXTERN PetscErrorCode QEPGetWhichEigenpairs(QEP,QEPWhich*);
134: PETSC_EXTERN PetscErrorCode QEPSetLeftVectorsWanted(QEP,PetscBool);
135: PETSC_EXTERN PetscErrorCode QEPGetLeftVectorsWanted(QEP,PetscBool*);
136: PETSC_EXTERN PetscErrorCode QEPSetEigenvalueComparison(QEP,PetscErrorCode (*func)(QEP,PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*),void*);
138: PETSC_EXTERN PetscErrorCode QEPMonitorAll(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
139: PETSC_EXTERN PetscErrorCode QEPMonitorFirst(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
140: PETSC_EXTERN PetscErrorCode QEPMonitorConverged(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
141: PETSC_EXTERN PetscErrorCode QEPMonitorLG(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
142: PETSC_EXTERN PetscErrorCode QEPMonitorLGAll(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
144: PETSC_EXTERN PetscErrorCode QEPSetTrackAll(QEP,PetscBool);
145: PETSC_EXTERN PetscErrorCode QEPGetTrackAll(QEP,PetscBool*);
147: PETSC_EXTERN PetscErrorCode QEPSetOptionsPrefix(QEP,const char*);
148: PETSC_EXTERN PetscErrorCode QEPAppendOptionsPrefix(QEP,const char*);
149: PETSC_EXTERN PetscErrorCode QEPGetOptionsPrefix(QEP,const char*[]);
151: /*E
152: QEPConvergedReason - Reason an eigensolver was said to
153: have converged or diverged
155: Level: beginner
157: .seealso: QEPSolve(), QEPGetConvergedReason(), QEPSetTolerances()
158: E*/
159: typedef enum {/* converged */
160: QEP_CONVERGED_TOL = 2,
161: /* diverged */
162: QEP_DIVERGED_ITS = -3,
163: QEP_DIVERGED_BREAKDOWN = -4,
164: QEP_CONVERGED_ITERATING = 0} QEPConvergedReason;
166: PETSC_EXTERN PetscErrorCode QEPGetConvergedReason(QEP,QEPConvergedReason *);
168: PETSC_EXTERN PetscErrorCode QEPSortEigenvalues(QEP,PetscInt,PetscScalar*,PetscScalar*,PetscInt*);
169: PETSC_EXTERN PetscErrorCode QEPCompareEigenvalues(QEP,PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*);
171: PETSC_EXTERN PetscFunctionList QEPList;
172: PETSC_EXTERN PetscBool QEPRegisterAllCalled;
173: PETSC_EXTERN PetscErrorCode QEPRegisterAll(void);
174: PETSC_EXTERN PetscErrorCode QEPRegister(const char[],PetscErrorCode(*)(QEP));
176: PETSC_EXTERN PetscErrorCode QEPSetWorkVecs(QEP,PetscInt);
178: /* --------- options specific to particular eigensolvers -------- */
180: PETSC_EXTERN PetscErrorCode QEPLinearSetCompanionForm(QEP,PetscInt);
181: PETSC_EXTERN PetscErrorCode QEPLinearGetCompanionForm(QEP,PetscInt*);
182: PETSC_EXTERN PetscErrorCode QEPLinearSetExplicitMatrix(QEP,PetscBool);
183: PETSC_EXTERN PetscErrorCode QEPLinearGetExplicitMatrix(QEP,PetscBool*);
184: PETSC_EXTERN PetscErrorCode QEPLinearSetEPS(QEP,EPS);
185: PETSC_EXTERN PetscErrorCode QEPLinearGetEPS(QEP,EPS*);
187: #endif