Actual source code: iterativ.c

  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>
  9: #include <petscdmshell.h>
 10: #include <petscdraw.h>

 12: /*@
 13:   KSPGetResidualNorm - Gets the last (possibly approximate and/or preconditioned) residual norm that has been computed.

 15:   Not Collective

 17:   Input Parameter:
 18: . ksp - the iterative context

 20:   Output Parameter:
 21: . rnorm - residual norm

 23:   Level: intermediate

 25:   Notes:
 26:   For some methods, such as `KSPGMRES`, the norm is not computed directly from the residual.

 28:   The type of norm used by the method can be controlled with `KSPSetNormType()`

 30:   Certain solvers, under certain conditions, may not compute the final residual norm in an iteration, in that case the previous norm is returned.

 32: .seealso: [](ch_ksp), `KSP`, `KSPSetNormType()`, `KSPBuildResidual()`, `KSPNormType`
 33: @*/
 34: PetscErrorCode KSPGetResidualNorm(KSP ksp, PetscReal *rnorm)
 35: {
 36:   PetscFunctionBegin;
 38:   PetscAssertPointer(rnorm, 2);
 39:   *rnorm = ksp->rnorm;
 40:   PetscFunctionReturn(PETSC_SUCCESS);
 41: }

 43: /*@
 44:   KSPGetIterationNumber - Gets the current iteration number; if the `KSPSolve()` is complete, returns the number of iterations used.

 46:   Not Collective

 48:   Input Parameter:
 49: . ksp - the iterative context

 51:   Output Parameter:
 52: . its - number of iterations

 54:   Level: intermediate

 56:   Note:
 57:   During the ith iteration this returns i-1

 59: .seealso: [](ch_ksp), `KSP`, `KSPGetResidualNorm()`, `KSPBuildResidual()`, `KSPGetTotalIterations()`
 60: @*/
 61: PetscErrorCode KSPGetIterationNumber(KSP ksp, PetscInt *its)
 62: {
 63:   PetscFunctionBegin;
 65:   PetscAssertPointer(its, 2);
 66:   *its = ksp->its;
 67:   PetscFunctionReturn(PETSC_SUCCESS);
 68: }

 70: /*@
 71:   KSPGetTotalIterations - Gets the total number of iterations this `KSP` object has performed since was created, counted over all linear solves

 73:   Not Collective

 75:   Input Parameter:
 76: . ksp - the iterative context

 78:   Output Parameter:
 79: . its - total number of iterations

 81:   Level: intermediate

 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

 87: .seealso: [](ch_ksp), `KSP`, `KSPBuildResidual()`, `KSPGetResidualNorm()`, `KSPGetIterationNumber()`
 88: @*/
 89: PetscErrorCode KSPGetTotalIterations(KSP ksp, PetscInt *its)
 90: {
 91:   PetscFunctionBegin;
 93:   PetscAssertPointer(its, 2);
 94:   *its = ksp->totalits;
 95:   PetscFunctionReturn(PETSC_SUCCESS);
 96: }

 98: /*@C
 99:   KSPMonitorResidual - Print the (possibly preconditioned) residual norm at each iteration of an iterative solver.

101:   Collective

103:   Input Parameters:
104: + ksp   - iterative context
105: . n     - iteration number
106: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
107: - vf    - The viewer context

109:   Options Database Key:
110: . -ksp_monitor - Activates `KSPMonitorResidual()`

112:   Level: intermediate

114:   Note:
115:   For some methods, such as `KSPGMRES`, the norm is not computed directly from the residual.

117:   The type of norm used by the method can be controlled with `KSPSetNormType()`

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.

122: .seealso: [](ch_ksp), `KSP`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `KSPMonitorResidualDraw()`, `KSPMonitorResidualDrawLG()`,
123:           `KSPMonitorResidualRange()`, `KSPMonitorTrueResidualDraw()`, `KSPMonitorTrueResidualDrawLG()`, `KSPMonitorTrueResidualMax()`,
124:           `KSPMonitorSingularValue()`, `KSPMonitorSolutionDrawLG()`, `KSPMonitorSolutionDraw()`, `KSPMonitorSolution()`,
125:           `KSPMonitorErrorDrawLG()`, `KSPMonitorErrorDraw()`, `KSPMonitorError()`
126: @*/
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;

134:   PetscFunctionBegin;
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: }

147: /*@C
148:   KSPMonitorResidualDraw - Plots the (possibly preconditioned) residual at each iteration of an iterative solver.

150:   Collective

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

158:   Options Database Key:
159: . -ksp_monitor draw - Activates `KSPMonitorResidualDraw()`

161:   Level: intermediate

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.

167: .seealso: [](ch_ksp), `KSP`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `KSPMonitorResidual()`, `KSPMonitorResidualDrawLG()`
168: @*/
169: PetscErrorCode KSPMonitorResidualDraw(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
170: {
171:   PetscViewer       viewer = vf->viewer;
172:   PetscViewerFormat format = vf->format;
173:   Vec               r;

175:   PetscFunctionBegin;
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: }

188: /*@C
189:   KSPMonitorResidualDrawLG - Plots the (possibly preconditioned) residual norm at each iteration of an iterative solver.

191:   Collective

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

199:   Options Database Key:
200: . -ksp_monitor draw::draw_lg - Activates `KSPMonitorResidualDrawLG()`

202:   Level: intermediate

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.

208:   Use `KSPMonitorResidualDrawLGCreate()` to create the context used with this monitor

210: .seealso: [](ch_ksp), `KSP`, `PETSCVIEWERDRAW`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `KSPMonitorResidualDraw()`, `KSPMonitorResidual()`
211: @*/
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     = vf->lg;
217:   KSPConvergedReason reason;
218:   PetscReal          x, y;

220:   PetscFunctionBegin;
223:   PetscCall(PetscViewerPushFormat(viewer, format));
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: }

238: /*@C
239:   KSPMonitorResidualDrawLGCreate - Creates the context for the (possibly preconditioned) residual norm monitor `KSPMonitorResidualDrawLG()`

241:   Collective

243:   Input Parameters:
244: + viewer - The `PetscViewer` of type `PETSCVIEWERDRAW`
245: . format - The viewer format
246: - ctx    - An optional user context

248:   Output Parameter:
249: . vf - The viewer context

251:   Level: intermediate

253: .seealso: [](ch_ksp), `KSP`, `PETSCVIEWERDRAW`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `KSPMonitorResidualDrawLG()`,
254:           `PetscViewerFormat`, `PetscViewer`, `PetscViewerAndFormat`
255: @*/
256: PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer viewer, PetscViewerFormat format, void *ctx, PetscViewerAndFormat **vf)
257: {
258:   PetscFunctionBegin;
259:   PetscCall(PetscViewerAndFormatCreate(viewer, format, vf));
260:   (*vf)->data = ctx;
261:   PetscCall(KSPMonitorLGCreate(PetscObjectComm((PetscObject)viewer), NULL, NULL, "Log Residual Norm", 1, NULL, PETSC_DECIDE, PETSC_DECIDE, 400, 300, &(*vf)->lg));
262:   PetscFunctionReturn(PETSC_SUCCESS);
263: }

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.

270:   Deprecated: Intentionally has no manual page
271: */
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;

279:   PetscFunctionBegin;
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: }

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;

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:   PetscCall(MPIU_Allreduce(&pwork, per, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)ksp)));
312:   *per = *per / N;
313:   PetscFunctionReturn(PETSC_SUCCESS);
314: }

316: /*@C
317:   KSPMonitorResidualRange - Prints the percentage of residual elements that are more than 10 percent of the maximum value.

319:   Collective

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

327:   Options Database Key:
328: . -ksp_monitor_range - Activates `KSPMonitorResidualRange()`

330:   Level: intermediate

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.

336: .seealso: [](ch_ksp), `KSP`, `KSPMonitorSet()`, `KSPMonitorResidual()`
337: @*/
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;

347:   PetscFunctionBegin;
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.0 * perc), (double)rel, (double)(rel / perc)));
359:   PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
360:   PetscCall(PetscViewerPopFormat(viewer));
361:   PetscFunctionReturn(PETSC_SUCCESS);
362: }

364: /*@C
365:   KSPMonitorTrueResidual - Prints the true residual norm, as well as the (possibly preconditioned) approximate residual norm, at each iteration of an iterative solver.

367:   Collective

369:   Input Parameters:
370: + ksp   - iterative context
371: . n     - iteration number
372: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
373: - vf    - The viewer context

375:   Options Database Key:
376: . -ksp_monitor_true_residual - Activates `KSPMonitorTrueResidual()`

378:   Level: intermediate

380:   Notes:
381:   When using right preconditioning, these values are equivalent.

383:   This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
384:   to be used during the `KSP` solve.

386: .seealso: [](ch_ksp), `KSP`, `KSPMonitorSet()`, `KSPMonitorResidual()`, `KSPMonitorTrueResidualMaxNorm()`, `PetscViewerAndFormat`
387: @*/
388: PetscErrorCode KSPMonitorTrueResidual(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
389: {
390:   PetscViewer       viewer = vf->viewer;
391:   PetscViewerFormat format = vf->format;
392:   Vec               r;
393:   PetscReal         truenorm, bnorm;
394:   char              normtype[256];
395:   PetscInt          tablevel;
396:   const char       *prefix;

398:   PetscFunctionBegin;
400:   PetscCall(PetscObjectGetTabLevel((PetscObject)ksp, &tablevel));
401:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ksp, &prefix));
402:   PetscCall(PetscStrncpy(normtype, KSPNormTypes[ksp->normtype], sizeof(normtype)));
403:   PetscCall(PetscStrtolower(normtype));
404:   PetscCall(KSPBuildResidual(ksp, NULL, NULL, &r));
405:   PetscCall(VecNorm(r, NORM_2, &truenorm));
406:   PetscCall(VecNorm(ksp->vec_rhs, NORM_2, &bnorm));
407:   PetscCall(VecDestroy(&r));

409:   PetscCall(PetscViewerPushFormat(viewer, format));
410:   PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
411:   if (n == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual norms for %s solve.\n", prefix));
412:   if (bnorm == 0) {
413:     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));
414:   } else {
415:     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)));
416:   }
417:   PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
418:   PetscCall(PetscViewerPopFormat(viewer));
419:   PetscFunctionReturn(PETSC_SUCCESS);
420: }

