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-2012, Universitat 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

 29: PETSC_EXTERN PetscErrorCode QEPInitializePackage(const char[]);

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

 35:    Level: beginner

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

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

 44:    Level: beginner

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

 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 } QEPWhich;

 81: PETSC_EXTERN PetscErrorCode QEPCreate(MPI_Comm,QEP*);
 82: PETSC_EXTERN PetscErrorCode QEPDestroy(QEP*);
 83: PETSC_EXTERN PetscErrorCode QEPReset(QEP);
 84: PETSC_EXTERN PetscErrorCode QEPSetType(QEP,const QEPType);
 85: PETSC_EXTERN PetscErrorCode QEPGetType(QEP,const QEPType*);
 86: PETSC_EXTERN PetscErrorCode QEPSetProblemType(QEP,QEPProblemType);
 87: PETSC_EXTERN PetscErrorCode QEPGetProblemType(QEP,QEPProblemType*);
 88: PETSC_EXTERN PetscErrorCode QEPSetOperators(QEP,Mat,Mat,Mat);
 89: PETSC_EXTERN PetscErrorCode QEPGetOperators(QEP,Mat*,Mat*,Mat*);
 90: PETSC_EXTERN PetscErrorCode QEPSetFromOptions(QEP);
 91: PETSC_EXTERN PetscErrorCode QEPSetUp(QEP);
 92: PETSC_EXTERN PetscErrorCode QEPSolve(QEP);
 93: PETSC_EXTERN PetscErrorCode QEPView(QEP,PetscViewer);
 94: PETSC_EXTERN PetscErrorCode QEPPrintSolution(QEP,PetscViewer);

 96: PETSC_EXTERN PetscErrorCode QEPSetIP(QEP,IP);
 97: PETSC_EXTERN PetscErrorCode QEPGetIP(QEP,IP*);
 98: PETSC_EXTERN PetscErrorCode QEPSetDS(QEP,DS);
 99: PETSC_EXTERN PetscErrorCode QEPGetDS(QEP,DS*);
