Line | Branch | Exec | Source |
---|---|---|---|
1 | /* | ||
2 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
3 | SLEPc - Scalable Library for Eigenvalue Problem Computations | ||
4 | Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain | ||
5 | |||
6 | This file is part of SLEPc. | ||
7 | SLEPc is distributed under a 2-clause BSD license (see LICENSE). | ||
8 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
9 | */ | ||
10 | /* | ||
11 | SVD routines related to monitors | ||
12 | */ | ||
13 | |||
14 | #include <slepc/private/svdimpl.h> /*I "slepcsvd.h" I*/ | ||
15 | #include <petscdraw.h> | ||
16 | |||
17 | /* | ||
18 | Runs the user provided monitor routines, if any. | ||
19 | */ | ||
20 | 49717 | PetscErrorCode SVDMonitor(SVD svd,PetscInt it,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest) | |
21 | { | ||
22 | 49717 | PetscInt i,n = svd->numbermonitors; | |
23 | |||
24 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
49717 | PetscFunctionBegin; |
25 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
49937 | for (i=0;i<n;i++) PetscCall((*svd->monitor[i])(svd,it,nconv,sigma,errest,nest,svd->monitorcontext[i])); |
26 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
9245 | PetscFunctionReturn(PETSC_SUCCESS); |
27 | } | ||
28 | |||
29 | /*@C | ||
30 | SVDMonitorSet - Sets an ADDITIONAL function to be called at every | ||
31 | iteration to monitor the error estimates for each requested singular triplet. | ||
32 | |||
33 | Logically Collective | ||
34 | |||
35 | Input Parameters: | ||
36 | + svd - singular value solver context obtained from SVDCreate() | ||
37 | . monitor - pointer to function (if this is NULL, it turns off monitoring), see SVDMonitorFn | ||
38 | . mctx - [optional] context for private data for the | ||
39 | monitor routine (use NULL if no context is desired) | ||
40 | - monitordestroy - [optional] routine that frees monitor context (may be NULL), | ||
41 | see PetscCtxDestroyFn for the calling sequence | ||
42 | |||
43 | Options Database Keys: | ||
44 | + -svd_monitor - print only the first error estimate | ||
45 | . -svd_monitor_all - print error estimates at each iteration | ||
46 | . -svd_monitor_conv - print the singular value approximations only when | ||
47 | convergence has been reached | ||
48 | . -svd_monitor_conditioning - print the condition number when available | ||
49 | . -svd_monitor draw::draw_lg - sets line graph monitor for the first unconverged | ||
50 | approximate singular value | ||
51 | . -svd_monitor_all draw::draw_lg - sets line graph monitor for all unconverged | ||
52 | approximate singular values | ||
53 | . -svd_monitor_conv draw::draw_lg - sets line graph monitor for convergence history | ||
54 | - -svd_monitor_cancel - cancels all monitors that have been hardwired into | ||
55 | a code by calls to SVDMonitorSet(), but does not cancel those set via | ||
56 | the options database. | ||
57 | |||
58 | Notes: | ||
59 | The options database option -svd_monitor and related options are the easiest way | ||
60 | to turn on SVD iteration monitoring. | ||
61 | |||
62 | SVDMonitorRegister() provides a way to associate an options database key with SVD | ||
63 | monitor function. | ||
64 | |||
65 | Several different monitoring routines may be set by calling SVDMonitorSet() multiple | ||
66 | times; all will be called in the order in which they were set. | ||
67 | |||
68 | Level: intermediate | ||
69 | |||
70 | .seealso: SVDMonitorFirst(), SVDMonitorAll(), SVDMonitorConditioning(), SVDMonitorCancel() | ||
71 | @*/ | ||
72 | 275 | PetscErrorCode SVDMonitorSet(SVD svd,SVDMonitorFn *monitor,void *mctx,PetscCtxDestroyFn *monitordestroy) | |
73 | { | ||
74 | 275 | PetscInt i; | |
75 | 275 | PetscBool identical; | |
76 | |||
77 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
275 | PetscFunctionBegin; |
78 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
275 | PetscValidHeaderSpecific(svd,SVD_CLASSID,1); |
79 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
333 | for (i=0;i<svd->numbermonitors;i++) { |
80 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
58 | PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))(PetscVoidFn*)monitor,mctx,monitordestroy,(PetscErrorCode (*)(void))(PetscVoidFn*)svd->monitor[i],svd->monitorcontext[i],svd->monitordestroy[i],&identical)); |
81 |
2/14✓ Branch 0 taken 8 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
|
58 | if (identical) PetscFunctionReturn(PETSC_SUCCESS); |
82 | } | ||
83 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
275 | PetscCheck(svd->numbermonitors<MAXSVDMONITORS,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Too many SVD monitors set"); |
84 | 275 | svd->monitor[svd->numbermonitors] = monitor; | |
85 | 275 | svd->monitorcontext[svd->numbermonitors] = mctx; | |
86 | 275 | svd->monitordestroy[svd->numbermonitors++] = monitordestroy; | |
87 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
275 | PetscFunctionReturn(PETSC_SUCCESS); |
88 | } | ||
89 | |||
90 | /*@ | ||
91 | SVDMonitorCancel - Clears all monitors for an SVD object. | ||
92 | |||
93 | Logically Collective | ||
94 | |||
95 | Input Parameters: | ||
96 | . svd - singular value solver context obtained from SVDCreate() | ||
97 | |||
98 | Options Database Key: | ||
99 | . -svd_monitor_cancel - Cancels all monitors that have been hardwired | ||
100 | into a code by calls to SVDMonitorSet(), | ||
101 | but does not cancel those set via the options database. | ||
102 | |||
103 | Level: intermediate | ||
104 | |||
105 | .seealso: SVDMonitorSet() | ||
106 | @*/ | ||
107 | 2549 | PetscErrorCode SVDMonitorCancel(SVD svd) | |
108 | { | ||
109 | 2549 | PetscInt i; | |
110 | |||
111 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
2549 | PetscFunctionBegin; |
112 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
2549 | PetscValidHeaderSpecific(svd,SVD_CLASSID,1); |
113 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2824 | for (i=0; i<svd->numbermonitors; i++) { |
114 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
275 | if (svd->monitordestroy[i]) PetscCall((*svd->monitordestroy[i])(&svd->monitorcontext[i])); |
115 | } | ||
116 | 2549 | svd->numbermonitors = 0; | |
117 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
2549 | PetscFunctionReturn(PETSC_SUCCESS); |
118 | } | ||
119 | |||
120 | /*@C | ||
121 | SVDGetMonitorContext - Gets the monitor context, as set by | ||
122 | SVDMonitorSet() for the FIRST monitor only. | ||
123 | |||
124 | Not Collective | ||
125 | |||
126 | Input Parameter: | ||
127 | . svd - singular value solver context obtained from SVDCreate() | ||
128 | |||
129 | Output Parameter: | ||
130 | . ctx - monitor context | ||
131 | |||
132 | Level: intermediate | ||
133 | |||
134 | .seealso: SVDMonitorSet() | ||
135 | @*/ | ||
136 | ✗ | PetscErrorCode SVDGetMonitorContext(SVD svd,void *ctx) | |
137 | { | ||
138 | ✗ | PetscFunctionBegin; | |
139 | ✗ | PetscValidHeaderSpecific(svd,SVD_CLASSID,1); | |
140 | ✗ | *(void**)ctx = svd->monitorcontext[0]; | |
141 | ✗ | PetscFunctionReturn(PETSC_SUCCESS); | |
142 | } | ||
143 | |||
144 | /*@C | ||
145 | SVDMonitorFirst - Print the first unconverged approximate value and | ||
146 | error estimate at each iteration of the singular value solver. | ||
147 | |||
148 | Collective | ||
149 | |||
150 | Input Parameters: | ||
151 | + svd - singular value solver context | ||
152 | . its - iteration number | ||
153 | . nconv - number of converged singular triplets so far | ||
154 | . sigma - singular values | ||
155 | . errest - error estimates | ||
156 | . nest - number of error estimates to display | ||
157 | - vf - viewer and format for monitoring | ||
158 | |||
159 | Options Database Key: | ||
160 | . -svd_monitor - activates SVDMonitorFirst() | ||
161 | |||
162 | Level: intermediate | ||
163 | |||
164 | .seealso: SVDMonitorSet(), SVDMonitorAll(), SVDMonitorConditioning(), SVDMonitorConverged() | ||
165 | @*/ | ||
166 | 70 | PetscErrorCode SVDMonitorFirst(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf) | |
167 | { | ||
168 | 70 | PetscViewer viewer = vf->viewer; | |
169 | |||
170 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
70 | PetscFunctionBegin; |
171 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
70 | PetscValidHeaderSpecific(svd,SVD_CLASSID,1); |
172 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
70 | PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,7); |
173 |
7/10✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
|
70 | if (its==1 && ((PetscObject)svd)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Singular value approximations and residual norms for %s solve.\n",((PetscObject)svd)->prefix)); |
174 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
70 | if (nconv<nest) { |
175 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
70 | PetscCall(PetscViewerPushFormat(viewer,vf->format)); |
176 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
70 | PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel)); |
177 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
70 | PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " SVD nconv=%" PetscInt_FMT " first unconverged value (error)",its,nconv)); |
178 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
70 | PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE)); |
179 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
70 | PetscCall(PetscViewerASCIIPrintf(viewer," %g (%10.8e)\n",(double)sigma[nconv],(double)errest[nconv])); |
180 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
70 | PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE)); |
181 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
70 | PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel)); |
182 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
70 | PetscCall(PetscViewerPopFormat(viewer)); |
183 | } | ||
184 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
14 | PetscFunctionReturn(PETSC_SUCCESS); |
185 | } | ||
186 | |||
187 | /*@C | ||
188 | SVDMonitorAll - Print the current approximate values and | ||
189 | error estimates at each iteration of the singular value solver. | ||
190 | |||
191 | Collective | ||
192 | |||
193 | Input Parameters: | ||
194 | + svd - singular value solver context | ||
195 | . its - iteration number | ||
196 | . nconv - number of converged singular triplets so far | ||
197 | . sigma - singular values | ||
198 | . errest - error estimates | ||
199 | . nest - number of error estimates to display | ||
200 | - vf - viewer and format for monitoring | ||
201 | |||
202 | Options Database Key: | ||
203 | . -svd_monitor_all - activates SVDMonitorAll() | ||
204 | |||
205 | Level: intermediate | ||
206 | |||
207 | .seealso: SVDMonitorSet(), SVDMonitorFirst(), SVDMonitorConditioning(), SVDMonitorConverged() | ||
208 | @*/ | ||
209 | 20 | PetscErrorCode SVDMonitorAll(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf) | |
210 | { | ||
211 | 20 | PetscInt i; | |
212 | 20 | PetscViewer viewer = vf->viewer; | |
213 | |||
214 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
20 | PetscFunctionBegin; |
215 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
20 | PetscValidHeaderSpecific(svd,SVD_CLASSID,1); |
216 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
20 | PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,7); |
217 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(PetscViewerPushFormat(viewer,vf->format)); |
218 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel)); |
219 |
2/10✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
20 | if (its==1 && ((PetscObject)svd)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Singular value approximations and residual norms for %s solve.\n",((PetscObject)svd)->prefix)); |
220 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " SVD nconv=%" PetscInt_FMT " Values (Errors)",its,nconv)); |
221 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE)); |
222 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
140 | for (i=0;i<nest;i++) PetscCall(PetscViewerASCIIPrintf(viewer," %g (%10.8e)",(double)sigma[i],(double)errest[i])); |
223 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(PetscViewerASCIIPrintf(viewer,"\n")); |
224 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE)); |
225 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel)); |
226 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(PetscViewerPopFormat(viewer)); |
227 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
4 | PetscFunctionReturn(PETSC_SUCCESS); |
228 | } | ||
229 | |||
230 | /*@C | ||
231 | SVDMonitorConverged - Print the approximate values and | ||
232 | error estimates as they converge. | ||
233 | |||
234 | Collective | ||
235 | |||
236 | Input Parameters: | ||
237 | + svd - singular value solver context | ||
238 | . its - iteration number | ||
239 | . nconv - number of converged singular triplets so far | ||
240 | . sigma - singular values | ||
241 | . errest - error estimates | ||
242 | . nest - number of error estimates to display | ||
243 | - vf - viewer and format for monitoring | ||
244 | |||
245 | Options Database Key: | ||
246 | . -svd_monitor_conv - activates SVDMonitorConverged() | ||
247 | |||
248 | Level: intermediate | ||
249 | |||
250 | .seealso: SVDMonitorSet(), SVDMonitorFirst(), SVDMonitorConditioning(), SVDMonitorAll() | ||
251 | @*/ | ||
252 | 20 | PetscErrorCode SVDMonitorConverged(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf) | |
253 | { | ||
254 | 20 | PetscInt i; | |
255 | 20 | PetscViewer viewer = vf->viewer; | |
256 | 20 | SlepcConvMon ctx; | |
257 | |||
258 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
20 | PetscFunctionBegin; |
259 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
20 | PetscValidHeaderSpecific(svd,SVD_CLASSID,1); |
260 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
20 | PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,7); |
261 | 20 | ctx = (SlepcConvMon)vf->data; | |
262 |
2/10✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
20 | if (its==1 && ((PetscObject)svd)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Convergence history for %s solve.\n",((PetscObject)svd)->prefix)); |
263 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
20 | if (its==1) ctx->oldnconv = 0; |
264 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
20 | if (ctx->oldnconv!=nconv) { |
265 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(PetscViewerPushFormat(viewer,vf->format)); |
266 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel)); |
267 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
140 | for (i=ctx->oldnconv;i<nconv;i++) { |
268 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
120 | PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " SVD converged value (error) #%" PetscInt_FMT,its,i)); |
269 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
120 | PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE)); |
270 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
120 | PetscCall(PetscViewerASCIIPrintf(viewer," %g (%10.8e)\n",(double)sigma[i],(double)errest[i])); |
271 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
120 | PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE)); |
272 | } | ||
273 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel)); |
274 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(PetscViewerPopFormat(viewer)); |
275 | 20 | ctx->oldnconv = nconv; | |
276 | } | ||
277 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
4 | PetscFunctionReturn(PETSC_SUCCESS); |
278 | } | ||
279 | |||
280 | 50 | PetscErrorCode SVDMonitorConvergedCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf) | |
281 | { | ||
282 | 50 | SlepcConvMon mctx; | |
283 | |||
284 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
50 | PetscFunctionBegin; |
285 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(PetscViewerAndFormatCreate(viewer,format,vf)); |
286 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(PetscNew(&mctx)); |
287 | 50 | mctx->ctx = ctx; | |
288 | 50 | (*vf)->data = (void*)mctx; | |
289 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
50 | PetscFunctionReturn(PETSC_SUCCESS); |
290 | } | ||
291 | |||
292 | 50 | PetscErrorCode SVDMonitorConvergedDestroy(PetscViewerAndFormat **vf) | |
293 | { | ||
294 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
50 | PetscFunctionBegin; |
295 |
2/14✓ Branch 0 taken 8 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
|
50 | if (!*vf) PetscFunctionReturn(PETSC_SUCCESS); |
296 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
50 | PetscCall(PetscFree((*vf)->data)); |
297 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(PetscViewerDestroy(&(*vf)->viewer)); |
298 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
50 | PetscCall(PetscFree(*vf)); |
299 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
4 | PetscFunctionReturn(PETSC_SUCCESS); |
300 | } | ||
301 | |||
302 | /*@C | ||
303 | SVDMonitorFirstDrawLG - Plots the error estimate of the first unconverged | ||
304 | approximation at each iteration of the singular value solver. | ||
305 | |||
306 | Collective | ||
307 | |||
308 | Input Parameters: | ||
309 | + svd - singular value solver context | ||
310 | . its - iteration number | ||
311 | . nconv - number of converged singular triplets so far | ||
312 | . sigma - singular values | ||
313 | . errest - error estimates | ||
314 | . nest - number of error estimates to display | ||
315 | - vf - viewer and format for monitoring | ||
316 | |||
317 | Options Database Key: | ||
318 | . -svd_monitor draw::draw_lg - activates SVDMonitorFirstDrawLG() | ||
319 | |||
320 | Level: intermediate | ||
321 | |||
322 | .seealso: SVDMonitorSet() | ||
323 | @*/ | ||
324 | 32 | PetscErrorCode SVDMonitorFirstDrawLG(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf) | |
325 | { | ||
326 | 32 | PetscViewer viewer = vf->viewer; | |
327 | 32 | PetscDrawLG lg; | |
328 | 32 | PetscReal x,y; | |
329 | |||
330 |
1/2✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
|
32 | PetscFunctionBegin; |
331 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 1 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
32 | PetscValidHeaderSpecific(svd,SVD_CLASSID,1); |
332 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 1 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
32 | PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,7); |
333 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
32 | PetscCall(PetscViewerPushFormat(viewer,vf->format)); |
334 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
32 | PetscCall(PetscViewerDrawGetDrawLG(viewer,0,&lg)); |
335 |
2/2✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
|
32 | if (its==1) { |
336 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
16 | PetscCall(PetscDrawLGReset(lg)); |
337 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
16 | PetscCall(PetscDrawLGSetDimension(lg,1)); |
338 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
16 | PetscCall(PetscDrawLGSetLimits(lg,1,1,PetscLog10Real(svd->tol)-2,0.0)); |
339 | } | ||
340 |
2/2✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
|
32 | if (nconv<nest) { |
341 | 24 | x = (PetscReal)its; | |
342 |
1/2✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
|
24 | if (errest[nconv] > 0.0) y = PetscLog10Real(errest[nconv]); |
343 | ✗ | else y = 0.0; | |
344 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
24 | PetscCall(PetscDrawLGAddPoint(lg,&x,&y)); |
345 |
1/6✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
24 | if (its <= 20 || !(its % 5) || svd->reason) { |
346 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
24 | PetscCall(PetscDrawLGDraw(lg)); |
347 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
24 | PetscCall(PetscDrawLGSave(lg)); |
348 | } | ||
349 | } | ||
350 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
32 | PetscCall(PetscViewerPopFormat(viewer)); |
351 |
6/12✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
|
4 | PetscFunctionReturn(PETSC_SUCCESS); |
352 | } | ||
353 | |||
354 | /*@C | ||
355 | SVDMonitorFirstDrawLGCreate - Creates the plotter for the first error estimate. | ||
356 | |||
357 | Collective | ||
358 | |||
359 | Input Parameters: | ||
360 | + viewer - the viewer | ||
361 | . format - the viewer format | ||
362 | - ctx - an optional user context | ||
363 | |||
364 | Output Parameter: | ||
365 | . vf - the viewer and format context | ||
366 | |||
367 | Level: intermediate | ||
368 | |||
369 | .seealso: SVDMonitorSet() | ||
370 | @*/ | ||
371 | 16 | PetscErrorCode SVDMonitorFirstDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf) | |
372 | { | ||
373 |
1/2✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
|
16 | PetscFunctionBegin; |
374 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
16 | PetscCall(PetscViewerAndFormatCreate(viewer,format,vf)); |
375 | 16 | (*vf)->data = ctx; | |
376 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
16 | PetscCall(PetscViewerMonitorLGSetUp(viewer,NULL,"First Error Estimate","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300)); |
377 |
6/12✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
|
2 | PetscFunctionReturn(PETSC_SUCCESS); |
378 | } | ||
379 | |||
380 | /*@C | ||
381 | SVDMonitorAllDrawLG - Plots the error estimate of all unconverged | ||
382 | approximations at each iteration of the singular value solver. | ||
383 | |||
384 | Collective | ||
385 | |||
386 | Input Parameters: | ||
387 | + svd - singular value solver context | ||
388 | . its - iteration number | ||
389 | . nconv - number of converged singular triplets so far | ||
390 | . sigma - singular values | ||
391 | . errest - error estimates | ||
392 | . nest - number of error estimates to display | ||
393 | - vf - viewer and format for monitoring | ||
394 | |||
395 | Options Database Key: | ||
396 | . -svd_monitor_all draw::draw_lg - activates SVDMonitorAllDrawLG() | ||
397 | |||
398 | Level: intermediate | ||
399 | |||
400 | .seealso: SVDMonitorSet() | ||
401 | @*/ | ||
402 | 8 | PetscErrorCode SVDMonitorAllDrawLG(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf) | |
403 | { | ||
404 | 8 | PetscViewer viewer = vf->viewer; | |
405 | 8 | PetscDrawLG lg; | |
406 | 8 | PetscInt i,n = PetscMin(svd->nsv,255); | |
407 | 8 | PetscReal *x,*y; | |
408 | |||
409 |
1/2✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
|
8 | PetscFunctionBegin; |
410 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 1 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
8 | PetscValidHeaderSpecific(svd,SVD_CLASSID,1); |
411 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 1 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
8 | PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,7); |
412 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
8 | PetscCall(PetscViewerPushFormat(viewer,vf->format)); |
413 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
8 | PetscCall(PetscViewerDrawGetDrawLG(viewer,0,&lg)); |
414 |
1/2✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
|
8 | if (its==1) { |
415 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
8 | PetscCall(PetscDrawLGReset(lg)); |
416 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
8 | PetscCall(PetscDrawLGSetDimension(lg,n)); |
417 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
8 | PetscCall(PetscDrawLGSetLimits(lg,1,1,PetscLog10Real(svd->tol)-2,0.0)); |
418 | } | ||
419 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
8 | PetscCall(PetscMalloc2(n,&x,n,&y)); |
420 |
2/2✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
|
16 | for (i=0;i<n;i++) { |
421 | 8 | x[i] = (PetscReal)its; | |
422 |
2/4✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 8 times.
|
8 | if (i < nest && errest[i] > 0.0) y[i] = PetscLog10Real(errest[i]); |
423 | 8 | else y[i] = 0.0; | |
424 | } | ||
425 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
8 | PetscCall(PetscDrawLGAddPoint(lg,x,y)); |
426 |
1/6✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
8 | if (its <= 20 || !(its % 5) || svd->reason) { |
427 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
8 | PetscCall(PetscDrawLGDraw(lg)); |
428 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
8 | PetscCall(PetscDrawLGSave(lg)); |
429 | } | ||
430 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
8 | PetscCall(PetscFree2(x,y)); |
431 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
8 | PetscCall(PetscViewerPopFormat(viewer)); |
432 |
6/12✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
|
1 | PetscFunctionReturn(PETSC_SUCCESS); |
433 | } | ||
434 | |||
435 | /*@C | ||
436 | SVDMonitorAllDrawLGCreate - Creates the plotter for all the error estimates. | ||
437 | |||
438 | Collective | ||
439 | |||
440 | Input Parameters: | ||
441 | + viewer - the viewer | ||
442 | . format - the viewer format | ||
443 | - ctx - an optional user context | ||
444 | |||
445 | Output Parameter: | ||
446 | . vf - the viewer and format context | ||
447 | |||
448 | Level: intermediate | ||
449 | |||
450 | .seealso: SVDMonitorSet() | ||
451 | @*/ | ||
452 | 8 | PetscErrorCode SVDMonitorAllDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf) | |
453 | { | ||
454 |
1/2✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
|
8 | PetscFunctionBegin; |
455 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
8 | PetscCall(PetscViewerAndFormatCreate(viewer,format,vf)); |
456 | 8 | (*vf)->data = ctx; | |
457 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 7 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
8 | PetscCall(PetscViewerMonitorLGSetUp(viewer,NULL,"All Error Estimates","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300)); |
458 |
6/12✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
|
1 | PetscFunctionReturn(PETSC_SUCCESS); |
459 | } | ||
460 | |||
461 | /*@C | ||
462 | SVDMonitorConvergedDrawLG - Plots the number of converged eigenvalues | ||
463 | at each iteration of the singular value solver. | ||
464 | |||
465 | Collective | ||
466 | |||
467 | Input Parameters: | ||
468 | + svd - singular value solver context | ||
469 | . its - iteration number | ||
470 | . nconv - number of converged singular triplets so far | ||
471 | . sigma - singular values | ||
472 | . errest - error estimates | ||
473 | . nest - number of error estimates to display | ||
474 | - vf - viewer and format for monitoring | ||
475 | |||
476 | Options Database Key: | ||
477 | . -svd_monitor_conv draw::draw_lg - activates SVDMonitorConvergedDrawLG() | ||
478 | |||
479 | Level: intermediate | ||
480 | |||
481 | .seealso: SVDMonitorSet() | ||
482 | @*/ | ||
483 | ✗ | PetscErrorCode SVDMonitorConvergedDrawLG(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf) | |
484 | { | ||
485 | ✗ | PetscViewer viewer = vf->viewer; | |
486 | ✗ | PetscDrawLG lg; | |
487 | ✗ | PetscReal x,y; | |
488 | |||
489 | ✗ | PetscFunctionBegin; | |
490 | ✗ | PetscValidHeaderSpecific(svd,SVD_CLASSID,1); | |
491 | ✗ | PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,7); | |
492 | ✗ | PetscCall(PetscViewerPushFormat(viewer,vf->format)); | |
493 | ✗ | PetscCall(PetscViewerDrawGetDrawLG(viewer,0,&lg)); | |
494 | ✗ | if (its==1) { | |
495 | ✗ | PetscCall(PetscDrawLGReset(lg)); | |
496 | ✗ | PetscCall(PetscDrawLGSetDimension(lg,1)); | |
497 | ✗ | PetscCall(PetscDrawLGSetLimits(lg,1,1,0,svd->nsv)); | |
498 | } | ||
499 | ✗ | x = (PetscReal)its; | |
500 | ✗ | y = (PetscReal)svd->nconv; | |
501 | ✗ | PetscCall(PetscDrawLGAddPoint(lg,&x,&y)); | |
502 | ✗ | if (its <= 20 || !(its % 5) || svd->reason) { | |
503 | ✗ | PetscCall(PetscDrawLGDraw(lg)); | |
504 | ✗ | PetscCall(PetscDrawLGSave(lg)); | |
505 | } | ||
506 | ✗ | PetscCall(PetscViewerPopFormat(viewer)); | |
507 | ✗ | PetscFunctionReturn(PETSC_SUCCESS); | |
508 | } | ||
509 | |||
510 | /*@C | ||
511 | SVDMonitorConvergedDrawLGCreate - Creates the plotter for the convergence history. | ||
512 | |||
513 | Collective | ||
514 | |||
515 | Input Parameters: | ||
516 | + viewer - the viewer | ||
517 | . format - the viewer format | ||
518 | - ctx - an optional user context | ||
519 | |||
520 | Output Parameter: | ||
521 | . vf - the viewer and format context | ||
522 | |||
523 | Level: intermediate | ||
524 | |||
525 | .seealso: SVDMonitorSet() | ||
526 | @*/ | ||
527 | ✗ | PetscErrorCode SVDMonitorConvergedDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf) | |
528 | { | ||
529 | ✗ | SlepcConvMon mctx; | |
530 | |||
531 | ✗ | PetscFunctionBegin; | |
532 | ✗ | PetscCall(PetscViewerAndFormatCreate(viewer,format,vf)); | |
533 | ✗ | PetscCall(PetscNew(&mctx)); | |
534 | ✗ | mctx->ctx = ctx; | |
535 | ✗ | (*vf)->data = (void*)mctx; | |
536 | ✗ | PetscCall(PetscViewerMonitorLGSetUp(viewer,NULL,"Convergence History","Number Converged",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300)); | |
537 | ✗ | PetscFunctionReturn(PETSC_SUCCESS); | |
538 | } | ||
539 | |||
540 | /*@C | ||
541 | SVDMonitorConditioning - Print the condition number at each iteration of the singular | ||
542 | value solver. | ||
543 | |||
544 | Collective | ||
545 | |||
546 | Input Parameters: | ||
547 | + svd - singular value solver context | ||
548 | . its - iteration number | ||
549 | . nconv - (unused) number of converged singular triplets so far | ||
550 | . sigma - (unused) singular values | ||
551 | . errest - (unused) error estimates | ||
552 | . nest - (unused) number of error estimates to display | ||
553 | - vf - viewer and format for monitoring | ||
554 | |||
555 | Options Database Key: | ||
556 | . -svd_monitor_conditioning - activates SVDMonitorConditioning() | ||
557 | |||
558 | Note: | ||
559 | Works only for solvers that use a DS of type GSVD. The printed information corresponds | ||
560 | to the maximum of the condition number of the two generated bidiagonal matrices. | ||
561 | |||
562 | Level: intermediate | ||
563 | |||
564 | .seealso: SVDMonitorSet(), SVDMonitorAll(), SVDMonitorFirst(), SVDMonitorConverged() | ||
565 | @*/ | ||
566 | 70 | PetscErrorCode SVDMonitorConditioning(SVD svd,PetscInt its,PetscInt nconv,PetscReal *sigma,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf) | |
567 | { | ||
568 | 70 | PetscViewer viewer = vf->viewer; | |
569 | 70 | PetscBool isgsvd; | |
570 | 70 | PetscReal cond; | |
571 | |||
572 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
70 | PetscFunctionBegin; |
573 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
70 | PetscValidHeaderSpecific(svd,SVD_CLASSID,1); |
574 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
70 | PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,7); |
575 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
70 | PetscCall(PetscObjectTypeCompare((PetscObject)svd->ds,DSGSVD,&isgsvd)); |
576 |
2/14✓ Branch 0 taken 8 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
|
70 | if (!isgsvd) PetscFunctionReturn(PETSC_SUCCESS); |
577 |
3/10✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
70 | if (its==1 && ((PetscObject)svd)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Condition number of bidiagonal matrices for %s solve.\n",((PetscObject)svd)->prefix)); |
578 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
70 | PetscCall(PetscViewerPushFormat(viewer,vf->format)); |
579 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
70 | PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)svd)->tablevel)); |
580 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
70 | PetscCall(DSCond(svd->ds,&cond)); |
581 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
70 | PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " SVD condition number = %g\n",its,(double)cond)); |
582 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
70 | PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)svd)->tablevel)); |
583 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
70 | PetscCall(PetscViewerPopFormat(viewer)); |
584 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
14 | PetscFunctionReturn(PETSC_SUCCESS); |
585 | } | ||
586 |