422: /*@C
423:   KSPMonitorTrueResidualDraw - Plots the true residual at each iteration of an iterative solver.

425:   Collective

427:   Input Parameters:
428: + ksp   - iterative context
429: . n     - iteration number
430: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
431: - vf    - The viewer context of type `PETSCVIEWERDRAW`

433:   Options Database Key:
434: . -ksp_monitor_true_residual draw - Activates `KSPMonitorResidualDraw()`

436:   Level: intermediate

438:   Note:
439:   This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
440:   to be used during the `KSP` solve.

442: .seealso: [](ch_ksp), `PETSCVIEWERDRAW`, `KSP`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `KSPMonitorResidual()`,
443:           `KSPMonitorTrueResidualDrawLG()`, `PetscViewerAndFormat`
444: @*/
445: PetscErrorCode KSPMonitorTrueResidualDraw(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
446: {
447:   PetscViewer       viewer = vf->viewer;
448:   PetscViewerFormat format = vf->format;
449:   Vec               r;

451:   PetscFunctionBegin;
453:   PetscCall(PetscViewerPushFormat(viewer, format));
454:   PetscCall(KSPBuildResidual(ksp, NULL, NULL, &r));
455:   PetscCall(PetscObjectSetName((PetscObject)r, "Residual"));
456:   PetscCall(PetscObjectCompose((PetscObject)r, "__Vec_bc_zero__", (PetscObject)ksp));
457:   PetscCall(VecView(r, viewer));
458:   PetscCall(PetscObjectCompose((PetscObject)r, "__Vec_bc_zero__", NULL));
459:   PetscCall(VecDestroy(&r));
460:   PetscCall(PetscViewerPopFormat(viewer));
461:   PetscFunctionReturn(PETSC_SUCCESS);
462: }

464: /*@C
465:   KSPMonitorTrueResidualDrawLG - Plots the true residual norm at each iteration of an iterative solver.

467:   Collective

469:   Input Parameters:
470: + ksp   - iterative context
471: . n     - iteration number
472: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
473: - vf    - The viewer context

475:   Options Database Key:
476: . -ksp_monitor_true_residual draw::draw_lg - Activates `KSPMonitorTrueResidualDrawLG()`

478:   Level: intermediate

480:   Notes:
481:   This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
482:   to be used during the `KSP` solve.

484:   Call `KSPMonitorTrueResidualDrawLGCreate()` to create the context needed for this monitor

486: .seealso: [](ch_ksp), `PETSCVIEWERDRAW`, `KSP`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `KSPMonitorTrueResidualDraw()`, `KSPMonitorResidual`,
487:           `KSPMonitorTrueResidualDrawLGCreate()`
488: @*/
489: PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
490: {
491:   PetscViewer        viewer = vf->viewer;
492:   PetscViewerFormat  format = vf->format;
493:   PetscDrawLG        lg     = vf->lg;
494:   Vec                r;
495:   KSPConvergedReason reason;
496:   PetscReal          truenorm, x[2], y[2];

498:   PetscFunctionBegin;
501:   PetscCall(KSPBuildResidual(ksp, NULL, NULL, &r));
502:   PetscCall(VecNorm(r, NORM_2, &truenorm));
503:   PetscCall(VecDestroy(&r));
504:   PetscCall(PetscViewerPushFormat(viewer, format));
505:   if (!n) PetscCall(PetscDrawLGReset(lg));
506:   x[0] = (PetscReal)n;
507:   if (rnorm > 0.0) y[0] = PetscLog10Real(rnorm);
508:   else y[0] = -15.0;
509:   x[1] = (PetscReal)n;
510:   if (truenorm > 0.0) y[1] = PetscLog10Real(truenorm);
511:   else y[1] = -15.0;
512:   PetscCall(PetscDrawLGAddPoint(lg, x, y));
513:   PetscCall(KSPGetConvergedReason(ksp, &reason));
514:   if (n <= 20 || !(n % 5) || reason) {
515:     PetscCall(PetscDrawLGDraw(lg));
516:     PetscCall(PetscDrawLGSave(lg));
517:   }
518:   PetscCall(PetscViewerPopFormat(viewer));
519:   PetscFunctionReturn(PETSC_SUCCESS);
520: }

522: /*@C
523:   KSPMonitorTrueResidualDrawLGCreate - Creates the context for the true residual monitor `KSPMonitorTrueResidualDrawLG()`

525:   Collective

527:   Input Parameters:
528: + viewer - The `PetscViewer` of type `PETSCVIEWERDRAW`
529: . format - The viewer format
530: - ctx    - An optional user context

532:   Output Parameter:
533: . vf - The viewer context

535:   Level: intermediate

537: .seealso: [](ch_ksp), `PETSCVIEWERDRAW`, `KSP`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `PetscViewerAndFormat`
538: @*/
539: PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer viewer, PetscViewerFormat format, void *ctx, PetscViewerAndFormat **vf)
540: {
541:   const char *names[] = {"preconditioned", "true"};

543:   PetscFunctionBegin;
544:   PetscCall(PetscViewerAndFormatCreate(viewer, format, vf));
545:   (*vf)->data = ctx;
546:   PetscCall(KSPMonitorLGCreate(PetscObjectComm((PetscObject)viewer), NULL, NULL, "Log Residual Norm", 2, names, PETSC_DECIDE, PETSC_DECIDE, 400, 300, &(*vf)->lg));
547:   PetscFunctionReturn(PETSC_SUCCESS);
548: }

550: /*@C
551:   KSPMonitorTrueResidualMax - Prints the true residual max norm at each iteration of an iterative solver.

553:   Collective

555:   Input Parameters:
556: + ksp   - iterative context
557: . n     - iteration number
558: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
559: - vf    - The viewer context

561:   Options Database Key:
562: . -ksp_monitor_true_residual_max - Activates `KSPMonitorTrueResidualMax()`

564:   Level: intermediate

566:   Note:
567:   This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
568:   to be used during the `KSP` solve.

570: .seealso: [](ch_ksp), `KSP`, `KSPMonitorSet()`, `KSPMonitorResidual()`, `KSPMonitorTrueResidualMaxNorm()`
571: @*/
572: PetscErrorCode KSPMonitorTrueResidualMax(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
573: {
574:   PetscViewer       viewer = vf->viewer;
575:   PetscViewerFormat format = vf->format;
576:   Vec               r;
577:   PetscReal         truenorm, bnorm;
578:   char              normtype[256];
579:   PetscInt          tablevel;
580:   const char       *prefix;

582:   PetscFunctionBegin;
584:   PetscCall(PetscObjectGetTabLevel((PetscObject)ksp, &tablevel));
585:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ksp, &prefix));
586:   PetscCall(PetscStrncpy(normtype, KSPNormTypes[ksp->normtype], sizeof(normtype)));
587:   PetscCall(PetscStrtolower(normtype));
588:   PetscCall(KSPBuildResidual(ksp, NULL, NULL, &r));
589:   PetscCall(VecNorm(r, NORM_INFINITY, &truenorm));
590:   PetscCall(VecNorm(ksp->vec_rhs, NORM_INFINITY, &bnorm));
591:   PetscCall(VecDestroy(&r));

593:   PetscCall(PetscViewerPushFormat(viewer, format));
594:   PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
595:   if (n == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual norms for %s solve.\n", prefix));
596:   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)));
597:   PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
598:   PetscCall(PetscViewerPopFormat(viewer));
599:   PetscFunctionReturn(PETSC_SUCCESS);
600: }

602: /*@C
603:   KSPMonitorError - Prints the error norm, as well as the (possibly preconditioned) residual norm, at each iteration of an iterative solver.

605:   Collective

607:   Input Parameters:
608: + ksp   - iterative context
609: . n     - iteration number
610: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
611: - vf    - The viewer context

613:   Options Database Key:
614: . -ksp_monitor_error - Activates `KSPMonitorError()`

616:   Level: intermediate

618:   Note:
619:   This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
620:   to be used during the `KSP` solve.

622: .seealso: [](ch_ksp), `KSP`, `KSPMonitorSet()`, `KSPMonitorResidual()`, `KSPMonitorTrueResidualMaxNorm()`
623: @*/
624: PetscErrorCode KSPMonitorError(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
625: {
626:   PetscViewer       viewer = vf->viewer;
627:   PetscViewerFormat format = vf->format;
628:   DM                dm;
629:   Vec               sol;
630:   PetscReal        *errors;
631:   PetscInt          Nf, f;
632:   PetscInt          tablevel;
633:   const char       *prefix;

635:   PetscFunctionBegin;
637:   PetscCall(PetscObjectGetTabLevel((PetscObject)ksp, &tablevel));
638:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ksp, &prefix));
639:   PetscCall(KSPGetDM(ksp, &dm));
640:   PetscCall(DMGetNumFields(dm, &Nf));
641:   PetscCall(DMGetGlobalVector(dm, &sol));
642:   PetscCall(KSPBuildSolution(ksp, sol, NULL));
643:   /* TODO: Make a different monitor that flips sign for SNES, Newton system is A dx = -b, so we need to negate the solution */
644:   PetscCall(VecScale(sol, -1.0));
645:   PetscCall(PetscCalloc1(Nf, &errors));
646:   PetscCall(DMComputeError(dm, sol, errors, NULL));

648:   PetscCall(PetscViewerPushFormat(viewer, format));
649:   PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
650:   if (n == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, "  Error norms for %s solve.\n", prefix));
651:   PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " KSP Error norm %s", n, Nf > 1 ? "[" : ""));
652:   PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
653:   for (f = 0; f < Nf; ++f) {
654:     if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
655:     PetscCall(PetscViewerASCIIPrintf(viewer, "%14.12e", (double)errors[f]));
656:   }
657:   PetscCall(PetscViewerASCIIPrintf(viewer, "%s resid norm %14.12e\n", Nf > 1 ? "]" : "", (double)rnorm));
658:   PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
659:   PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
660:   PetscCall(PetscViewerPopFormat(viewer));
661:   PetscCall(DMRestoreGlobalVector(dm, &sol));
662:   PetscFunctionReturn(PETSC_SUCCESS);
663: }

