/usr/lib/slepcdir/3.1/include/slepcqep.h is in libslepc3.1-dev 3.1-p6.dfsg-1.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 | /*
User interface for SLEPc's quadratic eigenvalue solvers.
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
SLEPc - Scalable Library for Eigenvalue Problem Computations
Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
This file is part of SLEPc.
SLEPc is free software: you can redistribute it and/or modify it under the
terms of version 3 of the GNU Lesser General Public License as published by
the Free Software Foundation.
SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
more details.
You should have received a copy of the GNU Lesser General Public License
along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
*/
#if !defined(__SLEPCQEP_H)
#define __SLEPCQEP_H
#include "slepcsys.h"
#include "slepceps.h"
PETSC_EXTERN_CXX_BEGIN
extern PetscCookie QEP_COOKIE;
/*S
QEP - Abstract SLEPc object that manages all the quadratic eigenvalue
problem solvers.
Level: beginner
.seealso: QEPCreate()
S*/
typedef struct _p_QEP* QEP;
/*E
QEPType - String with the name of a quadratic eigensolver
Level: beginner
.seealso: QEPSetType(), QEP
E*/
#define QEPType char*
#define QEPLINEAR "linear"
#define QEPQARNOLDI "qarnoldi"
/*E
QEPProblemType - determines the type of the quadratic eigenproblem
Level: intermediate
.seealso: QEPSetProblemType(), QEPGetProblemType()
E*/
typedef enum { QEP_GENERAL=1,
QEP_HERMITIAN, /* M, C, K Hermitian */
QEP_GYROSCOPIC /* M, K Hermitian, M>0, C skew-Hermitian */
} QEPProblemType;
/*E
QEPWhich - determines which part of the spectrum is requested
Level: intermediate
.seealso: QEPSetWhichEigenpairs(), QEPGetWhichEigenpairs()
E*/
typedef enum { QEP_LARGEST_MAGNITUDE=1,
QEP_SMALLEST_MAGNITUDE,
QEP_LARGEST_REAL,
QEP_SMALLEST_REAL,
QEP_LARGEST_IMAGINARY,
QEP_SMALLEST_IMAGINARY } QEPWhich;
EXTERN PetscErrorCode QEPCreate(MPI_Comm,QEP*);
EXTERN PetscErrorCode QEPDestroy(QEP);
EXTERN PetscErrorCode QEPSetType(QEP,const QEPType);
EXTERN PetscErrorCode QEPGetType(QEP,const QEPType*);
EXTERN PetscErrorCode QEPSetProblemType(QEP,QEPProblemType);
EXTERN PetscErrorCode QEPGetProblemType(QEP,QEPProblemType*);
EXTERN PetscErrorCode QEPSetOperators(QEP,Mat,Mat,Mat);
EXTERN PetscErrorCode QEPGetOperators(QEP,Mat*,Mat*,Mat*);
EXTERN PetscErrorCode QEPSetFromOptions(QEP);
EXTERN PetscErrorCode QEPSetUp(QEP);
EXTERN PetscErrorCode QEPSolve(QEP);
EXTERN PetscErrorCode QEPView(QEP,PetscViewer);
EXTERN PetscErrorCode QEPSetIP(QEP,IP);
EXTERN PetscErrorCode QEPGetIP(QEP,IP*);
EXTERN PetscErrorCode QEPSetTolerances(QEP,PetscReal,PetscInt);
EXTERN PetscErrorCode QEPGetTolerances(QEP,PetscReal*,PetscInt*);
EXTERN PetscErrorCode QEPSetConvergenceTest(QEP,PetscErrorCode (*)(QEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void*);
EXTERN PetscErrorCode QEPDefaultConverged(QEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
EXTERN PetscErrorCode QEPAbsoluteConverged(QEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
EXTERN PetscErrorCode QEPSetDimensions(QEP,PetscInt,PetscInt,PetscInt);
EXTERN PetscErrorCode QEPGetDimensions(QEP,PetscInt*,PetscInt*,PetscInt*);
EXTERN PetscErrorCode QEPSetScaleFactor(QEP,PetscReal);
EXTERN PetscErrorCode QEPGetScaleFactor(QEP,PetscReal*);
EXTERN PetscErrorCode QEPGetConverged(QEP,PetscInt*);
EXTERN PetscErrorCode QEPGetEigenpair(QEP,PetscInt,PetscScalar*,PetscScalar*,Vec,Vec);
EXTERN PetscErrorCode QEPComputeRelativeError(QEP,PetscInt,PetscReal*);
EXTERN PetscErrorCode QEPComputeResidualNorm(QEP,PetscInt,PetscReal*);
EXTERN PetscErrorCode QEPGetErrorEstimate(QEP,PetscInt,PetscReal*);
EXTERN PetscErrorCode QEPMonitorSet(QEP,PetscErrorCode (*)(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*),
void*,PetscErrorCode (*monitordestroy)(void*));
EXTERN PetscErrorCode QEPMonitorCancel(QEP);
EXTERN PetscErrorCode QEPGetMonitorContext(QEP,void **);
EXTERN PetscErrorCode QEPGetIterationNumber(QEP,PetscInt*);
EXTERN PetscErrorCode QEPGetOperationCounters(QEP,PetscInt*,PetscInt*,PetscInt*);
EXTERN PetscErrorCode QEPSetInitialSpace(QEP,PetscInt,Vec*);
EXTERN PetscErrorCode QEPSetInitialSpaceLeft(QEP,PetscInt,Vec*);
EXTERN PetscErrorCode QEPSetWhichEigenpairs(QEP,QEPWhich);
EXTERN PetscErrorCode QEPGetWhichEigenpairs(QEP,QEPWhich*);
EXTERN PetscErrorCode QEPSetLeftVectorsWanted(QEP,PetscTruth);
EXTERN PetscErrorCode QEPGetLeftVectorsWanted(QEP,PetscTruth*);
EXTERN PetscErrorCode QEPSetEigenvalueComparison(QEP,PetscErrorCode (*func)(QEP,PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*),void*);
EXTERN PetscErrorCode QEPMonitorAll(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
EXTERN PetscErrorCode QEPMonitorFirst(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
EXTERN PetscErrorCode QEPMonitorConverged(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
EXTERN PetscErrorCode QEPMonitorLG(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
EXTERN PetscErrorCode QEPMonitorLGAll(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
EXTERN PetscErrorCode QEPSetTrackAll(QEP,PetscTruth);
EXTERN PetscErrorCode QEPGetTrackAll(QEP,PetscTruth*);
EXTERN PetscErrorCode QEPSetOptionsPrefix(QEP,const char*);
EXTERN PetscErrorCode QEPAppendOptionsPrefix(QEP,const char*);
EXTERN PetscErrorCode QEPGetOptionsPrefix(QEP,const char*[]);
/*E
QEPConvergedReason - reason an eigensolver was said to
have converged or diverged
Level: beginner
.seealso: QEPSolve(), QEPGetConvergedReason(), QEPSetTolerances()
E*/
typedef enum {/* converged */
QEP_CONVERGED_TOL = 2,
/* diverged */
QEP_DIVERGED_ITS = -3,
QEP_DIVERGED_BREAKDOWN = -4,
QEP_CONVERGED_ITERATING = 0} QEPConvergedReason;
EXTERN PetscErrorCode QEPGetConvergedReason(QEP,QEPConvergedReason *);
EXTERN PetscErrorCode QEPSortEigenvalues(QEP,PetscInt,PetscScalar*,PetscScalar*,PetscInt*);
EXTERN PetscErrorCode QEPSortEigenvaluesReal(QEP,PetscInt,PetscReal*,PetscInt*);
EXTERN PetscErrorCode QEPCompareEigenvalues(QEP,PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*);
EXTERN PetscErrorCode QEPSortDenseSchur(QEP,PetscInt,PetscInt,PetscScalar*,PetscInt,PetscScalar*,PetscScalar*,PetscScalar*);
EXTERN PetscErrorCode QEPRegister(const char*,const char*,const char*,PetscErrorCode(*)(QEP));
#if defined(PETSC_USE_DYNAMIC_LIBRARIES)
#define QEPRegisterDynamic(a,b,c,d) QEPRegister(a,b,c,0)
#else
#define QEPRegisterDynamic(a,b,c,d) QEPRegister(a,b,c,d)
#endif
EXTERN PetscErrorCode QEPRegisterDestroy(void);
/* --------- options specific to particular eigensolvers -------- */
EXTERN PetscErrorCode QEPLinearSetCompanionForm(QEP,PetscInt);
EXTERN PetscErrorCode QEPLinearGetCompanionForm(QEP,PetscInt*);
EXTERN PetscErrorCode QEPLinearSetExplicitMatrix(QEP,PetscTruth);
EXTERN PetscErrorCode QEPLinearGetExplicitMatrix(QEP,PetscTruth*);
EXTERN PetscErrorCode QEPLinearSetEPS(QEP,EPS);
EXTERN PetscErrorCode QEPLinearGetEPS(QEP,EPS*);
PETSC_EXTERN_CXX_END
#endif
|