#include "slepcqep.h" PetscErrorCode QEPMonitorSet(QEP qep,PetscErrorCode (*monitor)(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*), void *mctx,PetscErrorCode (*monitordestroy)(void*))Collective on QEP
qep | - eigensolver context obtained from QEPCreate() | |
monitor | - pointer to function (if this is PETSC_NULL, it turns off monitoring) | |
mctx | - [optional] context for private data for the monitor routine (use PETSC_NULL if no context is desired) | |
monitordestroy | - [optional] routine that frees monitor context (may be PETSC_NULL) |
monitor (QEP qep, int its, int nconv, PetscScalar *eigr, PetscScalar *eigi, PetscReal* errest, int nest, void *mctx)
qep | - quadratic eigensolver context obtained from QEPCreate() | |
its | - iteration number | |
nconv | - number of converged eigenpairs | |
eigr | - real part of the eigenvalues | |
eigi | - imaginary part of the eigenvalues | |
errest | - relative error estimates for each eigenpair | |
nest | - number of error estimates | |
mctx | - optional monitoring context, as set by QEPMonitorSet() |
-qep_monitor | - print only the first error estimate | |
-qep_monitor_all | - print error estimates at each iteration | |
-qep_monitor_conv | - print the eigenvalue approximations only when convergence has been reached | |
-qep_monitor_draw | - sets line graph monitor for the first unconverged approximate eigenvalue | |
-qep_monitor_draw_all | - sets line graph monitor for all unconverged approximate eigenvalue | |
-qep_monitor_cancel | - cancels all monitors that have been hardwired into a code by calls to QEPMonitorSet(), but does not cancel those set via the options database. |
Location: src/qep/interface/qepmon.c
Index of all QEP routines
Table of Contents for all manual pages
Index of all manual pages