665: /*@C
666:   KSPMonitorErrorDraw - Plots the error at each iteration of an iterative solver.

668:   Collective

670:   Input Parameters:
671: + ksp   - iterative context
672: . n     - iteration number
673: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
674: - vf    - The viewer context

676:   Options Database Key:
677: . -ksp_monitor_error draw - Activates `KSPMonitorErrorDraw()`

679:   Level: intermediate

681:   Note:
682:   This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
683:   to be used during the `KSP` solve.

685: .seealso: [](ch_ksp), `PETSCVIEWERDRAW`, `KSP`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `KSPMonitorErrorDrawLG()`
686: @*/
687: PetscErrorCode KSPMonitorErrorDraw(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
688: {
689:   PetscViewer       viewer = vf->viewer;
690:   PetscViewerFormat format = vf->format;
691:   DM                dm;
692:   Vec               sol, e;

694:   PetscFunctionBegin;
696:   PetscCall(PetscViewerPushFormat(viewer, format));
697:   PetscCall(KSPGetDM(ksp, &dm));
698:   PetscCall(DMGetGlobalVector(dm, &sol));
699:   PetscCall(KSPBuildSolution(ksp, sol, NULL));
700:   PetscCall(DMComputeError(dm, sol, NULL, &e));
701:   PetscCall(PetscObjectSetName((PetscObject)e, "Error"));
702:   PetscCall(PetscObjectCompose((PetscObject)e, "__Vec_bc_zero__", (PetscObject)ksp));
703:   PetscCall(VecView(e, viewer));
704:   PetscCall(PetscObjectCompose((PetscObject)e, "__Vec_bc_zero__", NULL));
705:   PetscCall(VecDestroy(&e));
706:   PetscCall(DMRestoreGlobalVector(dm, &sol));
707:   PetscCall(PetscViewerPopFormat(viewer));
708:   PetscFunctionReturn(PETSC_SUCCESS);
709: }

711: /*@C
712:   KSPMonitorErrorDrawLG - Plots the error and residual norm at each iteration of an iterative solver.

714:   Collective

716:   Input Parameters:
717: + ksp   - iterative context
718: . n     - iteration number
719: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
720: - vf    - The viewer context

722:   Options Database Key:
723: . -ksp_monitor_error draw::draw_lg - Activates `KSPMonitorTrueResidualDrawLG()`

725:   Level: intermediate

727:   Notes:
728:   This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
729:   to be used during the `KSP` solve.

731:   Call `KSPMonitorErrorDrawLGCreate()` to create the context used with this monitor

733: .seealso: [](ch_ksp), `PETSCVIEWERDRAW`, `KSP`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `KSPMonitorErrorDraw()`
734: @*/
735: PetscErrorCode KSPMonitorErrorDrawLG(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
736: {
737:   PetscViewer        viewer = vf->viewer;
738:   PetscViewerFormat  format = vf->format;
739:   PetscDrawLG        lg     = vf->lg;
740:   DM                 dm;
741:   Vec                sol;
742:   KSPConvergedReason reason;
743:   PetscReal         *x, *errors;
744:   PetscInt           Nf, f;

746:   PetscFunctionBegin;
749:   PetscCall(KSPGetDM(ksp, &dm));
750:   PetscCall(DMGetNumFields(dm, &Nf));
751:   PetscCall(DMGetGlobalVector(dm, &sol));
752:   PetscCall(KSPBuildSolution(ksp, sol, NULL));
753:   /* TODO: Make a different monitor that flips sign for SNES, Newton system is A dx = -b, so we need to negate the solution */
754:   PetscCall(VecScale(sol, -1.0));
755:   PetscCall(PetscCalloc2(Nf + 1, &x, Nf + 1, &errors));
756:   PetscCall(DMComputeError(dm, sol, errors, NULL));

758:   PetscCall(PetscViewerPushFormat(viewer, format));
759:   if (!n) PetscCall(PetscDrawLGReset(lg));
760:   for (f = 0; f < Nf; ++f) {
761:     x[f]      = (PetscReal)n;
762:     errors[f] = errors[f] > 0.0 ? PetscLog10Real(errors[f]) : -15.;
763:   }
764:   x[Nf]      = (PetscReal)n;
765:   errors[Nf] = rnorm > 0.0 ? PetscLog10Real(rnorm) : -15.;
766:   PetscCall(PetscDrawLGAddPoint(lg, x, errors));
767:   PetscCall(KSPGetConvergedReason(ksp, &reason));
768:   if (n <= 20 || !(n % 5) || reason) {
769:     PetscCall(PetscDrawLGDraw(lg));
770:     PetscCall(PetscDrawLGSave(lg));
771:   }
772:   PetscCall(PetscViewerPopFormat(viewer));
773:   PetscFunctionReturn(PETSC_SUCCESS);
774: }

776: /*@C
777:   KSPMonitorErrorDrawLGCreate - Creates the context for the error and preconditioned residual plotter `KSPMonitorErrorDrawLG()`

779:   Collective

781:   Input Parameters:
782: + viewer - The `PetscViewer`
783: . format - The viewer format
784: - ctx    - An optional user context

786:   Output Parameter:
787: . vf - The viewer context

789:   Level: intermediate

791: .seealso: [](ch_ksp), `PETSCVIEWERDRAW`, `KSP`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `KSPMonitorErrorDrawLG()`
792: @*/
793: PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer viewer, PetscViewerFormat format, void *ctx, PetscViewerAndFormat **vf)
794: {
795:   KSP      ksp = (KSP)ctx;
796:   DM       dm;
797:   char   **names;
798:   PetscInt Nf, f;

800:   PetscFunctionBegin;
801:   PetscCall(KSPGetDM(ksp, &dm));
802:   PetscCall(DMGetNumFields(dm, &Nf));
803:   PetscCall(PetscMalloc1(Nf + 1, &names));
804:   for (f = 0; f < Nf; ++f) {
805:     PetscObject disc;
806:     const char *fname;
807:     char        lname[PETSC_MAX_PATH_LEN];

809:     PetscCall(DMGetField(dm, f, NULL, &disc));
810:     PetscCall(PetscObjectGetName(disc, &fname));
811:     PetscCall(PetscStrncpy(lname, fname, PETSC_MAX_PATH_LEN));
812:     PetscCall(PetscStrlcat(lname, " Error", PETSC_MAX_PATH_LEN));
813:     PetscCall(PetscStrallocpy(lname, &names[f]));
814:   }
815:   PetscCall(PetscStrallocpy("residual", &names[Nf]));
816:   PetscCall(PetscViewerAndFormatCreate(viewer, format, vf));
817:   (*vf)->data = ctx;
818:   PetscCall(KSPMonitorLGCreate(PetscObjectComm((PetscObject)viewer), NULL, NULL, "Log Error Norm", Nf + 1, (const char **)names, PETSC_DECIDE, PETSC_DECIDE, 400, 300, &(*vf)->lg));
819:   for (f = 0; f <= Nf; ++f) PetscCall(PetscFree(names[f]));
820:   PetscCall(PetscFree(names));
821:   PetscFunctionReturn(PETSC_SUCCESS);
822: }

824: /*@C
825:   KSPMonitorSolution - Print the solution norm at each iteration of an iterative solver.

827:   Collective

829:   Input Parameters:
830: + ksp   - iterative context
831: . n     - iteration number
832: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
833: - vf    - The viewer context

835:   Options Database Key:
836: . -ksp_monitor_solution - Activates `KSPMonitorSolution()`

838:   Level: intermediate

840:   Note:
841:   This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
842:   to be used during the `KSP` solve.

844: .seealso: [](ch_ksp), `KSPMonitorSet()`, `KSPMonitorTrueResidual()`
845: @*/
846: PetscErrorCode KSPMonitorSolution(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
847: {
848:   PetscViewer       viewer = vf->viewer;
849:   PetscViewerFormat format = vf->format;
850:   Vec               x;
851:   PetscReal         snorm;
852:   PetscInt          tablevel;
853:   const char       *prefix;

855:   PetscFunctionBegin;
857:   PetscCall(KSPBuildSolution(ksp, NULL, &x));
858:   PetscCall(VecNorm(x, NORM_2, &snorm));
859:   PetscCall(PetscObjectGetTabLevel((PetscObject)ksp, &tablevel));
860:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ksp, &prefix));
861:   PetscCall(PetscViewerPushFormat(viewer, format));
862:   PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
863:   if (n == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, "  Solution norms for %s solve.\n", prefix));
864:   PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " KSP Solution norm %14.12e\n", n, (double)snorm));
865:   PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
866:   PetscCall(PetscViewerPopFormat(viewer));
867:   PetscFunctionReturn(PETSC_SUCCESS);
868: }

