1 /*
2 This file contains some simple default routines.
3 These routines should be SHORT, since they will be included in every
4 executable image that uses the iterative routines (note that, through
5 the registry system, we provide a way to load only the truly necessary
6 files)
7 */
8 #include <petsc/private/kspimpl.h> /*I "petscksp.h" I*/
9 #include <petscdmshell.h>
10 #include <petscdraw.h>
11
12 /*@
13 KSPGetResidualNorm - Gets the last (possibly approximate and/or preconditioned) residual norm that has been computed.
14
15 Not Collective
16
17 Input Parameter:
18 . ksp - the iterative context
19
20 Output Parameter:
21 . rnorm - residual norm
22
23 Level: intermediate
24
25 Notes:
26 For some methods, such as `KSPGMRES`, the norm is not computed directly from the residual.
27
28 The type of norm used by the method can be controlled with `KSPSetNormType()`
29
30 Certain solvers, under certain conditions, may not compute the final residual norm in an iteration, in that case the previous norm is returned.
31
32 .seealso: [](ch_ksp), `KSP`, `KSPSetNormType()`, `KSPBuildResidual()`, `KSPNormType`
33 @*/
KSPGetResidualNorm(KSP ksp,PetscReal * rnorm)34 PetscErrorCode KSPGetResidualNorm(KSP ksp, PetscReal *rnorm)
35 {
36 PetscFunctionBegin;
37 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
38 PetscAssertPointer(rnorm, 2);
39 *rnorm = ksp->rnorm;
40 PetscFunctionReturn(PETSC_SUCCESS);
41 }
42
43 /*@
44 KSPGetIterationNumber - Gets the current iteration number; if the `KSPSolve()` is complete, returns the number of iterations used.
45
46 Not Collective
47
48 Input Parameter:
49 . ksp - the iterative context
50
51 Output Parameter:
52 . its - number of iterations
53
54 Level: intermediate
55
56 Note:
57 During the ith iteration this returns i-1
58
59 .seealso: [](ch_ksp), `KSP`, `KSPGetResidualNorm()`, `KSPBuildResidual()`, `KSPGetTotalIterations()`
60 @*/
KSPGetIterationNumber(KSP ksp,PetscInt * its)61 PetscErrorCode KSPGetIterationNumber(KSP ksp, PetscInt *its)
62 {
63 PetscFunctionBegin;
64 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
65 PetscAssertPointer(its, 2);
66 *its = ksp->its;
67 PetscFunctionReturn(PETSC_SUCCESS);
68 }
69
70 /*@
71 KSPGetTotalIterations - Gets the total number of iterations this `KSP` object has performed since was created, counted over all linear solves
72
73 Not Collective
74
75 Input Parameter:
76 . ksp - the iterative context
77
78 Output Parameter:
79 . its - total number of iterations
80
81 Level: intermediate
82
83 Note:
84 Use `KSPGetIterationNumber()` to get the count for the most recent solve only
85 If this is called within a `KSPSolve()` (such as in a `KSPMonitor` routine) then it does not include iterations within that current solve
86
87 .seealso: [](ch_ksp), `KSP`, `KSPBuildResidual()`, `KSPGetResidualNorm()`, `KSPGetIterationNumber()`
88 @*/
KSPGetTotalIterations(KSP ksp,PetscInt * its)89 PetscErrorCode KSPGetTotalIterations(KSP ksp, PetscInt *its)
90 {
91 PetscFunctionBegin;
92 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
93 PetscAssertPointer(its, 2);
94 *its = ksp->totalits;
95 PetscFunctionReturn(PETSC_SUCCESS);
96 }
97
98 /*@C
99 KSPMonitorResidual - Print the (possibly preconditioned, possibly approximate) residual norm at each iteration of an iterative solver.
100
101 Collective
102
103 Input Parameters:
104 + ksp - iterative context
105 . n - iteration number
106 . rnorm - (preconditioned) residual norm value (may be estimated).
107 - vf - The viewer context
108
109 Options Database Key:
110 . -ksp_monitor - Activates `KSPMonitorResidual()` to print the norm value at each iteration
111
112 Level: intermediate
113
114 Note:
115 For some methods, such as `KSPGMRES`, the norm is not computed directly from the residual.
116
117 The type of norm used by the method can be controlled with `KSPSetNormType()`
118
119 This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
120 to be used during the `KSP` solve.
121
122 .seealso: [](ch_ksp), `KSP`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `KSPMonitorResidualView()`, `KSPMonitorResidualDrawLG()`,
123 `KSPMonitorResidualRange()`, `KSPMonitorTrueResidualDraw()`, `KSPMonitorTrueResidualDrawLG()`, `KSPMonitorTrueResidualMax()`,
124 `KSPMonitorSingularValue()`, `KSPMonitorSolutionDrawLG()`, `KSPMonitorSolutionDraw()`, `KSPMonitorSolution()`,
125 `KSPMonitorErrorDrawLG()`, `KSPMonitorErrorDraw()`, `KSPMonitorError()`
126 @*/
KSPMonitorResidual(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat * vf)127 PetscErrorCode KSPMonitorResidual(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
128 {
129 PetscViewer viewer = vf->viewer;
130 PetscViewerFormat format = vf->format;
131 PetscInt tablevel;
132 const char *prefix;
133
134 PetscFunctionBegin;
135 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4);
136 PetscCall(PetscObjectGetTabLevel((PetscObject)ksp, &tablevel));
137 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ksp, &prefix));
138 PetscCall(PetscViewerPushFormat(viewer, format));
139 PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
140 if (n == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, " Residual norms for %s solve.\n", prefix));
141 PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " KSP Residual norm %14.12e\n", n, (double)rnorm));
142 PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
143 PetscCall(PetscViewerPopFormat(viewer));
144 PetscFunctionReturn(PETSC_SUCCESS);
145 }
146
147 /*@C
148 KSPMonitorResidualView - Plots the (possibly preconditioned) residual at each iteration of an iterative solver.
149
150 Collective
151
152 Input Parameters:
153 + ksp - iterative context
154 . n - iteration number
155 . rnorm - 2-norm (preconditioned) residual value (may be estimated).
156 - vf - The viewer context
157
158 Options Database Key:
159 . -ksp_monitor viewertype - Activates `KSPMonitorResidualView()`
160
161 Level: intermediate
162
163 Note:
164 This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
165 to be used during the `KSP` solve.
166
167 .seealso: [](ch_ksp), `KSP`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `KSPMonitorResidual()`, `KSPMonitorResidualDrawLG()`
168 @*/
KSPMonitorResidualView(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat * vf)169 PetscErrorCode KSPMonitorResidualView(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
170 {
171 PetscViewer viewer = vf->viewer;
172 PetscViewerFormat format = vf->format;
173 Vec r;
174
175 PetscFunctionBegin;
176 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4);
177 PetscCall(PetscViewerPushFormat(viewer, format));
178 PetscCall(KSPBuildResidual(ksp, NULL, NULL, &r));
179 PetscCall(PetscObjectSetName((PetscObject)r, "Residual"));
180 PetscCall(PetscObjectCompose((PetscObject)r, "__Vec_bc_zero__", (PetscObject)ksp));
181 PetscCall(VecView(r, viewer));
182 PetscCall(PetscObjectCompose((PetscObject)r, "__Vec_bc_zero__", NULL));
183 PetscCall(VecDestroy(&r));
184 PetscCall(PetscViewerPopFormat(viewer));
185 PetscFunctionReturn(PETSC_SUCCESS);
186 }
187
188 /*@C
189 KSPMonitorResidualDrawLG - Plots the (possibly preconditioned) residual norm at each iteration of an iterative solver.
190
191 Collective
192
193 Input Parameters:
194 + ksp - iterative context
195 . n - iteration number
196 . rnorm - 2-norm (preconditioned) residual value (may be estimated).
197 - vf - The viewer context
198
199 Options Database Key:
200 . -ksp_monitor draw::draw_lg - Activates `KSPMonitorResidualDrawLG()`
201
202 Level: intermediate
203
204 Notes:
205 This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
206 to be used during the `KSP` solve.
207
208 Use `KSPMonitorResidualDrawLGCreate()` to create the context used with this monitor
209
210 .seealso: [](ch_ksp), `KSP`, `PETSCVIEWERDRAW`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `KSPMonitorResidualView()`, `KSPMonitorResidual()`
211 @*/
KSPMonitorResidualDrawLG(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat * vf)212 PetscErrorCode KSPMonitorResidualDrawLG(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
213 {
214 PetscViewer viewer = vf->viewer;
215 PetscViewerFormat format = vf->format;
216 PetscDrawLG lg;
217 KSPConvergedReason reason;
218 PetscReal x, y;
219
220 PetscFunctionBegin;
221 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4);
222 PetscCall(PetscViewerPushFormat(viewer, format));
223 PetscCall(PetscViewerDrawGetDrawLG(viewer, 0, &lg));
224 if (!n) PetscCall(PetscDrawLGReset(lg));
225 x = (PetscReal)n;
226 if (rnorm > 0.0) y = PetscLog10Real(rnorm);
227 else y = -15.0;
228 PetscCall(PetscDrawLGAddPoint(lg, &x, &y));
229 PetscCall(KSPGetConvergedReason(ksp, &reason));
230 if (n <= 20 || !(n % 5) || reason) {
231 PetscCall(PetscDrawLGDraw(lg));
232 PetscCall(PetscDrawLGSave(lg));
233 }
234 PetscCall(PetscViewerPopFormat(viewer));
235 PetscFunctionReturn(PETSC_SUCCESS);
236 }
237
238 /*@C
239 KSPMonitorResidualDrawLGCreate - Creates the context for the (possibly preconditioned) residual norm monitor `KSPMonitorResidualDrawLG()`
240
241 Collective
242
243 Input Parameters:
244 + viewer - The `PetscViewer` of type `PETSCVIEWERDRAW`
245 . format - The viewer format
246 - ctx - An optional user context
247
248 Output Parameter:
249 . vf - The viewer context
250
251 Level: intermediate
252
253 .seealso: [](ch_ksp), `KSP`, `PETSCVIEWERDRAW`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `KSPMonitorResidualDrawLG()`,
254 `PetscViewerFormat`, `PetscViewer`, `PetscViewerAndFormat`
255 @*/
KSPMonitorResidualDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,PetscCtx ctx,PetscViewerAndFormat ** vf)256 PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer viewer, PetscViewerFormat format, PetscCtx ctx, PetscViewerAndFormat **vf)
257 {
258 PetscFunctionBegin;
259 PetscCall(PetscViewerAndFormatCreate(viewer, format, vf));
260 (*vf)->data = ctx;
261 PetscCall(PetscViewerMonitorLGSetUp(viewer, NULL, NULL, "Log Residual Norm", 1, NULL, PETSC_DECIDE, PETSC_DECIDE, 400, 300));
262 PetscFunctionReturn(PETSC_SUCCESS);
263 }
264
265 /*
266 This is the same as KSPMonitorResidual() except it prints fewer digits of the residual as the residual gets smaller.
267 This is because the later digits are meaningless and are often different on different machines; by using this routine different
268 machines will usually generate the same output.
269
270 Deprecated: Intentionally has no manual page
271 */
KSPMonitorResidualShort(KSP ksp,PetscInt its,PetscReal fnorm,PetscViewerAndFormat * vf)272 PetscErrorCode KSPMonitorResidualShort(KSP ksp, PetscInt its, PetscReal fnorm, PetscViewerAndFormat *vf)
273 {
274 PetscViewer viewer = vf->viewer;
275 PetscViewerFormat format = vf->format;
276 PetscInt tablevel;
277 const char *prefix;
278
279 PetscFunctionBegin;
280 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4);
281 PetscCall(PetscObjectGetTabLevel((PetscObject)ksp, &tablevel));
282 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ksp, &prefix));
283 PetscCall(PetscViewerPushFormat(viewer, format));
284 PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
285 if (its == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, " Residual norms for %s solve.\n", prefix));
286 if (fnorm > 1.e-9) PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " KSP Residual norm %g\n", its, (double)fnorm));
287 else if (fnorm > 1.e-11) PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " KSP Residual norm %5.3e\n", its, (double)fnorm));
288 else PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " KSP Residual norm < 1.e-11\n", its));
289 PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
290 PetscCall(PetscViewerPopFormat(viewer));
291 PetscFunctionReturn(PETSC_SUCCESS);
292 }
293
KSPMonitorRange_Private(KSP ksp,PetscInt it,PetscReal * per)294 PetscErrorCode KSPMonitorRange_Private(KSP ksp, PetscInt it, PetscReal *per)
295 {
296 Vec resid;
297 const PetscScalar *r;
298 PetscReal rmax, pwork;
299 PetscInt i, n, N;
300
301 PetscFunctionBegin;
302 PetscCall(KSPBuildResidual(ksp, NULL, NULL, &resid));
303 PetscCall(VecNorm(resid, NORM_INFINITY, &rmax));
304 PetscCall(VecGetLocalSize(resid, &n));
305 PetscCall(VecGetSize(resid, &N));
306 PetscCall(VecGetArrayRead(resid, &r));
307 pwork = 0.0;
308 for (i = 0; i < n; ++i) pwork += (PetscAbsScalar(r[i]) > .20 * rmax);
309 PetscCall(VecRestoreArrayRead(resid, &r));
310 PetscCall(VecDestroy(&resid));
311 PetscCallMPI(MPIU_Allreduce(&pwork, per, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)ksp)));
312 *per = *per / N;
313 PetscFunctionReturn(PETSC_SUCCESS);
314 }
315
316 /*@C
317 KSPMonitorResidualRange - Prints the percentage of residual elements that are more than 10 percent of the maximum value.
318
319 Collective
320
321 Input Parameters:
322 + ksp - iterative context
323 . it - iteration number
324 . rnorm - 2-norm (preconditioned) residual value (may be estimated).
325 - vf - The viewer context
326
327 Options Database Key:
328 . -ksp_monitor_range - Activates `KSPMonitorResidualRange()`
329
330 Level: intermediate
331
332 Note:
333 This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
334 to be used during the `KSP` solve.
335
336 .seealso: [](ch_ksp), `KSP`, `KSPMonitorSet()`, `KSPMonitorResidual()`
337 @*/
KSPMonitorResidualRange(KSP ksp,PetscInt it,PetscReal rnorm,PetscViewerAndFormat * vf)338 PetscErrorCode KSPMonitorResidualRange(KSP ksp, PetscInt it, PetscReal rnorm, PetscViewerAndFormat *vf)
339 {
340 static PetscReal prev;
341 PetscViewer viewer = vf->viewer;
342 PetscViewerFormat format = vf->format;
343 PetscInt tablevel;
344 const char *prefix;
345 PetscReal perc, rel;
346
347 PetscFunctionBegin;
348 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4);
349 PetscCall(PetscObjectGetTabLevel((PetscObject)ksp, &tablevel));
350 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ksp, &prefix));
351 PetscCall(PetscViewerPushFormat(viewer, format));
352 PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
353 if (!it) prev = rnorm;
354 if (it == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, " Residual norms for %s solve.\n", prefix));
355 PetscCall(KSPMonitorRange_Private(ksp, it, &perc));
356 rel = (prev - rnorm) / prev;
357 prev = rnorm;
358 PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " KSP preconditioned resid norm %14.12e Percent values above 20 percent of maximum %5.2f relative decrease %5.2e ratio %5.2e\n", it, (double)rnorm, (double)(100 * perc), (double)rel, (double)(rel / perc)));
359 PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
360 PetscCall(PetscViewerPopFormat(viewer));
361 PetscFunctionReturn(PETSC_SUCCESS);
362 }
363
364 /*@C
365 KSPMonitorTrueResidual - Prints the true residual norm, as well as the (possibly preconditioned, possibly approximate) residual norm,
366 at each iteration of a `KSPSolve()` iterative solver.
367
368 Collective
369
370 Input Parameters:
371 + ksp - iterative context
372 . n - iteration number
373 . rnorm - 2-norm (preconditioned) residual value (may be estimated).
374 - vf - The viewer context
375
376 Options Database Key:
377 . -ksp_monitor_true_residual - Activates `KSPMonitorTrueResidual()` to print both norm values at each iteration
378
379 Level: intermediate
380
381 Notes:
382 When using right preconditioning, the two norm values are equivalent.
383
384 This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
385 to be used during the `KSP` solve.
386
387 .seealso: [](ch_ksp), `KSP`, `KSPMonitorSet()`, `KSPMonitorResidual()`, `KSPMonitorTrueResidualMaxNorm()`, `PetscViewerAndFormat`
388 @*/
KSPMonitorTrueResidual(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat * vf)389 PetscErrorCode KSPMonitorTrueResidual(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
390 {
391 PetscViewer viewer = vf->viewer;
392 PetscViewerFormat format = vf->format;
393 Vec r;
394 PetscReal truenorm, bnorm;
395 char normtype[256];
396 PetscInt tablevel;
397 const char *prefix;
398
399 PetscFunctionBegin;
400 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4);
401 PetscCall(PetscObjectGetTabLevel((PetscObject)ksp, &tablevel));
402 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ksp, &prefix));
403 PetscCall(PetscStrncpy(normtype, KSPNormTypes[ksp->normtype], sizeof(normtype)));
404 PetscCall(PetscStrtolower(normtype));
405 PetscCall(KSPBuildResidual(ksp, NULL, NULL, &r));
406 PetscCall(VecNorm(r, NORM_2, &truenorm));
407 PetscCall(VecNorm(ksp->vec_rhs, NORM_2, &bnorm));
408 PetscCall(VecDestroy(&r));
409
410 PetscCall(PetscViewerPushFormat(viewer, format));
411 PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
412 if (n == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, " Residual norms for %s solve.\n", prefix));
413 if (bnorm == 0) {
414 PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " KSP %s resid norm %14.12e true resid norm %14.12e ||r(i)||/||b|| inf\n", n, normtype, (double)rnorm, (double)truenorm));
415 } else {
416 PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " KSP %s resid norm %14.12e true resid norm %14.12e ||r(i)||/||b|| %14.12e\n", n, normtype, (double)rnorm, (double)truenorm, (double)(truenorm / bnorm)));
417 }
418 PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
419 PetscCall(PetscViewerPopFormat(viewer));
420 PetscFunctionReturn(PETSC_SUCCESS);
421 }
422
423 /*@C
424 KSPMonitorTrueResidualView - Plots the true residual at each iteration of an iterative solver.
425
426 Collective
427
428 Input Parameters:
429 + ksp - iterative context
430 . n - iteration number
431 . rnorm - 2-norm (preconditioned) residual value (may be estimated).
432 - vf - The viewer context of type `PETSCVIEWERDRAW`
433
434 Options Database Key:
435 . -ksp_monitor_true_residual viewertype - Activates `KSPMonitorTrueResidualView()`
436
437 Level: intermediate
438
439 Note:
440 This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
441 to be used during the `KSP` solve.
442
443 .seealso: [](ch_ksp), `PETSCVIEWERDRAW`, `KSP`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `KSPMonitorResidual()`,
444 `KSPMonitorTrueResidualDrawLG()`, `PetscViewerAndFormat`
445 @*/
KSPMonitorTrueResidualView(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat * vf)446 PetscErrorCode KSPMonitorTrueResidualView(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
447 {
448 PetscViewer viewer = vf->viewer;
449 PetscViewerFormat format = vf->format;
450 Vec r;
451
452 PetscFunctionBegin;
453 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4);
454 PetscCall(PetscViewerPushFormat(viewer, format));
455 PetscCall(KSPBuildResidual(ksp, NULL, NULL, &r));
456 PetscCall(PetscObjectSetName((PetscObject)r, "True Residual"));
457 PetscCall(PetscObjectCompose((PetscObject)r, "__Vec_bc_zero__", (PetscObject)ksp));
458 PetscCall(VecView(r, viewer));
459 PetscCall(PetscObjectCompose((PetscObject)r, "__Vec_bc_zero__", NULL));
460 PetscCall(VecDestroy(&r));
461 PetscCall(PetscViewerPopFormat(viewer));
462 PetscFunctionReturn(PETSC_SUCCESS);
463 }
464
465 /*@C
466 KSPMonitorTrueResidualDrawLG - Plots the true residual norm at each iteration of an iterative solver.
467
468 Collective
469
470 Input Parameters:
471 + ksp - iterative context
472 . n - iteration number
473 . rnorm - 2-norm (preconditioned) residual value (may be estimated).
474 - vf - The viewer context
475
476 Options Database Key:
477 . -ksp_monitor_true_residual draw::draw_lg - Activates `KSPMonitorTrueResidualDrawLG()`
478
479 Level: intermediate
480
481 Notes:
482 This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
483 to be used during the `KSP` solve.
484
485 Call `KSPMonitorTrueResidualDrawLGCreate()` to create the context needed for this monitor
486
487 .seealso: [](ch_ksp), `PETSCVIEWERDRAW`, `KSP`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `KSPMonitorTrueResidualDraw()`, `KSPMonitorResidual`,
488 `KSPMonitorTrueResidualDrawLGCreate()`
489 @*/
KSPMonitorTrueResidualDrawLG(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat * vf)490 PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
491 {
492 PetscViewer viewer = vf->viewer;
493 PetscViewerFormat format = vf->format;
494 Vec r;
495 KSPConvergedReason reason;
496 PetscReal truenorm, x[2], y[2];
497 PetscDrawLG lg;
498
499 PetscFunctionBegin;
500 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4);
501 PetscCall(KSPBuildResidual(ksp, NULL, NULL, &r));
502 PetscCall(VecNorm(r, NORM_2, &truenorm));
503 PetscCall(VecDestroy(&r));
504 PetscCall(PetscViewerPushFormat(viewer, format));
505 PetscCall(PetscViewerDrawGetDrawLG(viewer, 0, &lg));
506 if (!n) PetscCall(PetscDrawLGReset(lg));
507 x[0] = (PetscReal)n;
508 if (rnorm > 0.0) y[0] = PetscLog10Real(rnorm);
509 else y[0] = -15.0;
510 x[1] = (PetscReal)n;
511 if (truenorm > 0.0) y[1] = PetscLog10Real(truenorm);
512 else y[1] = -15.0;
513 PetscCall(PetscDrawLGAddPoint(lg, x, y));
514 PetscCall(KSPGetConvergedReason(ksp, &reason));
515 if (n <= 20 || !(n % 5) || reason) {
516 PetscCall(PetscDrawLGDraw(lg));
517 PetscCall(PetscDrawLGSave(lg));
518 }
519 PetscCall(PetscViewerPopFormat(viewer));
520 PetscFunctionReturn(PETSC_SUCCESS);
521 }
522
523 /*@C
524 KSPMonitorTrueResidualDrawLGCreate - Creates the context for the true residual monitor `KSPMonitorTrueResidualDrawLG()`
525
526 Collective
527
528 Input Parameters:
529 + viewer - The `PetscViewer` of type `PETSCVIEWERDRAW`
530 . format - The viewer format
531 - ctx - An optional user context
532
533 Output Parameter:
534 . vf - The viewer context
535
536 Level: intermediate
537
538 .seealso: [](ch_ksp), `PETSCVIEWERDRAW`, `KSP`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `PetscViewerAndFormat`
539 @*/
KSPMonitorTrueResidualDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,PetscCtx ctx,PetscViewerAndFormat ** vf)540 PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer viewer, PetscViewerFormat format, PetscCtx ctx, PetscViewerAndFormat **vf)
541 {
542 const char *names[] = {"preconditioned", "true"};
543
544 PetscFunctionBegin;
545 PetscCall(PetscViewerAndFormatCreate(viewer, format, vf));
546 (*vf)->data = ctx;
547 PetscCall(PetscViewerMonitorLGSetUp(viewer, NULL, NULL, "Log Residual Norm", 2, names, PETSC_DECIDE, PETSC_DECIDE, 400, 300));
548 PetscFunctionReturn(PETSC_SUCCESS);
549 }
550
551 /*@C
552 KSPMonitorTrueResidualMax - Prints the true residual max norm at each iteration of an iterative solver.
553
554 Collective
555
556 Input Parameters:
557 + ksp - iterative context
558 . n - iteration number
559 . rnorm - 2-norm (preconditioned) residual value (may be estimated).
560 - vf - The viewer context
561
562 Options Database Key:
563 . -ksp_monitor_true_residual_max - Activates `KSPMonitorTrueResidualMax()`
564
565 Level: intermediate
566
567 Note:
568 This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
569 to be used during the `KSP` solve.
570
571 .seealso: [](ch_ksp), `KSP`, `KSPMonitorSet()`, `KSPMonitorResidual()`, `KSPMonitorTrueResidualMaxNorm()`
572 @*/
KSPMonitorTrueResidualMax(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat * vf)573 PetscErrorCode KSPMonitorTrueResidualMax(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
574 {
575 PetscViewer viewer = vf->viewer;
576 PetscViewerFormat format = vf->format;
577 Vec r;
578 PetscReal truenorm, bnorm;
579 char normtype[256];
580 PetscInt tablevel;
581 const char *prefix;
582
583 PetscFunctionBegin;
584 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4);
585 PetscCall(PetscObjectGetTabLevel((PetscObject)ksp, &tablevel));
586 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ksp, &prefix));
587 PetscCall(PetscStrncpy(normtype, KSPNormTypes[ksp->normtype], sizeof(normtype)));
588 PetscCall(PetscStrtolower(normtype));
589 PetscCall(KSPBuildResidual(ksp, NULL, NULL, &r));
590 PetscCall(VecNorm(r, NORM_INFINITY, &truenorm));
591 PetscCall(VecNorm(ksp->vec_rhs, NORM_INFINITY, &bnorm));
592 PetscCall(VecDestroy(&r));
593
594 PetscCall(PetscViewerPushFormat(viewer, format));
595 PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
596 if (n == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, " Residual norms for %s solve.\n", prefix));
597 PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " KSP %s true resid max norm %14.12e ||r(i)||/||b|| %14.12e\n", n, normtype, (double)truenorm, (double)(truenorm / bnorm)));
598 PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
599 PetscCall(PetscViewerPopFormat(viewer));
600 PetscFunctionReturn(PETSC_SUCCESS);
601 }
602
603 /*@C
604 KSPMonitorError - Prints the error norm, as well as the (possibly preconditioned) residual norm, at each iteration of an iterative solver.
605
606 Collective
607
608 Input Parameters:
609 + ksp - iterative context
610 . n - iteration number
611 . rnorm - 2-norm (preconditioned) residual value (may be estimated).
612 - vf - The viewer context
613
614 Options Database Key:
615 . -ksp_monitor_error - Activates `KSPMonitorError()`
616
617 Level: intermediate
618
619 Note:
620 This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
621 to be used during the `KSP` solve.
622
623 .seealso: [](ch_ksp), `KSP`, `KSPMonitorSet()`, `KSPMonitorResidual()`, `KSPMonitorTrueResidualMaxNorm()`
624 @*/
KSPMonitorError(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat * vf)625 PetscErrorCode KSPMonitorError(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
626 {
627 PetscViewer viewer = vf->viewer;
628 PetscViewerFormat format = vf->format;
629 DM dm;
630 Vec sol;
631 PetscReal *errors;
632 PetscInt Nf, f;
633 PetscInt tablevel;
634 const char *prefix;
635
636 PetscFunctionBegin;
637 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4);
638 PetscCall(PetscObjectGetTabLevel((PetscObject)ksp, &tablevel));
639 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ksp, &prefix));
640 PetscCall(KSPGetDM(ksp, &dm));
641 PetscCall(DMGetNumFields(dm, &Nf));
642 PetscCall(DMGetGlobalVector(dm, &sol));
643 PetscCall(KSPBuildSolution(ksp, sol, NULL));
644 /* TODO: Make a different monitor that flips sign for SNES, Newton system is A dx = -b, so we need to negate the solution */
645 PetscCall(VecScale(sol, -1.0));
646 PetscCall(PetscCalloc1(Nf, &errors));
647 PetscCall(DMComputeError(dm, sol, errors, NULL));
648
649 PetscCall(PetscViewerPushFormat(viewer, format));
650 PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
651 if (n == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, " Error norms for %s solve.\n", prefix));
652 PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " KSP Error norm %s", n, Nf > 1 ? "[" : ""));
653 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
654 for (f = 0; f < Nf; ++f) {
655 if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
656 PetscCall(PetscViewerASCIIPrintf(viewer, "%14.12e", (double)errors[f]));
657 }
658 PetscCall(PetscViewerASCIIPrintf(viewer, "%s resid norm %14.12e\n", Nf > 1 ? "]" : "", (double)rnorm));
659 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
660 PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
661 PetscCall(PetscViewerPopFormat(viewer));
662 PetscCall(DMRestoreGlobalVector(dm, &sol));
663 PetscFunctionReturn(PETSC_SUCCESS);
664 }
665
666 /*@C
667 KSPMonitorErrorDraw - Plots the error at each iteration of an iterative solver.
668
669 Collective
670
671 Input Parameters:
672 + ksp - iterative context
673 . n - iteration number
674 . rnorm - 2-norm (preconditioned) residual value (may be estimated).
675 - vf - The viewer context
676
677 Options Database Key:
678 . -ksp_monitor_error draw - Activates `KSPMonitorErrorDraw()`
679
680 Level: intermediate
681
682 Note:
683 This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
684 to be used during the `KSP` solve.
685
686 .seealso: [](ch_ksp), `PETSCVIEWERDRAW`, `KSP`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `KSPMonitorErrorDrawLG()`
687 @*/
KSPMonitorErrorDraw(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat * vf)688 PetscErrorCode KSPMonitorErrorDraw(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
689 {
690 PetscViewer viewer = vf->viewer;
691 PetscViewerFormat format = vf->format;
692 DM dm;
693 Vec sol, e;
694
695 PetscFunctionBegin;
696 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4);
697 PetscCall(PetscViewerPushFormat(viewer, format));
698 PetscCall(KSPGetDM(ksp, &dm));
699 PetscCall(DMGetGlobalVector(dm, &sol));
700 PetscCall(KSPBuildSolution(ksp, sol, NULL));
701 PetscCall(DMComputeError(dm, sol, NULL, &e));
702 PetscCall(PetscObjectSetName((PetscObject)e, "Error"));
703 PetscCall(PetscObjectCompose((PetscObject)e, "__Vec_bc_zero__", (PetscObject)ksp));
704 PetscCall(VecView(e, viewer));
705 PetscCall(PetscObjectCompose((PetscObject)e, "__Vec_bc_zero__", NULL));
706 PetscCall(VecDestroy(&e));
707 PetscCall(DMRestoreGlobalVector(dm, &sol));
708 PetscCall(PetscViewerPopFormat(viewer));
709 PetscFunctionReturn(PETSC_SUCCESS);
710 }
711
712 /*@C
713 KSPMonitorErrorDrawLG - Plots the error and residual norm at each iteration of an iterative solver.
714
715 Collective
716
717 Input Parameters:
718 + ksp - iterative context
719 . n - iteration number
720 . rnorm - 2-norm (preconditioned) residual value (may be estimated).
721 - vf - The viewer context
722
723 Options Database Key:
724 . -ksp_monitor_error draw::draw_lg - Activates `KSPMonitorTrueResidualDrawLG()`
725
726 Level: intermediate
727
728 Notes:
729 This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
730 to be used during the `KSP` solve.
731
732 Call `KSPMonitorErrorDrawLGCreate()` to create the context used with this monitor
733
734 .seealso: [](ch_ksp), `PETSCVIEWERDRAW`, `KSP`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `KSPMonitorErrorDraw()`
735 @*/
KSPMonitorErrorDrawLG(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat * vf)736 PetscErrorCode KSPMonitorErrorDrawLG(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
737 {
738 PetscViewer viewer = vf->viewer;
739 PetscViewerFormat format = vf->format;
740 PetscDrawLG lg;
741 DM dm;
742 Vec sol;
743 KSPConvergedReason reason;
744 PetscReal *x, *errors;
745 PetscInt Nf, f;
746
747 PetscFunctionBegin;
748 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4);
749 PetscCall(PetscViewerDrawGetDrawLG(viewer, 0, &lg));
750 PetscCall(KSPGetDM(ksp, &dm));
751 PetscCall(DMGetNumFields(dm, &Nf));
752 PetscCall(DMGetGlobalVector(dm, &sol));
753 PetscCall(KSPBuildSolution(ksp, sol, NULL));
754 /* TODO: Make a different monitor that flips sign for SNES, Newton system is A dx = -b, so we need to negate the solution */
755 PetscCall(VecScale(sol, -1.0));
756 PetscCall(PetscCalloc2(Nf + 1, &x, Nf + 1, &errors));
757 PetscCall(DMComputeError(dm, sol, errors, NULL));
758
759 PetscCall(PetscViewerPushFormat(viewer, format));
760 if (!n) PetscCall(PetscDrawLGReset(lg));
761 for (f = 0; f < Nf; ++f) {
762 x[f] = (PetscReal)n;
763 errors[f] = errors[f] > 0.0 ? PetscLog10Real(errors[f]) : -15.;
764 }
765 x[Nf] = (PetscReal)n;
766 errors[Nf] = rnorm > 0.0 ? PetscLog10Real(rnorm) : -15.;
767 PetscCall(PetscDrawLGAddPoint(lg, x, errors));
768 PetscCall(KSPGetConvergedReason(ksp, &reason));
769 if (n <= 20 || !(n % 5) || reason) {
770 PetscCall(PetscDrawLGDraw(lg));
771 PetscCall(PetscDrawLGSave(lg));
772 }
773 PetscCall(PetscViewerPopFormat(viewer));
774 PetscFunctionReturn(PETSC_SUCCESS);
775 }
776
777 /*@C
778 KSPMonitorErrorDrawLGCreate - Creates the context for the error and preconditioned residual plotter `KSPMonitorErrorDrawLG()`
779
780 Collective
781
782 Input Parameters:
783 + viewer - The `PetscViewer`
784 . format - The viewer format
785 - ctx - An optional user context
786
787 Output Parameter:
788 . vf - The viewer context
789
790 Level: intermediate
791
792 .seealso: [](ch_ksp), `PETSCVIEWERDRAW`, `KSP`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `KSPMonitorErrorDrawLG()`
793 @*/
KSPMonitorErrorDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,PetscCtx ctx,PetscViewerAndFormat ** vf)794 PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer viewer, PetscViewerFormat format, PetscCtx ctx, PetscViewerAndFormat **vf)
795 {
796 KSP ksp = (KSP)ctx;
797 DM dm;
798 char **names;
799 PetscInt Nf, f;
800
801 PetscFunctionBegin;
802 PetscCall(KSPGetDM(ksp, &dm));
803 PetscCall(DMGetNumFields(dm, &Nf));
804 PetscCall(PetscMalloc1(Nf + 1, &names));
805 for (f = 0; f < Nf; ++f) {
806 PetscObject disc;
807 const char *fname;
808 char lname[PETSC_MAX_PATH_LEN];
809
810 PetscCall(DMGetField(dm, f, NULL, &disc));
811 PetscCall(PetscObjectGetName(disc, &fname));
812 PetscCall(PetscStrncpy(lname, fname, PETSC_MAX_PATH_LEN));
813 PetscCall(PetscStrlcat(lname, " Error", PETSC_MAX_PATH_LEN));
814 PetscCall(PetscStrallocpy(lname, &names[f]));
815 }
816 PetscCall(PetscStrallocpy("residual", &names[Nf]));
817 PetscCall(PetscViewerAndFormatCreate(viewer, format, vf));
818 (*vf)->data = ctx;
819 PetscCall(PetscViewerMonitorLGSetUp(viewer, NULL, NULL, "Log Error Norm", Nf + 1, (const char **)names, PETSC_DECIDE, PETSC_DECIDE, 400, 300));
820 for (f = 0; f <= Nf; ++f) PetscCall(PetscFree(names[f]));
821 PetscCall(PetscFree(names));
822 PetscFunctionReturn(PETSC_SUCCESS);
823 }
824
825 /*@C
826 KSPMonitorSolution - Print the solution norm at each iteration of an iterative solver.
827
828 Collective
829
830 Input Parameters:
831 + ksp - iterative context
832 . n - iteration number
833 . rnorm - 2-norm (preconditioned) residual value (may be estimated).
834 - vf - The viewer context
835
836 Options Database Key:
837 . -ksp_monitor_solution - Activates `KSPMonitorSolution()`
838
839 Level: intermediate
840
841 Note:
842 This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
843 to be used during the `KSP` solve.
844
845 .seealso: [](ch_ksp), `KSPMonitorSet()`, `KSPMonitorTrueResidual()`
846 @*/
KSPMonitorSolution(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat * vf)847 PetscErrorCode KSPMonitorSolution(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
848 {
849 PetscViewer viewer = vf->viewer;
850 PetscViewerFormat format = vf->format;
851 Vec x;
852 PetscReal snorm;
853 PetscInt tablevel;
854 const char *prefix;
855
856 PetscFunctionBegin;
857 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4);
858 PetscCall(KSPBuildSolution(ksp, NULL, &x));
859 PetscCall(VecNorm(x, NORM_2, &snorm));
860 PetscCall(PetscObjectGetTabLevel((PetscObject)ksp, &tablevel));
861 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ksp, &prefix));
862 PetscCall(PetscViewerPushFormat(viewer, format));
863 PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
864 if (n == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, " Solution norms for %s solve.\n", prefix));
865 PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " KSP Solution norm %14.12e\n", n, (double)snorm));
866 PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
867 PetscCall(PetscViewerPopFormat(viewer));
868 PetscFunctionReturn(PETSC_SUCCESS);
869 }
870
871 /*@C
872 KSPMonitorSolutionDraw - Plots the solution at each iteration of an iterative solver.
873
874 Collective
875
876 Input Parameters:
877 + ksp - iterative context
878 . n - iteration number
879 . rnorm - 2-norm (preconditioned) residual value (may be estimated).
880 - vf - The viewer context
881
882 Options Database Key:
883 . -ksp_monitor_solution draw - Activates `KSPMonitorSolutionDraw()`
884
885 Level: intermediate
886
887 Note:
888 This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
889 to be used during the `KSP` solve.
890
891 .seealso: [](ch_ksp), `KSPMonitorSet()`, `KSPMonitorTrueResidual()`
892 @*/
KSPMonitorSolutionDraw(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat * vf)893 PetscErrorCode KSPMonitorSolutionDraw(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
894 {
895 PetscViewer viewer = vf->viewer;
896 PetscViewerFormat format = vf->format;
897 Vec x;
898
899 PetscFunctionBegin;
900 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4);
901 PetscCall(KSPBuildSolution(ksp, NULL, &x));
902 PetscCall(PetscViewerPushFormat(viewer, format));
903 PetscCall(PetscObjectSetName((PetscObject)x, "Solution"));
904 PetscCall(PetscObjectCompose((PetscObject)x, "__Vec_bc_zero__", (PetscObject)ksp));
905 PetscCall(VecView(x, viewer));
906 PetscCall(PetscObjectCompose((PetscObject)x, "__Vec_bc_zero__", NULL));
907 PetscCall(PetscViewerPopFormat(viewer));
908 PetscFunctionReturn(PETSC_SUCCESS);
909 }
910
911 /*@C
912 KSPMonitorSolutionDrawLG - Plots the solution norm at each iteration of an iterative solver.
913
914 Collective
915
916 Input Parameters:
917 + ksp - iterative context
918 . n - iteration number
919 . rnorm - 2-norm (preconditioned) residual value (may be estimated).
920 - vf - The viewer context
921
922 Options Database Key:
923 . -ksp_monitor_solution draw::draw_lg - Activates `KSPMonitorSolutionDrawLG()`
924
925 Level: intermediate
926
927 Notes:
928 This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
929 to be used during the `KSP` solve.
930
931 Call `KSPMonitorSolutionDrawLGCreate()` to create the context needed with this monitor
932
933 .seealso: [](ch_ksp), `KSP`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `KSPMonitorSolutionDrawLGCreate()`
934 @*/
KSPMonitorSolutionDrawLG(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat * vf)935 PetscErrorCode KSPMonitorSolutionDrawLG(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
936 {
937 PetscViewer viewer = vf->viewer;
938 PetscViewerFormat format = vf->format;
939 PetscDrawLG lg;
940 Vec u;
941 KSPConvergedReason reason;
942 PetscReal snorm, x, y;
943
944 PetscFunctionBegin;
945 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4);
946 PetscCall(PetscViewerDrawGetDrawLG(viewer, 0, &lg));
947 PetscCall(KSPBuildSolution(ksp, NULL, &u));
948 PetscCall(VecNorm(u, NORM_2, &snorm));
949 PetscCall(PetscViewerPushFormat(viewer, format));
950 if (!n) PetscCall(PetscDrawLGReset(lg));
951 x = (PetscReal)n;
952 if (snorm > 0.0) y = PetscLog10Real(snorm);
953 else y = -15.0;
954 PetscCall(PetscDrawLGAddPoint(lg, &x, &y));
955 PetscCall(KSPGetConvergedReason(ksp, &reason));
956 if (n <= 20 || !(n % 5) || reason) {
957 PetscCall(PetscDrawLGDraw(lg));
958 PetscCall(PetscDrawLGSave(lg));
959 }
960 PetscCall(PetscViewerPopFormat(viewer));
961 PetscFunctionReturn(PETSC_SUCCESS);
962 }
963
964 /*@C
965 KSPMonitorSolutionDrawLGCreate - Creates the context for the `KSP` monitor `KSPMonitorSolutionDrawLG()`
966
967 Collective
968
969 Input Parameters:
970 + viewer - The `PetscViewer`
971 . format - The viewer format
972 - ctx - An optional user context
973
974 Output Parameter:
975 . vf - The viewer context
976
977 Level: intermediate
978
979 Note:
980 This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
981 to be used during the `KSP` solve.
982
983 .seealso: [](ch_ksp), `KSPMonitorSet()`, `KSPMonitorTrueResidual()`
984 @*/
KSPMonitorSolutionDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,PetscCtx ctx,PetscViewerAndFormat ** vf)985 PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer viewer, PetscViewerFormat format, PetscCtx ctx, PetscViewerAndFormat **vf)
986 {
987 PetscFunctionBegin;
988 PetscCall(PetscViewerAndFormatCreate(viewer, format, vf));
989 (*vf)->data = ctx;
990 PetscCall(PetscViewerMonitorLGSetUp(viewer, NULL, NULL, "Log Solution Norm", 1, NULL, PETSC_DECIDE, PETSC_DECIDE, 400, 300));
991 PetscFunctionReturn(PETSC_SUCCESS);
992 }
993
994 /*@C
995 KSPMonitorSingularValue - Prints the two norm of the true residual and estimation of the extreme singular values of the preconditioned problem at each iteration.
996
997 Logically Collective
998
999 Input Parameters:
1000 + ksp - the iterative context
1001 . n - the iteration
1002 . rnorm - the two norm of the residual
1003 - vf - The viewer context
1004
1005 Options Database Key:
1006 . -ksp_monitor_singular_value - Activates `KSPMonitorSingularValue()`
1007
1008 Level: intermediate
1009
1010 Notes:
1011 The `KSPCG` solver uses the Lanczos technique for eigenvalue computation,
1012 while `KSPGMRES` uses the Arnoldi technique; other iterative methods do
1013 not currently compute singular values.
1014
1015 This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
1016 to be used during the `KSP` solve.
1017
1018 Call `KSPMonitorSingularValueCreate()` to create the context needed by this monitor
1019
1020 .seealso: [](ch_ksp), `KSP`, `KSPMonitorSet()`, `KSPComputeExtremeSingularValues()`, `KSPMonitorSingularValueCreate()`
1021 @*/
KSPMonitorSingularValue(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat * vf)1022 PetscErrorCode KSPMonitorSingularValue(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
1023 {
1024 PetscViewer viewer = vf->viewer;
1025 PetscViewerFormat format = vf->format;
1026 PetscReal emin, emax;
1027 PetscInt tablevel;
1028 const char *prefix;
1029
1030 PetscFunctionBegin;
1031 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
1032 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4);
1033 PetscCall(PetscObjectGetTabLevel((PetscObject)ksp, &tablevel));
1034 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ksp, &prefix));
1035 PetscCall(PetscViewerPushFormat(viewer, format));
1036 PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
1037 if (n == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, " Residual norms for %s solve.\n", prefix));
1038 if (!ksp->calc_sings) {
1039 PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " KSP Residual norm %14.12e\n", n, (double)rnorm));
1040 } else {
1041 PetscCall(KSPComputeExtremeSingularValues(ksp, &emax, &emin));
1042 PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " KSP Residual norm %14.12e %% max %14.12e min %14.12e max/min %14.12e\n", n, (double)rnorm, (double)emax, (double)emin, (double)(emax / emin)));
1043 }
1044 PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
1045 PetscCall(PetscViewerPopFormat(viewer));
1046 PetscFunctionReturn(PETSC_SUCCESS);
1047 }
1048
1049 /*@C
1050 KSPMonitorSingularValueCreate - Creates the singular value monitor context needed by `KSPMonitorSingularValue()`
1051
1052 Collective
1053
1054 Input Parameters:
1055 + viewer - The PetscViewer
1056 . format - The viewer format
1057 - ctx - An optional user context
1058
1059 Output Parameter:
1060 . vf - The viewer context
1061
1062 Level: intermediate
1063
1064 .seealso: [](ch_ksp), `KSP`, `KSPMonitorSet()`, `KSPMonitorSingularValue()`, `PetscViewer`
1065 @*/
KSPMonitorSingularValueCreate(PetscViewer viewer,PetscViewerFormat format,PetscCtx ctx,PetscViewerAndFormat ** vf)1066 PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer viewer, PetscViewerFormat format, PetscCtx ctx, PetscViewerAndFormat **vf)
1067 {
1068 KSP ksp = (KSP)ctx;
1069
1070 PetscFunctionBegin;
1071 PetscCall(PetscViewerAndFormatCreate(viewer, format, vf));
1072 (*vf)->data = ctx;
1073 PetscCall(KSPSetComputeSingularValues(ksp, PETSC_TRUE));
1074 PetscFunctionReturn(PETSC_SUCCESS);
1075 }
1076
1077 /*@C
1078 KSPMonitorDynamicToleranceCreate - Creates the context used by `KSPMonitorDynamicTolerance()`
1079
1080 Logically Collective
1081
1082 Output Parameter:
1083 . ctx - a void pointer
1084
1085 Options Database Key:
1086 . -sub_ksp_dynamic_tolerance <coef> - coefficient of dynamic tolerance for inner solver, default is 1.0
1087
1088 Level: advanced
1089
1090 Note:
1091 Use before calling `KSPMonitorSet()` with `KSPMonitorDynamicTolerance()`
1092
1093 The default coefficient for the tolerance can be changed with `KSPMonitorDynamicToleranceSetCoefficient()`
1094
1095 .seealso: [](sec_flexibleksp), `KSP`, `KSPMonitorDynamicTolerance()`, `KSPMonitorDynamicToleranceDestroy()`, `KSPMonitorDynamicToleranceSetCoefficient()`
1096 @*/
KSPMonitorDynamicToleranceCreate(PetscCtx ctx)1097 PetscErrorCode KSPMonitorDynamicToleranceCreate(PetscCtx ctx)
1098 {
1099 KSPDynTolCtx *scale;
1100
1101 PetscFunctionBegin;
1102 PetscCall(PetscMalloc1(1, &scale));
1103 scale->bnrm = -1.0;
1104 scale->coef = 1.0;
1105 *(void **)ctx = scale;
1106 PetscFunctionReturn(PETSC_SUCCESS);
1107 }
1108
1109 /*@C
1110 KSPMonitorDynamicToleranceSetCoefficient - Sets the coefficient in the context used by `KSPMonitorDynamicTolerance()`
1111
1112 Logically Collective
1113
1114 Output Parameters:
1115 + ctx - the context for `KSPMonitorDynamicTolerance()`
1116 - coeff - the coefficient, default is 1.0
1117
1118 Options Database Key:
1119 . -sub_ksp_dynamic_tolerance <coef> - coefficient of dynamic tolerance for inner solver, default is 1.0
1120
1121 Level: advanced
1122
1123 Note:
1124 Use before calling `KSPMonitorSet()` and after `KSPMonitorDynamicToleranceCreate()`
1125
1126 .seealso: [](sec_flexibleksp), `KSP`, `KSPMonitorDynamicTolerance()`, `KSPMonitorDynamicToleranceDestroy()`, `KSPMonitorDynamicToleranceCreate()`
1127 @*/
KSPMonitorDynamicToleranceSetCoefficient(PetscCtx ctx,PetscReal coeff)1128 PetscErrorCode KSPMonitorDynamicToleranceSetCoefficient(PetscCtx ctx, PetscReal coeff)
1129 {
1130 KSPDynTolCtx *scale = (KSPDynTolCtx *)ctx;
1131
1132 PetscFunctionBegin;
1133 scale->coef = coeff;
1134 PetscFunctionReturn(PETSC_SUCCESS);
1135 }
1136
1137 /*@C
1138 KSPMonitorDynamicTolerance - A monitor that changes the inner tolerance of nested preconditioners in every outer iteration in an adaptive way.
1139
1140 Collective
1141
1142 Input Parameters:
1143 + ksp - iterative context
1144 . its - iteration number (not used)
1145 . fnorm - the current residual norm
1146 - ctx - context used by monitor
1147
1148 Options Database Key:
1149 . -sub_ksp_dynamic_tolerance <coef> - coefficient of dynamic tolerance for inner solver, default is 1.0
1150
1151 Level: advanced
1152
1153 Notes:
1154 Applies for `PCKSP`, `PCBJACOBI`, and `PCDEFLATION` preconditioners
1155
1156 This may be useful for a flexible preconditioned Krylov method, such as `KSPFGMRES`, [](sec_flexibleksp) to
1157 control the accuracy of the inner solves needed to guarantee convergence of the outer iterations.
1158
1159 This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
1160 to be used during the `KSP` solve.
1161
1162 Use `KSPMonitorDynamicToleranceCreate()` and `KSPMonitorDynamicToleranceSetCoefficient()` to create the context needed by this
1163 monitor function.
1164
1165 Pass the context and `KSPMonitorDynamicToleranceDestroy()` to `KSPMonitorSet()`
1166
1167 .seealso: [](sec_flexibleksp), `KSP`, `KSPMonitorDynamicToleranceCreate()`, `KSPMonitorDynamicToleranceDestroy()`, `KSPMonitorDynamicToleranceSetCoefficient()`
1168 @*/
KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,PetscCtx ctx)1169 PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp, PetscInt its, PetscReal fnorm, PetscCtx ctx)
1170 {
1171 PC pc;
1172 PetscReal outer_rtol, outer_abstol, outer_dtol, inner_rtol;
1173 PetscInt outer_maxits, nksp, first, i;
1174 KSPDynTolCtx *scale = (KSPDynTolCtx *)ctx;
1175 KSP *subksp = NULL;
1176 KSP kspinner;
1177 PetscBool flg;
1178
1179 PetscFunctionBegin;
1180 PetscCall(KSPGetPC(ksp, &pc));
1181
1182 /* compute inner_rtol */
1183 if (scale->bnrm < 0.0) {
1184 Vec b;
1185 PetscCall(KSPGetRhs(ksp, &b));
1186 PetscCall(VecNorm(b, NORM_2, &scale->bnrm));
1187 }
1188 PetscCall(KSPGetTolerances(ksp, &outer_rtol, &outer_abstol, &outer_dtol, &outer_maxits));
1189 inner_rtol = PetscMin(scale->coef * scale->bnrm * outer_rtol / fnorm, 0.999);
1190
1191 /* if pc is ksp */
1192 PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCKSP, &flg));
1193 if (flg) {
1194 PetscCall(PCKSPGetKSP(pc, &kspinner));
1195 PetscCall(KSPSetTolerances(kspinner, inner_rtol, outer_abstol, outer_dtol, outer_maxits));
1196 PetscFunctionReturn(PETSC_SUCCESS);
1197 }
1198
1199 /* if pc is bjacobi */
1200 PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCBJACOBI, &flg));
1201 if (flg) {
1202 PetscCall(PCBJacobiGetSubKSP(pc, &nksp, &first, &subksp));
1203 if (subksp) {
1204 for (i = 0; i < nksp; i++) PetscCall(KSPSetTolerances(subksp[i], inner_rtol, outer_abstol, outer_dtol, outer_maxits));
1205 PetscFunctionReturn(PETSC_SUCCESS);
1206 }
1207 }
1208
1209 /* if pc is deflation*/
1210 PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCDEFLATION, &flg));
1211 if (flg) {
1212 PetscCall(PCDeflationGetCoarseKSP(pc, &kspinner));
1213 PetscCall(KSPSetTolerances(kspinner, inner_rtol, outer_abstol, outer_dtol, PETSC_CURRENT));
1214 PetscFunctionReturn(PETSC_SUCCESS);
1215 }
1216
1217 /* TODO: dynamic tolerance may apply to other types of pc */
1218 PetscFunctionReturn(PETSC_SUCCESS);
1219 }
1220
1221 /*@C
1222 KSPMonitorDynamicToleranceDestroy - Destroy the monitor context used in `KSPMonitorDynamicTolerance()`
1223
1224 Input Parameter:
1225 . ctx - the monitor context
1226
1227 Level: advanced
1228
1229 Note:
1230 This is not called directly but is passed to `KSPMonitorSet()` along with `KSPMonitorDynamicTolerance()`
1231
1232 .seealso: [](ch_ksp), `KSP`, `KSPMonitorDynamicTolerance()`, `KSPMonitorSet()`, `KSPMonitorDynamicToleranceCreate()`
1233 @*/
KSPMonitorDynamicToleranceDestroy(PetscCtxRt ctx)1234 PetscErrorCode KSPMonitorDynamicToleranceDestroy(PetscCtxRt ctx)
1235 {
1236 PetscFunctionBegin;
1237 PetscCall(PetscFree(*(void **)ctx));
1238 PetscFunctionReturn(PETSC_SUCCESS);
1239 }
1240
1241 /*@C
1242 KSPConvergedSkip - Convergence test that do not return as converged
1243 until the maximum number of iterations is reached.
1244
1245 Collective
1246
1247 Input Parameters:
1248 + ksp - iterative context
1249 . n - iteration number
1250 . rnorm - 2-norm residual value (may be estimated)
1251 - dtx - unused convergence context
1252
1253 Output Parameter:
1254 . reason - `KSP_CONVERGED_ITERATING` or `KSP_CONVERGED_ITS`
1255
1256 Options Database Key:
1257 . -ksp_convergence_test skip - skips the test
1258
1259 Level: advanced
1260
1261 Note:
1262 This should be used as the convergence test with the option
1263 `KSPSetNormType`(ksp,`KSP_NORM_NONE`), since norms of the residual are
1264 not computed. Convergence is then declared after the maximum number
1265 of iterations have been reached. Useful when one is using `KSPCG` or
1266 `KSPBCGS`. [](sec_flexibleksp)
1267
1268 .seealso: [](ch_ksp), `KSP`, `KSPCG`, `KSPBCGS`, `KSPConvergenceTestFn`, `KSPSetConvergenceTest()`, `KSPSetTolerances()`, `KSPSetNormType()`, [](sec_flexibleksp),
1269 `KSPConvergedReason`
1270 @*/
KSPConvergedSkip(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason * reason,PetscCtx dtx)1271 PetscErrorCode KSPConvergedSkip(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, PetscCtx dtx)
1272 {
1273 PetscFunctionBegin;
1274 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
1275 PetscAssertPointer(reason, 4);
1276 *reason = KSP_CONVERGED_ITERATING;
1277 if (n >= ksp->max_it) *reason = KSP_CONVERGED_ITS;
1278 PetscFunctionReturn(PETSC_SUCCESS);
1279 }
1280
1281 /*@
1282 KSPSetConvergedNegativeCurvature - Allows to declare convergence and return `KSP_CONVERGED_NEG_CURVE` when negative curvature is detected
1283
1284 Collective
1285
1286 Input Parameters:
1287 + ksp - iterative context
1288 - flg - the Boolean value
1289
1290 Options Database Key:
1291 . -ksp_converged_neg_curve <bool> - Declare convergence if negative curvature is detected
1292
1293 Level: advanced
1294
1295 Note:
1296 This is currently used only by a subset of the Krylov solvers, namely `KSPCG`, `KSPSTCG`, `KSPQCG`, `KSPGLTR`, `KSPNASH`, and `KSPMINRES`.
1297
1298 .seealso: [](ch_ksp), `KSP`, `KSPConvergedReason`, `KSPGetConvergedNegativeCurvature()`
1299 @*/
KSPSetConvergedNegativeCurvature(KSP ksp,PetscBool flg)1300 PetscErrorCode KSPSetConvergedNegativeCurvature(KSP ksp, PetscBool flg)
1301 {
1302 PetscFunctionBegin;
1303 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
1304 PetscValidLogicalCollectiveBool(ksp, flg, 2);
1305 ksp->converged_neg_curve = flg;
1306 PetscFunctionReturn(PETSC_SUCCESS);
1307 }
1308
1309 /*@
1310 KSPGetConvergedNegativeCurvature - Get the flag to declare convergence if negative curvature is detected
1311
1312 Collective
1313
1314 Input Parameter:
1315 . ksp - iterative context
1316
1317 Output Parameter:
1318 . flg - the Boolean value
1319
1320 Level: advanced
1321
1322 .seealso: [](ch_ksp), `KSP`, `KSPConvergedReason`, `KSPSetConvergedNegativeCurvature()`
1323 @*/
KSPGetConvergedNegativeCurvature(KSP ksp,PetscBool * flg)1324 PetscErrorCode KSPGetConvergedNegativeCurvature(KSP ksp, PetscBool *flg)
1325 {
1326 PetscFunctionBegin;
1327 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
1328 PetscAssertPointer(flg, 2);
1329 *flg = ksp->converged_neg_curve;
1330 PetscFunctionReturn(PETSC_SUCCESS);
1331 }
1332
1333 /*@C
1334 KSPConvergedDefaultCreate - Creates and initializes the context used by the `KSPConvergedDefault()` function
1335
1336 Not Collective
1337
1338 Output Parameter:
1339 . ctx - convergence context
1340
1341 Level: intermediate
1342
1343 .seealso: [](ch_ksp), `KSP`, `KSPConvergedDefault()`, `KSPConvergedDefaultDestroy()`, `KSPSetConvergenceTest()`, `KSPSetTolerances()`,
1344 `KSPConvergedSkip()`, `KSPConvergedReason`, `KSPGetConvergedReason()`, `KSPConvergedDefaultSetUIRNorm()`, `KSPConvergedDefaultSetUMIRNorm()`,
1345 `KSPConvergedDefaultSetConvergedMaxits()`
1346 @*/
KSPConvergedDefaultCreate(void ** ctx)1347 PetscErrorCode KSPConvergedDefaultCreate(void **ctx) PeNS
1348 {
1349 KSPConvergedDefaultCtx *cctx;
1350
1351 PetscFunctionBegin;
1352 PetscCall(PetscNew(&cctx));
1353 *ctx = cctx;
1354 PetscFunctionReturn(PETSC_SUCCESS);
1355 }
1356
1357 /*@
1358 KSPConvergedDefaultSetUIRNorm - makes the default convergence test use $ || B*(b - A*(initial guess))||$
1359 instead of $ || B*b ||$. In the case of right preconditioner or if `KSPSetNormType`(ksp,`KSP_NORM_UNPRECONDITIONED`)
1360 is used there is no B in the above formula.
1361
1362 Collective
1363
1364 Input Parameters:
1365 . ksp - iterative context
1366
1367 Options Database Key:
1368 . -ksp_converged_use_initial_residual_norm <bool> - Use initial residual norm for computing relative convergence
1369
1370 Level: intermediate
1371
1372 Notes:
1373 UIRNorm is short for Use Initial Residual Norm.
1374
1375 Use `KSPSetTolerances()` to alter the defaults for rtol, abstol, dtol.
1376
1377 The precise values of reason are macros such as `KSP_CONVERGED_RTOL`, which
1378 are defined in petscksp.h.
1379
1380 If the convergence test is not `KSPConvergedDefault()` then this is ignored.
1381
1382 If right preconditioning is being used then B does not appear in the above formula.
1383
1384 .seealso: [](ch_ksp), `KSP`, `KSPSetConvergenceTest()`, `KSPSetTolerances()`, `KSPConvergedSkip()`, `KSPConvergedReason`, `KSPGetConvergedReason()`, `KSPConvergedDefaultSetUMIRNorm()`, `KSPConvergedDefaultSetConvergedMaxits()`
1385 @*/
KSPConvergedDefaultSetUIRNorm(KSP ksp)1386 PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP ksp)
1387 {
1388 KSPConvergedDefaultCtx *ctx = (KSPConvergedDefaultCtx *)ksp->cnvP;
1389
1390 PetscFunctionBegin;
1391 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
1392 if (ksp->converged != KSPConvergedDefault) PetscFunctionReturn(PETSC_SUCCESS);
1393 PetscCheck(!ctx->mininitialrtol, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_WRONGSTATE, "Cannot use KSPConvergedDefaultSetUIRNorm() and KSPConvergedDefaultSetUMIRNorm() together");
1394 ctx->initialrtol = PETSC_TRUE;
1395 PetscFunctionReturn(PETSC_SUCCESS);
1396 }
1397
1398 /*@
1399 KSPConvergedDefaultSetUMIRNorm - makes the default convergence test use min(|| B*(b - A*(initial guess))||,|| B*b ||)
1400 In the case of right preconditioner or if `KSPSetNormType`(ksp,`KSP_NORM_UNPRECONDITIONED`)
1401 is used there is no B in the above formula.
1402
1403 Collective
1404
1405 Input Parameters:
1406 . ksp - iterative context
1407
1408 Options Database Key:
1409 . -ksp_converged_use_min_initial_residual_norm <bool> - Use minimum of initial residual norm and b for computing relative convergence
1410
1411 Level: intermediate
1412
1413 Notes:
1414 UMIRNorm is short for Use Minimum Initial Residual Norm.
1415
1416 Use `KSPSetTolerances()` to alter the defaults for rtol, abstol, dtol.
1417
1418 .seealso: [](ch_ksp), `KSP`, `KSPSetConvergenceTest()`, `KSPSetTolerances()`, `KSPConvergedSkip()`, `KSPConvergedReason`, `KSPGetConvergedReason()`, `KSPConvergedDefaultSetUIRNorm()`, `KSPConvergedDefaultSetConvergedMaxits()`
1419 @*/
KSPConvergedDefaultSetUMIRNorm(KSP ksp)1420 PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP ksp)
1421 {
1422 KSPConvergedDefaultCtx *ctx = (KSPConvergedDefaultCtx *)ksp->cnvP;
1423
1424 PetscFunctionBegin;
1425 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
1426 if (ksp->converged != KSPConvergedDefault) PetscFunctionReturn(PETSC_SUCCESS);
1427 PetscCheck(!ctx->initialrtol, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_WRONGSTATE, "Cannot use KSPConvergedDefaultSetUIRNorm() and KSPConvergedDefaultSetUMIRNorm() together");
1428 ctx->mininitialrtol = PETSC_TRUE;
1429 PetscFunctionReturn(PETSC_SUCCESS);
1430 }
1431
1432 /*@
1433 KSPConvergedDefaultSetConvergedMaxits - allows the default convergence test to declare convergence and return `KSP_CONVERGED_ITS` if the maximum number of iterations is reached
1434
1435 Collective
1436
1437 Input Parameters:
1438 + ksp - iterative context
1439 - flg - boolean flag
1440
1441 Options Database Key:
1442 . -ksp_converged_maxits <bool> - Declare convergence if the maximum number of iterations is reached
1443
1444 Level: intermediate
1445
1446 .seealso: [](ch_ksp), `KSP`, `KSPSetConvergenceTest()`, `KSPSetTolerances()`, `KSPConvergedSkip()`, `KSPConvergedReason`, `KSPGetConvergedReason()`, `KSPConvergedDefaultSetUMIRNorm()`, `KSPConvergedDefaultSetUIRNorm()`
1447 @*/
KSPConvergedDefaultSetConvergedMaxits(KSP ksp,PetscBool flg)1448 PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP ksp, PetscBool flg)
1449 {
1450 KSPConvergedDefaultCtx *ctx = (KSPConvergedDefaultCtx *)ksp->cnvP;
1451
1452 PetscFunctionBegin;
1453 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
1454 PetscValidLogicalCollectiveBool(ksp, flg, 2);
1455 if (ksp->converged != KSPConvergedDefault) PetscFunctionReturn(PETSC_SUCCESS);
1456 ctx->convmaxits = flg;
1457 PetscFunctionReturn(PETSC_SUCCESS);
1458 }
1459
1460 /*@C
1461 KSPConvergedDefault - Default code to determine convergence of the linear iterative solvers
1462
1463 Collective
1464
1465 Input Parameters:
1466 + ksp - iterative context
1467 . n - iteration number
1468 . rnorm - residual norm (may be estimated, depending on the method may be the preconditioned residual norm)
1469 - ctx - convergence context which must be created by `KSPConvergedDefaultCreate()`
1470
1471 Output Parameter:
1472 . reason - the convergence reason; it is positive if the iteration has converged,
1473 negative if the iteration has diverged, and `KSP_CONVERGED_ITERATING` otherwise
1474
1475 Options Database Keys:
1476 + -ksp_max_it - maximum number of linear iterations
1477 . -ksp_min_it - minimum number of linear iterations, defaults to 0
1478 . -ksp_rtol rtol - relative tolerance used in default determination of convergence, i.e. if residual norm decreases by this factor than convergence is declared
1479 . -ksp_atol abstol - absolute tolerance used in default convergence test, i.e. if residual norm is less than this then convergence is declared
1480 . -ksp_divtol tol - if residual norm increases by this factor than divergence is declared
1481 . -ksp_converged_use_initial_residual_norm - see `KSPConvergedDefaultSetUIRNorm()`
1482 . -ksp_converged_use_min_initial_residual_norm - see `KSPConvergedDefaultSetUMIRNorm()`
1483 - -ksp_converged_maxits - see `KSPConvergedDefaultSetConvergedMaxits()`
1484
1485 Level: advanced
1486
1487 Notes:
1488 `KSPConvergedDefault()` reaches convergence when rnorm < MAX (rtol * rnorm_0, abstol);
1489 Divergence is detected if rnorm > dtol * rnorm_0, or when failures are detected throughout the iteration.
1490 By default, reaching the maximum number of iterations is considered divergence (i.e. `KSP_DIVERGED_ITS`).
1491 In order to have PETSc declaring convergence in such a case (i.e. `KSP_CONVERGED_ITS`), users can use `KSPConvergedDefaultSetConvergedMaxits()`
1492
1493 where\:
1494 + `rtol` - relative tolerance,
1495 . `abstol` - absolute tolerance.
1496 . `dtol` - divergence tolerance,
1497 - `rnorm_0` - the two norm of the right-hand side (or the preconditioned norm, depending on what was set with
1498 `KSPSetNormType()`). When initial guess is non-zero you
1499 can call `KSPConvergedDefaultSetUIRNorm()` to use the norm of (b - A*(initial guess))
1500 as the starting point for relative norm convergence testing, that is as `rnorm_0`.
1501 Call `KSPConvergedDefaultSetUMIRNorm()` to use the minimum of the norm of (b - A*(initial guess)) and the norm of b as the starting point.
1502
1503 Use `KSPSetTolerances()` to alter the defaults for `rtol`, `abstol`, `dtol`.
1504
1505 Use `KSPSetNormType()` (or `-ksp_norm_type <none,preconditioned,unpreconditioned,natural>`) to change the norm used for computing rnorm
1506
1507 The precise values of reason are available in `KSPConvergedReason`
1508
1509 This routine is used by `KSP` by default so the user generally never needs call it directly.
1510
1511 Use `KSPSetConvergenceTest()` to provide your own test instead of using this one.
1512
1513 Call `KSPSetConvergenceTest()` with the `ctx`, as created above and the destruction function `KSPConvergedDefaultDestroy()`
1514
1515 .seealso: [](ch_ksp), `KSP`, `KSPSetConvergenceTest()`, `KSPSetTolerances()`, `KSPConvergedSkip()`, `KSPConvergedReason`, `KSPGetConvergedReason()`,
1516 `KSPSetMinimumIterations()`, `KSPConvergenceTestFn`,
1517 `KSPConvergedDefaultSetUIRNorm()`, `KSPConvergedDefaultSetUMIRNorm()`, `KSPConvergedDefaultSetConvergedMaxits()`, `KSPConvergedDefaultCreate()`, `KSPConvergedDefaultDestroy()`
1518 @*/
KSPConvergedDefault(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason * reason,PetscCtx ctx)1519 PetscErrorCode KSPConvergedDefault(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, PetscCtx ctx)
1520 {
1521 KSPConvergedDefaultCtx *cctx = (KSPConvergedDefaultCtx *)ctx;
1522 KSPNormType normtype;
1523
1524 PetscFunctionBegin;
1525 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
1526 PetscValidLogicalCollectiveInt(ksp, n, 2);
1527 PetscAssertPointer(reason, 4);
1528 PetscCheck(cctx, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_NULL, "Convergence context must have been created with KSPConvergedDefaultCreate()");
1529 *reason = KSP_CONVERGED_ITERATING;
1530
1531 if (cctx->convmaxits && n >= ksp->max_it) {
1532 *reason = KSP_CONVERGED_ITS;
1533 PetscCall(PetscInfo(ksp, "Linear solver has converged. Maximum number of iterations reached %" PetscInt_FMT "\n", n));
1534 PetscFunctionReturn(PETSC_SUCCESS);
1535 }
1536 PetscCall(KSPGetNormType(ksp, &normtype));
1537 if (normtype == KSP_NORM_NONE) PetscFunctionReturn(PETSC_SUCCESS);
1538
1539 if (!n) {
1540 /* if user gives initial guess need to compute norm of b */
1541 if (!ksp->guess_zero && !cctx->initialrtol) {
1542 PetscReal snorm = 0.0;
1543 if (ksp->normtype == KSP_NORM_UNPRECONDITIONED || ksp->pc_side == PC_RIGHT) {
1544 PetscCall(PetscInfo(ksp, "user has provided nonzero initial guess, computing 2-norm of RHS\n"));
1545 PetscCall(VecNorm(ksp->vec_rhs, NORM_2, &snorm)); /* <- b'*b */
1546 } else {
1547 Vec z;
1548 /* Should avoid allocating the z vector each time but cannot stash it in cctx because if KSPReset() is called the vector size might change */
1549 PetscCall(VecDuplicate(ksp->vec_rhs, &z));
1550 PetscCall(KSP_PCApply(ksp, ksp->vec_rhs, z));
1551 if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
1552 PetscCall(PetscInfo(ksp, "user has provided nonzero initial guess, computing 2-norm of preconditioned RHS\n"));
1553 PetscCall(VecNorm(z, NORM_2, &snorm)); /* dp <- b'*B'*B*b */
1554 } else if (ksp->normtype == KSP_NORM_NATURAL) {
1555 PetscScalar norm;
1556 PetscCall(PetscInfo(ksp, "user has provided nonzero initial guess, computing natural norm of RHS\n"));
1557 PetscCall(VecDot(ksp->vec_rhs, z, &norm));
1558 snorm = PetscSqrtReal(PetscAbsScalar(norm)); /* dp <- b'*B*b */
1559 }
1560 PetscCall(VecDestroy(&z));
1561 }
1562 /* handle special case of zero RHS and nonzero guess */
1563 if (!snorm) {
1564 PetscCall(PetscInfo(ksp, "Special case, user has provided nonzero initial guess and zero RHS\n"));
1565 snorm = rnorm;
1566 }
1567 if (cctx->mininitialrtol) ksp->rnorm0 = PetscMin(snorm, rnorm);
1568 else ksp->rnorm0 = snorm;
1569 } else {
1570 ksp->rnorm0 = rnorm;
1571 }
1572 ksp->ttol = PetscMax(ksp->rtol * ksp->rnorm0, ksp->abstol);
1573 }
1574
1575 if (n <= ksp->chknorm) PetscFunctionReturn(PETSC_SUCCESS);
1576
1577 if (PetscIsInfOrNanReal(rnorm)) {
1578 PCFailedReason pcreason;
1579 PetscCall(PCReduceFailedReason(ksp->pc));
1580 PetscCall(PCGetFailedReason(ksp->pc, &pcreason));
1581 if (pcreason) {
1582 *reason = KSP_DIVERGED_PC_FAILED;
1583 PetscCall(PetscInfo(ksp, "Linear solver pcsetup fails, declaring divergence \n"));
1584 } else {
1585 *reason = KSP_DIVERGED_NANORINF;
1586 PetscCall(PetscInfo(ksp, "Linear solver has created a not a number (NaN) as the residual norm, declaring divergence \n"));
1587 }
1588 PetscFunctionReturn(PETSC_SUCCESS);
1589 }
1590
1591 if (n < ksp->min_it) PetscFunctionReturn(PETSC_SUCCESS);
1592
1593 if (rnorm <= ksp->ttol) {
1594 if (rnorm < ksp->abstol) {
1595 PetscCall(PetscInfo(ksp, "Linear solver has converged. Residual norm %14.12e is less than absolute tolerance %14.12e at iteration %" PetscInt_FMT "\n", (double)rnorm, (double)ksp->abstol, n));
1596 *reason = KSP_CONVERGED_ATOL;
1597 } else {
1598 if (cctx->initialrtol) {
1599 PetscCall(PetscInfo(ksp, "Linear solver has converged. Residual norm %14.12e is less than relative tolerance %14.12e times initial residual norm %14.12e at iteration %" PetscInt_FMT "\n", (double)rnorm, (double)ksp->rtol, (double)ksp->rnorm0, n));
1600 } else {
1601 PetscCall(PetscInfo(ksp, "Linear solver has converged. Residual norm %14.12e is less than relative tolerance %14.12e times initial right-hand side norm %14.12e at iteration %" PetscInt_FMT "\n", (double)rnorm, (double)ksp->rtol, (double)ksp->rnorm0, n));
1602 }
1603 *reason = KSP_CONVERGED_RTOL;
1604 }
1605 } else if (rnorm >= ksp->divtol * ksp->rnorm0) {
1606 PetscCall(PetscInfo(ksp, "Linear solver is diverging. Initial right hand size norm %14.12e, current residual norm %14.12e at iteration %" PetscInt_FMT "\n", (double)ksp->rnorm0, (double)rnorm, n));
1607 *reason = KSP_DIVERGED_DTOL;
1608 }
1609 PetscFunctionReturn(PETSC_SUCCESS);
1610 }
1611
1612 /*@C
1613 KSPConvergedDefaultDestroy - Frees the space used by the `KSPConvergedDefault()` function context
1614
1615 Not Collective
1616
1617 Input Parameter:
1618 . ctx - convergence context
1619
1620 Level: intermediate
1621
1622 Note:
1623 Pass this function name into `KSPSetConvergenceTest()` along with the context obtained with `KSPConvergedDefaultCreate()` and `KSPConvergedDefault()`
1624
1625 .seealso: [](ch_ksp), `KSP`, `KSPConvergedDefault()`, `KSPConvergedDefaultCreate()`, `KSPSetConvergenceTest()`, `KSPSetTolerances()`, `KSPConvergedSkip()`,
1626 `KSPConvergedReason`, `KSPGetConvergedReason()`, `KSPConvergedDefaultSetUIRNorm()`, `KSPConvergedDefaultSetUMIRNorm()`
1627 @*/
KSPConvergedDefaultDestroy(PetscCtxRt ctx)1628 PetscErrorCode KSPConvergedDefaultDestroy(PetscCtxRt ctx)
1629 {
1630 KSPConvergedDefaultCtx *cctx = *(KSPConvergedDefaultCtx **)ctx;
1631
1632 PetscFunctionBegin;
1633 PetscCall(VecDestroy(&cctx->work));
1634 PetscCall(PetscFree(cctx));
1635 PetscFunctionReturn(PETSC_SUCCESS);
1636 }
1637
1638 // PetscClangLinter pragma disable: -fdoc-sowing-chars
1639 /*
1640 KSPBuildSolutionDefault - Default code to build/move the solution.
1641
1642 Collective
1643
1644 Input Parameters:
1645 + ksp - iterative context
1646 - v - pointer to the user's vector
1647
1648 Output Parameter:
1649 . V - pointer to a vector containing the solution
1650
1651 Level: advanced
1652
1653 Note:
1654 Some `KSP` methods such as `KSPGMRES` do not compute the explicit solution at each iteration, this routine takes the information
1655 they have computed during the previous iterations and uses it to compute the explicit solution
1656
1657 Developer Note:
1658 This is `PETSC_EXTERN` because it may be used by user written plugin `KSPType` implementations
1659
1660 .seealso: [](ch_ksp), `KSP`, `KSPGetSolution()`, `KSPBuildResidualDefault()`
1661 */
KSPBuildSolutionDefault(KSP ksp,Vec v,Vec * V)1662 PetscErrorCode KSPBuildSolutionDefault(KSP ksp, Vec v, Vec *V)
1663 {
1664 PetscFunctionBegin;
1665 if (ksp->pc_side == PC_RIGHT) {
1666 if (ksp->pc) {
1667 PetscCheck(v, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Not working with right preconditioner");
1668 PetscCall(KSP_PCApply(ksp, ksp->vec_sol, v));
1669 *V = v;
1670 } else {
1671 if (v) {
1672 PetscCall(VecCopy(ksp->vec_sol, v));
1673 *V = v;
1674 } else *V = ksp->vec_sol;
1675 }
1676 } else if (ksp->pc_side == PC_SYMMETRIC) {
1677 if (ksp->pc) {
1678 PetscCheck(!ksp->transpose_solve, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Not working with symmetric preconditioner and transpose solve");
1679 PetscCheck(v, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Not working with symmetric preconditioner");
1680 PetscCall(PCApplySymmetricRight(ksp->pc, ksp->vec_sol, v));
1681 *V = v;
1682 } else {
1683 if (v) {
1684 PetscCall(VecCopy(ksp->vec_sol, v));
1685 *V = v;
1686 } else *V = ksp->vec_sol;
1687 }
1688 } else {
1689 if (v) {
1690 PetscCall(VecCopy(ksp->vec_sol, v));
1691 *V = v;
1692 } else *V = ksp->vec_sol;
1693 }
1694 PetscFunctionReturn(PETSC_SUCCESS);
1695 }
1696
1697 /*@
1698 KSPBuildResidualDefault - Default code to compute the residual.
1699
1700 Collecive on ksp
1701
1702 Input Parameters:
1703 + ksp - iterative context
1704 . t - pointer to temporary vector
1705 - v - pointer to user vector
1706
1707 Output Parameter:
1708 . V - pointer to a vector containing the residual
1709
1710 Level: advanced
1711
1712 Note:
1713 Some `KSP` methods such as `KSPGMRES` do not compute the explicit residual at each iteration, this routine takes the information
1714 they have computed during the previous iterations and uses it to compute the explicit residual via the formula r = b - A*x.
1715
1716 Developer Note:
1717 This is `PETSC_EXTERN` because it may be used by user written plugin `KSPType` implementations
1718
1719 .seealso: [](ch_ksp), `KSP`, `KSPBuildSolutionDefault()`
1720 @*/
KSPBuildResidualDefault(KSP ksp,Vec t,Vec v,Vec * V)1721 PetscErrorCode KSPBuildResidualDefault(KSP ksp, Vec t, Vec v, Vec *V)
1722 {
1723 Mat Amat, Pmat;
1724
1725 PetscFunctionBegin;
1726 if (!ksp->pc) PetscCall(KSPGetPC(ksp, &ksp->pc));
1727 PetscCall(PCGetOperators(ksp->pc, &Amat, &Pmat));
1728 PetscCall(KSPBuildSolution(ksp, t, NULL));
1729 PetscCall(KSP_MatMult(ksp, Amat, t, v));
1730 PetscCall(VecAYPX(v, -1.0, ksp->vec_rhs));
1731 *V = v;
1732 PetscFunctionReturn(PETSC_SUCCESS);
1733 }
1734
1735 /*@C
1736 KSPCreateVecs - Gets a number of work vectors suitably sized for the operator in the `KSP`
1737
1738 Collective
1739
1740 Input Parameters:
1741 + ksp - iterative context
1742 . rightn - number of right work vectors to allocate
1743 - leftn - number of left work vectors to allocate
1744
1745 Output Parameters:
1746 + right - the array of vectors created
1747 - left - the array of left vectors
1748
1749 Level: advanced
1750
1751 Notes:
1752 The right vector has as many elements as the matrix has columns. The left
1753 vector has as many elements as the matrix has rows, see `MatSetSizes()` for details on the layout of the vectors.
1754
1755 The vectors are new vectors that are not owned by the `KSP`, they should be destroyed with calls to `VecDestroyVecs()` when no longer needed.
1756
1757 Developer Note:
1758 First tries to duplicate the rhs and solution vectors of the `KSP`, if they do not exist tries to get them from the matrix with `MatCreateVecs()`, if
1759 that does not exist tries to get them from the `DM` (if it is provided) with `DMCreateGlobalVectors()`.
1760
1761 .seealso: [](ch_ksp), `MatCreateVecs()`, `VecDestroyVecs()`, `KSPSetWorkVecs()`
1762 @*/
KSPCreateVecs(KSP ksp,PetscInt rightn,Vec * right[],PetscInt leftn,Vec * left[])1763 PetscErrorCode KSPCreateVecs(KSP ksp, PetscInt rightn, Vec *right[], PetscInt leftn, Vec *left[])
1764 {
1765 Vec vecr = NULL, vecl = NULL;
1766 PetscBool matset, pmatset, isshell, preferdm = PETSC_FALSE;
1767 Mat mat = NULL;
1768
1769 PetscFunctionBegin;
1770 if (ksp->dm) {
1771 PetscCall(PetscObjectTypeCompare((PetscObject)ksp->dm, DMSHELL, &isshell));
1772 preferdm = isshell ? PETSC_FALSE : PETSC_TRUE;
1773 }
1774 if (rightn) {
1775 PetscCheck(right, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_INCOMP, "You asked for right vectors but did not pass a pointer to hold them");
1776 if (ksp->vec_sol) vecr = ksp->vec_sol;
1777 else {
1778 if (preferdm) {
1779 PetscCall(DMGetGlobalVector(ksp->dm, &vecr));
1780 } else if (ksp->pc) {
1781 PetscCall(PCGetOperatorsSet(ksp->pc, &matset, &pmatset));
1782 /* check for mat before pmat because for KSPLSQR pmat may be a different size than mat since pmat maybe mat'*mat */
1783 if (matset) {
1784 PetscCall(PCGetOperators(ksp->pc, &mat, NULL));
1785 PetscCall(MatCreateVecs(mat, &vecr, NULL));
1786 } else if (pmatset) {
1787 PetscCall(PCGetOperators(ksp->pc, NULL, &mat));
1788 PetscCall(MatCreateVecs(mat, &vecr, NULL));
1789 }
1790 }
1791 if (!vecr && ksp->dm) PetscCall(DMGetGlobalVector(ksp->dm, &vecr));
1792 PetscCheck(vecr, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_WRONGSTATE, "You requested a vector from a KSP that cannot provide one");
1793 }
1794 PetscCall(VecDuplicateVecs(vecr, rightn, right));
1795 if (!ksp->vec_sol) {
1796 if (preferdm) {
1797 PetscCall(DMRestoreGlobalVector(ksp->dm, &vecr));
1798 } else if (mat) {
1799 PetscCall(VecDestroy(&vecr));
1800 } else if (ksp->dm) {
1801 PetscCall(DMRestoreGlobalVector(ksp->dm, &vecr));
1802 }
1803 }
1804 }
1805 if (leftn) {
1806 PetscCheck(left, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_INCOMP, "You asked for left vectors but did not pass a pointer to hold them");
1807 if (ksp->vec_rhs) vecl = ksp->vec_rhs;
1808 else {
1809 if (preferdm) {
1810 PetscCall(DMGetGlobalVector(ksp->dm, &vecl));
1811 } else if (ksp->pc) {
1812 PetscCall(PCGetOperatorsSet(ksp->pc, &matset, &pmatset));
1813 /* check for mat before pmat because for KSPLSQR pmat may be a different size than mat since pmat maybe mat'*mat */
1814 if (matset) {
1815 PetscCall(PCGetOperators(ksp->pc, &mat, NULL));
1816 PetscCall(MatCreateVecs(mat, NULL, &vecl));
1817 } else if (pmatset) {
1818 PetscCall(PCGetOperators(ksp->pc, NULL, &mat));
1819 PetscCall(MatCreateVecs(mat, NULL, &vecl));
1820 }
1821 }
1822 if (!vecl && ksp->dm) PetscCall(DMGetGlobalVector(ksp->dm, &vecl));
1823 PetscCheck(vecl, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_WRONGSTATE, "You requested a vector from a KSP that cannot provide one");
1824 }
1825 PetscCall(VecDuplicateVecs(vecl, leftn, left));
1826 if (!ksp->vec_rhs) {
1827 if (preferdm) {
1828 PetscCall(DMRestoreGlobalVector(ksp->dm, &vecl));
1829 } else if (mat) {
1830 PetscCall(VecDestroy(&vecl));
1831 } else if (ksp->dm) {
1832 PetscCall(DMRestoreGlobalVector(ksp->dm, &vecl));
1833 }
1834 }
1835 }
1836 PetscFunctionReturn(PETSC_SUCCESS);
1837 }
1838
1839 /*@
1840 KSPSetWorkVecs - Sets a number of work vectors into a `KSP` object
1841
1842 Collective
1843
1844 Input Parameters:
1845 + ksp - iterative context
1846 - nw - number of work vectors to allocate
1847
1848 Level: developer
1849
1850 Developer Note:
1851 This is `PETSC_EXTERN` because it may be used by user written plugin `KSPType` implementations
1852
1853 .seealso: [](ch_ksp), `KSP`, `KSPCreateVecs()`
1854 @*/
KSPSetWorkVecs(KSP ksp,PetscInt nw)1855 PetscErrorCode KSPSetWorkVecs(KSP ksp, PetscInt nw)
1856 {
1857 PetscFunctionBegin;
1858 PetscCall(VecDestroyVecs(ksp->nwork, &ksp->work));
1859 ksp->nwork = nw;
1860 PetscCall(KSPCreateVecs(ksp, nw, &ksp->work, 0, NULL));
1861 PetscFunctionReturn(PETSC_SUCCESS);
1862 }
1863
1864 // PetscClangLinter pragma disable: -fdoc-sowing-chars
1865 /*
1866 KSPDestroyDefault - Destroys a iterative context variable for methods with no separate context. Preferred calling sequence `KSPDestroy()`.
1867
1868 Input Parameter:
1869 . ksp - the iterative context
1870
1871 Level: advanced
1872
1873 Developer Note:
1874 This is `PETSC_EXTERN` because it may be used by user written plugin `KSPType` implementations
1875
1876 .seealso: [](ch_ksp), `KSP`, `KSPDestroy()`
1877 */
KSPDestroyDefault(KSP ksp)1878 PetscErrorCode KSPDestroyDefault(KSP ksp)
1879 {
1880 PetscFunctionBegin;
1881 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
1882 PetscCall(PetscFree(ksp->data));
1883 PetscFunctionReturn(PETSC_SUCCESS);
1884 }
1885
1886 /*@
1887 KSPGetConvergedReason - Gets the reason the `KSP` iteration was stopped.
1888
1889 Not Collective
1890
1891 Input Parameter:
1892 . ksp - the `KSP` context
1893
1894 Output Parameter:
1895 . reason - negative value indicates diverged, positive value converged, see `KSPConvergedReason` for the possible values
1896
1897 Options Database Key:
1898 . -ksp_converged_reason - prints the reason to standard out when the solve ends
1899
1900 Level: intermediate
1901
1902 Note:
1903 If this routine is called before or doing the `KSPSolve()` the value of `KSP_CONVERGED_ITERATING` is returned
1904
1905 .seealso: [](ch_ksp), `KSPConvergedReason`, `KSP`, `KSPSetConvergenceTest()`, `KSPConvergedDefault()`, `KSPSetTolerances()`,
1906 `KSPConvergedReasonView()`, `KSPGetConvergedReasonString()`
1907 @*/
KSPGetConvergedReason(KSP ksp,KSPConvergedReason * reason)1908 PetscErrorCode KSPGetConvergedReason(KSP ksp, KSPConvergedReason *reason)
1909 {
1910 PetscFunctionBegin;
1911 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
1912 PetscAssertPointer(reason, 2);
1913 *reason = ksp->reason;
1914 PetscFunctionReturn(PETSC_SUCCESS);
1915 }
1916
1917 /*@C
1918 KSPGetConvergedReasonString - Return a human readable string for a `KSPConvergedReason`
1919
1920 Not Collective
1921
1922 Input Parameter:
1923 . ksp - the `KSP` context
1924
1925 Output Parameter:
1926 . strreason - a human readable string that describes ksp converged reason
1927
1928 Level: beginner
1929
1930 .seealso: [](ch_ksp), `KSP`, `KSPGetConvergedReason()`
1931 @*/
KSPGetConvergedReasonString(KSP ksp,const char * strreason[])1932 PetscErrorCode KSPGetConvergedReasonString(KSP ksp, const char *strreason[])
1933 {
1934 PetscFunctionBegin;
1935 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
1936 PetscAssertPointer(strreason, 2);
1937 *strreason = KSPConvergedReasons[ksp->reason];
1938 PetscFunctionReturn(PETSC_SUCCESS);
1939 }
1940
1941 #include <petsc/private/dmimpl.h>
1942 /*@
1943 KSPSetDM - Sets the `DM` that may be used by some preconditioners and that may be used to construct the linear system
1944
1945 Logically Collective
1946
1947 Input Parameters:
1948 + ksp - the `KSP`
1949 - dm - the `DM`, cannot be `NULL` to remove a previously set `DM`
1950
1951 Level: intermediate
1952
1953 Notes:
1954 If this is used then the `KSP` will attempt to use the `DM` to create the matrix and use the routine set with
1955 `DMKSPSetComputeOperators()`. Use `KSPSetDMActive`(ksp, `KSP_DMACTIVE_OPERATOR`, `PETSC_FALSE`) to instead use the matrix you've provided with
1956 `KSPSetOperators()`.
1957
1958 A `DM` can only be used for solving one problem at a time because information about the problem is stored on the `DM`,
1959 even when not using interfaces like `DMKSPSetComputeOperators()`. Use `DMClone()` to get a distinct `DM` when solving
1960 different problems using the same function space.
1961
1962 .seealso: [](ch_ksp), `KSP`, `DM`, `KSPGetDM()`, `KSPSetDMActive()`, `KSPSetComputeOperators()`, `KSPSetComputeRHS()`, `KSPSetComputeInitialGuess()`, `DMKSPSetComputeOperators()`, `DMKSPSetComputeRHS()`, `DMKSPSetComputeInitialGuess()`
1963 @*/
KSPSetDM(KSP ksp,DM dm)1964 PetscErrorCode KSPSetDM(KSP ksp, DM dm)
1965 {
1966 PetscFunctionBegin;
1967 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
1968 PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1969 PetscCall(PetscObjectReference((PetscObject)dm));
1970 if (ksp->dm) { /* Move the DMSNES context over to the new DM unless the new DM already has one */
1971 if (ksp->dm->dmksp && !dm->dmksp) {
1972 DMKSP kdm;
1973 PetscCall(DMCopyDMKSP(ksp->dm, dm));
1974 PetscCall(DMGetDMKSP(ksp->dm, &kdm));
1975 if (kdm->originaldm == ksp->dm) kdm->originaldm = dm; /* Grant write privileges to the replacement DM */
1976 }
1977 PetscCall(DMDestroy(&ksp->dm));
1978 }
1979 ksp->dm = dm;
1980 ksp->dmAuto = PETSC_FALSE;
1981 ksp->dmActive = KSP_DMACTIVE_ALL;
1982 if (ksp->pc) PetscCall(PCSetDM(ksp->pc, dm));
1983 PetscFunctionReturn(PETSC_SUCCESS);
1984 }
1985
1986 /*@
1987 KSPSetDMActive - Indicates the `DM` should be used to generate the linear system matrix, the right-hand side vector, and the initial guess
1988
1989 Logically Collective
1990
1991 Input Parameters:
1992 + ksp - the `KSP`
1993 . active - one of `KSP_DMACTIVE_OPERATOR`, `KSP_DMACTIVE_RHS`, or `KSP_DMACTIVE_INITIAL_GUESS`
1994 - flg - use the `DM`
1995
1996 Level: intermediate
1997
1998 Note:
1999 By default `KSPSetDM()` sets the `DM` as active.
2000
2001 .seealso: [](ch_ksp), `KSP`, `DM`, `KSPGetDM()`, `KSPSetDM()`, `SNESSetDM()`, `KSPSetComputeOperators()`, `KSPSetComputeRHS()`, `KSPSetComputeInitialGuess()`
2002 @*/
KSPSetDMActive(KSP ksp,KSPDMActive active,PetscBool flg)2003 PetscErrorCode KSPSetDMActive(KSP ksp, KSPDMActive active, PetscBool flg)
2004 {
2005 PetscFunctionBegin;
2006 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
2007 PetscValidLogicalCollectiveBool(ksp, flg, 3);
2008 if (flg) ksp->dmActive = (KSPDMActive)(ksp->dmActive | active);
2009 else ksp->dmActive = (KSPDMActive)(ksp->dmActive & ~active);
2010 PetscFunctionReturn(PETSC_SUCCESS);
2011 }
2012
2013 /*@
2014 KSPGetDM - Gets the `DM` that may be used by some preconditioners and that may be used to construct the linear system
2015
2016 Not Collective
2017
2018 Input Parameter:
2019 . ksp - the `KSP`
2020
2021 Output Parameter:
2022 . dm - the `DM`
2023
2024 Level: intermediate
2025
2026 .seealso: [](ch_ksp), `KSP`, `DM`, `KSPSetDM()`, `KSPSetDMActive()`
2027 @*/
KSPGetDM(KSP ksp,DM * dm)2028 PetscErrorCode KSPGetDM(KSP ksp, DM *dm)
2029 {
2030 PetscFunctionBegin;
2031 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
2032 if (!ksp->dm) {
2033 PetscCall(DMShellCreate(PetscObjectComm((PetscObject)ksp), &ksp->dm));
2034 ksp->dmAuto = PETSC_TRUE;
2035 }
2036 *dm = ksp->dm;
2037 PetscFunctionReturn(PETSC_SUCCESS);
2038 }
2039
2040 /*@
2041 KSPSetApplicationContext - Sets the optional user-defined context for the linear solver.
2042
2043 Logically Collective
2044
2045 Input Parameters:
2046 + ksp - the `KSP` context
2047 - ctx - user context
2048
2049 Level: intermediate
2050
2051 Notes:
2052 The user context is a way for users to attach any information to the `KSP` that they may need later when interacting with the `KSP`
2053
2054 Use `KSPGetApplicationContext()` to get access to the context at a later time.
2055
2056 Fortran Note:
2057 This only works when `ctx` is a Fortran derived type (it cannot be a `PetscObject`), we recommend writing a Fortran interface definition for this
2058 function that tells the Fortran compiler the derived data type that is passed in as the `ctx` argument. See `KSPGetApplicationContext()` for
2059 an example.
2060
2061 .seealso: [](ch_ksp), `KSP`, `KSPGetApplicationContext()`
2062 @*/
KSPSetApplicationContext(KSP ksp,PetscCtx ctx)2063 PetscErrorCode KSPSetApplicationContext(KSP ksp, PetscCtx ctx)
2064 {
2065 PC pc;
2066
2067 PetscFunctionBegin;
2068 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
2069 ksp->ctx = ctx;
2070 PetscCall(KSPGetPC(ksp, &pc));
2071 PetscCall(PCSetApplicationContext(pc, ctx));
2072 PetscFunctionReturn(PETSC_SUCCESS);
2073 }
2074
2075 /*@
2076 KSPGetApplicationContext - Gets the user-defined context for the linear solver set with `KSPSetApplicationContext()`
2077
2078 Not Collective
2079
2080 Input Parameter:
2081 . ksp - `KSP` context
2082
2083 Output Parameter:
2084 . ctx - a pointer to the application context
2085
2086 Level: intermediate
2087
2088 Fortran Notes:
2089 This only works when the context is a Fortran derived type or a `PetscObject`. Define `ctx` with
2090 .vb
2091 type(tUsertype), pointer :: ctx
2092 .ve
2093
2094 .seealso: [](ch_ksp), `KSP`, `KSPSetApplicationContext()`
2095 @*/
KSPGetApplicationContext(KSP ksp,PetscCtxRt ctx)2096 PetscErrorCode KSPGetApplicationContext(KSP ksp, PetscCtxRt ctx)
2097 {
2098 PetscFunctionBegin;
2099 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
2100 *(void **)ctx = ksp->ctx;
2101 PetscFunctionReturn(PETSC_SUCCESS);
2102 }
2103
2104 #include <petsc/private/pcimpl.h>
2105
2106 /*@
2107 KSPCheckSolve - Checks if the `PCSetUp()` or `KSPSolve()` failed and set the error flag for the outer `PC`. A `KSP_DIVERGED_ITS` is
2108 not considered a failure in this context
2109
2110 Collective
2111
2112 Input Parameters:
2113 + ksp - the linear solver `KSP` context.
2114 . pc - the preconditioner context
2115 - vec - a vector that will be initialized with infinity to indicate lack of convergence
2116
2117 Level: developer
2118
2119 Note:
2120 This is called within `PCApply()` implementations to check if an error has been detected on any particular MPI processes. By initializing the vector
2121 with infinity the next call to `KSPCheckNorm()` or `KSPCheckDot()` will provide the same information to all the MPI processes that an error occurred on
2122 at least one of the processes.
2123
2124 This may be called by a subset of the processes in the `PC`.
2125
2126 Developer Note:
2127 This is used to manage returning with appropriate information from preconditioners whose inner `KSP` solvers have failed in some way
2128
2129 .seealso: [](ch_ksp), `KSP`, `KSPCreate()`, `KSPSetType()`, `KSPCheckNorm()`, `KSPCheckDot()`
2130 @*/
KSPCheckSolve(KSP ksp,PC pc,Vec vec)2131 PetscErrorCode KSPCheckSolve(KSP ksp, PC pc, Vec vec)
2132 {
2133 PCFailedReason pcreason;
2134 PC subpc;
2135
2136 PetscFunctionBegin;
2137 PetscCall(KSPGetPC(ksp, &subpc));
2138 PetscCall(PCGetFailedReason(subpc, &pcreason));
2139 PetscCall(VecFlag(vec, pcreason || (ksp->reason < 0 && ksp->reason != KSP_DIVERGED_ITS)));
2140 if (pcreason || (ksp->reason < 0 && ksp->reason != KSP_DIVERGED_ITS)) {
2141 PetscCheck(!pc->erroriffailure, PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "Detected not converged in KSP inner solve: KSP reason %s PC reason %s", KSPConvergedReasons[ksp->reason], PCFailedReasons[pcreason]);
2142 PetscCall(PetscInfo(ksp, "Detected not converged in KSP inner solve: KSP reason %s PC reason %s\n", KSPConvergedReasons[ksp->reason], PCFailedReasons[pcreason]));
2143 pc->failedreason = PC_SUBPC_ERROR;
2144 }
2145 PetscFunctionReturn(PETSC_SUCCESS);
2146 }
2147