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 | LME routines related to monitors | ||
12 | */ | ||
13 | |||
14 | #include <slepc/private/lmeimpl.h> /*I "slepclme.h" I*/ | ||
15 | #include <petscdraw.h> | ||
16 | |||
17 | /* | ||
18 | Runs the user provided monitor routines, if any. | ||
19 | */ | ||
20 | 789 | PetscErrorCode LMEMonitor(LME lme,PetscInt it,PetscReal errest) | |
21 | { | ||
22 | 789 | PetscInt i,n = lme->numbermonitors; | |
23 | |||
24 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
789 | 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.
|
861 | for (i=0;i<n;i++) PetscCall((*lme->monitor[i])(lme,it,errest,lme->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.
|
126 | PetscFunctionReturn(PETSC_SUCCESS); |
27 | } | ||
28 | |||
29 | /*@C | ||
30 | LMEMonitorSet - Sets an ADDITIONAL function to be called at every | ||
31 | iteration to monitor convergence. | ||
32 | |||
33 | Logically Collective | ||
34 | |||
35 | Input Parameters: | ||
36 | + lme - linear matrix equation solver context obtained from LMECreate() | ||
37 | . monitor - pointer to function (if this is NULL, it turns off monitoring), see LMEMonitorFn | ||
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 | + -lme_monitor - print the error estimate | ||
45 | . -lme_monitor draw::draw_lg - sets line graph monitor for the error estimate | ||
46 | - -lme_monitor_cancel - cancels all monitors that have been hardwired into | ||
47 | a code by calls to LMEMonitorSet(), but does not cancel those set via | ||
48 | the options database. | ||
49 | |||
50 | Notes: | ||
51 | The options database option -lme_monitor and related options are the easiest way | ||
52 | to turn on LME iteration monitoring. | ||
53 | |||
54 | LMEMonitorRegister() provides a way to associate an options database key with LME | ||
55 | monitor function. | ||
56 | |||
57 | Several different monitoring routines may be set by calling LMEMonitorSet() multiple | ||
58 | times; all will be called in the order in which they were set. | ||
59 | |||
60 | Level: intermediate | ||
61 | |||
62 | .seealso: LMEMonitorCancel() | ||
63 | @*/ | ||
64 | 46 | PetscErrorCode LMEMonitorSet(LME lme,LMEMonitorFn *monitor,void *mctx,PetscCtxDestroyFn *monitordestroy) | |
65 | { | ||
66 | 46 | PetscInt i; | |
67 | 46 | PetscBool identical; | |
68 | |||
69 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
46 | PetscFunctionBegin; |
70 |
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.
|
46 | PetscValidHeaderSpecific(lme,LME_CLASSID,1); |
71 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
46 | for (i=0;i<lme->numbermonitors;i++) { |
72 |
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.
|
10 | PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))(PetscVoidFn*)monitor,mctx,monitordestroy,(PetscErrorCode (*)(void))(PetscVoidFn*)lme->monitor[i],lme->monitorcontext[i],lme->monitordestroy[i],&identical)); |
73 |
8/14✓ 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.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
|
10 | if (identical) PetscFunctionReturn(PETSC_SUCCESS); |
74 | } | ||
75 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
36 | PetscCheck(lme->numbermonitors<MAXLMEMONITORS,PetscObjectComm((PetscObject)lme),PETSC_ERR_ARG_OUTOFRANGE,"Too many LME monitors set"); |
76 | 36 | lme->monitor[lme->numbermonitors] = monitor; | |
77 | 36 | lme->monitorcontext[lme->numbermonitors] = mctx; | |
78 | 36 | lme->monitordestroy[lme->numbermonitors++] = monitordestroy; | |
79 |
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.
|
36 | PetscFunctionReturn(PETSC_SUCCESS); |
80 | } | ||
81 | |||
82 | /*@ | ||
83 | LMEMonitorCancel - Clears all monitors for an LME object. | ||
84 | |||
85 | Logically Collective | ||
86 | |||
87 | Input Parameters: | ||
88 | . lme - linear matrix equation solver context obtained from LMECreate() | ||
89 | |||
90 | Options Database Key: | ||
91 | . -lme_monitor_cancel - cancels all monitors that have been hardwired | ||
92 | into a code by calls to LMEMonitorSet(), | ||
93 | but does not cancel those set via the options database. | ||
94 | |||
95 | Level: intermediate | ||
96 | |||
97 | .seealso: LMEMonitorSet() | ||
98 | @*/ | ||
99 | 218 | PetscErrorCode LMEMonitorCancel(LME lme) | |
100 | { | ||
101 | 218 | PetscInt i; | |
102 | |||
103 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
218 | PetscFunctionBegin; |
104 |
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.
|
218 | PetscValidHeaderSpecific(lme,LME_CLASSID,1); |
105 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
254 | for (i=0; i<lme->numbermonitors; i++) { |
106 |
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.
|
36 | if (lme->monitordestroy[i]) PetscCall((*lme->monitordestroy[i])(&lme->monitorcontext[i])); |
107 | } | ||
108 | 218 | lme->numbermonitors = 0; | |
109 |
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.
|
218 | PetscFunctionReturn(PETSC_SUCCESS); |
110 | } | ||
111 | |||
112 | /*@C | ||
113 | LMEGetMonitorContext - Gets the monitor context, as set by | ||
114 | LMEMonitorSet() for the FIRST monitor only. | ||
115 | |||
116 | Not Collective | ||
117 | |||
118 | Input Parameter: | ||
119 | . lme - linear matrix equation solver context obtained from LMECreate() | ||
120 | |||
121 | Output Parameter: | ||
122 | . ctx - monitor context | ||
123 | |||
124 | Level: intermediate | ||
125 | |||
126 | .seealso: LMEMonitorSet() | ||
127 | @*/ | ||
128 | ✗ | PetscErrorCode LMEGetMonitorContext(LME lme,void *ctx) | |
129 | { | ||
130 | ✗ | PetscFunctionBegin; | |
131 | ✗ | PetscValidHeaderSpecific(lme,LME_CLASSID,1); | |
132 | ✗ | *(void**)ctx = lme->monitorcontext[0]; | |
133 | ✗ | PetscFunctionReturn(PETSC_SUCCESS); | |
134 | } | ||
135 | |||
136 | /*@C | ||
137 | LMEMonitorDefault - Print the error estimate of the current approximation at each | ||
138 | iteration of the linear matrix equation solver. | ||
139 | |||
140 | Collective | ||
141 | |||
142 | Input Parameters: | ||
143 | + lme - linear matrix equation solver context | ||
144 | . its - iteration number | ||
145 | . errest - error estimate | ||
146 | - vf - viewer and format for monitoring | ||
147 | |||
148 | Options Database Key: | ||
149 | . -lme_monitor - activates LMEMonitorDefault() | ||
150 | |||
151 | Level: intermediate | ||
152 | |||
153 | .seealso: LMEMonitorSet() | ||
154 | @*/ | ||
155 | 40 | PetscErrorCode LMEMonitorDefault(LME lme,PetscInt its,PetscReal errest,PetscViewerAndFormat *vf) | |
156 | { | ||
157 | 40 | PetscViewer viewer = vf->viewer; | |
158 | |||
159 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
40 | PetscFunctionBegin; |
160 |
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.
|
40 | PetscValidHeaderSpecific(lme,LME_CLASSID,1); |
161 |
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.
|
40 | PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); |
162 |
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.
|
40 | PetscCall(PetscViewerPushFormat(viewer,vf->format)); |
163 |
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.
|
40 | PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)lme)->tablevel)); |
164 |
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.
|
40 | if (its == 1 && ((PetscObject)lme)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Error estimates for %s solve.\n",((PetscObject)lme)->prefix)); |
165 |
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.
|
40 | PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " LME Error estimate %14.12e\n",its,(double)errest)); |
166 |
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.
|
40 | PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)lme)->tablevel)); |
167 |
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.
|
40 | PetscCall(PetscViewerPopFormat(viewer)); |
168 |
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.
|
8 | PetscFunctionReturn(PETSC_SUCCESS); |
169 | } | ||
170 | |||
171 | /*@C | ||
172 | LMEMonitorDefaultDrawLG - Plots the error estimate of the current approximation at each | ||
173 | iteration of the linear matrix equation solver. | ||
174 | |||
175 | Collective | ||
176 | |||
177 | Input Parameters: | ||
178 | + lme - linear matrix equation solver context | ||
179 | . its - iteration number | ||
180 | . errest - error estimate | ||
181 | - vf - viewer and format for monitoring | ||
182 | |||
183 | Options Database Key: | ||
184 | . -lme_monitor draw::draw_lg - activates LMEMonitorDefaultDrawLG() | ||
185 | |||
186 | Level: intermediate | ||
187 | |||
188 | .seealso: LMEMonitorSet() | ||
189 | @*/ | ||
190 | 32 | PetscErrorCode LMEMonitorDefaultDrawLG(LME lme,PetscInt its,PetscReal errest,PetscViewerAndFormat *vf) | |
191 | { | ||
192 | 32 | PetscViewer viewer = vf->viewer; | |
193 | 32 | PetscDrawLG lg; | |
194 | 32 | PetscReal x,y; | |
195 | |||
196 |
1/2✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
|
32 | PetscFunctionBegin; |
197 |
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(lme,LME_CLASSID,1); |
198 |
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,4); |
199 |
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)); |
200 |
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)); |
201 |
2/2✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
|
32 | if (its==1) { |
202 |
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)); |
203 |
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,1)); |
204 |
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.0,PetscLog10Real(lme->tol)-2,0.0)); |
205 | } | ||
206 | 32 | x = (PetscReal)its; | |
207 |
1/2✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
|
32 | if (errest > 0.0) y = PetscLog10Real(errest); |
208 | ✗ | else y = 0.0; | |
209 |
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(PetscDrawLGAddPoint(lg,&x,&y)); |
210 |
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.
|
32 | if (its <= 20 || !(its % 5) || lme->reason) { |
211 |
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(PetscDrawLGDraw(lg)); |
212 |
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(PetscDrawLGSave(lg)); |
213 | } | ||
214 |
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)); |
215 |
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); |
216 | } | ||
217 | |||
218 | /*@C | ||
219 | LMEMonitorDefaultDrawLGCreate - Creates the plotter for the error estimate. | ||
220 | |||
221 | Collective | ||
222 | |||
223 | Input Parameters: | ||
224 | + viewer - the viewer | ||
225 | . format - the viewer format | ||
226 | - ctx - an optional user context | ||
227 | |||
228 | Output Parameter: | ||
229 | . vf - the viewer and format context | ||
230 | |||
231 | Level: intermediate | ||
232 | |||
233 | .seealso: LMEMonitorSet() | ||
234 | @*/ | ||
235 | 8 | PetscErrorCode LMEMonitorDefaultDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf) | |
236 | { | ||
237 |
1/2✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
|
8 | PetscFunctionBegin; |
238 |
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)); |
239 | 8 | (*vf)->data = ctx; | |
240 |
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,"Error Estimate","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300)); |
241 |
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); |
242 | } | ||
243 |