870: /*@C
871:   KSPMonitorSolutionDraw - Plots the solution at each iteration of an iterative solver.

873:   Collective

875:   Input Parameters:
876: + ksp   - iterative context
877: . n     - iteration number
878: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
879: - vf    - The viewer context

881:   Options Database Key:
882: . -ksp_monitor_solution draw - Activates `KSPMonitorSolutionDraw()`

884:   Level: intermediate

886:   Note:
887:   This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
888:   to be used during the `KSP` solve.

890: .seealso: [](ch_ksp), `KSPMonitorSet()`, `KSPMonitorTrueResidual()`
891: @*/
892: PetscErrorCode KSPMonitorSolutionDraw(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
893: {
894:   PetscViewer       viewer = vf->viewer;
895:   PetscViewerFormat format = vf->format;
896:   Vec               x;

898:   PetscFunctionBegin;
900:   PetscCall(KSPBuildSolution(ksp, NULL, &x));
901:   PetscCall(PetscViewerPushFormat(viewer, format));
902:   PetscCall(PetscObjectSetName((PetscObject)x, "Solution"));
903:   PetscCall(PetscObjectCompose((PetscObject)x, "__Vec_bc_zero__", (PetscObject)ksp));
904:   PetscCall(VecView(x, viewer));
905:   PetscCall(PetscObjectCompose((PetscObject)x, "__Vec_bc_zero__", NULL));
906:   PetscCall(PetscViewerPopFormat(viewer));
907:   PetscFunctionReturn(PETSC_SUCCESS);
908: }

910: /*@C
911:   KSPMonitorSolutionDrawLG - Plots the solution norm at each iteration of an iterative solver.

913:   Collective

915:   Input Parameters:
916: + ksp   - iterative context
917: . n     - iteration number
918: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
919: - vf    - The viewer context

921:   Options Database Key:
922: . -ksp_monitor_solution draw::draw_lg - Activates `KSPMonitorSolutionDrawLG()`

924:   Level: intermediate

926:   Notes:
927:   This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
928:   to be used during the `KSP` solve.

930:   Call `KSPMonitorSolutionDrawLGCreate()` to create the context needed with this monitor

932: .seealso: [](ch_ksp), `KSP`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `KSPMonitorSolutionDrawLGCreate()`
933: @*/
934: PetscErrorCode KSPMonitorSolutionDrawLG(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
935: {
936:   PetscViewer        viewer = vf->viewer;
937:   PetscViewerFormat  format = vf->format;
938:   PetscDrawLG        lg     = vf->lg;
939:   Vec                u;
940:   KSPConvergedReason reason;
941:   PetscReal          snorm, x, y;

943:   PetscFunctionBegin;
946:   PetscCall(KSPBuildSolution(ksp, NULL, &u));
947:   PetscCall(VecNorm(u, NORM_2, &snorm));
948:   PetscCall(PetscViewerPushFormat(viewer, format));
949:   if (!n) PetscCall(PetscDrawLGReset(lg));
950:   x = (PetscReal)n;
951:   if (snorm > 0.0) y = PetscLog10Real(snorm);
952:   else y = -15.0;
953:   PetscCall(PetscDrawLGAddPoint(lg, &x, &y));
954:   PetscCall(KSPGetConvergedReason(ksp, &reason));
955:   if (n <= 20 || !(n % 5) || reason) {
956:     PetscCall(PetscDrawLGDraw(lg));
957:     PetscCall(PetscDrawLGSave(lg));
958:   }
959:   PetscCall(PetscViewerPopFormat(viewer));
960:   PetscFunctionReturn(PETSC_SUCCESS);
961: }

963: /*@C
964:   KSPMonitorSolutionDrawLGCreate - Creates the context for the `KSP` monitor `KSPMonitorSolutionDrawLG()`

966:   Collective

968:   Input Parameters:
969: + viewer - The `PetscViewer`
970: . format - The viewer format
971: - ctx    - An optional user context

973:   Output Parameter:
974: . vf - The viewer context

976:   Level: intermediate

978:   Note:
979:   This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
980:   to be used during the `KSP` solve.

982: .seealso: [](ch_ksp), `KSPMonitorSet()`, `KSPMonitorTrueResidual()`
983: @*/
984: PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer viewer, PetscViewerFormat format, void *ctx, PetscViewerAndFormat **vf)
985: {
986:   PetscFunctionBegin;
987:   PetscCall(PetscViewerAndFormatCreate(viewer, format, vf));
988:   (*vf)->data = ctx;
989:   PetscCall(KSPMonitorLGCreate(PetscObjectComm((PetscObject)viewer), NULL, NULL, "Log Solution Norm", 1, NULL, PETSC_DECIDE, PETSC_DECIDE, 400, 300, &(*vf)->lg));
990:   PetscFunctionReturn(PETSC_SUCCESS);
991: }

993: /*@C
994:   KSPMonitorSingularValue - Prints the two norm of the true residual and estimation of the extreme singular values of the preconditioned problem at each iteration.

996:   Logically Collective

998:   Input Parameters:
999: + ksp   - the iterative context
1000: . n     - the iteration
1001: . rnorm - the two norm of the residual
1002: - vf    - The viewer context

1004:   Options Database Key:
1005: . -ksp_monitor_singular_value - Activates `KSPMonitorSingularValue()`

1007:   Level: intermediate

1009:   Notes:
1010:   The `KSPCG` solver uses the Lanczos technique for eigenvalue computation,
1011:   while `KSPGMRES` uses the Arnoldi technique; other iterative methods do
1012:   not currently compute singular values.

1014:   This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
1015:   to be used during the `KSP` solve.

1017:   Call `KSPMonitorSingularValueCreate()` to create the context needed by this monitor

1019: .seealso: [](ch_ksp), `KSP`, `KSPMonitorSet()`, `KSPComputeExtremeSingularValues()`, `KSPMonitorSingularValueCreate()`
1020: @*/
1021: PetscErrorCode KSPMonitorSingularValue(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
1022: {
1023:   PetscViewer       viewer = vf->viewer;
1024:   PetscViewerFormat format = vf->format;
1025:   PetscReal         emin, emax;
1026:   PetscInt          tablevel;
1027:   const char       *prefix;

1029:   PetscFunctionBegin;
1032:   PetscCall(PetscObjectGetTabLevel((PetscObject)ksp, &tablevel));
1033:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ksp, &prefix));
1034:   PetscCall(PetscViewerPushFormat(viewer, format));
1035:   PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
1036:   if (n == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual norms for %s solve.\n", prefix));
1037:   if (!ksp->calc_sings) {
1038:     PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " KSP Residual norm %14.12e\n", n, (double)rnorm));
1039:   } else {
1040:     PetscCall(KSPComputeExtremeSingularValues(ksp, &emax, &emin));
1041:     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)));
1042:   }
1043:   PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
1044:   PetscCall(PetscViewerPopFormat(viewer));
1045:   PetscFunctionReturn(PETSC_SUCCESS);
1046: }

1048: /*@C
1049:   KSPMonitorSingularValueCreate - Creates the singular value monitor context needed by `KSPMonitorSingularValue()`

1051:   Collective

1053:   Input Parameters:
1054: + viewer - The PetscViewer
1055: . format - The viewer format
1056: - ctx    - An optional user context

1058:   Output Parameter:
1059: . vf - The viewer context

1061:   Level: intermediate

1063: .seealso: [](ch_ksp), `KSP`, `KSPMonitorSet()`, `KSPMonitorSingularValue()`, `PetscViewer`
1064: @*/
1065: PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer viewer, PetscViewerFormat format, void *ctx, PetscViewerAndFormat **vf)
1066: {
1067:   KSP ksp = (KSP)ctx;

1069:   PetscFunctionBegin;
1070:   PetscCall(PetscViewerAndFormatCreate(viewer, format, vf));
1071:   (*vf)->data = ctx;
1072:   PetscCall(KSPSetComputeSingularValues(ksp, PETSC_TRUE));
1073:   PetscFunctionReturn(PETSC_SUCCESS);
1074: }

1076: /*@C
1077:   KSPMonitorDynamicToleranceCreate - Creates the context used by `KSPMonitorDynamicTolerance()`

1079:   Logically Collective

1081:   Output Parameter:
1082: . ctx - a void pointer

1084:   Options Database Key:
1085: . -sub_ksp_dynamic_tolerance <coef> - coefficient of dynamic tolerance for inner solver, default is 1.0

1087:   Level: advanced

1089:   Note:
1090:   Use before calling `KSPMonitorSet()` with `KSPMonitorDynamicTolerance()`

1092:   The default coefficient for the tolerance can be changed with `KSPMonitorDynamicToleranceSetCoefficient()`

1094: .seealso: [](sec_flexibleksp), `KSP`, `KSPMonitorDynamicTolerance()`, `KSPMonitorDynamicToleranceDestroy()`, `KSPMonitorDynamicToleranceSetCoefficient()`
1095: @*/
1096: PetscErrorCode KSPMonitorDynamicToleranceCreate(void *ctx)
1097: {
1098:   KSPDynTolCtx *scale;

1100:   PetscFunctionBegin;
1101:   PetscCall(PetscMalloc1(1, &scale));
1102:   scale->bnrm   = -1.0;
1103:   scale->coef   = 1.0;
1104:   *(void **)ctx = scale;
1105:   PetscFunctionReturn(PETSC_SUCCESS);
1106: }

1108: /*@C
1109:   KSPMonitorDynamicToleranceSetCoefficient - Sets the coefficient in the context used by `KSPMonitorDynamicTolerance()`

1111:   Logically Collective

1113:   Output Parameters:
1114: + ctx   - the context for `KSPMonitorDynamicTolerance()`
1115: - coeff - the coefficient, default is 1.0

1117:   Options Database Key:
1118: . -sub_ksp_dynamic_tolerance <coef> - coefficient of dynamic tolerance for inner solver, default is 1.0

1120:   Level: advanced

1122:   Note:
1123:   Use before calling `KSPMonitorSet()` and after `KSPMonitorDynamicToleranceCreate()`

1125: .seealso: [](sec_flexibleksp), `KSP`, `KSPMonitorDynamicTolerance()`, `KSPMonitorDynamicToleranceDestroy()`, `KSPMonitorDynamicToleranceCreate()`
1126: @*/
1127: PetscErrorCode KSPMonitorDynamicToleranceSetCoefficient(void *ctx, PetscReal coeff)
1128: {
1129:   KSPDynTolCtx *scale = (KSPDynTolCtx *)ctx;

1131:   PetscFunctionBegin;
1132:   scale->coef = coeff;
1133:   PetscFunctionReturn(PETSC_SUCCESS);
1134: }

1136: /*@C
1137:   KSPMonitorDynamicTolerance - A monitor that changes the inner tolerance of nested preconditioners in every outer iteration in an adaptive way.

1139:   Collective

1141:   Input Parameters:
1142: + ksp   - iterative context
1143: . its   - iteration number (not used)
1144: . fnorm - the current residual norm
1145: - ctx   - context used by monitor

1147:   Options Database Key:
1148: . -sub_ksp_dynamic_tolerance <coef> - coefficient of dynamic tolerance for inner solver, default is 1.0

1150:   Level: advanced

1152:   Notes:
1153:   Applies for `PCKSP`, `PCBJACOBI`, and `PCDEFLATION` preconditioners

1155:   This may be useful for a flexible preconditioned Krylov method, such as `KSPFGMRES`, [](sec_flexibleksp) to
1156:   control the accuracy of the inner solves needed to guarantee convergence of the outer iterations.

1158:   This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
1159:   to be used during the `KSP` solve.

1161:   Use `KSPMonitorDynamicToleranceCreate()` and `KSPMonitorDynamicToleranceSetCoefficient()` to create the context needed by this
1162:   monitor function.

1164:   Pass the context and `KSPMonitorDynamicToleranceDestroy()` to `KSPMonitorSet()`

1166: .seealso: [](sec_flexibleksp), `KSP`, `KSPMonitorDynamicToleranceCreate()`, `KSPMonitorDynamicToleranceDestroy()`, `KSPMonitorDynamicToleranceSetCoefficient()`
1167: @*/
1168: PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp, PetscInt its, PetscReal fnorm, void *ctx)
1169: {
1170:   PC            pc;
1171:   PetscReal     outer_rtol, outer_abstol, outer_dtol, inner_rtol;
1172:   PetscInt      outer_maxits, nksp, first, i;
1173:   KSPDynTolCtx *scale  = (KSPDynTolCtx *)ctx;
1174:   KSP          *subksp = NULL;
1175:   KSP           kspinner;
1176:   PetscBool     flg;

1178:   PetscFunctionBegin;
1179:   PetscCall(KSPGetPC(ksp, &pc));

1181:   /* compute inner_rtol */
1182:   if (scale->bnrm < 0.0) {
1183:     Vec b;
1184:     PetscCall(KSPGetRhs(ksp, &b));
1185:     PetscCall(VecNorm(b, NORM_2, &scale->bnrm));
1186:   }
1187:   PetscCall(KSPGetTolerances(ksp, &outer_rtol, &outer_abstol, &outer_dtol, &outer_maxits));
1188:   inner_rtol = PetscMin(scale->coef * scale->bnrm * outer_rtol / fnorm, 0.999);

1190:   /* if pc is ksp */
1191:   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCKSP, &flg));
1192:   if (flg) {
1193:     PetscCall(PCKSPGetKSP(pc, &kspinner));
1194:     PetscCall(KSPSetTolerances(kspinner, inner_rtol, outer_abstol, outer_dtol, outer_maxits));
1195:     PetscFunctionReturn(PETSC_SUCCESS);
1196:   }

1198:   /* if pc is bjacobi */
1199:   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCBJACOBI, &flg));
1200:   if (flg) {
1201:     PetscCall(PCBJacobiGetSubKSP(pc, &nksp, &first, &subksp));
1202:     if (subksp) {
1203:       for (i = 0; i < nksp; i++) PetscCall(KSPSetTolerances(subksp[i], inner_rtol, outer_abstol, outer_dtol, outer_maxits));
1204:       PetscFunctionReturn(PETSC_SUCCESS);
1205:     }
1206:   }

