SVDMonitorFirstDrawLG#
Plots the error estimate of the first unconverged approximation at each iteration of the singular value solver.
Synopsis#
#include "slepcsvd.h"
PetscErrorCode SVDMonitorFirstDrawLG(SVD svd,PetscInt its,PetscInt nconv,PetscReal sigma[],PetscReal errest[],PetscInt nest,PetscViewerAndFormat *vf)
Collective
Input Parameters#
svd - the singular value solver context
its - iteration number
nconv - number of converged singular triplets so far
sigma - singular values
errest - error estimates
nest - number of error estimates to display
vf - viewer and format for monitoring
Options Database Key#
-svd_monitor draw::draw_lg - activates
SVDMonitorFirstDrawLG()
Notes#
This is not called directly by users, rather one calls SVDMonitorSet(), with this
function as an argument, to cause the monitor to be used during the SVD solve.
Call SVDMonitorFirstDrawLGCreate() to create the context used with this monitor.
See Also#
SVD: Singular Value Decomposition, SVDMonitorSet(), SVDMonitorFirstDrawLGCreate()
Level#
intermediate
Location#
Index of all SVD routines Table of Contents for all manual pages Index of all manual pages