File: SVDMonitorSet.html

package info (click to toggle)
slepc 3.4.2.dfsg-2
  • links: PTS, VCS
  • area: main
  • in suites: jessie, jessie-kfreebsd
  • size: 15,584 kB
  • ctags: 66,037
  • sloc: ansic: 51,275; makefile: 2,762; python: 1,577; fortran: 785; f90: 186; sh: 9
file content (83 lines) | stat: -rw-r--r-- 4,685 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD>
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<link rel="stylesheet" href="/slepc/slepc.css" type="text/css">
<TITLE>SVDMonitorSet</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">

<H1>SVDMonitorSet</H1>
Sets an ADDITIONAL function to be called at every iteration to monitor the error estimates for each requested singular triplet. 
<H3><FONT COLOR="#883300">Synopsis</FONT></H3>
<PRE>
#include "slepcsvd.h" 
PetscErrorCode SVDMonitorSet(SVD svd,PetscErrorCode (*monitor)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
</PRE>
Collective on <A HREF="../SVD/SVD.html#SVD">SVD</A>
<P>
<H3><FONT COLOR="#883300">Input Parameters</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>svd     </B></TD><TD> - singular value solver context obtained from <A HREF="../SVD/SVDCreate.html#SVDCreate">SVDCreate</A>()
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>monitor </B></TD><TD> - pointer to function (if this is NULL, it turns off monitoring)
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>mctx    </B></TD><TD> - [optional] context for private data for the
monitor routine (use NULL if no context is desired)
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#883300">Calling Sequence of monitor</FONT></H3>
<pre>
    monitor (<A HREF="../SVD/SVD.html#SVD">SVD</A> svd, PetscInt its, PetscInt nconv, PetscReal *sigma, PetscReal* errest, PetscInt nest, void *mctx)
</pre>
<P>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>svd    </B></TD><TD> - singular value solver context obtained from <A HREF="../SVD/SVDCreate.html#SVDCreate">SVDCreate</A>()
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>its    </B></TD><TD> - iteration number
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>nconv  </B></TD><TD> - number of converged singular triplets
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>sigma  </B></TD><TD> - singular values
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>errest </B></TD><TD> - relative error estimates for each singular triplet
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>nest   </B></TD><TD> - number of error estimates
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>mctx   </B></TD><TD> - optional monitoring context, as set by <A HREF="../SVD/SVDMonitorSet.html#SVDMonitorSet">SVDMonitorSet</A>()
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#883300">Options Database Keys</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-svd_monitor        </B></TD><TD> - print only the first error estimate
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-svd_monitor_all    </B></TD><TD> - print error estimates at each iteration
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-svd_monitor_conv   </B></TD><TD> - print the singular value approximations only when
convergence has been reached
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-svd_monitor_lg     </B></TD><TD> - sets line graph monitor for the first unconverged
approximate singular value
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-svd_monitor_lg_all </B></TD><TD> - sets line graph monitor for all unconverged
approximate singular values
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-svd_monitor_cancel </B></TD><TD> - cancels all monitors that have been hardwired into
a code by calls to <A HREF="../SVD/SVDMonitorSet.html#SVDMonitorSet">SVDMonitorSet</A>(), but does not cancel those set via
the options database.
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#883300">Notes</FONT></H3>
Several different monitoring routines may be set by calling
<A HREF="../SVD/SVDMonitorSet.html#SVDMonitorSet">SVDMonitorSet</A>() multiple times; all will be called in the
order in which they were set.
<P>

<P>
<H3><FONT COLOR="#883300">See Also</FONT></H3>
 <A HREF="../SVD/SVDMonitorFirst.html#SVDMonitorFirst">SVDMonitorFirst</A>(), <A HREF="../SVD/SVDMonitorAll.html#SVDMonitorAll">SVDMonitorAll</A>(), <A HREF="../SVD/SVDMonitorCancel.html#SVDMonitorCancel">SVDMonitorCancel</A>()
<BR><P><B><FONT COLOR="#883300">Location: </FONT></B><A HREF="../../../src/svd/interface/svdmon.c.html#SVDMonitorSet">src/svd/interface/svdmon.c</A>
<BR><A HREF="./index.html">Index of all SVD routines</A>
<BR><A HREF="../../index.html">Table of Contents for all manual pages</A>
<BR><A HREF="../singleindex.html">Index of all manual pages</A>
</BODY></HTML>