1208:   /* if pc is deflation*/
1209:   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCDEFLATION, &flg));
1210:   if (flg) {
1211:     PetscCall(PCDeflationGetCoarseKSP(pc, &kspinner));
1212:     PetscCall(KSPSetTolerances(kspinner, inner_rtol, outer_abstol, outer_dtol, PETSC_DEFAULT));
1213:     PetscFunctionReturn(PETSC_SUCCESS);
1214:   }

1216:   /* todo: dynamic tolerance may apply to other types of pc */
1217:   PetscFunctionReturn(PETSC_SUCCESS);
1218: }

1220: /*@C
1221:   KSPMonitorDynamicToleranceDestroy - Destroy the monitor context used in `KSPMonitorDynamicTolerance()`

1223:   Input Parameter:
1224: . ctx - the monitor context

1226:   Level: advanced

1228:   Note:
1229:   This is not called directly but is passed to `KSPMonitorSet()` along with `KSPMonitorDynamicTolerance()`

1231: .seealso: [](ch_ksp), `KSP`, `KSPMonitorDynamicTolerance()`, `KSPMonitorSet()`, `KSPMonitorDynamicToleranceCreate()`
1232: @*/
1233: PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **ctx)
1234: {
1235:   PetscFunctionBegin;
1236:   PetscCall(PetscFree(*ctx));
1237:   PetscFunctionReturn(PETSC_SUCCESS);
1238: }

1240: /*@C
1241:   KSPConvergedSkip - Convergence test that do not return as converged
1242:   until the maximum number of iterations is reached.

1244:   Collective

1246:   Input Parameters:
1247: + ksp   - iterative context
1248: . n     - iteration number
1249: . rnorm - 2-norm residual value (may be estimated)
1250: - dtx   - unused convergence context

1252:   Output Parameter:
1253: . reason - `KSP_CONVERGED_ITERATING` or `KSP_CONVERGED_ITS`

1255:   Options Database Key:
1256: . -ksp_convergence_test skip - skips the test

1258:   Level: advanced

1260:   Note:
1261:   This should be used as the convergence test with the option
1262:   `KSPSetNormType`(ksp,`KSP_NORM_NONE`), since norms of the residual are
1263:   not computed. Convergence is then declared after the maximum number
1264:   of iterations have been reached. Useful when one is using `KSPCG` or
1265:   `KSPBCGS`. [](sec_flexibleksp)

1267: .seealso: [](ch_ksp), `KSP`, `KSPCG`, `KSPBCGS`, `KSPSetConvergenceTest()`, `KSPSetTolerances()`, `KSPSetNormType()`, [](sec_flexibleksp),
1268:           `KSPConvergedReason`
1269: @*/
1270: PetscErrorCode KSPConvergedSkip(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, void *dtx)
1271: {
1272:   PetscFunctionBegin;
1274:   PetscAssertPointer(reason, 4);
1275:   *reason = KSP_CONVERGED_ITERATING;
1276:   if (n >= ksp->max_it) *reason = KSP_CONVERGED_ITS;
1277:   PetscFunctionReturn(PETSC_SUCCESS);
1278: }

1280: /*@
1281:   KSPSetConvergedNegativeCurvature - Allows to declare convergence and return `KSP_CONVERGED_NEG_CURVE` when negative curvature is detected

1283:   Collective

1285:   Input Parameters:
1286: + ksp - iterative context
1287: - flg - the Boolean value

1289:   Options Database Key:
1290: . -ksp_converged_neg_curve <bool> - Declare convergence if negative curvature is detected

1292:   Level: advanced

1294:   Note:
1295:   This is currently used only by a subset of the Krylov solvers, namely `KSPCG`, `KSPSTCG`, `KSPQCG`, `KSPGLTR`, `KSPNASH`, and `KSPMINRES`.

1297: .seealso: [](ch_ksp), `KSP`, `KSPConvergedReason`, `KSPGetConvergedNegativeCurvature()`
1298: @*/
1299: PetscErrorCode KSPSetConvergedNegativeCurvature(KSP ksp, PetscBool flg)
1300: {
1301:   PetscFunctionBegin;
1304:   ksp->converged_neg_curve = flg;
1305:   PetscFunctionReturn(PETSC_SUCCESS);
1306: }

1308: /*@
1309:   KSPGetConvergedNegativeCurvature - Get the flag to declare convergence if negative curvature is detected

1311:   Collective

1313:   Input Parameter:
1314: . ksp - iterative context

1316:   Output Parameter:
1317: . flg - the Boolean value

1319:   Level: advanced

1321: .seealso: [](ch_ksp), `KSP`, `KSPConvergedReason`, `KSPSetConvergedNegativeCurvature()`
1322: @*/
1323: PetscErrorCode KSPGetConvergedNegativeCurvature(KSP ksp, PetscBool *flg)
1324: {
1325:   PetscFunctionBegin;
1327:   PetscAssertPointer(flg, 2);
1328:   *flg = ksp->converged_neg_curve;
1329:   PetscFunctionReturn(PETSC_SUCCESS);
1330: }

1332: /*@C
1333:   KSPConvergedDefaultCreate - Creates and initializes the context used by the `KSPConvergedDefault()` function

1335:   Not Collective

1337:   Output Parameter:
1338: . ctx - convergence context

1340:   Level: intermediate

1342: .seealso: [](ch_ksp), `KSP`, `KSPConvergedDefault()`, `KSPConvergedDefaultDestroy()`, `KSPSetConvergenceTest()`, `KSPSetTolerances()`,
1343:           `KSPConvergedSkip()`, `KSPConvergedReason`, `KSPGetConvergedReason()`, `KSPConvergedDefaultSetUIRNorm()`, `KSPConvergedDefaultSetUMIRNorm()`,
1344:           `KSPConvergedDefaultSetConvergedMaxits()`
1345: @*/
1346: PetscErrorCode KSPConvergedDefaultCreate(void **ctx)
1347: {
1348:   KSPConvergedDefaultCtx *cctx;

1350:   PetscFunctionBegin;
1351:   PetscCall(PetscNew(&cctx));
1352:   *ctx = cctx;
1353:   PetscFunctionReturn(PETSC_SUCCESS);
1354: }

1356: /*@
1357:   KSPConvergedDefaultSetUIRNorm - makes the default convergence test use $ || B*(b - A*(initial guess))||$
1358:   instead of $ || B*b ||$. In the case of right preconditioner or if `KSPSetNormType`(ksp,`KSP_NORM_UNPRECONDITIONED`)
1359:   is used there is no B in the above formula.

1361:   Collective

1363:   Input Parameters:
1364: . ksp - iterative context

1366:   Options Database Key:
1367: . -ksp_converged_use_initial_residual_norm <bool> - Use initial residual norm for computing relative convergence

1369:   Level: intermediate

1371:   Notes:
1372:   UIRNorm is short for Use Initial Residual Norm.

1374:   Use `KSPSetTolerances()` to alter the defaults for rtol, abstol, dtol.

1376:   The precise values of reason are macros such as `KSP_CONVERGED_RTOL`, which
1377:   are defined in petscksp.h.

1379:   If the convergence test is not `KSPConvergedDefault()` then this is ignored.

1381:   If right preconditioning is being used then B does not appear in the above formula.

1383: .seealso: [](ch_ksp), `KSP`, `KSPSetConvergenceTest()`, `KSPSetTolerances()`, `KSPConvergedSkip()`, `KSPConvergedReason`, `KSPGetConvergedReason()`, `KSPConvergedDefaultSetUMIRNorm()`, `KSPConvergedDefaultSetConvergedMaxits()`
1384: @*/
1385: PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP ksp)
1386: {
1387:   KSPConvergedDefaultCtx *ctx = (KSPConvergedDefaultCtx *)ksp->cnvP;

1389:   PetscFunctionBegin;
1391:   if (ksp->converged != KSPConvergedDefault) PetscFunctionReturn(PETSC_SUCCESS);
1392:   PetscCheck(!ctx->mininitialrtol, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_WRONGSTATE, "Cannot use KSPConvergedDefaultSetUIRNorm() and KSPConvergedDefaultSetUMIRNorm() together");
1393:   ctx->initialrtol = PETSC_TRUE;
1394:   PetscFunctionReturn(PETSC_SUCCESS);
1395: }

1397: /*@
1398:   KSPConvergedDefaultSetUMIRNorm - makes the default convergence test use min(|| B*(b - A*(initial guess))||,|| B*b ||)
1399:   In the case of right preconditioner or if `KSPSetNormType`(ksp,`KSP_NORM_UNPRECONDITIONED`)
1400:   is used there is no B in the above formula.

1402:   Collective

1404:   Input Parameters:
1405: . ksp - iterative context

1407:   Options Database Key:
1408: . -ksp_converged_use_min_initial_residual_norm <bool> - Use minimum of initial residual norm and b for computing relative convergence

1410:   Level: intermediate

1412:   Notes:
1413:   UMIRNorm is short for Use Minimum Initial Residual Norm.

1415:   Use `KSPSetTolerances()` to alter the defaults for rtol, abstol, dtol.

1417: .seealso: [](ch_ksp), `KSP`, `KSPSetConvergenceTest()`, `KSPSetTolerances()`, `KSPConvergedSkip()`, `KSPConvergedReason`, `KSPGetConvergedReason()`, `KSPConvergedDefaultSetUIRNorm()`, `KSPConvergedDefaultSetConvergedMaxits()`
1418: @*/
1419: PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP ksp)
1420: {
1421:   KSPConvergedDefaultCtx *ctx = (KSPConvergedDefaultCtx *)ksp->cnvP;

1423:   PetscFunctionBegin;
1425:   if (ksp->converged != KSPConvergedDefault) PetscFunctionReturn(PETSC_SUCCESS);
1426:   PetscCheck(!ctx->initialrtol, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_WRONGSTATE, "Cannot use KSPConvergedDefaultSetUIRNorm() and KSPConvergedDefaultSetUMIRNorm() together");
1427:   ctx->mininitialrtol = PETSC_TRUE;
1428:   PetscFunctionReturn(PETSC_SUCCESS);
1429: }

