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-2010, Universidad Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.
  9:       
 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 slepcsys.h
 27:  #include slepceps.h


 32: /*S
 33:      QEP - Abstract SLEPc object that manages all the quadratic eigenvalue 
 34:      problem solvers.

 36:    Level: beginner

 38: .seealso:  QEPCreate()
 39: S*/
 40: typedef struct _p_QEP* QEP;

 42: /*E
 43:     QEPType - String with the name of a quadratic eigensolver

 45:    Level: beginner

 47: .seealso: QEPSetType(), QEP
 48: E*/
 49: #define QEPType      char*
 50: #define QEPLINEAR    "linear"
 51: #define QEPQARNOLDI  "qarnoldi"

 53: /*E
 54:     QEPProblemType - determines the type of the quadratic eigenproblem

 56:     Level: intermediate

 58: .seealso: QEPSetProblemType(), QEPGetProblemType()
 59: E*/
 60: typedef enum { QEP_GENERAL=1,
 61:                QEP_HERMITIAN,   /* M, C, K  Hermitian */
 62:                QEP_GYROSCOPIC   /* M, K  Hermitian, M>0, C skew-Hermitian */
 63:              } QEPProblemType;

 65: /*E
 66:     QEPWhich - determines which part of the spectrum is requested

 68:     Level: intermediate

 70: .seealso: QEPSetWhichEigenpairs(), QEPGetWhichEigenpairs()
 71: E*/
 72: typedef enum { QEP_LARGEST_MAGNITUDE=1,
 73:                QEP_SMALLEST_MAGNITUDE,
 74:                QEP_LARGEST_REAL,
 75:                QEP_SMALLEST_REAL,
 76:                QEP_LARGEST_IMAGINARY,
 77:                QEP_SMALLEST_IMAGINARY } QEPWhich;

 79: EXTERN PetscErrorCode QEPCreate(MPI_Comm,QEP*);
 80: EXTERN PetscErrorCode QEPDestroy(QEP);
 81: EXTERN PetscErrorCode QEPSetType(QEP,const QEPType);
 82: EXTERN PetscErrorCode QEPGetType(QEP,const QEPType*);
 83: EXTERN PetscErrorCode QEPSetProblemType(QEP,QEPProblemType);
 84: EXTERN PetscErrorCode QEPGetProblemType(QEP,QEPProblemType*);
 85: EXTERN PetscErrorCode QEPSetOperators(QEP,Mat,Mat,Mat);
 86: EXTERN PetscErrorCode QEPGetOperators(QEP,Mat*,Mat*,Mat*);
 87: EXTERN PetscErrorCode QEPSetFromOptions(QEP);
 88: EXTERN PetscErrorCode QEPSetUp(QEP);
 89: EXTERN PetscErrorCode QEPSolve(QEP);
 90: EXTERN PetscErrorCode QEPView(QEP,PetscViewer);

 92: EXTERN PetscErrorCode QEPSetIP(QEP,IP);
 93: EXTERN PetscErrorCode QEPGetIP(QEP,IP*);
 94: EXTERN PetscErrorCode QEPSetTolerances(QEP,PetscReal,PetscInt);
 95: EXTERN PetscErrorCode QEPGetTolerances(QEP,PetscReal*,PetscInt*);
 96: EXTERN PetscErrorCode QEPSetConvergenceTest(QEP,PetscErrorCode (*)(QEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void*);
 97: EXTERN PetscErrorCode QEPDefaultConverged(QEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
 98: EXTERN PetscErrorCode QEPAbsoluteConverged(QEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
 99: EXTERN PetscErrorCode QEPSetDimensions(QEP,PetscInt,PetscInt,PetscInt);
100: EXTERN PetscErrorCode QEPGetDimensions(QEP,PetscInt*,PetscInt*,PetscInt*);
101: EXTERN PetscErrorCode QEPSetScaleFactor(QEP,PetscReal);
102: EXTERN PetscErrorCode QEPGetScaleFactor(QEP,PetscReal*);

104: EXTERN PetscErrorCode QEPGetConverged(QEP,PetscInt*);
105: EXTERN PetscErrorCode QEPGetEigenpair(QEP,PetscInt,PetscScalar*,PetscScalar*,Vec,Vec);
106: EXTERN PetscErrorCode QEPComputeRelativeError(QEP,PetscInt,PetscReal*);
107: EXTERN PetscErrorCode QEPComputeResidualNorm(QEP,PetscInt,PetscReal*);
108: EXTERN PetscErrorCode QEPGetErrorEstimate(QEP,PetscInt,PetscReal*);

110: EXTERN PetscErrorCode QEPMonitorSet(QEP,PetscErrorCode (*)(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*),
111:                                     void*,PetscErrorCode (*monitordestroy)(void*));
112: EXTERN PetscErrorCode QEPMonitorCancel(QEP);
113: EXTERN PetscErrorCode QEPGetMonitorContext(QEP,void **);
114: EXTERN PetscErrorCode QEPGetIterationNumber(QEP,PetscInt*);
115: EXTERN PetscErrorCode QEPGetOperationCounters(QEP,PetscInt*,PetscInt*,PetscInt*);

117: EXTERN PetscErrorCode QEPSetInitialSpace(QEP,PetscInt,Vec*);
118: EXTERN PetscErrorCode QEPSetInitialSpaceLeft(QEP,PetscInt,Vec*);
119: EXTERN PetscErrorCode QEPSetWhichEigenpairs(QEP,QEPWhich);
120: EXTERN PetscErrorCode QEPGetWhichEigenpairs(QEP,QEPWhich*);
121: EXTERN PetscErrorCode QEPSetLeftVectorsWanted(QEP,PetscTruth);
122: EXTERN PetscErrorCode QEPGetLeftVectorsWanted(QEP,PetscTruth*);
123: EXTERN PetscErrorCode QEPSetEigenvalueComparison(QEP,PetscErrorCode (*func)(QEP,PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*),void*);

125: EXTERN PetscErrorCode QEPMonitorAll(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
126: EXTERN PetscErrorCode QEPMonitorFirst(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
127: EXTERN PetscErrorCode QEPMonitorConverged(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
128: EXTERN PetscErrorCode QEPMonitorLG(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
129: EXTERN PetscErrorCode QEPMonitorLGAll(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);

131: EXTERN PetscErrorCode QEPSetTrackAll(QEP,PetscTruth);
132: EXTERN PetscErrorCode QEPGetTrackAll(QEP,PetscTruth*);

134: EXTERN PetscErrorCode QEPSetOptionsPrefix(QEP,const char*);
135: EXTERN PetscErrorCode QEPAppendOptionsPrefix(QEP,const char*);
136: EXTERN PetscErrorCode QEPGetOptionsPrefix(QEP,const char*[]);

138: /*E
139:     QEPConvergedReason - reason an eigensolver was said to 
140:          have converged or diverged

142:     Level: beginner

144: .seealso: QEPSolve(), QEPGetConvergedReason(), QEPSetTolerances()
145: E*/
146: typedef enum {/* converged */
147:               QEP_CONVERGED_TOL                =  2,
148:               /* diverged */
149:               QEP_DIVERGED_ITS                 = -3,
150:               QEP_DIVERGED_BREAKDOWN           = -4,
151:               QEP_CONVERGED_ITERATING          =  0} QEPConvergedReason;

153: EXTERN PetscErrorCode QEPGetConvergedReason(QEP,QEPConvergedReason *);

155: EXTERN PetscErrorCode QEPSortEigenvalues(QEP,PetscInt,PetscScalar*,PetscScalar*,PetscInt*);
156: EXTERN PetscErrorCode QEPSortEigenvaluesReal(QEP,PetscInt,PetscReal*,PetscInt*);
157: EXTERN PetscErrorCode QEPCompareEigenvalues(QEP,PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*);
158: EXTERN PetscErrorCode QEPSortDenseSchur(QEP,PetscInt,PetscInt,PetscScalar*,PetscInt,PetscScalar*,PetscScalar*,PetscScalar*);

160: EXTERN PetscErrorCode QEPRegister(const char*,const char*,const char*,PetscErrorCode(*)(QEP));
161: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
162: #define QEPRegisterDynamic(a,b,c,d) QEPRegister(a,b,c,0)
163: #else
164: #define QEPRegisterDynamic(a,b,c,d) QEPRegister(a,b,c,d)
165: #endif
166: EXTERN PetscErrorCode QEPRegisterDestroy(void);

168: /* --------- options specific to particular eigensolvers -------- */

170: EXTERN PetscErrorCode QEPLinearSetCompanionForm(QEP,PetscInt);
171: EXTERN PetscErrorCode QEPLinearGetCompanionForm(QEP,PetscInt*);
172: EXTERN PetscErrorCode QEPLinearSetExplicitMatrix(QEP,PetscTruth);
173: EXTERN PetscErrorCode QEPLinearGetExplicitMatrix(QEP,PetscTruth*);
174: EXTERN PetscErrorCode QEPLinearSetEPS(QEP,EPS);
175: EXTERN PetscErrorCode QEPLinearGetEPS(QEP,EPS*);

178: #endif