100: PETSC_EXTERN PetscErrorCode QEPSetTolerances(QEP,PetscReal,PetscInt);
101: PETSC_EXTERN PetscErrorCode QEPGetTolerances(QEP,PetscReal*,PetscInt*);
102: PETSC_EXTERN PetscErrorCode QEPSetConvergenceTest(QEP,PetscErrorCode (*)(QEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void*);
103: PETSC_EXTERN PetscErrorCode QEPConvergedDefault(QEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
104: PETSC_EXTERN PetscErrorCode QEPConvergedAbsolute(QEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
105: PETSC_EXTERN PetscErrorCode QEPSetDimensions(QEP,PetscInt,PetscInt,PetscInt);
106: PETSC_EXTERN PetscErrorCode QEPGetDimensions(QEP,PetscInt*,PetscInt*,PetscInt*);
107: PETSC_EXTERN PetscErrorCode QEPSetScaleFactor(QEP,PetscReal);
108: PETSC_EXTERN PetscErrorCode QEPGetScaleFactor(QEP,PetscReal*);

110: PETSC_EXTERN PetscErrorCode QEPGetConverged(QEP,PetscInt*);
111: PETSC_EXTERN PetscErrorCode QEPGetEigenpair(QEP,PetscInt,PetscScalar*,PetscScalar*,Vec,Vec);
112: PETSC_EXTERN PetscErrorCode QEPComputeRelativeError(QEP,PetscInt,PetscReal*);
113: PETSC_EXTERN PetscErrorCode QEPComputeResidualNorm(QEP,PetscInt,PetscReal*);
114: PETSC_EXTERN PetscErrorCode QEPGetErrorEstimate(QEP,PetscInt,PetscReal*);

116: PETSC_EXTERN PetscErrorCode QEPMonitorSet(QEP,PetscErrorCode (*)(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*),void*,PetscErrorCode (*)(void**));
117: PETSC_EXTERN PetscErrorCode QEPMonitorCancel(QEP);
118: PETSC_EXTERN PetscErrorCode QEPGetMonitorContext(QEP,void **);
119: PETSC_EXTERN PetscErrorCode QEPGetIterationNumber(QEP,PetscInt*);
120: PETSC_EXTERN PetscErrorCode QEPGetOperationCounters(QEP,PetscInt*,PetscInt*,PetscInt*);

122: PETSC_EXTERN PetscErrorCode QEPSetInitialSpace(QEP,PetscInt,Vec*);
123: PETSC_EXTERN PetscErrorCode QEPSetInitialSpaceLeft(QEP,PetscInt,Vec*);
124: PETSC_EXTERN PetscErrorCode QEPSetWhichEigenpairs(QEP,QEPWhich);
125: PETSC_EXTERN PetscErrorCode QEPGetWhichEigenpairs(QEP,QEPWhich*);
126: PETSC_EXTERN PetscErrorCode QEPSetLeftVectorsWanted(QEP,PetscBool);
127: PETSC_EXTERN PetscErrorCode QEPGetLeftVectorsWanted(QEP,PetscBool*);
128: PETSC_EXTERN PetscErrorCode QEPSetEigenvalueComparison(QEP,PetscErrorCode (*func)(QEP,PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*),void*);

130: PETSC_EXTERN PetscErrorCode QEPMonitorAll(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
131: PETSC_EXTERN PetscErrorCode QEPMonitorFirst(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
132: PETSC_EXTERN PetscErrorCode QEPMonitorConverged(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
133: PETSC_EXTERN PetscErrorCode QEPMonitorLG(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
134: PETSC_EXTERN PetscErrorCode QEPMonitorLGAll(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);

136: PETSC_EXTERN PetscErrorCode QEPSetTrackAll(QEP,PetscBool);
137: PETSC_EXTERN PetscErrorCode QEPGetTrackAll(QEP,PetscBool*);

139: PETSC_EXTERN PetscErrorCode QEPSetOptionsPrefix(QEP,const char*);
140: PETSC_EXTERN PetscErrorCode QEPAppendOptionsPrefix(QEP,const char*);
141: PETSC_EXTERN PetscErrorCode QEPGetOptionsPrefix(QEP,const char*[]);

143: /*E
144:     QEPConvergedReason - reason an eigensolver was said to 
145:          have converged or diverged

147:     Level: beginner

149: .seealso: QEPSolve(), QEPGetConvergedReason(), QEPSetTolerances()
150: E*/
151: typedef enum {/* converged */
152:               QEP_CONVERGED_TOL                =  2,
153:               /* diverged */
154:               QEP_DIVERGED_ITS                 = -3,
155:               QEP_DIVERGED_BREAKDOWN           = -4,
156:               QEP_CONVERGED_ITERATING          =  0} QEPConvergedReason;

158: PETSC_EXTERN PetscErrorCode QEPGetConvergedReason(QEP,QEPConvergedReason *);

160: PETSC_EXTERN PetscErrorCode QEPSortEigenvalues(QEP,PetscInt,PetscScalar*,PetscScalar*,PetscInt*);
161: PETSC_EXTERN PetscErrorCode QEPCompareEigenvalues(QEP,PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*);

163: PETSC_EXTERN PetscFList QEPList;
164: PETSC_EXTERN PetscBool  QEPRegisterAllCalled;
165: PETSC_EXTERN PetscErrorCode QEPRegisterAll(const char[]);
166: PETSC_EXTERN PetscErrorCode QEPRegisterDestroy(void);
167: PETSC_EXTERN PetscErrorCode QEPRegister(const char[],const char[],const char[],PetscErrorCode(*)(QEP));

169: /*MC
170:    QEPRegisterDynamic - Adds a method to the quadratic eigenproblem solver package.

172:    Synopsis:
173:    PetscErrorCode QEPRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(QEP))

175:    Not Collective

177:    Input Parameters:
178: +  name_solver - name of a new user-defined solver
179: .  path - path (either absolute or relative) the library containing this solver
180: .  name_create - name of routine to create the solver context
181: -  routine_create - routine to create the solver context

183:    Notes:
184:    QEPRegisterDynamic() may be called multiple times to add several user-defined solvers.

186:    If dynamic libraries are used, then the fourth input argument (routine_create)
187:    is ignored.

189:    Sample usage:
190: .vb
191:    QEPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
192:                "MySolverCreate",MySolverCreate);
193: .ve

195:    Then, your solver can be chosen with the procedural interface via
196: $     QEPSetType(qep,"my_solver")
197:    or at runtime via the option
198: $     -qep_type my_solver

200:    Level: advanced

202: .seealso: QEPRegisterDestroy(), QEPRegisterAll()
203: M*/
204: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
205: #define QEPRegisterDynamic(a,b,c,d) QEPRegister(a,b,c,0)
206: #else
207: #define QEPRegisterDynamic(a,b,c,d) QEPRegister(a,b,c,d)
208: #endif

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

212: PETSC_EXTERN PetscErrorCode QEPLinearSetCompanionForm(QEP,PetscInt);
213: PETSC_EXTERN PetscErrorCode QEPLinearGetCompanionForm(QEP,PetscInt*);
214: PETSC_EXTERN PetscErrorCode QEPLinearSetExplicitMatrix(QEP,PetscBool);
215: PETSC_EXTERN PetscErrorCode QEPLinearGetExplicitMatrix(QEP,PetscBool*);
216: PETSC_EXTERN PetscErrorCode QEPLinearSetEPS(QEP,EPS);
217: PETSC_EXTERN PetscErrorCode QEPLinearGetEPS(QEP,EPS*);

219: #endif