#include "slepceps.h" PetscErrorCode EPSMonitorSet(EPS eps,PetscErrorCode (*monitor)(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *mctx),void *mctx,PetscErrorCode (*monitordestroy)(void**))Logically Collective
eps | - eigensolver context obtained from EPSCreate() | |
monitor | - pointer to function (if this is NULL, it turns off monitoring) | |
mctx | - [optional] context for private data for the monitor routine (use NULL if no context is desired) | |
monitordestroy | - [optional] routine that frees monitor context (may be NULL) |
PetscErrorCode monitor(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *mctx)
eps | - eigensolver context obtained from EPSCreate() | |
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 EPSMonitorSet() |
-eps_monitor | - print only the first error estimate | |
-eps_monitor_all | - print error estimates at each iteration | |
-eps_monitor_conv | - print the eigenvalue approximations only when convergence has been reached | |
-eps_monitor draw::draw_lg | - sets line graph monitor for the first unconverged approximate eigenvalue | |
-eps_monitor_all draw::draw_lg | - sets line graph monitor for all unconverged approximate eigenvalues | |
-eps_monitor_conv draw::draw_lg | - sets line graph monitor for convergence history | |
-eps_monitor_cancel | - cancels all monitors that have been hardwired into a code by calls to EPSMonitorSet(), but does not cancel those set via the options database. |