1431: /*@
1432:   KSPConvergedDefaultSetConvergedMaxits - allows the default convergence test to declare convergence and return `KSP_CONVERGED_ITS` if the maximum number of iterations is reached

1434:   Collective

1436:   Input Parameters:
1437: + ksp - iterative context
1438: - flg - boolean flag

1440:   Options Database Key:
1441: . -ksp_converged_maxits <bool> - Declare convergence if the maximum number of iterations is reached

1443:   Level: intermediate

1445: .seealso: [](ch_ksp), `KSP`, `KSPSetConvergenceTest()`, `KSPSetTolerances()`, `KSPConvergedSkip()`, `KSPConvergedReason`, `KSPGetConvergedReason()`, `KSPConvergedDefaultSetUMIRNorm()`, `KSPConvergedDefaultSetUIRNorm()`
1446: @*/
1447: PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP ksp, PetscBool flg)
1448: {
1449:   KSPConvergedDefaultCtx *ctx = (KSPConvergedDefaultCtx *)ksp->cnvP;

1451:   PetscFunctionBegin;
1454:   if (ksp->converged != KSPConvergedDefault) PetscFunctionReturn(PETSC_SUCCESS);
1455:   ctx->convmaxits = flg;
1456:   PetscFunctionReturn(PETSC_SUCCESS);
1457: }

1459: /*@C
1460:   KSPConvergedDefault - Default code to determine convergence of the linear iterative solvers

1462:   Collective

1464:   Input Parameters:
1465: + ksp   - iterative context
1466: . n     - iteration number
1467: . rnorm - residual norm (may be estimated, depending on the method may be the preconditioned residual norm)
1468: - ctx   - convergence context which must be created by `KSPConvergedDefaultCreate()`

1470:   Output Parameter:
1471: . reason - the convergence reason; it is positive if the iteration has converged,
1472:            negative if the iteration has diverged, and `KSP_CONVERGED_ITERATING` otherwise

1474:   Options Database Keys:
1475: + -ksp_max_it                                  - maximum number of linear iterations
1476: . -ksp_min_it                                  - minimum number of linear iterations, defaults to 0
1477: . -ksp_rtol rtol                               - relative tolerance used in default determination of convergence, i.e. if residual norm decreases by this factor than convergence is declared
1478: . -ksp_atol abstol                             - absolute tolerance used in default convergence test, i.e. if residual norm is less than this then convergence is declared
1479: . -ksp_divtol tol                              - if residual norm increases by this factor than divergence is declared
1480: . -ksp_converged_use_initial_residual_norm     - see `KSPConvergedDefaultSetUIRNorm()`
1481: . -ksp_converged_use_min_initial_residual_norm - see `KSPConvergedDefaultSetUMIRNorm()`
1482: - -ksp_converged_maxits                        - see `KSPConvergedDefaultSetConvergedMaxits()`

1484:   Level: advanced

1486:   Notes:
1487:   `KSPConvergedDefault()` reaches convergence when   rnorm < MAX (rtol * rnorm_0, abstol);
1488:   Divergence is detected if rnorm > dtol * rnorm_0, or when failures are detected throughout the iteration.
1489:   By default, reaching the maximum number of iterations is considered divergence (i.e. KSP_DIVERGED_ITS).
1490:   In order to have PETSc declaring convergence in such a case (i.e. `KSP_CONVERGED_ITS`), users can use `KSPConvergedDefaultSetConvergedMaxits()`

1492:   where\:
1493: +     `rtol` - relative tolerance,
1494: .     `abstol` - absolute tolerance.
1495: .     `dtol` - divergence tolerance,
1496: -     `rnorm_0` - the two norm of the right-hand side (or the preconditioned norm, depending on what was set with
1497:   `KSPSetNormType()`. When initial guess is non-zero you
1498:   can call `KSPConvergedDefaultSetUIRNorm()` to use the norm of (b - A*(initial guess))
1499:   as the starting point for relative norm convergence testing, that is as `rnorm_0`.
1500:   Call `KSPConvergedDefaultSetUMIRNorm()` to use the minimum of the norm of (b - A*(initial guess)) and the norm of b as the starting point.

1502:   Use `KSPSetTolerances()` to alter the defaults for `rtol`, `abstol`, `dtol`.

1504:   Use `KSPSetNormType()` (or `-ksp_norm_type <none,preconditioned,unpreconditioned,natural>`) to change the norm used for computing rnorm

1506:   The precise values of reason are available in `KSPConvergedReason`

1508:   This routine is used by `KSP` by default so the user generally never needs call it directly.

1510:   Use `KSPSetConvergenceTest()` to provide your own test instead of using this one.

1512:   Call `KSPSetConvergenceTest()` with the `ctx`, as created above and the destruction function `KSPConvergedDefaultDestroy()`

1514: .seealso: [](ch_ksp), `KSP`, `KSPSetConvergenceTest()`, `KSPSetTolerances()`, `KSPConvergedSkip()`, `KSPConvergedReason`, `KSPGetConvergedReason()`, `KSPSetMinimumIterations()`,
1515:           `KSPConvergedDefaultSetUIRNorm()`, `KSPConvergedDefaultSetUMIRNorm()`, `KSPConvergedDefaultSetConvergedMaxits()`, `KSPConvergedDefaultCreate()`, `KSPConvergedDefaultDestroy()`
1516: @*/
1517: PetscErrorCode KSPConvergedDefault(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, void *ctx)
1518: {
1519:   KSPConvergedDefaultCtx *cctx = (KSPConvergedDefaultCtx *)ctx;
1520:   KSPNormType             normtype;

1522:   PetscFunctionBegin;
1525:   PetscAssertPointer(reason, 4);
1526:   PetscCheck(cctx, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_NULL, "Convergence context must have been created with KSPConvergedDefaultCreate()");
1527:   *reason = KSP_CONVERGED_ITERATING;

1529:   if (cctx->convmaxits && n >= ksp->max_it) {
1530:     *reason = KSP_CONVERGED_ITS;
1531:     PetscCall(PetscInfo(ksp, "Linear solver has converged. Maximum number of iterations reached %" PetscInt_FMT "\n", n));
1532:     PetscFunctionReturn(PETSC_SUCCESS);
1533:   }
1534:   PetscCall(KSPGetNormType(ksp, &normtype));
1535:   if (normtype == KSP_NORM_NONE) PetscFunctionReturn(PETSC_SUCCESS);

1537:   if (!n) {
1538:     /* if user gives initial guess need to compute norm of b */
1539:     if (!ksp->guess_zero && !cctx->initialrtol) {
1540:       PetscReal snorm = 0.0;
1541:       if (ksp->normtype == KSP_NORM_UNPRECONDITIONED || ksp->pc_side == PC_RIGHT) {
1542:         PetscCall(PetscInfo(ksp, "user has provided nonzero initial guess, computing 2-norm of RHS\n"));
1543:         PetscCall(VecNorm(ksp->vec_rhs, NORM_2, &snorm)); /*     <- b'*b */
1544:       } else {
1545:         Vec z;
1546:         /* Should avoid allocating the z vector each time but cannot stash it in cctx because if KSPReset() is called the vector size might change */
1547:         PetscCall(VecDuplicate(ksp->vec_rhs, &z));
1548:         PetscCall(KSP_PCApply(ksp, ksp->vec_rhs, z));
1549:         if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
1550:           PetscCall(PetscInfo(ksp, "user has provided nonzero initial guess, computing 2-norm of preconditioned RHS\n"));
1551:           PetscCall(VecNorm(z, NORM_2, &snorm)); /*    dp <- b'*B'*B*b */
1552:         } else if (ksp->normtype == KSP_NORM_NATURAL) {
1553:           PetscScalar norm;
1554:           PetscCall(PetscInfo(ksp, "user has provided nonzero initial guess, computing natural norm of RHS\n"));
1555:           PetscCall(VecDot(ksp->vec_rhs, z, &norm));
1556:           snorm = PetscSqrtReal(PetscAbsScalar(norm)); /*    dp <- b'*B*b */
1557:         }
1558:         PetscCall(VecDestroy(&z));
1559:       }
1560:       /* handle special case of zero RHS and nonzero guess */
1561:       if (!snorm) {
1562:         PetscCall(PetscInfo(ksp, "Special case, user has provided nonzero initial guess and zero RHS\n"));
1563:         snorm = rnorm;
1564:       }
1565:       if (cctx->mininitialrtol) ksp->rnorm0 = PetscMin(snorm, rnorm);
1566:       else ksp->rnorm0 = snorm;
1567:     } else {
1568:       ksp->rnorm0 = rnorm;
1569:     }
1570:     ksp->ttol = PetscMax(ksp->rtol * ksp->rnorm0, ksp->abstol);
1571:   }

1573:   if (n <= ksp->chknorm) PetscFunctionReturn(PETSC_SUCCESS);

1575:   if (PetscIsInfOrNanReal(rnorm)) {
1576:     PCFailedReason pcreason;
1577:     PetscCall(PCReduceFailedReason(ksp->pc));
1578:     PetscCall(PCGetFailedReason(ksp->pc, &pcreason));
1579:     if (pcreason) {
1580:       *reason = KSP_DIVERGED_PC_FAILED;
1581:       PetscCall(PetscInfo(ksp, "Linear solver pcsetup fails, declaring divergence \n"));
1582:     } else {
1583:       *reason = KSP_DIVERGED_NANORINF;
1584:       PetscCall(PetscInfo(ksp, "Linear solver has created a not a number (NaN) as the residual norm, declaring divergence \n"));
1585:     }
1586:     PetscFunctionReturn(PETSC_SUCCESS);
1587:   }

1589:   if (n < ksp->min_it) PetscFunctionReturn(PETSC_SUCCESS);

1591:   if (rnorm <= ksp->ttol) {
1592:     if (rnorm < ksp->abstol) {
1593:       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));
1594:       *reason = KSP_CONVERGED_ATOL;
1595:     } else {
1596:       if (cctx->initialrtol) {
1597:         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));
1598:       } else {
1599:         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));
1600:       }
1601:       *reason = KSP_CONVERGED_RTOL;
1602:     }
1603:   } else if (rnorm >= ksp->divtol * ksp->rnorm0) {
1604:     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));
1605:     *reason = KSP_DIVERGED_DTOL;
1606:   }
1607:   PetscFunctionReturn(PETSC_SUCCESS);
1608: }

1610: /*@C
1611:   KSPConvergedDefaultDestroy - Frees the space used by the `KSPConvergedDefault()` function context

1613:   Not Collective

1615:   Input Parameter:
1616: . ctx - convergence context

1618:   Level: intermediate

1620:   Note:
1621:   Pass this function name into `KSPSetConvergenceTest()` along with the context obtained with `KSPConvergedDefaultCreate()` and `KSPConvergedDefault()`

1623: .seealso: [](ch_ksp), `KSP`, `KSPConvergedDefault()`, `KSPConvergedDefaultCreate()`, `KSPSetConvergenceTest()`, `KSPSetTolerances()`, `KSPConvergedSkip()`,
1624:           `KSPConvergedReason`, `KSPGetConvergedReason()`, `KSPConvergedDefaultSetUIRNorm()`, `KSPConvergedDefaultSetUMIRNorm()`
1625: @*/
1626: PetscErrorCode KSPConvergedDefaultDestroy(void *ctx)
1627: {
1628:   KSPConvergedDefaultCtx *cctx = (KSPConvergedDefaultCtx *)ctx;

1630:   PetscFunctionBegin;
1631:   PetscCall(VecDestroy(&cctx->work));
1632:   PetscCall(PetscFree(ctx));
1633:   PetscFunctionReturn(PETSC_SUCCESS);
1634: }

1636: // PetscClangLinter pragma disable: -fdoc-sowing-chars
1637: /*
1638:   KSPBuildSolutionDefault - Default code to build/move the solution.

1640:   Collective

1642:   Input Parameters:
1643: + ksp - iterative context
1644: - v   - pointer to the user's vector

1646:   Output Parameter:
1647: . V - pointer to a vector containing the solution

1649:   Level: advanced

1651:   Note:
1652:   Some `KSP` methods such as `KSPGMRES` do not compute the explicit solution at each iteration, this routine takes the information
1653:   they have computed during the previous iterations and uses it to compute the explicit solution

1655:   Developer Note:
1656:   This is `PETSC_EXTERN` because it may be used by user written plugin `KSPType` implementations

1658: .seealso: [](ch_ksp), `KSP`, `KSPGetSolution()`, `KSPBuildResidualDefault()`
1659: */
1660: PetscErrorCode KSPBuildSolutionDefault(KSP ksp, Vec v, Vec *V)
1661: {
1662:   PetscFunctionBegin;
1663:   if (ksp->pc_side == PC_RIGHT) {
1664:     if (ksp->pc) {
1665:       PetscCheck(v, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Not working with right preconditioner");
1666:       PetscCall(KSP_PCApply(ksp, ksp->vec_sol, v));
1667:       *V = v;
1668:     } else {
1669:       if (v) {
1670:         PetscCall(VecCopy(ksp->vec_sol, v));
1671:         *V = v;
1672:       } else *V = ksp->vec_sol;
1673:     }
1674:   } else if (ksp->pc_side == PC_SYMMETRIC) {
1675:     if (ksp->pc) {
1676:       PetscCheck(!ksp->transpose_solve, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Not working with symmetric preconditioner and transpose solve");
1677:       PetscCheck(v, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Not working with symmetric preconditioner");
1678:       PetscCall(PCApplySymmetricRight(ksp->pc, ksp->vec_sol, v));
1679:       *V = v;
1680:     } else {
1681:       if (v) {
1682:         PetscCall(VecCopy(ksp->vec_sol, v));
1683:         *V = v;
1684:       } else *V = ksp->vec_sol;
1685:     }
1686:   } else {
1687:     if (v) {
1688:       PetscCall(VecCopy(ksp->vec_sol, v));
1689:       *V = v;
1690:     } else *V = ksp->vec_sol;
1691:   }
1692:   PetscFunctionReturn(PETSC_SUCCESS);
1693: }

1695: /*@
1696:   KSPBuildResidualDefault - Default code to compute the residual.

1698:   Collecive on ksp

1700:   Input Parameters:
1701: + ksp - iterative context
1702: . t   - pointer to temporary vector
1703: - v   - pointer to user vector

1705:   Output Parameter:
1706: . V - pointer to a vector containing the residual

1708:   Level: advanced

1710:   Note:
1711:   Some `KSP` methods such as `KSPGMRES` do not compute the explicit residual at each iteration, this routine takes the information
1712:   they have computed during the previous iterations and uses it to compute the explicit residual via the formula r = b - A*x.

1714:   Developer Note:
1715:   This is `PETSC_EXTERN` because it may be used by user written plugin `KSPType` implementations

1717: .seealso: [](ch_ksp), `KSP`, `KSPBuildSolutionDefault()`
1718: @*/
1719: PetscErrorCode KSPBuildResidualDefault(KSP ksp, Vec t, Vec v, Vec *V)
1720: {
1721:   Mat Amat, Pmat;

1723:   PetscFunctionBegin;
1724:   if (!ksp->pc) PetscCall(KSPGetPC(ksp, &ksp->pc));
1725:   PetscCall(PCGetOperators(ksp->pc, &Amat, &Pmat));
1726:   PetscCall(KSPBuildSolution(ksp, t, NULL));
1727:   PetscCall(KSP_MatMult(ksp, Amat, t, v));
1728:   PetscCall(VecAYPX(v, -1.0, ksp->vec_rhs));
1729:   *V = v;
1730:   PetscFunctionReturn(PETSC_SUCCESS);
1731: }

1733: /*@C
1734:   KSPCreateVecs - Gets a number of work vectors suitably sized for the operator in the `KSP`

1736:   Collective

1738:   Input Parameters:
1739: + ksp    - iterative context
1740: . rightn - number of right work vectors
1741: - leftn  - number of left work vectors to allocate

1743:   Output Parameters:
1744: + right - the array of vectors created
1745: - left  - the array of left vectors

1747:   Level: advanced

1749:   Notes:
1750:   The right vector has as many elements as the matrix has columns. The left
1751:   vector has as many elements as the matrix has rows, see `MatSetSizes()` for details on the layout of the vectors.

1753:   The vectors are new vectors that are not owned by the `KSP`, they should be destroyed with calls to `VecDestroyVecs()` when no longer needed.

1755:   Developer Note:
1756:   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
1757:   that does not exist tries to get them from the `DM` (if it is provided) with `DMCreateGlobalVectors()`.

1759: .seealso: [](ch_ksp), `MatCreateVecs()`, `VecDestroyVecs()`, `KSPSetWorkVecs()`
1760: @*/
1761: PetscErrorCode KSPCreateVecs(KSP ksp, PetscInt rightn, Vec **right, PetscInt leftn, Vec **left)
1762: {
1763:   Vec       vecr = NULL, vecl = NULL;
1764:   PetscBool matset, pmatset, isshell, preferdm = PETSC_FALSE;
1765:   Mat       mat = NULL;

1767:   PetscFunctionBegin;
1768:   if (ksp->dm) {
1769:     PetscCall(PetscObjectTypeCompare((PetscObject)ksp->dm, DMSHELL, &isshell));
1770:     preferdm = isshell ? PETSC_FALSE : PETSC_TRUE;
1771:   }
1772:   if (rightn) {
1773:     PetscCheck(right, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_INCOMP, "You asked for right vectors but did not pass a pointer to hold them");
1774:     if (ksp->vec_sol) vecr = ksp->vec_sol;
1775:     else {
1776:       if (preferdm) {
1777:         PetscCall(DMGetGlobalVector(ksp->dm, &vecr));
1778:       } else if (ksp->pc) {
1779:         PetscCall(PCGetOperatorsSet(ksp->pc, &matset, &pmatset));
1780:         /* check for mat before pmat because for KSPLSQR pmat may be a different size than mat since pmat maybe mat'*mat */
1781:         if (matset) {
1782:           PetscCall(PCGetOperators(ksp->pc, &mat, NULL));
1783:           PetscCall(MatCreateVecs(mat, &vecr, NULL));
1784:         } else if (pmatset) {
1785:           PetscCall(PCGetOperators(ksp->pc, NULL, &mat));
1786:           PetscCall(MatCreateVecs(mat, &vecr, NULL));
1787:         }
1788:       }
1789:       if (!vecr && ksp->dm) PetscCall(DMGetGlobalVector(ksp->dm, &vecr));
1790:       PetscCheck(vecr, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_WRONGSTATE, "You requested a vector from a KSP that cannot provide one");
1791:     }
1792:     PetscCall(VecDuplicateVecs(vecr, rightn, right));
1793:     if (!ksp->vec_sol) {
1794:       if (preferdm) {
1795:         PetscCall(DMRestoreGlobalVector(ksp->dm, &vecr));
1796:       } else if (mat) {
1797:         PetscCall(VecDestroy(&vecr));
1798:       } else if (ksp->dm) {
1799:         PetscCall(DMRestoreGlobalVector(ksp->dm, &vecr));
1800:       }
1801:     }
1802:   }
1803:   if (leftn) {
1804:     PetscCheck(left, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_INCOMP, "You asked for left vectors but did not pass a pointer to hold them");
1805:     if (ksp->vec_rhs) vecl = ksp->vec_rhs;
1806:     else {
1807:       if (preferdm) {
1808:         PetscCall(DMGetGlobalVector(ksp->dm, &vecl));
1809:       } else if (ksp->pc) {
1810:         PetscCall(PCGetOperatorsSet(ksp->pc, &matset, &pmatset));
1811:         /* check for mat before pmat because for KSPLSQR pmat may be a different size than mat since pmat maybe mat'*mat */
1812:         if (matset) {
1813:           PetscCall(PCGetOperators(ksp->pc, &mat, NULL));
1814:           PetscCall(MatCreateVecs(mat, NULL, &vecl));
1815:         } else if (pmatset) {
1816:           PetscCall(PCGetOperators(ksp->pc, NULL, &mat));
1817:           PetscCall(MatCreateVecs(mat, NULL, &vecl));
1818:         }
1819:       }
1820:       if (!vecl && ksp->dm) PetscCall(DMGetGlobalVector(ksp->dm, &vecl));
1821:       PetscCheck(vecl, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_WRONGSTATE, "You requested a vector from a KSP that cannot provide one");
1822:     }
1823:     PetscCall(VecDuplicateVecs(vecl, leftn, left));
1824:     if (!ksp->vec_rhs) {
1825:       if (preferdm) {
1826:         PetscCall(DMRestoreGlobalVector(ksp->dm, &vecl));
1827:       } else if (mat) {
1828:         PetscCall(VecDestroy(&vecl));
1829:       } else if (ksp->dm) {
1830:         PetscCall(DMRestoreGlobalVector(ksp->dm, &vecl));
1831:       }
1832:     }
1833:   }
1834:   PetscFunctionReturn(PETSC_SUCCESS);
1835: }

1837: /*@C
1838:   KSPSetWorkVecs - Sets a number of work vectors into a `KSP` object

1840:   Collective

1842:   Input Parameters:
1843: + ksp - iterative context
1844: - nw  - number of work vectors to allocate

1846:   Level: developer

1848:   Developer Note:
1849:   This is `PETSC_EXTERN` because it may be used by user written plugin `KSPType` implementations

1851: .seealso: [](ch_ksp), `KSP`, `KSPCreateVecs()`
1852: @*/
1853: PetscErrorCode KSPSetWorkVecs(KSP ksp, PetscInt nw)
1854: {
1855:   PetscFunctionBegin;
1856:   PetscCall(VecDestroyVecs(ksp->nwork, &ksp->work));
1857:   ksp->nwork = nw;
1858:   PetscCall(KSPCreateVecs(ksp, nw, &ksp->work, 0, NULL));
1859:   PetscFunctionReturn(PETSC_SUCCESS);
1860: }

1862: // PetscClangLinter pragma disable: -fdoc-sowing-chars
1863: /*
1864:   KSPDestroyDefault - Destroys a iterative context variable for methods with no separate context.  Preferred calling sequence `KSPDestroy()`.

1866:   Input Parameter:
1867: . ksp - the iterative context

1869:   Level: advanced

1871:   Developer Note:
1872:   This is `PETSC_EXTERN` because it may be used by user written plugin `KSPType` implementations

1874: .seealso: [](ch_ksp), `KSP`, `KSPDestroy()`
1875: */
1876: PetscErrorCode KSPDestroyDefault(KSP ksp)
1877: {
1878:   PetscFunctionBegin;
1880:   PetscCall(PetscFree(ksp->data));
1881:   PetscFunctionReturn(PETSC_SUCCESS);
1882: }

1884: /*@
1885:   KSPGetConvergedReason - Gets the reason the `KSP` iteration was stopped.

1887:   Not Collective

1889:   Input Parameter:
1890: . ksp - the `KSP` context

1892:   Output Parameter:
1893: . reason - negative value indicates diverged, positive value converged, see `KSPConvergedReason` for the possible values

1895:   Options Database Key:
1896: . -ksp_converged_reason - prints the reason to standard out when the solve ends

1898:   Level: intermediate

1900:   Note:
1901:   If this routine is called before or doing the `KSPSolve()` the value of `KSP_CONVERGED_ITERATING` is returned

1903: .seealso: [](ch_ksp), `KSPConvergedReason`, `KSP`, `KSPSetConvergenceTest()`, `KSPConvergedDefault()`, `KSPSetTolerances()`,
1904:           `KSPConvergedReasonView()`, `KSPGetConvergedReasonString()`
1905: @*/
1906: PetscErrorCode KSPGetConvergedReason(KSP ksp, KSPConvergedReason *reason)
1907: {
1908:   PetscFunctionBegin;
1910:   PetscAssertPointer(reason, 2);
1911:   *reason = ksp->reason;
1912:   PetscFunctionReturn(PETSC_SUCCESS);
1913: }

1915: /*@C
1916:   KSPGetConvergedReasonString - Return a human readable string for a `KSPConvergedReason`

1918:   Not Collective

1920:   Input Parameter:
1921: . ksp - the `KSP` context

1923:   Output Parameter:
1924: . strreason - a human readable string that describes ksp converged reason

1926:   Level: beginner

1928: .seealso: [](ch_ksp), `KSP`, `KSPGetConvergedReason()`
1929: @*/
1930: PetscErrorCode KSPGetConvergedReasonString(KSP ksp, const char **strreason)
1931: {
1932:   PetscFunctionBegin;
1934:   PetscAssertPointer(strreason, 2);
1935:   *strreason = KSPConvergedReasons[ksp->reason];
1936:   PetscFunctionReturn(PETSC_SUCCESS);
1937: }

1939: #include <petsc/private/dmimpl.h>
1940: /*@
1941:   KSPSetDM - Sets the `DM` that may be used by some preconditioners and that may be used to construct the linear system

1943:   Logically Collective

1945:   Input Parameters:
1946: + ksp - the `KSP`
1947: - dm  - the `DM`, cannot be `NULL` to remove a previously set `DM`

1949:   Level: intermediate

1951:   Notes:
1952:   If this is used then the `KSP` will attempt to use the `DM` to create the matrix and use the routine set with
1953:   `DMKSPSetComputeOperators()`. Use `KSPSetDMActive`(ksp,`PETSC_FALSE`) to instead use the matrix you've provided with
1954:   `KSPSetOperators()`.

1956:   A `DM` can only be used for solving one problem at a time because information about the problem is stored on the `DM`,
1957:   even when not using interfaces like `DMKSPSetComputeOperators()`.  Use `DMClone()` to get a distinct `DM` when solving
1958:   different problems using the same function space.

1960: .seealso: [](ch_ksp), `KSP`, `DM`, `KSPGetDM()`, `KSPSetDMActive()`, `KSPSetComputeOperators()`, `KSPSetComputeRHS()`, `KSPSetComputeInitialGuess()`, `DMKSPSetComputeOperators()`, `DMKSPSetComputeRHS()`, `DMKSPSetComputeInitialGuess()`
1961: @*/
1962: PetscErrorCode KSPSetDM(KSP ksp, DM dm)
1963: {
1964:   PC pc;

1966:   PetscFunctionBegin;
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:   PetscCall(KSPGetPC(ksp, &pc));
1982:   PetscCall(PCSetDM(pc, dm));
1983:   ksp->dmActive = PETSC_TRUE;
1984:   PetscFunctionReturn(PETSC_SUCCESS);
1985: }

1987: /*@
1988:   KSPSetDMActive - Indicates the `DM` should be used to generate the linear system matrix and right-hand side vector

1990:   Logically Collective

1992:   Input Parameters:
1993: + ksp - the `KSP`
1994: - flg - use the `DM`

1996:   Level: intermediate

1998:   Note:
1999:   By default `KSPSetDM()` sets the `DM` as active, call `KSPSetDMActive`(ksp,`PETSC_FALSE`); after `KSPSetDM`(ksp,dm) to not have the `KSP` object use the `DM` to generate the matrices.

2001: .seealso: [](ch_ksp), `KSP`, `DM`, `KSPGetDM()`, `KSPSetDM()`, `SNESSetDM()`, `KSPSetComputeOperators()`, `KSPSetComputeRHS()`, `KSPSetComputeInitialGuess()`
2002: @*/
2003: PetscErrorCode KSPSetDMActive(KSP ksp, PetscBool flg)
2004: {
2005:   PetscFunctionBegin;
2008:   ksp->dmActive = flg;
2009:   PetscFunctionReturn(PETSC_SUCCESS);
2010: }

2012: /*@
2013:   KSPGetDM - Gets the `DM` that may be used by some preconditioners and that may be used to construct the linear system

2015:   Not Collective

2017:   Input Parameter:
2018: . ksp - the `KSP`

2020:   Output Parameter:
2021: . dm - the `DM`

2023:   Level: intermediate

2025: .seealso: [](ch_ksp), `KSP`, `DM`, `KSPSetDM()`, `KSPSetDMActive()`
2026: @*/
2027: PetscErrorCode KSPGetDM(KSP ksp, DM *dm)
2028: {
2029:   PetscFunctionBegin;
2031:   if (!ksp->dm) {
2032:     PetscCall(DMShellCreate(PetscObjectComm((PetscObject)ksp), &ksp->dm));
2033:     ksp->dmAuto = PETSC_TRUE;
2034:   }
2035:   *dm = ksp->dm;
2036:   PetscFunctionReturn(PETSC_SUCCESS);
2037: }

2039: /*@
2040:   KSPSetApplicationContext - Sets the optional user-defined context for the linear solver.

2042:   Logically Collective

2044:   Input Parameters:
2045: + ksp - the `KSP` context
2046: - ctx - optional user context

2048:   Level: intermediate

2050:   Notes:
2051:   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:   Use `KSPGetApplicationContext()` to get access to the context at a later time.

2055:   Fortran Note:
2056:   To use this from Fortran you must write a Fortran interface definition for this
2057:   function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.

2059: .seealso: [](ch_ksp), `KSP`, `KSPGetApplicationContext()`
2060: @*/
2061: PetscErrorCode KSPSetApplicationContext(KSP ksp, void *ctx)
2062: {
2063:   PC pc;

2065:   PetscFunctionBegin;
2067:   ksp->user = ctx;
2068:   PetscCall(KSPGetPC(ksp, &pc));
2069:   PetscCall(PCSetApplicationContext(pc, ctx));
2070:   PetscFunctionReturn(PETSC_SUCCESS);
2071: }

2073: /*@
2074:   KSPGetApplicationContext - Gets the user-defined context for the linear solver set with `KSPSetApplicationContext()`

2076:   Not Collective

2078:   Input Parameter:
2079: . ksp - `KSP` context

2081:   Output Parameter:
2082: . ctx - user context

2084:   Level: intermediate

2086:   Fortran Note:
2087:   You may need to write a Fortran interface definition for this
2088:   function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.

2090: .seealso: [](ch_ksp), `KSP`, `KSPSetApplicationContext()`
2091: @*/
2092: PetscErrorCode KSPGetApplicationContext(KSP ksp, void *ctx)
2093: {
2094:   PetscFunctionBegin;
2096:   *(void **)ctx = ksp->user;
2097:   PetscFunctionReturn(PETSC_SUCCESS);
2098: }

2100: #include <petsc/private/pcimpl.h>

2102: /*@
2103:   KSPCheckSolve - Checks if the `PCSetUp()` or `KSPSolve()` failed and set the error flag for the outer `PC`. A `KSP_DIVERGED_ITS` is
2104:   not considered a failure in this context

2106:   Collective

2108:   Input Parameters:
2109: + ksp - the linear solver `KSP` context.
2110: . pc  - the preconditioner context
2111: - vec - a vector that will be initialized with Inf to indicate lack of convergence

2113:   Level: developer

2115:   Note:
2116:   This is called within `PCApply()` implementations to check if an error has been detected on any particular MPI processes. By initializing the vector
2117:   with Inf the next call to `KSPCheckNorm()` or `KSPCheckDot()` will provide the same information to all the MPI processes that an error occurred on
2118:   at least one of the processes.

2120:   This may be called by a subset of the processes in the `PC`.

2122:   Developer Note:
2123:   This is used to manage returning with appropriate information from preconditioners whose inner `KSP` solvers have failed in some way

2125: .seealso: [](ch_ksp), `KSP`, `KSPCreate()`, `KSPSetType()`, `KSPCheckNorm()`, `KSPCheckDot()`
2126: @*/
2127: PetscErrorCode KSPCheckSolve(KSP ksp, PC pc, Vec vec)
2128: {
2129:   PCFailedReason pcreason;
2130:   PC             subpc;

2132:   PetscFunctionBegin;
2133:   PetscCall(KSPGetPC(ksp, &subpc));
2134:   PetscCall(PCGetFailedReason(subpc, &pcreason));
2135:   if (pcreason || (ksp->reason < 0 && ksp->reason != KSP_DIVERGED_ITS)) {
2136:     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]);
2137:     PetscCall(PetscInfo(ksp, "Detected not converged in KSP inner solve: KSP reason %s PC reason %s\n", KSPConvergedReasons[ksp->reason], PCFailedReasons[pcreason]));
2138:     pc->failedreason = PC_SUBPC_ERROR;
2139:     if (vec) PetscCall(VecSetInf(vec));
2140:   }
2141:   PetscFunctionReturn(PETSC_SUCCESS);
2142: }