Actual source code: itfunc.c

  1: /*
  2:       Interface KSP routines that the user calls.
  3: */

  5: #include <petsc/private/kspimpl.h>
  6: #include <petsc/private/matimpl.h>
  7: #include <petscdm.h>

  9: /* number of nested levels of KSPSetUp/Solve(). This is used to determine if KSP_DIVERGED_ITS should be fatal. */
 10: static PetscInt level = 0;

 12: static inline PetscErrorCode ObjectView(PetscObject obj, PetscViewer viewer, PetscViewerFormat format)
 13: {
 14:   PetscCall(PetscViewerPushFormat(viewer, format));
 15:   PetscCall(PetscObjectView(obj, viewer));
 16:   PetscCall(PetscViewerPopFormat(viewer));
 17:   return PETSC_SUCCESS;
 18: }

 20: /*@
 21:   KSPComputeExtremeSingularValues - Computes the extreme singular values
 22:   for the preconditioned operator. Called after or during `KSPSolve()`.

 24:   Not Collective

 26:   Input Parameter:
 27: . ksp - iterative solver obtained from `KSPCreate()`

 29:   Output Parameters:
 30: + emax - maximum estimated singular value
 31: - emin - minimum estimated singular value

 33:   Options Database Key:
 34: . -ksp_view_singularvalues - compute extreme singular values and print when `KSPSolve()` completes.

 36:   Level: advanced

 38:   Notes:
 39:   One must call `KSPSetComputeSingularValues()` before calling `KSPSetUp()`
 40:   (or use the option `-ksp_view_singularvalues`) in order for this routine to work correctly.

 42:   Many users may just want to use the monitoring routine
 43:   `KSPMonitorSingularValue()` (which can be set with option `-ksp_monitor_singular_value`)
 44:   to print the extreme singular values at each iteration of the linear solve.

 46:   Estimates of the smallest singular value may be very inaccurate, especially if the Krylov method has not converged.
 47:   The largest singular value is usually accurate to within a few percent if the method has converged, but is still not
 48:   intended for eigenanalysis. Consider the excellent package SLEPc if accurate values are required.

 50:   Disable restarts if using `KSPGMRES`, otherwise this estimate will only be using those iterations after the last
 51:   restart. See `KSPGMRESSetRestart()` for more details.

 53: .seealso: [](ch_ksp), `KSPSetComputeSingularValues()`, `KSPMonitorSingularValue()`, `KSPComputeEigenvalues()`, `KSP`, `KSPComputeRitz()`
 54: @*/
 55: PetscErrorCode KSPComputeExtremeSingularValues(KSP ksp, PetscReal *emax, PetscReal *emin)
 56: {
 57:   PetscFunctionBegin;
 59:   PetscAssertPointer(emax, 2);
 60:   PetscAssertPointer(emin, 3);
 61:   PetscCheck(ksp->calc_sings, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_WRONGSTATE, "Singular values not requested before KSPSetUp()");

 63:   if (ksp->ops->computeextremesingularvalues) PetscUseTypeMethod(ksp, computeextremesingularvalues, emax, emin);
 64:   else {
 65:     *emin = -1.0;
 66:     *emax = -1.0;
 67:   }
 68:   PetscFunctionReturn(PETSC_SUCCESS);
 69: }

 71: /*@
 72:   KSPComputeEigenvalues - Computes the extreme eigenvalues for the
 73:   preconditioned operator. Called after or during `KSPSolve()`.

 75:   Not Collective

 77:   Input Parameters:
 78: + ksp - iterative solver obtained from `KSPCreate()`
 79: - n   - size of arrays `r` and `c`. The number of eigenvalues computed `neig` will, in
 80:        general, be less than this.

 82:   Output Parameters:
 83: + r    - real part of computed eigenvalues, provided by user with a dimension of at least `n`
 84: . c    - complex part of computed eigenvalues, provided by user with a dimension of at least `n`
 85: - neig - actual number of eigenvalues computed (will be less than or equal to `n`)

 87:   Options Database Key:
 88: . -ksp_view_eigenvalues - Prints eigenvalues to stdout

 90:   Level: advanced

 92:   Notes:
 93:   The number of eigenvalues estimated depends on the size of the Krylov space
 94:   generated during the `KSPSolve()` ; for example, with
 95:   `KSPCG` it corresponds to the number of CG iterations, for `KSPGMRES` it is the number
 96:   of GMRES iterations SINCE the last restart. Any extra space in `r` and `c`
 97:   will be ignored.

 99:   `KSPComputeEigenvalues()` does not usually provide accurate estimates; it is
100:   intended only for assistance in understanding the convergence of iterative
101:   methods, not for eigenanalysis. For accurate computation of eigenvalues we recommend using
102:   the excellent package SLEPc.

104:   One must call `KSPSetComputeEigenvalues()` before calling `KSPSetUp()`
105:   in order for this routine to work correctly.

107:   Many users may just want to use the monitoring routine
108:   `KSPMonitorSingularValue()` (which can be set with option `-ksp_monitor_singular_value`)
109:   to print the singular values at each iteration of the linear solve.

111:   `KSPComputeRitz()` provides estimates for both the eigenvalues and their corresponding eigenvectors.

113: .seealso: [](ch_ksp), `KSPSetComputeEigenvalues()`, `KSPSetComputeSingularValues()`, `KSPMonitorSingularValue()`, `KSPComputeExtremeSingularValues()`, `KSP`, `KSPComputeRitz()`
114: @*/
115: PetscErrorCode KSPComputeEigenvalues(KSP ksp, PetscInt n, PetscReal r[], PetscReal c[], PetscInt *neig)
116: {
117:   PetscFunctionBegin;
119:   if (n) PetscAssertPointer(r, 3);
120:   if (n) PetscAssertPointer(c, 4);
121:   PetscCheck(n >= 0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Requested < 0 Eigenvalues");
122:   PetscAssertPointer(neig, 5);
123:   PetscCheck(ksp->calc_sings, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_WRONGSTATE, "Eigenvalues not requested before KSPSetUp()");

125:   if (n && ksp->ops->computeeigenvalues) PetscUseTypeMethod(ksp, computeeigenvalues, n, r, c, neig);
126:   else *neig = 0;
127:   PetscFunctionReturn(PETSC_SUCCESS);
128: }

130: /*@
131:   KSPComputeRitz - Computes the Ritz or harmonic Ritz pairs associated with the
132:   smallest or largest in modulus, for the preconditioned operator.

134:   Not Collective

136:   Input Parameters:
137: + ksp   - iterative solver obtained from `KSPCreate()`
138: . ritz  - `PETSC_TRUE` or `PETSC_FALSE` for Ritz pairs or harmonic Ritz pairs, respectively
139: - small - `PETSC_TRUE` or `PETSC_FALSE` for smallest or largest (harmonic) Ritz values, respectively

141:   Output Parameters:
142: + nrit  - On input number of (harmonic) Ritz pairs to compute; on output, actual number of computed (harmonic) Ritz pairs
143: . S     - an array of the Ritz vectors, pass in an array of vectors of size `nrit`
144: . tetar - real part of the Ritz values, pass in an array of size `nrit`
145: - tetai - imaginary part of the Ritz values, pass in an array of size `nrit`

147:   Level: advanced

149:   Notes:
150:   This only works with a `KSPType` of `KSPGMRES`.

152:   One must call `KSPSetComputeRitz()` before calling `KSPSetUp()` in order for this routine to work correctly.

154:   This routine must be called after `KSPSolve()`.

156:   In `KSPGMRES`, the (harmonic) Ritz pairs are computed from the Hessenberg matrix obtained during
157:   the last complete cycle of the GMRES solve, or during the partial cycle if the solve ended before
158:   a restart (that is a complete GMRES cycle was never achieved).

160:   The number of actual (harmonic) Ritz pairs computed is less than or equal to the restart
161:   parameter for GMRES if a complete cycle has been performed or less or equal to the number of GMRES
162:   iterations.

164:   `KSPComputeEigenvalues()` provides estimates for only the eigenvalues (Ritz values).

166:   For real matrices, the (harmonic) Ritz pairs can be complex-valued. In such a case,
167:   the routine selects the complex (harmonic) Ritz value and its conjugate, and two successive entries of the
168:   vectors `S` are equal to the real and the imaginary parts of the associated vectors.
169:   When PETSc has been built with complex scalars, the real and imaginary parts of the Ritz
170:   values are still returned in `tetar` and `tetai`, as is done in `KSPComputeEigenvalues()`, but
171:   the Ritz vectors S are complex.

173:   The (harmonic) Ritz pairs are given in order of increasing (harmonic) Ritz values in modulus.

175:   The Ritz pairs do not necessarily accurately reflect the eigenvalues and eigenvectors of the operator, consider the
176:   excellent package SLEPc if accurate values are required.

178: .seealso: [](ch_ksp), `KSPSetComputeRitz()`, `KSP`, `KSPGMRES`, `KSPComputeEigenvalues()`, `KSPSetComputeSingularValues()`, `KSPMonitorSingularValue()`
179: @*/
180: PetscErrorCode KSPComputeRitz(KSP ksp, PetscBool ritz, PetscBool small, PetscInt *nrit, Vec S[], PetscReal tetar[], PetscReal tetai[])
181: {
182:   PetscFunctionBegin;
184:   PetscCheck(ksp->calc_ritz, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_WRONGSTATE, "Ritz pairs not requested before KSPSetUp()");
185:   PetscTryTypeMethod(ksp, computeritz, ritz, small, nrit, S, tetar, tetai);
186:   PetscFunctionReturn(PETSC_SUCCESS);
187: }

189: /*@
190:   KSPSetUpOnBlocks - Sets up the preconditioner for each block in
191:   the block Jacobi `PCJACOBI`, overlapping Schwarz `PCASM`, and fieldsplit `PCFIELDSPLIT` preconditioners

193:   Collective

195:   Input Parameter:
196: . ksp - the `KSP` context

198:   Level: advanced

200:   Notes:
201:   `KSPSetUpOnBlocks()` is a routine that the user can optionally call for
202:   more precise profiling (via `-log_view`) of the setup phase for these
203:   block preconditioners.  If the user does not call `KSPSetUpOnBlocks()`,
204:   it will automatically be called from within `KSPSolve()`.

206:   Calling `KSPSetUpOnBlocks()` is the same as calling `PCSetUpOnBlocks()`
207:   on the `PC` context within the `KSP` context.

209: .seealso: [](ch_ksp), `PCSetUpOnBlocks()`, `KSPSetUp()`, `PCSetUp()`, `KSP`
210: @*/
211: PetscErrorCode KSPSetUpOnBlocks(KSP ksp)
212: {
213:   PC             pc;
214:   PCFailedReason pcreason;

216:   PetscFunctionBegin;
218:   level++;
219:   PetscCall(KSPGetPC(ksp, &pc));
220:   PetscCall(PCSetUpOnBlocks(pc));
221:   PetscCall(PCGetFailedReason(pc, &pcreason));
222:   level--;
223:   /*
224:      This is tricky since only a subset of MPI ranks may set this; each KSPSolve_*() is responsible for checking
225:      this flag and initializing an appropriate vector with VecFlag() so that the first norm computation can
226:      produce a result at KSPCheckNorm() thus communicating the known problem to all MPI ranks so they may
227:      terminate the Krylov solve. For many KSP implementations this is handled within KSPInitialResidual()
228:   */
229:   if (pcreason) ksp->reason = KSP_DIVERGED_PC_FAILED;
230:   PetscFunctionReturn(PETSC_SUCCESS);
231: }

233: /*@
234:   KSPSetReusePreconditioner - reuse the current preconditioner for future `KSPSolve()`, do not construct a new preconditioner even if the `Mat` operator
235:   in the `KSP` has different values

237:   Collective

239:   Input Parameters:
240: + ksp  - iterative solver obtained from `KSPCreate()`
241: - flag - `PETSC_TRUE` to reuse the current preconditioner, or `PETSC_FALSE` to construct a new preconditioner

243:   Options Database Key:
244: . -ksp_reuse_preconditioner <true,false> - reuse the previously computed preconditioner

246:   Level: intermediate

248:   Note:
249:   When using `SNES` one can use `SNESSetLagPreconditioner()` to determine when preconditioners are reused.

251: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGetReusePreconditioner()`,
252:           `SNESSetLagPreconditioner()`, `SNES`
253: @*/
254: PetscErrorCode KSPSetReusePreconditioner(KSP ksp, PetscBool flag)
255: {
256:   PC pc;

258:   PetscFunctionBegin;
260:   PetscCall(KSPGetPC(ksp, &pc));
261:   PetscCall(PCSetReusePreconditioner(pc, flag));
262:   PetscFunctionReturn(PETSC_SUCCESS);
263: }

265: /*@
266:   KSPGetReusePreconditioner - Determines if the `KSP` reuses the current preconditioner even if the `Mat` operator in the `KSP` has changed.

268:   Collective

270:   Input Parameter:
271: . ksp - iterative solver obtained from `KSPCreate()`

273:   Output Parameter:
274: . flag - the boolean flag indicating if the current preconditioner should be reused

276:   Level: intermediate

278: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSolve()`, `KSPDestroy()`, `KSPSetReusePreconditioner()`, `KSP`
279: @*/
280: PetscErrorCode KSPGetReusePreconditioner(KSP ksp, PetscBool *flag)
281: {
282:   PetscFunctionBegin;
284:   PetscAssertPointer(flag, 2);
285:   *flag = PETSC_FALSE;
286:   if (ksp->pc) PetscCall(PCGetReusePreconditioner(ksp->pc, flag));
287:   PetscFunctionReturn(PETSC_SUCCESS);
288: }

290: /*@
291:   KSPSetSkipPCSetFromOptions - prevents `KSPSetFromOptions()` from calling `PCSetFromOptions()`.
292:   This is used if the same `PC` is shared by more than one `KSP` so its options are not reset for each `KSP`

294:   Collective

296:   Input Parameters:
297: + ksp  - iterative solver obtained from `KSPCreate()`
298: - flag - `PETSC_TRUE` to skip calling the `PCSetFromOptions()`

300:   Level: intermediate

302: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSolve()`, `KSPDestroy()`, `PCSetReusePreconditioner()`, `KSP`
303: @*/
304: PetscErrorCode KSPSetSkipPCSetFromOptions(KSP ksp, PetscBool flag)
305: {
306:   PetscFunctionBegin;
308:   ksp->skippcsetfromoptions = flag;
309:   PetscFunctionReturn(PETSC_SUCCESS);
310: }

312: /*@
313:   KSPSetUp - Sets up the internal data structures for the
314:   later use `KSPSolve()` the `KSP` linear iterative solver.

316:   Collective

318:   Input Parameter:
319: . ksp - iterative solver, `KSP`, obtained from `KSPCreate()`

321:   Level: developer

323:   Note:
324:   This is called automatically by `KSPSolve()` so usually does not need to be called directly.

326: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPSetUpOnBlocks()`
327: @*/
328: PetscErrorCode KSPSetUp(KSP ksp)
329: {
330:   Mat            A, B;
331:   Mat            mat, pmat;
332:   MatNullSpace   nullsp;
333:   PCFailedReason pcreason;
334:   PC             pc;
335:   PetscBool      pcmpi;

337:   PetscFunctionBegin;
339:   PetscCall(KSPGetPC(ksp, &pc));
340:   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCMPI, &pcmpi));
341:   if (pcmpi) {
342:     PetscBool ksppreonly;
343:     PetscCall(PetscObjectTypeCompare((PetscObject)ksp, KSPPREONLY, &ksppreonly));
344:     if (!ksppreonly) PetscCall(KSPSetType(ksp, KSPPREONLY));
345:   }
346:   level++;

348:   /* reset the convergence flag from the previous solves */
349:   ksp->reason = KSP_CONVERGED_ITERATING;

351:   if (!((PetscObject)ksp)->type_name) PetscCall(KSPSetType(ksp, KSPGMRES));
352:   PetscCall(KSPSetUpNorms_Private(ksp, PETSC_TRUE, &ksp->normtype, &ksp->pc_side));

354:   if (ksp->dmActive && !ksp->setupstage) {
355:     /* first time in so build matrix and vector data structures using DM */
356:     if (!ksp->vec_rhs) PetscCall(DMCreateGlobalVector(ksp->dm, &ksp->vec_rhs));
357:     if (!ksp->vec_sol) PetscCall(DMCreateGlobalVector(ksp->dm, &ksp->vec_sol));
358:     PetscCall(DMCreateMatrix(ksp->dm, &A));
359:     PetscCall(KSPSetOperators(ksp, A, A));
360:     PetscCall(PetscObjectDereference((PetscObject)A));
361:   }

363:   if (ksp->dmActive) {
364:     DMKSP kdm;
365:     PetscCall(DMGetDMKSP(ksp->dm, &kdm));

367:     if (kdm->ops->computeinitialguess && ksp->setupstage != KSP_SETUP_NEWRHS) {
368:       /* only computes initial guess the first time through */
369:       PetscCallBack("KSP callback initial guess", (*kdm->ops->computeinitialguess)(ksp, ksp->vec_sol, kdm->initialguessctx));
370:       PetscCall(KSPSetInitialGuessNonzero(ksp, PETSC_TRUE));
371:     }
372:     if (kdm->ops->computerhs) PetscCallBack("KSP callback rhs", (*kdm->ops->computerhs)(ksp, ksp->vec_rhs, kdm->rhsctx));

374:     if (ksp->setupstage != KSP_SETUP_NEWRHS) {
375:       PetscCheck(kdm->ops->computeoperators, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_WRONGSTATE, "You called KSPSetDM() but did not use DMKSPSetComputeOperators() or KSPSetDMActive(ksp,PETSC_FALSE);");
376:       PetscCall(KSPGetOperators(ksp, &A, &B));
377:       PetscCallBack("KSP callback operators", (*kdm->ops->computeoperators)(ksp, A, B, kdm->operatorsctx));
378:     }
379:   }

381:   if (ksp->setupstage == KSP_SETUP_NEWRHS) {
382:     level--;
383:     PetscFunctionReturn(PETSC_SUCCESS);
384:   }
385:   PetscCall(PetscLogEventBegin(KSP_SetUp, ksp, ksp->vec_rhs, ksp->vec_sol, 0));

387:   switch (ksp->setupstage) {
388:   case KSP_SETUP_NEW:
389:     PetscUseTypeMethod(ksp, setup);
390:     break;
391:   case KSP_SETUP_NEWMATRIX: /* This should be replaced with a more general mechanism */
392:     if (ksp->setupnewmatrix) PetscUseTypeMethod(ksp, setup);
393:     break;
394:   default:
395:     break;
396:   }

398:   if (!ksp->pc) PetscCall(KSPGetPC(ksp, &ksp->pc));
399:   PetscCall(PCGetOperators(ksp->pc, &mat, &pmat));
400:   /* scale the matrix if requested */
401:   if (ksp->dscale) {
402:     PetscScalar *xx;
403:     PetscInt     i, n;
404:     PetscBool    zeroflag = PETSC_FALSE;

406:     if (!ksp->diagonal) { /* allocate vector to hold diagonal */
407:       PetscCall(MatCreateVecs(pmat, &ksp->diagonal, NULL));
408:     }
409:     PetscCall(MatGetDiagonal(pmat, ksp->diagonal));
410:     PetscCall(VecGetLocalSize(ksp->diagonal, &n));
411:     PetscCall(VecGetArray(ksp->diagonal, &xx));
412:     for (i = 0; i < n; i++) {
413:       if (xx[i] != 0.0) xx[i] = 1.0 / PetscSqrtReal(PetscAbsScalar(xx[i]));
414:       else {
415:         xx[i]    = 1.0;
416:         zeroflag = PETSC_TRUE;
417:       }
418:     }
419:     PetscCall(VecRestoreArray(ksp->diagonal, &xx));
420:     if (zeroflag) PetscCall(PetscInfo(ksp, "Zero detected in diagonal of matrix, using 1 at those locations\n"));
421:     PetscCall(MatDiagonalScale(pmat, ksp->diagonal, ksp->diagonal));
422:     if (mat != pmat) PetscCall(MatDiagonalScale(mat, ksp->diagonal, ksp->diagonal));
423:     ksp->dscalefix2 = PETSC_FALSE;
424:   }
425:   PetscCall(PetscLogEventEnd(KSP_SetUp, ksp, ksp->vec_rhs, ksp->vec_sol, 0));
426:   PetscCall(PCSetErrorIfFailure(ksp->pc, ksp->errorifnotconverged));
427:   PetscCall(PCSetUp(ksp->pc));
428:   PetscCall(PCGetFailedReason(ksp->pc, &pcreason));
429:   /* TODO: this code was wrong and is still wrong, there is no way to propagate the failure to all processes; their is no code to handle a ksp->reason on only some ranks */
430:   if (pcreason) ksp->reason = KSP_DIVERGED_PC_FAILED;

432:   PetscCall(MatGetNullSpace(mat, &nullsp));
433:   if (nullsp) {
434:     PetscBool test = PETSC_FALSE;
435:     PetscCall(PetscOptionsGetBool(((PetscObject)ksp)->options, ((PetscObject)ksp)->prefix, "-ksp_test_null_space", &test, NULL));
436:     if (test) PetscCall(MatNullSpaceTest(nullsp, mat, NULL));
437:   }
438:   ksp->setupstage = KSP_SETUP_NEWRHS;
439:   level--;
440:   PetscFunctionReturn(PETSC_SUCCESS);
441: }

443: /*@
444:   KSPConvergedReasonView - Displays the reason a `KSP` solve converged or diverged, `KSPConvergedReason` to a `PetscViewer`

446:   Collective

448:   Input Parameters:
449: + ksp    - iterative solver obtained from `KSPCreate()`
450: - viewer - the `PetscViewer` on which to display the reason

452:   Options Database Keys:
453: + -ksp_converged_reason          - print reason for converged or diverged, also prints number of iterations
454: - -ksp_converged_reason ::failed - only print reason and number of iterations when diverged

456:   Level: beginner

458:   Note:
459:   Use `KSPConvergedReasonViewFromOptions()` to display the reason based on values in the PETSc options database.

461:   To change the format of the output call `PetscViewerPushFormat`(`viewer`,`format`) before this call. Use `PETSC_VIEWER_DEFAULT` for the default,
462:   use `PETSC_VIEWER_FAILED` to only display a reason if it fails.

464: .seealso: [](ch_ksp), `KSPConvergedReasonViewFromOptions()`, `KSPCreate()`, `KSPSetUp()`, `KSPDestroy()`, `KSPSetTolerances()`, `KSPConvergedDefault()`,
465:           `KSPSolveTranspose()`, `KSPGetIterationNumber()`, `KSP`, `KSPGetConvergedReason()`, `PetscViewerPushFormat()`, `PetscViewerPopFormat()`
466: @*/
467: PetscErrorCode KSPConvergedReasonView(KSP ksp, PetscViewer viewer)
468: {
469:   PetscBool         isAscii;
470:   PetscViewerFormat format;

472:   PetscFunctionBegin;
473:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)ksp));
474:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isAscii));
475:   if (isAscii) {
476:     PetscCall(PetscViewerGetFormat(viewer, &format));
477:     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ksp)->tablevel + 1));
478:     if (ksp->reason > 0 && format != PETSC_VIEWER_FAILED) {
479:       if (((PetscObject)ksp)->prefix) {
480:         PetscCall(PetscViewerASCIIPrintf(viewer, "Linear %s solve converged due to %s iterations %" PetscInt_FMT "\n", ((PetscObject)ksp)->prefix, KSPConvergedReasons[ksp->reason], ksp->its));
481:       } else {
482:         PetscCall(PetscViewerASCIIPrintf(viewer, "Linear solve converged due to %s iterations %" PetscInt_FMT "\n", KSPConvergedReasons[ksp->reason], ksp->its));
483:       }
484:     } else if (ksp->reason <= 0) {
485:       if (((PetscObject)ksp)->prefix) {
486:         PetscCall(PetscViewerASCIIPrintf(viewer, "Linear %s solve did not converge due to %s iterations %" PetscInt_FMT "\n", ((PetscObject)ksp)->prefix, KSPConvergedReasons[ksp->reason], ksp->its));
487:       } else {
488:         PetscCall(PetscViewerASCIIPrintf(viewer, "Linear solve did not converge due to %s iterations %" PetscInt_FMT "\n", KSPConvergedReasons[ksp->reason], ksp->its));
489:       }
490:       if (ksp->reason == KSP_DIVERGED_PC_FAILED) {
491:         PCFailedReason reason;
492:         PetscCall(PCGetFailedReason(ksp->pc, &reason));
493:         PetscCall(PetscViewerASCIIPrintf(viewer, "               PC failed due to %s\n", PCFailedReasons[reason]));
494:       }
495:     }
496:     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ksp)->tablevel + 1));
497:   }
498:   PetscFunctionReturn(PETSC_SUCCESS);
499: }

501: /*@C
502:   KSPConvergedReasonViewSet - Sets an ADDITIONAL function that is to be used at the
503:   end of the linear solver to display the convergence reason of the linear solver.

505:   Logically Collective

507:   Input Parameters:
508: + ksp               - the `KSP` context
509: . f                 - the ksp converged reason view function
510: . vctx              - [optional] user-defined context for private data for the
511:                       `KSPConvergedReason` view routine (use `NULL` if no context is desired)
512: - reasonviewdestroy - [optional] routine that frees `vctx` (may be `NULL`), see `PetscCtxDestroyFn` for the calling sequence

514:   Options Database Keys:
515: + -ksp_converged_reason             - sets a default `KSPConvergedReasonView()`
516: - -ksp_converged_reason_view_cancel - cancels all converged reason viewers that have been hardwired into a code by
517:                                       calls to `KSPConvergedReasonViewSet()`, but does not cancel those set via the options database.

519:   Level: intermediate

521:   Note:
522:   Several different converged reason view routines may be set by calling
523:   `KSPConvergedReasonViewSet()` multiple times; all will be called in the
524:   order in which they were set.

526:   Developer Note:
527:   Should be named KSPConvergedReasonViewAdd().

529: .seealso: [](ch_ksp), `KSPConvergedReasonView()`, `KSPConvergedReasonViewCancel()`, `PetscCtxDestroyFn`
530: @*/
531: PetscErrorCode KSPConvergedReasonViewSet(KSP ksp, PetscErrorCode (*f)(KSP, void *), void *vctx, PetscCtxDestroyFn *reasonviewdestroy)
532: {
533:   PetscInt  i;
534:   PetscBool identical;

536:   PetscFunctionBegin;
538:   for (i = 0; i < ksp->numberreasonviews; i++) {
539:     PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))f, vctx, reasonviewdestroy, (PetscErrorCode (*)(void))ksp->reasonview[i], ksp->reasonviewcontext[i], ksp->reasonviewdestroy[i], &identical));
540:     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
541:   }
542:   PetscCheck(ksp->numberreasonviews < MAXKSPREASONVIEWS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many KSP reasonview set");
543:   ksp->reasonview[ksp->numberreasonviews]          = f;
544:   ksp->reasonviewdestroy[ksp->numberreasonviews]   = reasonviewdestroy;
545:   ksp->reasonviewcontext[ksp->numberreasonviews++] = vctx;
546:   PetscFunctionReturn(PETSC_SUCCESS);
547: }

549: /*@
550:   KSPConvergedReasonViewCancel - Clears all the reasonview functions for a `KSP` object set with `KSPConvergedReasonViewSet()`
551:   as well as the default viewer.

553:   Collective

555:   Input Parameter:
556: . ksp - iterative solver obtained from `KSPCreate()`

558:   Level: intermediate

560: .seealso: [](ch_ksp), `KSPCreate()`, `KSPDestroy()`, `KSPReset()`, `KSPConvergedReasonViewSet()`
561: @*/
562: PetscErrorCode KSPConvergedReasonViewCancel(KSP ksp)
563: {
564:   PetscInt i;

566:   PetscFunctionBegin;
568:   for (i = 0; i < ksp->numberreasonviews; i++) {
569:     if (ksp->reasonviewdestroy[i]) PetscCall((*ksp->reasonviewdestroy[i])(&ksp->reasonviewcontext[i]));
570:   }
571:   ksp->numberreasonviews = 0;
572:   PetscCall(PetscViewerDestroy(&ksp->convergedreasonviewer));
573:   PetscFunctionReturn(PETSC_SUCCESS);
574: }

576: /*@
577:   KSPConvergedReasonViewFromOptions - Processes command line options to determine if/how a `KSPReason` is to be viewed.

579:   Collective

581:   Input Parameter:
582: . ksp - the `KSP` object

584:   Level: intermediate

586:   Note:
587:   This is called automatically at the conclusion of `KSPSolve()` so is rarely called directly by user code.

589: .seealso: [](ch_ksp), `KSPConvergedReasonView()`, `KSPConvergedReasonViewSet()`
590: @*/
591: PetscErrorCode KSPConvergedReasonViewFromOptions(KSP ksp)
592: {
593:   PetscFunctionBegin;
594:   /* Call all user-provided reason review routines */
595:   for (PetscInt i = 0; i < ksp->numberreasonviews; i++) PetscCall((*ksp->reasonview[i])(ksp, ksp->reasonviewcontext[i]));

597:   /* Call the default PETSc routine */
598:   if (ksp->convergedreasonviewer) {
599:     PetscCall(PetscViewerPushFormat(ksp->convergedreasonviewer, ksp->convergedreasonformat));
600:     PetscCall(KSPConvergedReasonView(ksp, ksp->convergedreasonviewer));
601:     PetscCall(PetscViewerPopFormat(ksp->convergedreasonviewer));
602:   }
603:   PetscFunctionReturn(PETSC_SUCCESS);
604: }

606: /*@
607:   KSPConvergedRateView - Displays the convergence rate <https://en.wikipedia.org/wiki/Coefficient_of_determination> of `KSPSolve()` to a viewer

609:   Collective

611:   Input Parameters:
612: + ksp    - iterative solver obtained from `KSPCreate()`
613: - viewer - the `PetscViewer` to display the reason

615:   Options Database Key:
616: . -ksp_converged_rate - print reason for convergence or divergence and the convergence rate (or 0.0 for divergence)

618:   Level: intermediate

620:   Notes:
621:   To change the format of the output, call `PetscViewerPushFormat`(`viewer`,`format`) before this call.

623:   Suppose that the residual is reduced linearly, $r_k = c^k r_0$, which means $\log r_k = \log r_0 + k \log c$. After linear regression,
624:   the slope is $\log c$. The coefficient of determination is given by $1 - \frac{\sum_i (y_i - f(x_i))^2}{\sum_i (y_i - \bar y)}$,

626: .seealso: [](ch_ksp), `KSPConvergedReasonView()`, `KSPGetConvergedRate()`, `KSPSetTolerances()`, `KSPConvergedDefault()`
627: @*/
628: PetscErrorCode KSPConvergedRateView(KSP ksp, PetscViewer viewer)
629: {
630:   PetscViewerFormat format;
631:   PetscBool         isAscii;
632:   PetscReal         rrate, rRsq, erate = 0.0, eRsq = 0.0;
633:   PetscInt          its;
634:   const char       *prefix, *reason = KSPConvergedReasons[ksp->reason];

636:   PetscFunctionBegin;
637:   PetscCall(KSPGetIterationNumber(ksp, &its));
638:   PetscCall(KSPComputeConvergenceRate(ksp, &rrate, &rRsq, &erate, &eRsq));
639:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)ksp));
640:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isAscii));
641:   if (isAscii) {
642:     PetscCall(KSPGetOptionsPrefix(ksp, &prefix));
643:     PetscCall(PetscViewerGetFormat(viewer, &format));
644:     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ksp)->tablevel));
645:     if (ksp->reason > 0) {
646:       if (prefix) PetscCall(PetscViewerASCIIPrintf(viewer, "Linear %s solve converged due to %s iterations %" PetscInt_FMT, prefix, reason, its));
647:       else PetscCall(PetscViewerASCIIPrintf(viewer, "Linear solve converged due to %s iterations %" PetscInt_FMT, reason, its));
648:       PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
649:       if (rRsq >= 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, " res rate %g R^2 %g", (double)rrate, (double)rRsq));
650:       if (eRsq >= 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, " error rate %g R^2 %g", (double)erate, (double)eRsq));
651:       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
652:       PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
653:     } else if (ksp->reason <= 0) {
654:       if (prefix) PetscCall(PetscViewerASCIIPrintf(viewer, "Linear %s solve did not converge due to %s iterations %" PetscInt_FMT, prefix, reason, its));
655:       else PetscCall(PetscViewerASCIIPrintf(viewer, "Linear solve did not converge due to %s iterations %" PetscInt_FMT, reason, its));
656:       PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
657:       if (rRsq >= 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, " res rate %g R^2 %g", (double)rrate, (double)rRsq));
658:       if (eRsq >= 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, " error rate %g R^2 %g", (double)erate, (double)eRsq));
659:       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
660:       PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
661:       if (ksp->reason == KSP_DIVERGED_PC_FAILED) {
662:         PCFailedReason reason;
663:         PetscCall(PCGetFailedReason(ksp->pc, &reason));
664:         PetscCall(PetscViewerASCIIPrintf(viewer, "               PC failed due to %s\n", PCFailedReasons[reason]));
665:       }
666:     }
667:     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ksp)->tablevel));
668:   }
669:   PetscFunctionReturn(PETSC_SUCCESS);
670: }

672: #include <petscdraw.h>

674: static PetscErrorCode KSPViewEigenvalues_Internal(KSP ksp, PetscBool isExplicit, PetscViewer viewer, PetscViewerFormat format)
675: {
676:   PetscReal  *r, *c;
677:   PetscInt    n, i, neig;
678:   PetscBool   isascii, isdraw;
679:   PetscMPIInt rank;

681:   PetscFunctionBegin;
682:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)ksp), &rank));
683:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
684:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
685:   if (isExplicit) {
686:     PetscCall(VecGetSize(ksp->vec_sol, &n));
687:     PetscCall(PetscMalloc2(n, &r, n, &c));
688:     PetscCall(KSPComputeEigenvaluesExplicitly(ksp, n, r, c));
689:     neig = n;
690:   } else {
691:     PetscInt nits;

693:     PetscCall(KSPGetIterationNumber(ksp, &nits));
694:     n = nits + 2;
695:     if (!nits) {
696:       PetscCall(PetscViewerASCIIPrintf(viewer, "Zero iterations in solver, cannot approximate any eigenvalues\n"));
697:       PetscFunctionReturn(PETSC_SUCCESS);
698:     }
699:     PetscCall(PetscMalloc2(n, &r, n, &c));
700:     PetscCall(KSPComputeEigenvalues(ksp, n, r, c, &neig));
701:   }
702:   if (isascii) {
703:     PetscCall(PetscViewerASCIIPrintf(viewer, "%s computed eigenvalues\n", isExplicit ? "Explicitly" : "Iteratively"));
704:     for (i = 0; i < neig; ++i) {
705:       if (c[i] >= 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, "%g + %gi\n", (double)r[i], (double)c[i]));
706:       else PetscCall(PetscViewerASCIIPrintf(viewer, "%g - %gi\n", (double)r[i], -(double)c[i]));
707:     }
708:   } else if (isdraw && rank == 0) {
709:     PetscDraw   draw;
710:     PetscDrawSP drawsp;

712:     if (format == PETSC_VIEWER_DRAW_CONTOUR) {
713:       PetscCall(KSPPlotEigenContours_Private(ksp, neig, r, c));
714:     } else {
715:       PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
716:       PetscCall(PetscDrawSPCreate(draw, 1, &drawsp));
717:       PetscCall(PetscDrawSPReset(drawsp));
718:       for (i = 0; i < neig; ++i) PetscCall(PetscDrawSPAddPoint(drawsp, r + i, c + i));
719:       PetscCall(PetscDrawSPDraw(drawsp, PETSC_TRUE));
720:       PetscCall(PetscDrawSPSave(drawsp));
721:       PetscCall(PetscDrawSPDestroy(&drawsp));
722:     }
723:   }
724:   PetscCall(PetscFree2(r, c));
725:   PetscFunctionReturn(PETSC_SUCCESS);
726: }

728: static PetscErrorCode KSPViewSingularvalues_Internal(KSP ksp, PetscViewer viewer, PetscViewerFormat format)
729: {
730:   PetscReal smax, smin;
731:   PetscInt  nits;
732:   PetscBool isascii;

734:   PetscFunctionBegin;
735:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
736:   PetscCall(KSPGetIterationNumber(ksp, &nits));
737:   if (!nits) {
738:     PetscCall(PetscViewerASCIIPrintf(viewer, "Zero iterations in solver, cannot approximate any singular values\n"));
739:     PetscFunctionReturn(PETSC_SUCCESS);
740:   }
741:   PetscCall(KSPComputeExtremeSingularValues(ksp, &smax, &smin));
742:   if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, "Iteratively computed extreme %svalues: max %g min %g max/min %g\n", smin < 0 ? "eigen" : "singular ", (double)smax, (double)smin, (double)(smax / smin)));
743:   PetscFunctionReturn(PETSC_SUCCESS);
744: }

746: static PetscErrorCode KSPViewFinalResidual_Internal(KSP ksp, PetscViewer viewer, PetscViewerFormat format)
747: {
748:   PetscBool isascii;

750:   PetscFunctionBegin;
751:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
752:   PetscCheck(!ksp->dscale || ksp->dscalefix, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_WRONGSTATE, "Cannot compute final scale with -ksp_diagonal_scale except also with -ksp_diagonal_scale_fix");
753:   if (isascii) {
754:     Mat       A;
755:     Vec       t;
756:     PetscReal norm;

758:     PetscCall(PCGetOperators(ksp->pc, &A, NULL));
759:     PetscCall(VecDuplicate(ksp->vec_rhs, &t));
760:     PetscCall(KSP_MatMult(ksp, A, ksp->vec_sol, t));
761:     PetscCall(VecAYPX(t, -1.0, ksp->vec_rhs));
762:     PetscCall(VecViewFromOptions(t, (PetscObject)ksp, "-ksp_view_final_residual_vec"));
763:     PetscCall(VecNorm(t, NORM_2, &norm));
764:     PetscCall(VecDestroy(&t));
765:     PetscCall(PetscViewerASCIIPrintf(viewer, "KSP final norm of residual %g\n", (double)norm));
766:   }
767:   PetscFunctionReturn(PETSC_SUCCESS);
768: }

770: PETSC_EXTERN PetscErrorCode PetscMonitorPauseFinal_Internal(PetscInt n, void *ctx[])
771: {
772:   PetscFunctionBegin;
773:   for (PetscInt i = 0; i < n; ++i) {
774:     PetscViewerAndFormat *vf = (PetscViewerAndFormat *)ctx[i];
775:     PetscDraw             draw;
776:     PetscReal             lpause;
777:     PetscBool             isdraw;

779:     if (!vf) continue;
780:     if (!PetscCheckPointer(vf->viewer, PETSC_OBJECT)) continue;
781:     if (((PetscObject)vf->viewer)->classid != PETSC_VIEWER_CLASSID) continue;
782:     PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERDRAW, &isdraw));
783:     if (!isdraw) continue;

785:     PetscCall(PetscViewerDrawGetDraw(vf->viewer, 0, &draw));
786:     PetscCall(PetscDrawGetPause(draw, &lpause));
787:     PetscCall(PetscDrawSetPause(draw, -1.0));
788:     PetscCall(PetscDrawPause(draw));
789:     PetscCall(PetscDrawSetPause(draw, lpause));
790:   }
791:   PetscFunctionReturn(PETSC_SUCCESS);
792: }

794: static PetscErrorCode KSPMonitorPauseFinal_Internal(KSP ksp)
795: {
796:   PetscFunctionBegin;
797:   if (!ksp->pauseFinal) PetscFunctionReturn(PETSC_SUCCESS);
798:   PetscCall(PetscMonitorPauseFinal_Internal(ksp->numbermonitors, ksp->monitorcontext));
799:   PetscFunctionReturn(PETSC_SUCCESS);
800: }

802: static PetscErrorCode KSPSolve_Private(KSP ksp, Vec b, Vec x)
803: {
804:   PetscBool    flg = PETSC_FALSE, inXisinB = PETSC_FALSE, guess_zero;
805:   Mat          mat, pmat;
806:   MPI_Comm     comm;
807:   MatNullSpace nullsp;
808:   Vec          btmp, vec_rhs = NULL;

810:   PetscFunctionBegin;
811:   level++;
812:   comm = PetscObjectComm((PetscObject)ksp);
813:   if (x && x == b) {
814:     PetscCheck(ksp->guess_zero, comm, PETSC_ERR_ARG_INCOMP, "Cannot use x == b with nonzero initial guess");
815:     PetscCall(VecDuplicate(b, &x));
816:     inXisinB = PETSC_TRUE;
817:   }
818:   if (b) {
819:     PetscCall(PetscObjectReference((PetscObject)b));
820:     PetscCall(VecDestroy(&ksp->vec_rhs));
821:     ksp->vec_rhs = b;
822:   }
823:   if (x) {
824:     PetscCall(PetscObjectReference((PetscObject)x));
825:     PetscCall(VecDestroy(&ksp->vec_sol));
826:     ksp->vec_sol = x;
827:   }

829:   if (ksp->viewPre) PetscCall(ObjectView((PetscObject)ksp, ksp->viewerPre, ksp->formatPre));

831:   if (ksp->presolve) PetscCall((*ksp->presolve)(ksp, ksp->vec_rhs, ksp->vec_sol, ksp->prectx));

833:   /* reset the residual history list if requested */
834:   if (ksp->res_hist_reset) ksp->res_hist_len = 0;
835:   if (ksp->err_hist_reset) ksp->err_hist_len = 0;

837:   /* KSPSetUp() scales the matrix if needed */
838:   PetscCall(KSPSetUp(ksp));
839:   PetscCall(KSPSetUpOnBlocks(ksp));

841:   if (ksp->guess) {
842:     PetscObjectState ostate, state;

844:     PetscCall(KSPGuessSetUp(ksp->guess));
845:     PetscCall(PetscObjectStateGet((PetscObject)ksp->vec_sol, &ostate));
846:     PetscCall(KSPGuessFormGuess(ksp->guess, ksp->vec_rhs, ksp->vec_sol));
847:     PetscCall(PetscObjectStateGet((PetscObject)ksp->vec_sol, &state));
848:     if (state != ostate) {
849:       ksp->guess_zero = PETSC_FALSE;
850:     } else {
851:       PetscCall(PetscInfo(ksp, "Using zero initial guess since the KSPGuess object did not change the vector\n"));
852:       ksp->guess_zero = PETSC_TRUE;
853:     }
854:   }

856:   PetscCall(VecSetErrorIfLocked(ksp->vec_sol, 3));

858:   PetscCall(PetscLogEventBegin(!ksp->transpose_solve ? KSP_Solve : KSP_SolveTranspose, ksp, ksp->vec_rhs, ksp->vec_sol, 0));
859:   PetscCall(PCGetOperators(ksp->pc, &mat, &pmat));
860:   /* diagonal scale RHS if called for */
861:   if (ksp->dscale) {
862:     PetscCall(VecPointwiseMult(ksp->vec_rhs, ksp->vec_rhs, ksp->diagonal));
863:     /* second time in, but matrix was scaled back to original */
864:     if (ksp->dscalefix && ksp->dscalefix2) {
865:       Mat mat, pmat;

867:       PetscCall(PCGetOperators(ksp->pc, &mat, &pmat));
868:       PetscCall(MatDiagonalScale(pmat, ksp->diagonal, ksp->diagonal));
869:       if (mat != pmat) PetscCall(MatDiagonalScale(mat, ksp->diagonal, ksp->diagonal));
870:     }

872:     /* scale initial guess */
873:     if (!ksp->guess_zero) {
874:       if (!ksp->truediagonal) {
875:         PetscCall(VecDuplicate(ksp->diagonal, &ksp->truediagonal));
876:         PetscCall(VecCopy(ksp->diagonal, ksp->truediagonal));
877:         PetscCall(VecReciprocal(ksp->truediagonal));
878:       }
879:       PetscCall(VecPointwiseMult(ksp->vec_sol, ksp->vec_sol, ksp->truediagonal));
880:     }
881:   }
882:   PetscCall(PCPreSolve(ksp->pc, ksp));

884:   if (ksp->guess_zero && !ksp->guess_not_read) PetscCall(VecSet(ksp->vec_sol, 0.0));
885:   if (ksp->guess_knoll) { /* The Knoll trick is independent on the KSPGuess specified */
886:     PetscCall(PCApply(ksp->pc, ksp->vec_rhs, ksp->vec_sol));
887:     PetscCall(KSP_RemoveNullSpace(ksp, ksp->vec_sol));
888:     ksp->guess_zero = PETSC_FALSE;
889:   }

891:   /* can we mark the initial guess as zero for this solve? */
892:   guess_zero = ksp->guess_zero;
893:   if (!ksp->guess_zero) {
894:     PetscReal norm;

896:     PetscCall(VecNormAvailable(ksp->vec_sol, NORM_2, &flg, &norm));
897:     if (flg && !norm) ksp->guess_zero = PETSC_TRUE;
898:   }
899:   if (ksp->transpose_solve) {
900:     PetscCall(MatGetNullSpace(mat, &nullsp));
901:   } else {
902:     PetscCall(MatGetTransposeNullSpace(mat, &nullsp));
903:   }
904:   if (nullsp) {
905:     PetscCall(VecDuplicate(ksp->vec_rhs, &btmp));
906:     PetscCall(VecCopy(ksp->vec_rhs, btmp));
907:     PetscCall(MatNullSpaceRemove(nullsp, btmp));
908:     vec_rhs      = ksp->vec_rhs;
909:     ksp->vec_rhs = btmp;
910:   }
911:   PetscCall(VecLockReadPush(ksp->vec_rhs));
912:   PetscUseTypeMethod(ksp, solve);
913:   PetscCall(KSPMonitorPauseFinal_Internal(ksp));

915:   PetscCall(VecLockReadPop(ksp->vec_rhs));
916:   if (nullsp) {
917:     ksp->vec_rhs = vec_rhs;
918:     PetscCall(VecDestroy(&btmp));
919:   }

921:   ksp->guess_zero = guess_zero;

923:   PetscCheck(ksp->reason, comm, PETSC_ERR_PLIB, "Internal error, solver returned without setting converged reason");
924:   ksp->totalits += ksp->its;

926:   PetscCall(KSPConvergedReasonViewFromOptions(ksp));

928:   if (ksp->viewRate) {
929:     PetscCall(PetscViewerPushFormat(ksp->viewerRate, ksp->formatRate));
930:     PetscCall(KSPConvergedRateView(ksp, ksp->viewerRate));
931:     PetscCall(PetscViewerPopFormat(ksp->viewerRate));
932:   }
933:   PetscCall(PCPostSolve(ksp->pc, ksp));

935:   /* diagonal scale solution if called for */
936:   if (ksp->dscale) {
937:     PetscCall(VecPointwiseMult(ksp->vec_sol, ksp->vec_sol, ksp->diagonal));
938:     /* unscale right-hand side and matrix */
939:     if (ksp->dscalefix) {
940:       Mat mat, pmat;

942:       PetscCall(VecReciprocal(ksp->diagonal));
943:       PetscCall(VecPointwiseMult(ksp->vec_rhs, ksp->vec_rhs, ksp->diagonal));
944:       PetscCall(PCGetOperators(ksp->pc, &mat, &pmat));
945:       PetscCall(MatDiagonalScale(pmat, ksp->diagonal, ksp->diagonal));
946:       if (mat != pmat) PetscCall(MatDiagonalScale(mat, ksp->diagonal, ksp->diagonal));
947:       PetscCall(VecReciprocal(ksp->diagonal));
948:       ksp->dscalefix2 = PETSC_TRUE;
949:     }
950:   }
951:   PetscCall(PetscLogEventEnd(!ksp->transpose_solve ? KSP_Solve : KSP_SolveTranspose, ksp, ksp->vec_rhs, ksp->vec_sol, 0));
952:   if (ksp->guess) PetscCall(KSPGuessUpdate(ksp->guess, ksp->vec_rhs, ksp->vec_sol));
953:   if (ksp->postsolve) PetscCall((*ksp->postsolve)(ksp, ksp->vec_rhs, ksp->vec_sol, ksp->postctx));

955:   PetscCall(PCGetOperators(ksp->pc, &mat, &pmat));
956:   if (ksp->viewEV) PetscCall(KSPViewEigenvalues_Internal(ksp, PETSC_FALSE, ksp->viewerEV, ksp->formatEV));
957:   if (ksp->viewEVExp) PetscCall(KSPViewEigenvalues_Internal(ksp, PETSC_TRUE, ksp->viewerEVExp, ksp->formatEVExp));
958:   if (ksp->viewSV) PetscCall(KSPViewSingularvalues_Internal(ksp, ksp->viewerSV, ksp->formatSV));
959:   if (ksp->viewFinalRes) PetscCall(KSPViewFinalResidual_Internal(ksp, ksp->viewerFinalRes, ksp->formatFinalRes));
960:   if (ksp->viewMat) PetscCall(ObjectView((PetscObject)mat, ksp->viewerMat, ksp->formatMat));
961:   if (ksp->viewPMat) PetscCall(ObjectView((PetscObject)pmat, ksp->viewerPMat, ksp->formatPMat));
962:   if (ksp->viewRhs) PetscCall(ObjectView((PetscObject)ksp->vec_rhs, ksp->viewerRhs, ksp->formatRhs));
963:   if (ksp->viewSol) PetscCall(ObjectView((PetscObject)ksp->vec_sol, ksp->viewerSol, ksp->formatSol));
964:   if (ksp->view) PetscCall(ObjectView((PetscObject)ksp, ksp->viewer, ksp->format));
965:   if (ksp->viewDScale) PetscCall(ObjectView((PetscObject)ksp->diagonal, ksp->viewerDScale, ksp->formatDScale));
966:   if (ksp->viewMatExp) {
967:     Mat A, B;

969:     PetscCall(PCGetOperators(ksp->pc, &A, NULL));
970:     if (ksp->transpose_solve) {
971:       Mat AT;

973:       PetscCall(MatCreateTranspose(A, &AT));
974:       PetscCall(MatComputeOperator(AT, MATAIJ, &B));
975:       PetscCall(MatDestroy(&AT));
976:     } else {
977:       PetscCall(MatComputeOperator(A, MATAIJ, &B));
978:     }
979:     PetscCall(ObjectView((PetscObject)B, ksp->viewerMatExp, ksp->formatMatExp));
980:     PetscCall(MatDestroy(&B));
981:   }
982:   if (ksp->viewPOpExp) {
983:     Mat B;

985:     PetscCall(KSPComputeOperator(ksp, MATAIJ, &B));
986:     PetscCall(ObjectView((PetscObject)B, ksp->viewerPOpExp, ksp->formatPOpExp));
987:     PetscCall(MatDestroy(&B));
988:   }

990:   if (inXisinB) {
991:     PetscCall(VecCopy(x, b));
992:     PetscCall(VecDestroy(&x));
993:   }
994:   PetscCall(PetscObjectSAWsBlock((PetscObject)ksp));
995:   if (ksp->errorifnotconverged && ksp->reason < 0 && ((level == 1) || (ksp->reason != KSP_DIVERGED_ITS))) {
996:     PCFailedReason reason;

998:     PetscCheck(ksp->reason == KSP_DIVERGED_PC_FAILED, comm, PETSC_ERR_NOT_CONVERGED, "KSPSolve%s() has not converged, reason %s", !ksp->transpose_solve ? "" : "Transpose", KSPConvergedReasons[ksp->reason]);
999:     PetscCall(PCGetFailedReason(ksp->pc, &reason));
1000:     SETERRQ(comm, PETSC_ERR_NOT_CONVERGED, "KSPSolve%s() has not converged, reason %s PC failed due to %s", !ksp->transpose_solve ? "" : "Transpose", KSPConvergedReasons[ksp->reason], PCFailedReasons[reason]);
1001:   }
1002:   level--;
1003:   PetscFunctionReturn(PETSC_SUCCESS);
1004: }

1006: /*@
1007:   KSPSolve - Solves a linear system associated with `KSP` object

1009:   Collective

1011:   Input Parameters:
1012: + ksp - iterative solver obtained from `KSPCreate()`
1013: . b   - the right-hand side vector
1014: - x   - the solution (this may be the same vector as `b`, then `b` will be overwritten with the answer)

1016:   Options Database Keys:
1017: + -ksp_view_eigenvalues                      - compute preconditioned operators eigenvalues
1018: . -ksp_view_eigenvalues_explicit             - compute the eigenvalues by forming the dense operator and using LAPACK
1019: . -ksp_view_mat binary                       - save matrix to the default binary viewer
1020: . -ksp_view_pmat binary                      - save matrix used to build preconditioner to the default binary viewer
1021: . -ksp_view_rhs binary                       - save right-hand side vector to the default binary viewer
1022: . -ksp_view_solution binary                  - save computed solution vector to the default binary viewer
1023:                                                (can be read later with src/ksp/tutorials/ex10.c for testing solvers)
1024: . -ksp_view_mat_explicit                     - for matrix-free operators, computes the matrix entries and views them
1025: . -ksp_view_preconditioned_operator_explicit - computes the product of the preconditioner and matrix as an explicit matrix and views it
1026: . -ksp_converged_reason                      - print reason for converged or diverged, also prints number of iterations
1027: . -ksp_view_final_residual                   - print 2-norm of true linear system residual at the end of the solution process
1028: . -ksp_error_if_not_converged                - stop the program as soon as an error is detected in a `KSPSolve()`
1029: . -ksp_view_pre                              - print the ksp data structure before the system solution
1030: - -ksp_view                                  - print the ksp data structure at the end of the system solution

1032:   Level: beginner

1034:   Notes:
1035:   See `KSPSetFromOptions()` for additional options database keys that affect `KSPSolve()`

1037:   If one uses `KSPSetDM()` then `x` or `b` need not be passed. Use `KSPGetSolution()` to access the solution in this case.

1039:   The operator is specified with `KSPSetOperators()`.

1041:   `KSPSolve()` will normally return without generating an error regardless of whether the linear system was solved or if constructing the preconditioner failed.
1042:   Call `KSPGetConvergedReason()` to determine if the solver converged or failed and why. The option -ksp_error_if_not_converged or function `KSPSetErrorIfNotConverged()`
1043:   will cause `KSPSolve()` to error as soon as an error occurs in the linear solver.  In inner `KSPSolve()` `KSP_DIVERGED_ITS` is not treated as an error because when using nested solvers
1044:   it may be fine that inner solvers in the preconditioner do not converge during the solution process.

1046:   The number of iterations can be obtained from `KSPGetIterationNumber()`.

1048:   If you provide a matrix that has a `MatSetNullSpace()` and `MatSetTransposeNullSpace()` this will use that information to solve singular systems
1049:   in the least squares sense with a norm minimizing solution.

1051:   $A x = b $  where $b = b_p + b_t$ where $b_t$ is not in the range of $A$ (and hence by the fundamental theorem of linear algebra is in the nullspace(A'), see `MatSetNullSpace()`).

1053:   `KSP` first removes $b_t$ producing the linear system $ A x = b_p $ (which has multiple solutions) and solves this to find the $\|x\|$ minimizing solution (and hence
1054:   it finds the solution $x$ orthogonal to the nullspace(A). The algorithm is simply in each iteration of the Krylov method we remove the nullspace(A) from the search
1055:   direction thus the solution which is a linear combination of the search directions has no component in the nullspace(A).

1057:   We recommend always using `KSPGMRES` for such singular systems.
1058:   If nullspace(A) = nullspace(A') (note symmetric matrices always satisfy this property) then both left and right preconditioning will work
1059:   If nullspace(A) != nullspace(A') then left preconditioning will work but right preconditioning may not work (or it may).

1061:   Developer Notes:
1062:   The reason we cannot always solve  nullspace(A) != nullspace(A') systems with right preconditioning is because we need to remove at each iteration
1063:   the nullspace(AB) from the search direction. While we know the nullspace(A) the nullspace(AB) equals $B^-1$ times the nullspace(A) but except for trivial preconditioners
1064:   such as diagonal scaling we cannot apply the inverse of the preconditioner to a vector and thus cannot compute the nullspace(AB).

1066:   If using a direct method (e.g., via the `KSP` solver
1067:   `KSPPREONLY` and a preconditioner such as `PCLU` or `PCCHOLESKY` then usually one iteration of the `KSP` method will be needed for convergence.

1069:   To solve a linear system with the transpose of the matrix use `KSPSolveTranspose()`.

1071:   Understanding Convergence\:
1072:   The manual pages `KSPMonitorSet()`, `KSPComputeEigenvalues()`, and
1073:   `KSPComputeEigenvaluesExplicitly()` provide information on additional
1074:   options to monitor convergence and print eigenvalue information.

1076: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPDestroy()`, `KSPSetTolerances()`, `KSPConvergedDefault()`,
1077:           `KSPSolveTranspose()`, `KSPGetIterationNumber()`, `MatNullSpaceCreate()`, `MatSetNullSpace()`, `MatSetTransposeNullSpace()`, `KSP`,
1078:           `KSPConvergedReasonView()`, `KSPCheckSolve()`, `KSPSetErrorIfNotConverged()`
1079: @*/
1080: PetscErrorCode KSPSolve(KSP ksp, Vec b, Vec x)
1081: {
1082:   PetscBool isPCMPI;

1084:   PetscFunctionBegin;
1088:   ksp->transpose_solve = PETSC_FALSE;
1089:   PetscCall(KSPSolve_Private(ksp, b, x));
1090:   PetscCall(PetscObjectTypeCompare((PetscObject)ksp->pc, PCMPI, &isPCMPI));
1091:   if (PCMPIServerActive && isPCMPI) {
1092:     KSP subksp;

1094:     PetscCall(PCMPIGetKSP(ksp->pc, &subksp));
1095:     ksp->its    = subksp->its;
1096:     ksp->reason = subksp->reason;
1097:   }
1098:   PetscFunctionReturn(PETSC_SUCCESS);
1099: }

1101: /*@
1102:   KSPSolveTranspose - Solves a linear system with the transpose of the matrix associated with the `KSP` object, $ A^T x = b$.

1104:   Collective

1106:   Input Parameters:
1107: + ksp - iterative solver obtained from `KSPCreate()`
1108: . b   - right-hand side vector
1109: - x   - solution vector

1111:   Level: developer

1113:   Note:
1114:   For complex numbers this solve the non-Hermitian transpose system.

1116:   Developer Note:
1117:   We need to implement a `KSPSolveHermitianTranspose()`

1119: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPDestroy()`, `KSPSetTolerances()`, `KSPConvergedDefault()`,
1120:           `KSPSolve()`, `KSP`, `KSPSetOperators()`
1121: @*/
1122: PetscErrorCode KSPSolveTranspose(KSP ksp, Vec b, Vec x)
1123: {
1124:   PetscFunctionBegin;
1128:   if (ksp->transpose.use_explicittranspose) {
1129:     Mat J, Jpre;
1130:     PetscCall(KSPGetOperators(ksp, &J, &Jpre));
1131:     if (!ksp->transpose.reuse_transpose) {
1132:       PetscCall(MatTranspose(J, MAT_INITIAL_MATRIX, &ksp->transpose.AT));
1133:       if (J != Jpre) PetscCall(MatTranspose(Jpre, MAT_INITIAL_MATRIX, &ksp->transpose.BT));
1134:       ksp->transpose.reuse_transpose = PETSC_TRUE;
1135:     } else {
1136:       PetscCall(MatTranspose(J, MAT_REUSE_MATRIX, &ksp->transpose.AT));
1137:       if (J != Jpre) PetscCall(MatTranspose(Jpre, MAT_REUSE_MATRIX, &ksp->transpose.BT));
1138:     }
1139:     if (J == Jpre && ksp->transpose.BT != ksp->transpose.AT) {
1140:       PetscCall(PetscObjectReference((PetscObject)ksp->transpose.AT));
1141:       ksp->transpose.BT = ksp->transpose.AT;
1142:     }
1143:     PetscCall(KSPSetOperators(ksp, ksp->transpose.AT, ksp->transpose.BT));
1144:   } else {
1145:     ksp->transpose_solve = PETSC_TRUE;
1146:   }
1147:   PetscCall(KSPSolve_Private(ksp, b, x));
1148:   PetscFunctionReturn(PETSC_SUCCESS);
1149: }

1151: static PetscErrorCode KSPViewFinalMatResidual_Internal(KSP ksp, Mat B, Mat X, PetscViewer viewer, PetscViewerFormat format, PetscInt shift)
1152: {
1153:   Mat        A, R;
1154:   PetscReal *norms;
1155:   PetscInt   i, N;
1156:   PetscBool  flg;

1158:   PetscFunctionBegin;
1159:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &flg));
1160:   if (flg) {
1161:     PetscCall(PCGetOperators(ksp->pc, &A, NULL));
1162:     if (!ksp->transpose_solve) PetscCall(MatMatMult(A, X, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &R));
1163:     else PetscCall(MatTransposeMatMult(A, X, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &R));
1164:     PetscCall(MatAYPX(R, -1.0, B, SAME_NONZERO_PATTERN));
1165:     PetscCall(MatGetSize(R, NULL, &N));
1166:     PetscCall(PetscMalloc1(N, &norms));
1167:     PetscCall(MatGetColumnNorms(R, NORM_2, norms));
1168:     PetscCall(MatDestroy(&R));
1169:     for (i = 0; i < N; ++i) PetscCall(PetscViewerASCIIPrintf(viewer, "%s #%" PetscInt_FMT " %g\n", i == 0 ? "KSP final norm of residual" : "                          ", shift + i, (double)norms[i]));
1170:     PetscCall(PetscFree(norms));
1171:   }
1172:   PetscFunctionReturn(PETSC_SUCCESS);
1173: }

1175: static PetscErrorCode KSPMatSolve_Private(KSP ksp, Mat B, Mat X)
1176: {
1177:   Mat       A, P, vB, vX;
1178:   Vec       cb, cx;
1179:   PetscInt  n1, N1, n2, N2, Bbn = PETSC_DECIDE;
1180:   PetscBool match;

1182:   PetscFunctionBegin;
1186:   PetscCheckSameComm(ksp, 1, B, 2);
1187:   PetscCheckSameComm(ksp, 1, X, 3);
1188:   PetscCheckSameType(B, 2, X, 3);
1189:   PetscCheck(B->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
1190:   MatCheckPreallocated(X, 3);
1191:   if (!X->assembled) {
1192:     PetscCall(MatSetOption(X, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
1193:     PetscCall(MatAssemblyBegin(X, MAT_FINAL_ASSEMBLY));
1194:     PetscCall(MatAssemblyEnd(X, MAT_FINAL_ASSEMBLY));
1195:   }
1196:   PetscCheck(B != X, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_IDN, "B and X must be different matrices");
1197:   PetscCheck(!ksp->transpose_solve || !ksp->transpose.use_explicittranspose, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "KSPMatSolveTranspose() does not support -ksp_use_explicittranspose");
1198:   PetscCall(KSPGetOperators(ksp, &A, &P));
1199:   PetscCall(MatGetLocalSize(B, NULL, &n2));
1200:   PetscCall(MatGetLocalSize(X, NULL, &n1));
1201:   PetscCall(MatGetSize(B, NULL, &N2));
1202:   PetscCall(MatGetSize(X, NULL, &N1));
1203:   PetscCheck(n1 == n2 && N1 == N2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible number of columns between block of right-hand sides (n,N) = (%" PetscInt_FMT ",%" PetscInt_FMT ") and block of solutions (n,N) = (%" PetscInt_FMT ",%" PetscInt_FMT ")", n2, N2, n1, N1);
1204:   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)B, &match, MATSEQDENSE, MATMPIDENSE, ""));
1205:   PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of right-hand sides not stored in a dense Mat");
1206:   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)X, &match, MATSEQDENSE, MATMPIDENSE, ""));
1207:   PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of solutions not stored in a dense Mat");
1208:   PetscCall(KSPSetUp(ksp));
1209:   PetscCall(KSPSetUpOnBlocks(ksp));
1210:   if (ksp->ops->matsolve) {
1211:     level++;
1212:     if (ksp->guess_zero) PetscCall(MatZeroEntries(X));
1213:     PetscCall(PetscLogEventBegin(!ksp->transpose_solve ? KSP_MatSolve : KSP_MatSolveTranspose, ksp, B, X, 0));
1214:     PetscCall(KSPGetMatSolveBatchSize(ksp, &Bbn));
1215:     /* by default, do a single solve with all columns */
1216:     if (Bbn == PETSC_DECIDE) Bbn = N2;
1217:     else PetscCheck(Bbn >= 1, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "KSPMatSolve() batch size %" PetscInt_FMT " must be positive", Bbn);
1218:     PetscCall(PetscInfo(ksp, "KSP type %s%s solving using batches of width at most %" PetscInt_FMT "\n", ((PetscObject)ksp)->type_name, ksp->transpose_solve ? " transpose" : "", Bbn));
1219:     /* if -ksp_matsolve_batch_size is greater than the actual number of columns, do a single solve with all columns */
1220:     if (Bbn >= N2) {
1221:       PetscUseTypeMethod(ksp, matsolve, B, X);
1222:       if (ksp->viewFinalRes) PetscCall(KSPViewFinalMatResidual_Internal(ksp, B, X, ksp->viewerFinalRes, ksp->formatFinalRes, 0));

1224:       PetscCall(KSPConvergedReasonViewFromOptions(ksp));

1226:       if (ksp->viewRate) {
1227:         PetscCall(PetscViewerPushFormat(ksp->viewerRate, PETSC_VIEWER_DEFAULT));
1228:         PetscCall(KSPConvergedRateView(ksp, ksp->viewerRate));
1229:         PetscCall(PetscViewerPopFormat(ksp->viewerRate));
1230:       }
1231:     } else {
1232:       for (n2 = 0; n2 < N2; n2 += Bbn) {
1233:         PetscCall(MatDenseGetSubMatrix(B, PETSC_DECIDE, PETSC_DECIDE, n2, PetscMin(n2 + Bbn, N2), &vB));
1234:         PetscCall(MatDenseGetSubMatrix(X, PETSC_DECIDE, PETSC_DECIDE, n2, PetscMin(n2 + Bbn, N2), &vX));
1235:         PetscUseTypeMethod(ksp, matsolve, vB, vX);
1236:         if (ksp->viewFinalRes) PetscCall(KSPViewFinalMatResidual_Internal(ksp, vB, vX, ksp->viewerFinalRes, ksp->formatFinalRes, n2));

1238:         PetscCall(KSPConvergedReasonViewFromOptions(ksp));

1240:         if (ksp->viewRate) {
1241:           PetscCall(PetscViewerPushFormat(ksp->viewerRate, PETSC_VIEWER_DEFAULT));
1242:           PetscCall(KSPConvergedRateView(ksp, ksp->viewerRate));
1243:           PetscCall(PetscViewerPopFormat(ksp->viewerRate));
1244:         }
1245:         PetscCall(MatDenseRestoreSubMatrix(B, &vB));
1246:         PetscCall(MatDenseRestoreSubMatrix(X, &vX));
1247:       }
1248:     }
1249:     if (ksp->viewMat) PetscCall(ObjectView((PetscObject)A, ksp->viewerMat, ksp->formatMat));
1250:     if (ksp->viewPMat) PetscCall(ObjectView((PetscObject)P, ksp->viewerPMat, ksp->formatPMat));
1251:     if (ksp->viewRhs) PetscCall(ObjectView((PetscObject)B, ksp->viewerRhs, ksp->formatRhs));
1252:     if (ksp->viewSol) PetscCall(ObjectView((PetscObject)X, ksp->viewerSol, ksp->formatSol));
1253:     if (ksp->view) PetscCall(KSPView(ksp, ksp->viewer));
1254:     PetscCall(PetscLogEventEnd(!ksp->transpose_solve ? KSP_MatSolve : KSP_MatSolveTranspose, ksp, B, X, 0));
1255:     if (ksp->errorifnotconverged && ksp->reason < 0 && (level == 1 || ksp->reason != KSP_DIVERGED_ITS)) {
1256:       PCFailedReason reason;

1258:       PetscCheck(ksp->reason == KSP_DIVERGED_PC_FAILED, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "KSPMatSolve%s() has not converged, reason %s", !ksp->transpose_solve ? "" : "Transpose", KSPConvergedReasons[ksp->reason]);
1259:       PetscCall(PCGetFailedReason(ksp->pc, &reason));
1260:       SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "KSPMatSolve%s() has not converged, reason %s PC failed due to %s", !ksp->transpose_solve ? "" : "Transpose", KSPConvergedReasons[ksp->reason], PCFailedReasons[reason]);
1261:     }
1262:     level--;
1263:   } else {
1264:     PetscCall(PetscInfo(ksp, "KSP type %s solving column by column\n", ((PetscObject)ksp)->type_name));
1265:     for (n2 = 0; n2 < N2; ++n2) {
1266:       PetscCall(MatDenseGetColumnVecRead(B, n2, &cb));
1267:       PetscCall(MatDenseGetColumnVecWrite(X, n2, &cx));
1268:       PetscCall(KSPSolve_Private(ksp, cb, cx));
1269:       PetscCall(MatDenseRestoreColumnVecWrite(X, n2, &cx));
1270:       PetscCall(MatDenseRestoreColumnVecRead(B, n2, &cb));
1271:     }
1272:   }
1273:   PetscFunctionReturn(PETSC_SUCCESS);
1274: }

1276: /*@
1277:   KSPMatSolve - Solves a linear system with multiple right-hand sides stored as a `MATDENSE`.

1279:   Input Parameters:
1280: + ksp - iterative solver
1281: - B   - block of right-hand sides

1283:   Output Parameter:
1284: . X - block of solutions

1286:   Level: intermediate

1288:   Notes:
1289:   This is a stripped-down version of `KSPSolve()`, which only handles `-ksp_view`, `-ksp_converged_reason`, `-ksp_converged_rate`, and `-ksp_view_final_residual`.

1291:   Unlike with `KSPSolve()`, `B` and `X` must be different matrices.

1293: .seealso: [](ch_ksp), `KSPSolve()`, `MatMatSolve()`, `KSPMatSolveTranspose()`, `MATDENSE`, `KSPHPDDM`, `PCBJACOBI`, `PCASM`, `KSPSetMatSolveBatchSize()`
1294: @*/
1295: PetscErrorCode KSPMatSolve(KSP ksp, Mat B, Mat X)
1296: {
1297:   PetscFunctionBegin;
1298:   ksp->transpose_solve = PETSC_FALSE;
1299:   PetscCall(KSPMatSolve_Private(ksp, B, X));
1300:   PetscFunctionReturn(PETSC_SUCCESS);
1301: }

1303: /*@
1304:   KSPMatSolveTranspose - Solves a linear system with the transposed matrix with multiple right-hand sides stored as a `MATDENSE`.

1306:   Input Parameters:
1307: + ksp - iterative solver
1308: - B   - block of right-hand sides

1310:   Output Parameter:
1311: . X - block of solutions

1313:   Level: intermediate

1315:   Notes:
1316:   This is a stripped-down version of `KSPSolveTranspose()`, which only handles `-ksp_view`, `-ksp_converged_reason`, `-ksp_converged_rate`, and `-ksp_view_final_residual`.

1318:   Unlike `KSPSolveTranspose()`,
1319:   `B` and `X` must be different matrices and the transposed matrix cannot be assembled explicitly for the user.

1321: .seealso: [](ch_ksp), `KSPSolveTranspose()`, `MatMatTransposeSolve()`, `KSPMatSolve()`, `MATDENSE`, `KSPHPDDM`, `PCBJACOBI`, `PCASM`
1322: @*/
1323: PetscErrorCode KSPMatSolveTranspose(KSP ksp, Mat B, Mat X)
1324: {
1325:   PetscFunctionBegin;
1326:   ksp->transpose_solve = PETSC_TRUE;
1327:   PetscCall(KSPMatSolve_Private(ksp, B, X));
1328:   PetscFunctionReturn(PETSC_SUCCESS);
1329: }

1331: /*@
1332:   KSPSetMatSolveBatchSize - Sets the maximum number of columns treated simultaneously in `KSPMatSolve()`.

1334:   Logically Collective

1336:   Input Parameters:
1337: + ksp - the `KSP` iterative solver
1338: - bs  - batch size

1340:   Level: advanced

1342: .seealso: [](ch_ksp), `KSPMatSolve()`, `KSPGetMatSolveBatchSize()`, `-mat_mumps_icntl_27`, `-matmatmult_Bbn`
1343: @*/
1344: PetscErrorCode KSPSetMatSolveBatchSize(KSP ksp, PetscInt bs)
1345: {
1346:   PetscFunctionBegin;
1349:   ksp->nmax = bs;
1350:   PetscFunctionReturn(PETSC_SUCCESS);
1351: }

1353: /*@
1354:   KSPGetMatSolveBatchSize - Gets the maximum number of columns treated simultaneously in `KSPMatSolve()`.

1356:   Input Parameter:
1357: . ksp - iterative solver context

1359:   Output Parameter:
1360: . bs - batch size

1362:   Level: advanced

1364: .seealso: [](ch_ksp), `KSPMatSolve()`, `KSPSetMatSolveBatchSize()`, `-mat_mumps_icntl_27`, `-matmatmult_Bbn`
1365: @*/
1366: PetscErrorCode KSPGetMatSolveBatchSize(KSP ksp, PetscInt *bs)
1367: {
1368:   PetscFunctionBegin;
1370:   PetscAssertPointer(bs, 2);
1371:   *bs = ksp->nmax;
1372:   PetscFunctionReturn(PETSC_SUCCESS);
1373: }

1375: /*@
1376:   KSPResetViewers - Resets all the viewers set from the options database during `KSPSetFromOptions()`

1378:   Collective

1380:   Input Parameter:
1381: . ksp - the `KSP` iterative solver context obtained from `KSPCreate()`

1383:   Level: beginner

1385: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPSolve()`, `KSPSetFromOptions()`, `KSP`
1386: @*/
1387: PetscErrorCode KSPResetViewers(KSP ksp)
1388: {
1389:   PetscFunctionBegin;
1391:   if (!ksp) PetscFunctionReturn(PETSC_SUCCESS);
1392:   PetscCall(PetscViewerDestroy(&ksp->viewer));
1393:   PetscCall(PetscViewerDestroy(&ksp->viewerPre));
1394:   PetscCall(PetscViewerDestroy(&ksp->viewerRate));
1395:   PetscCall(PetscViewerDestroy(&ksp->viewerMat));
1396:   PetscCall(PetscViewerDestroy(&ksp->viewerPMat));
1397:   PetscCall(PetscViewerDestroy(&ksp->viewerRhs));
1398:   PetscCall(PetscViewerDestroy(&ksp->viewerSol));
1399:   PetscCall(PetscViewerDestroy(&ksp->viewerMatExp));
1400:   PetscCall(PetscViewerDestroy(&ksp->viewerEV));
1401:   PetscCall(PetscViewerDestroy(&ksp->viewerSV));
1402:   PetscCall(PetscViewerDestroy(&ksp->viewerEVExp));
1403:   PetscCall(PetscViewerDestroy(&ksp->viewerFinalRes));
1404:   PetscCall(PetscViewerDestroy(&ksp->viewerPOpExp));
1405:   PetscCall(PetscViewerDestroy(&ksp->viewerDScale));
1406:   ksp->view         = PETSC_FALSE;
1407:   ksp->viewPre      = PETSC_FALSE;
1408:   ksp->viewMat      = PETSC_FALSE;
1409:   ksp->viewPMat     = PETSC_FALSE;
1410:   ksp->viewRhs      = PETSC_FALSE;
1411:   ksp->viewSol      = PETSC_FALSE;
1412:   ksp->viewMatExp   = PETSC_FALSE;
1413:   ksp->viewEV       = PETSC_FALSE;
1414:   ksp->viewSV       = PETSC_FALSE;
1415:   ksp->viewEVExp    = PETSC_FALSE;
1416:   ksp->viewFinalRes = PETSC_FALSE;
1417:   ksp->viewPOpExp   = PETSC_FALSE;
1418:   ksp->viewDScale   = PETSC_FALSE;
1419:   PetscFunctionReturn(PETSC_SUCCESS);
1420: }

1422: /*@
1423:   KSPReset - Removes any allocated `Vec` and `Mat` from the `KSP` data structures.

1425:   Collective

1427:   Input Parameter:
1428: . ksp - iterative solver obtained from `KSPCreate()`

1430:   Level: beginner

1432:   Notes:
1433:   Any options set in the `KSP`, including those set with `KSPSetFromOptions()` remain.

1435:   Call `KSPReset()` only before you call `KSPSetOperators()` with a different sized matrix than the previous matrix used with the `KSP`.

1437: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPSolve()`, `KSP`
1438: @*/
1439: PetscErrorCode KSPReset(KSP ksp)
1440: {
1441:   PetscFunctionBegin;
1443:   if (!ksp) PetscFunctionReturn(PETSC_SUCCESS);
1444:   PetscTryTypeMethod(ksp, reset);
1445:   if (ksp->pc) PetscCall(PCReset(ksp->pc));
1446:   if (ksp->guess) {
1447:     KSPGuess guess = ksp->guess;
1448:     PetscTryTypeMethod(guess, reset);
1449:   }
1450:   PetscCall(VecDestroyVecs(ksp->nwork, &ksp->work));
1451:   PetscCall(VecDestroy(&ksp->vec_rhs));
1452:   PetscCall(VecDestroy(&ksp->vec_sol));
1453:   PetscCall(VecDestroy(&ksp->diagonal));
1454:   PetscCall(VecDestroy(&ksp->truediagonal));

1456:   ksp->setupstage = KSP_SETUP_NEW;
1457:   ksp->nmax       = PETSC_DECIDE;
1458:   PetscFunctionReturn(PETSC_SUCCESS);
1459: }

1461: /*@
1462:   KSPDestroy - Destroys a `KSP` context.

1464:   Collective

1466:   Input Parameter:
1467: . ksp - iterative solver obtained from `KSPCreate()`

1469:   Level: beginner

1471: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPSolve()`, `KSP`
1472: @*/
1473: PetscErrorCode KSPDestroy(KSP *ksp)
1474: {
1475:   PC pc;

1477:   PetscFunctionBegin;
1478:   if (!*ksp) PetscFunctionReturn(PETSC_SUCCESS);
1480:   if (--((PetscObject)*ksp)->refct > 0) {
1481:     *ksp = NULL;
1482:     PetscFunctionReturn(PETSC_SUCCESS);
1483:   }

1485:   PetscCall(PetscObjectSAWsViewOff((PetscObject)*ksp));

1487:   /*
1488:    Avoid a cascading call to PCReset(ksp->pc) from the following call:
1489:    PCReset() shouldn't be called from KSPDestroy() as it is unprotected by pc's
1490:    refcount (and may be shared, e.g., by other ksps).
1491:    */
1492:   pc         = (*ksp)->pc;
1493:   (*ksp)->pc = NULL;
1494:   PetscCall(KSPReset(*ksp));
1495:   PetscCall(KSPResetViewers(*ksp));
1496:   (*ksp)->pc = pc;
1497:   PetscTryTypeMethod(*ksp, destroy);

1499:   if ((*ksp)->transpose.use_explicittranspose) {
1500:     PetscCall(MatDestroy(&(*ksp)->transpose.AT));
1501:     PetscCall(MatDestroy(&(*ksp)->transpose.BT));
1502:     (*ksp)->transpose.reuse_transpose = PETSC_FALSE;
1503:   }

1505:   PetscCall(KSPGuessDestroy(&(*ksp)->guess));
1506:   PetscCall(DMDestroy(&(*ksp)->dm));
1507:   PetscCall(PCDestroy(&(*ksp)->pc));
1508:   PetscCall(PetscFree((*ksp)->res_hist_alloc));
1509:   PetscCall(PetscFree((*ksp)->err_hist_alloc));
1510:   if ((*ksp)->convergeddestroy) PetscCall((*(*ksp)->convergeddestroy)((*ksp)->cnvP));
1511:   PetscCall(KSPMonitorCancel(*ksp));
1512:   PetscCall(KSPConvergedReasonViewCancel(*ksp));
1513:   PetscCall(PetscHeaderDestroy(ksp));
1514:   PetscFunctionReturn(PETSC_SUCCESS);
1515: }

1517: /*@
1518:   KSPSetPCSide - Sets the preconditioning side.

1520:   Logically Collective

1522:   Input Parameter:
1523: . ksp - iterative solver obtained from `KSPCreate()`

1525:   Output Parameter:
1526: . side - the preconditioning side, where side is one of
1527: .vb
1528:   PC_LEFT      - left preconditioning (default)
1529:   PC_RIGHT     - right preconditioning
1530:   PC_SYMMETRIC - symmetric preconditioning
1531: .ve

1533:   Options Database Key:
1534: . -ksp_pc_side <right,left,symmetric> - `KSP` preconditioner side

1536:   Level: intermediate

1538:   Notes:
1539:   Left preconditioning is used by default for most Krylov methods except `KSPFGMRES` which only supports right preconditioning.

1541:   For methods changing the side of the preconditioner changes the norm type that is used, see `KSPSetNormType()`.

1543:   Symmetric preconditioning is currently available only for the `KSPQCG` method. However, note that
1544:   symmetric preconditioning can be emulated by using either right or left
1545:   preconditioning, modifying the application of the matrix (with a custom `Mat` argument to `KSPSetOperators()`,
1546:   and using a pre 'KSPSetPreSolve()` or post processing `KSPSetPostSolve()` step).

1548:   Setting the `PCSide` often affects the default norm type. See `KSPSetNormType()` for details.

1550: .seealso: [](ch_ksp), `KSPGetPCSide()`, `KSPSetNormType()`, `KSPGetNormType()`, `KSP`, `KSPSetPreSolve()`, `KSPSetPostSolve()`
1551: @*/
1552: PetscErrorCode KSPSetPCSide(KSP ksp, PCSide side)
1553: {
1554:   PetscFunctionBegin;
1557:   ksp->pc_side = ksp->pc_side_set = side;
1558:   PetscFunctionReturn(PETSC_SUCCESS);
1559: }

1561: /*@
1562:   KSPGetPCSide - Gets the preconditioning side.

1564:   Not Collective

1566:   Input Parameter:
1567: . ksp - iterative solver obtained from `KSPCreate()`

1569:   Output Parameter:
1570: . side - the preconditioning side, where side is one of
1571: .vb
1572:   PC_LEFT      - left preconditioning (default)
1573:   PC_RIGHT     - right preconditioning
1574:   PC_SYMMETRIC - symmetric preconditioning
1575: .ve

1577:   Level: intermediate

1579: .seealso: [](ch_ksp), `KSPSetPCSide()`, `KSP`
1580: @*/
1581: PetscErrorCode KSPGetPCSide(KSP ksp, PCSide *side)
1582: {
1583:   PetscFunctionBegin;
1585:   PetscAssertPointer(side, 2);
1586:   PetscCall(KSPSetUpNorms_Private(ksp, PETSC_TRUE, &ksp->normtype, &ksp->pc_side));
1587:   *side = ksp->pc_side;
1588:   PetscFunctionReturn(PETSC_SUCCESS);
1589: }

1591: /*@
1592:   KSPGetTolerances - Gets the relative, absolute, divergence, and maximum
1593:   iteration tolerances used by the default `KSP` convergence tests.

1595:   Not Collective

1597:   Input Parameter:
1598: . ksp - the Krylov subspace context

1600:   Output Parameters:
1601: + rtol   - the relative convergence tolerance
1602: . abstol - the absolute convergence tolerance
1603: . dtol   - the divergence tolerance
1604: - maxits - maximum number of iterations

1606:   Level: intermediate

1608:   Note:
1609:   The user can specify `NULL` for any parameter that is not needed.

1611: .seealso: [](ch_ksp), `KSPSetTolerances()`, `KSP`, `KSPSetMinimumIterations()`, `KSPGetMinimumIterations()`
1612: @*/
1613: PetscErrorCode KSPGetTolerances(KSP ksp, PeOp PetscReal *rtol, PeOp PetscReal *abstol, PeOp PetscReal *dtol, PeOp PetscInt *maxits)
1614: {
1615:   PetscFunctionBegin;
1617:   if (abstol) *abstol = ksp->abstol;
1618:   if (rtol) *rtol = ksp->rtol;
1619:   if (dtol) *dtol = ksp->divtol;
1620:   if (maxits) *maxits = ksp->max_it;
1621:   PetscFunctionReturn(PETSC_SUCCESS);
1622: }

1624: /*@
1625:   KSPSetTolerances - Sets the relative, absolute, divergence, and maximum
1626:   iteration tolerances used by the default `KSP` convergence testers.

1628:   Logically Collective

1630:   Input Parameters:
1631: + ksp    - the Krylov subspace context
1632: . rtol   - the relative convergence tolerance, relative decrease in the (possibly preconditioned) residual norm
1633: . abstol - the absolute convergence tolerance   absolute size of the (possibly preconditioned) residual norm
1634: . dtol   - the divergence tolerance,   amount (possibly preconditioned) residual norm can increase before `KSPConvergedDefault()` concludes that the method is diverging
1635: - maxits - maximum number of iterations to use

1637:   Options Database Keys:
1638: + -ksp_atol <abstol>   - Sets `abstol`
1639: . -ksp_rtol <rtol>     - Sets `rtol`
1640: . -ksp_divtol <dtol>   - Sets `dtol`
1641: - -ksp_max_it <maxits> - Sets `maxits`

1643:   Level: intermediate

1645:   Notes:
1646:   The tolerances are with respect to a norm of the residual of the equation $ \| b - A x^n \|$, they do not directly use the error of the equation.
1647:   The norm used depends on the `KSPNormType` that has been set with `KSPSetNormType()`, the default depends on the `KSPType` used.

1649:   All parameters must be non-negative.

1651:   Use `PETSC_CURRENT` to retain the current value of any of the parameters. The deprecated `PETSC_DEFAULT` also retains the current value (though the name is confusing).

1653:   Use `PETSC_DETERMINE` to use the default value for the given `KSP`. The default value is the value when the object's type is set.

1655:   For `dtol` and `maxits` use `PETSC_UMLIMITED` to indicate there is no upper bound on these values

1657:   See `KSPConvergedDefault()` for details how these parameters are used in the default convergence test.  See also `KSPSetConvergenceTest()`
1658:   for setting user-defined stopping criteria.

1660:   Fortran Note:
1661:   Use `PETSC_CURRENT_INTEGER`, `PETSC_CURRENT_REAL`, `PETSC_DETERMINE_INTEGER`, or `PETSC_DETERMINE_REAL`

1663: .seealso: [](ch_ksp), `KSPGetTolerances()`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSP`, `KSPSetMinimumIterations()`
1664: @*/
1665: PetscErrorCode KSPSetTolerances(KSP ksp, PetscReal rtol, PetscReal abstol, PetscReal dtol, PetscInt maxits)
1666: {
1667:   PetscFunctionBegin;

1674:   if (rtol == (PetscReal)PETSC_DETERMINE) {
1675:     ksp->rtol = ksp->default_rtol;
1676:   } else if (rtol != (PetscReal)PETSC_CURRENT) {
1677:     PetscCheck(rtol >= 0.0 && rtol < 1.0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Relative tolerance %g must be non-negative and less than 1.0", (double)rtol);
1678:     ksp->rtol = rtol;
1679:   }
1680:   if (abstol == (PetscReal)PETSC_DETERMINE) {
1681:     ksp->abstol = ksp->default_abstol;
1682:   } else if (abstol != (PetscReal)PETSC_CURRENT) {
1683:     PetscCheck(abstol >= 0.0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %g must be non-negative", (double)abstol);
1684:     ksp->abstol = abstol;
1685:   }
1686:   if (dtol == (PetscReal)PETSC_DETERMINE) {
1687:     ksp->divtol = ksp->default_divtol;
1688:   } else if (dtol == (PetscReal)PETSC_UNLIMITED) {
1689:     ksp->divtol = PETSC_MAX_REAL;
1690:   } else if (dtol != (PetscReal)PETSC_CURRENT) {
1691:     PetscCheck(dtol >= 0.0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Divergence tolerance %g must be larger than 1.0", (double)dtol);
1692:     ksp->divtol = dtol;
1693:   }
1694:   if (maxits == PETSC_DETERMINE) {
1695:     ksp->max_it = ksp->default_max_it;
1696:   } else if (maxits == PETSC_UNLIMITED) {
1697:     ksp->max_it = PETSC_INT_MAX;
1698:   } else if (maxits != PETSC_CURRENT) {
1699:     PetscCheck(maxits >= 0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", maxits);
1700:     ksp->max_it = maxits;
1701:   }
1702:   PetscFunctionReturn(PETSC_SUCCESS);
1703: }

1705: /*@
1706:   KSPSetMinimumIterations - Sets the minimum number of iterations to use, regardless of the tolerances

1708:   Logically Collective

1710:   Input Parameters:
1711: + ksp   - the Krylov subspace context
1712: - minit - minimum number of iterations to use

1714:   Options Database Key:
1715: . -ksp_min_it <minits> - Sets `minit`

1717:   Level: intermediate

1719:   Notes:
1720:   Use `KSPSetTolerances()` to set a variety of other tolerances

1722:   See `KSPConvergedDefault()` for details on how these parameters are used in the default convergence test. See also `KSPSetConvergenceTest()`
1723:   for setting user-defined stopping criteria.

1725: .seealso: [](ch_ksp), `KSPGetTolerances()`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSP`, `KSPSetTolerances()`, `KSPGetMinimumIterations()`
1726: @*/
1727: PetscErrorCode KSPSetMinimumIterations(KSP ksp, PetscInt minit)
1728: {
1729:   PetscFunctionBegin;

1733:   PetscCheck(minit >= 0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Minimum number of iterations %" PetscInt_FMT " must be non-negative", minit);
1734:   ksp->min_it = minit;
1735:   PetscFunctionReturn(PETSC_SUCCESS);
1736: }

1738: /*@
1739:   KSPGetMinimumIterations - Gets the minimum number of iterations to use, regardless of the tolerances, that was set with `KSPSetMinimumIterations()` or `-ksp_min_it`

1741:   Not Collective

1743:   Input Parameter:
1744: . ksp - the Krylov subspace context

1746:   Output Parameter:
1747: . minit - minimum number of iterations to use

1749:   Level: intermediate

1751: .seealso: [](ch_ksp), `KSPGetTolerances()`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSP`, `KSPSetTolerances()`, `KSPSetMinimumIterations()`
1752: @*/
1753: PetscErrorCode KSPGetMinimumIterations(KSP ksp, PetscInt *minit)
1754: {
1755:   PetscFunctionBegin;
1757:   PetscAssertPointer(minit, 2);

1759:   *minit = ksp->min_it;
1760:   PetscFunctionReturn(PETSC_SUCCESS);
1761: }

1763: /*@
1764:   KSPSetInitialGuessNonzero - Tells the iterative solver that the
1765:   initial guess is nonzero; otherwise `KSP` assumes the initial guess
1766:   is to be zero (and thus zeros it out before solving).

1768:   Logically Collective

1770:   Input Parameters:
1771: + ksp - iterative solver obtained from `KSPCreate()`
1772: - flg - ``PETSC_TRUE`` indicates the guess is non-zero, `PETSC_FALSE` indicates the guess is zero

1774:   Options Database Key:
1775: . -ksp_initial_guess_nonzero <true,false> - use nonzero initial guess

1777:   Level: beginner

1779:   Note:
1780:   If this is not called the X vector is zeroed in the call to `KSPSolve()`.

1782: .seealso: [](ch_ksp), `KSPGetInitialGuessNonzero()`, `KSPGuessSetType()`, `KSPGuessType`, `KSP`
1783: @*/
1784: PetscErrorCode KSPSetInitialGuessNonzero(KSP ksp, PetscBool flg)
1785: {
1786:   PetscFunctionBegin;
1789:   ksp->guess_zero = (PetscBool)!flg;
1790:   PetscFunctionReturn(PETSC_SUCCESS);
1791: }

1793: /*@
1794:   KSPGetInitialGuessNonzero - Determines whether the `KSP` solver is using
1795:   a zero initial guess.

1797:   Not Collective

1799:   Input Parameter:
1800: . ksp - iterative solver obtained from `KSPCreate()`

1802:   Output Parameter:
1803: . flag - `PETSC_TRUE` if guess is nonzero, else `PETSC_FALSE`

1805:   Level: intermediate

1807: .seealso: [](ch_ksp), `KSPSetInitialGuessNonzero()`, `KSP`
1808: @*/
1809: PetscErrorCode KSPGetInitialGuessNonzero(KSP ksp, PetscBool *flag)
1810: {
1811:   PetscFunctionBegin;
1813:   PetscAssertPointer(flag, 2);
1814:   if (ksp->guess_zero) *flag = PETSC_FALSE;
1815:   else *flag = PETSC_TRUE;
1816:   PetscFunctionReturn(PETSC_SUCCESS);
1817: }

1819: /*@
1820:   KSPSetErrorIfNotConverged - Causes `KSPSolve()` to generate an error if the solver has not converged as soon as the error is detected.

1822:   Logically Collective

1824:   Input Parameters:
1825: + ksp - iterative solver obtained from `KSPCreate()`
1826: - flg - `PETSC_TRUE` indicates you want the error generated

1828:   Options Database Key:
1829: . -ksp_error_if_not_converged <true,false> - generate an error and stop the program

1831:   Level: intermediate

1833:   Notes:
1834:   Normally PETSc continues if a linear solver fails to converge, you can call `KSPGetConvergedReason()` after a `KSPSolve()`
1835:   to determine if it has converged.

1837:   A `KSP_DIVERGED_ITS` will not generate an error in a `KSPSolve()` inside a nested linear solver

1839: .seealso: [](ch_ksp), `KSPGetErrorIfNotConverged()`, `KSP`
1840: @*/
1841: PetscErrorCode KSPSetErrorIfNotConverged(KSP ksp, PetscBool flg)
1842: {
1843:   PetscFunctionBegin;
1846:   ksp->errorifnotconverged = flg;
1847:   PetscFunctionReturn(PETSC_SUCCESS);
1848: }

1850: /*@
1851:   KSPGetErrorIfNotConverged - Will `KSPSolve()` generate an error if the solver does not converge?

1853:   Not Collective

1855:   Input Parameter:
1856: . ksp - iterative solver obtained from KSPCreate()

1858:   Output Parameter:
1859: . flag - `PETSC_TRUE` if it will generate an error, else `PETSC_FALSE`

1861:   Level: intermediate

1863: .seealso: [](ch_ksp), `KSPSetErrorIfNotConverged()`, `KSP`
1864: @*/
1865: PetscErrorCode KSPGetErrorIfNotConverged(KSP ksp, PetscBool *flag)
1866: {
1867:   PetscFunctionBegin;
1869:   PetscAssertPointer(flag, 2);
1870:   *flag = ksp->errorifnotconverged;
1871:   PetscFunctionReturn(PETSC_SUCCESS);
1872: }

1874: /*@
1875:   KSPSetInitialGuessKnoll - Tells the iterative solver to use `PCApply()` on the right hand side vector to compute the initial guess (The Knoll trick)

1877:   Logically Collective

1879:   Input Parameters:
1880: + ksp - iterative solver obtained from `KSPCreate()`
1881: - flg - `PETSC_TRUE` or `PETSC_FALSE`

1883:   Level: advanced

1885:   Developer Note:
1886:   The Knoll trick is not currently implemented using the `KSPGuess` class

1888: .seealso: [](ch_ksp), `KSPGetInitialGuessKnoll()`, `KSPSetInitialGuessNonzero()`, `KSPGetInitialGuessNonzero()`, `KSP`
1889: @*/
1890: PetscErrorCode KSPSetInitialGuessKnoll(KSP ksp, PetscBool flg)
1891: {
1892:   PetscFunctionBegin;
1895:   ksp->guess_knoll = flg;
1896:   PetscFunctionReturn(PETSC_SUCCESS);
1897: }

1899: /*@
1900:   KSPGetInitialGuessKnoll - Determines whether the `KSP` solver is using the Knoll trick (using PCApply(pc,b,...) to compute
1901:   the initial guess

1903:   Not Collective

1905:   Input Parameter:
1906: . ksp - iterative solver obtained from `KSPCreate()`

1908:   Output Parameter:
1909: . flag - `PETSC_TRUE` if using Knoll trick, else `PETSC_FALSE`

1911:   Level: advanced

1913: .seealso: [](ch_ksp), `KSPSetInitialGuessKnoll()`, `KSPSetInitialGuessNonzero()`, `KSPGetInitialGuessNonzero()`, `KSP`
1914: @*/
1915: PetscErrorCode KSPGetInitialGuessKnoll(KSP ksp, PetscBool *flag)
1916: {
1917:   PetscFunctionBegin;
1919:   PetscAssertPointer(flag, 2);
1920:   *flag = ksp->guess_knoll;
1921:   PetscFunctionReturn(PETSC_SUCCESS);
1922: }

1924: /*@
1925:   KSPGetComputeSingularValues - Gets the flag indicating whether the extreme singular
1926:   values will be calculated via a Lanczos or Arnoldi process as the linear
1927:   system is solved.

1929:   Not Collective

1931:   Input Parameter:
1932: . ksp - iterative solver obtained from `KSPCreate()`

1934:   Output Parameter:
1935: . flg - `PETSC_TRUE` or `PETSC_FALSE`

1937:   Options Database Key:
1938: . -ksp_monitor_singular_value - Activates `KSPSetComputeSingularValues()`

1940:   Level: advanced

1942:   Notes:
1943:   This option is not valid for `KSPType`.

1945:   Many users may just want to use the monitoring routine
1946:   `KSPMonitorSingularValue()` (which can be set with option `-ksp_monitor_singular_value`)
1947:   to print the singular values at each iteration of the linear solve.

1949: .seealso: [](ch_ksp), `KSPComputeExtremeSingularValues()`, `KSPMonitorSingularValue()`, `KSP`
1950: @*/
1951: PetscErrorCode KSPGetComputeSingularValues(KSP ksp, PetscBool *flg)
1952: {
1953:   PetscFunctionBegin;
1955:   PetscAssertPointer(flg, 2);
1956:   *flg = ksp->calc_sings;
1957:   PetscFunctionReturn(PETSC_SUCCESS);
1958: }

1960: /*@
1961:   KSPSetComputeSingularValues - Sets a flag so that the extreme singular
1962:   values will be calculated via a Lanczos or Arnoldi process as the linear
1963:   system is solved.

1965:   Logically Collective

1967:   Input Parameters:
1968: + ksp - iterative solver obtained from `KSPCreate()`
1969: - flg - `PETSC_TRUE` or `PETSC_FALSE`

1971:   Options Database Key:
1972: . -ksp_monitor_singular_value - Activates `KSPSetComputeSingularValues()`

1974:   Level: advanced

1976:   Notes:
1977:   This option is not valid for all iterative methods.

1979:   Many users may just want to use the monitoring routine
1980:   `KSPMonitorSingularValue()` (which can be set with option `-ksp_monitor_singular_value`)
1981:   to print the singular values at each iteration of the linear solve.

1983: .seealso: [](ch_ksp), `KSPComputeExtremeSingularValues()`, `KSPMonitorSingularValue()`, `KSP`, `KSPSetComputeRitz()`
1984: @*/
1985: PetscErrorCode KSPSetComputeSingularValues(KSP ksp, PetscBool flg)
1986: {
1987:   PetscFunctionBegin;
1990:   ksp->calc_sings = flg;
1991:   PetscFunctionReturn(PETSC_SUCCESS);
1992: }

1994: /*@
1995:   KSPGetComputeEigenvalues - Gets the flag indicating that the extreme eigenvalues
1996:   values will be calculated via a Lanczos or Arnoldi process as the linear
1997:   system is solved.

1999:   Not Collective

2001:   Input Parameter:
2002: . ksp - iterative solver obtained from `KSPCreate()`

2004:   Output Parameter:
2005: . flg - `PETSC_TRUE` or `PETSC_FALSE`

2007:   Level: advanced

2009:   Note:
2010:   Currently this option is not valid for all iterative methods.

2012: .seealso: [](ch_ksp), `KSPComputeEigenvalues()`, `KSPComputeEigenvaluesExplicitly()`, `KSP`, `KSPSetComputeRitz()`
2013: @*/
2014: PetscErrorCode KSPGetComputeEigenvalues(KSP ksp, PetscBool *flg)
2015: {
2016:   PetscFunctionBegin;
2018:   PetscAssertPointer(flg, 2);
2019:   *flg = ksp->calc_sings;
2020:   PetscFunctionReturn(PETSC_SUCCESS);
2021: }

2023: /*@
2024:   KSPSetComputeEigenvalues - Sets a flag so that the extreme eigenvalues
2025:   values will be calculated via a Lanczos or Arnoldi process as the linear
2026:   system is solved.

2028:   Logically Collective

2030:   Input Parameters:
2031: + ksp - iterative solver obtained from `KSPCreate()`
2032: - flg - `PETSC_TRUE` or `PETSC_FALSE`

2034:   Level: advanced

2036:   Note:
2037:   Currently this option is not valid for all iterative methods.

2039: .seealso: [](ch_ksp), `KSPComputeEigenvalues()`, `KSPComputeEigenvaluesExplicitly()`, `KSP`, `KSPSetComputeRitz()`
2040: @*/
2041: PetscErrorCode KSPSetComputeEigenvalues(KSP ksp, PetscBool flg)
2042: {
2043:   PetscFunctionBegin;
2046:   ksp->calc_sings = flg;
2047:   PetscFunctionReturn(PETSC_SUCCESS);
2048: }

2050: /*@
2051:   KSPSetComputeRitz - Sets a flag so that the Ritz or harmonic Ritz pairs
2052:   will be calculated via a Lanczos or Arnoldi process as the linear
2053:   system is solved.

2055:   Logically Collective

2057:   Input Parameters:
2058: + ksp - iterative solver obtained from `KSPCreate()`
2059: - flg - `PETSC_TRUE` or `PETSC_FALSE`

2061:   Level: advanced

2063:   Note:
2064:   Currently this option is only valid for the `KSPGMRES` method.

2066: .seealso: [](ch_ksp), `KSPComputeRitz()`, `KSP`, `KSPComputeEigenvalues()`, `KSPComputeExtremeSingularValues()`
2067: @*/
2068: PetscErrorCode KSPSetComputeRitz(KSP ksp, PetscBool flg)
2069: {
2070:   PetscFunctionBegin;
2073:   ksp->calc_ritz = flg;
2074:   PetscFunctionReturn(PETSC_SUCCESS);
2075: }

2077: /*@
2078:   KSPGetRhs - Gets the right-hand-side vector for the linear system to
2079:   be solved.

2081:   Not Collective

2083:   Input Parameter:
2084: . ksp - iterative solver obtained from `KSPCreate()`

2086:   Output Parameter:
2087: . r - right-hand-side vector

2089:   Level: developer

2091: .seealso: [](ch_ksp), `KSPGetSolution()`, `KSPSolve()`, `KSP`
2092: @*/
2093: PetscErrorCode KSPGetRhs(KSP ksp, Vec *r)
2094: {
2095:   PetscFunctionBegin;
2097:   PetscAssertPointer(r, 2);
2098:   *r = ksp->vec_rhs;
2099:   PetscFunctionReturn(PETSC_SUCCESS);
2100: }

2102: /*@
2103:   KSPGetSolution - Gets the location of the solution for the
2104:   linear system to be solved.

2106:   Not Collective

2108:   Input Parameter:
2109: . ksp - iterative solver obtained from `KSPCreate()`

2111:   Output Parameter:
2112: . v - solution vector

2114:   Level: developer

2116:   Note:
2117:   If this is called during a `KSPSolve()` the vector's values may not represent the solution
2118:   to the linear system.

2120: .seealso: [](ch_ksp), `KSPGetRhs()`, `KSPBuildSolution()`, `KSPSolve()`, `KSP`
2121: @*/
2122: PetscErrorCode KSPGetSolution(KSP ksp, Vec *v)
2123: {
2124:   PetscFunctionBegin;
2126:   PetscAssertPointer(v, 2);
2127:   *v = ksp->vec_sol;
2128:   PetscFunctionReturn(PETSC_SUCCESS);
2129: }

2131: /*@
2132:   KSPSetPC - Sets the preconditioner to be used to calculate the
2133:   application of the preconditioner on a vector into a `KSP`.

2135:   Collective

2137:   Input Parameters:
2138: + ksp - the `KSP` iterative solver obtained from `KSPCreate()`
2139: - pc  - the preconditioner object (if `NULL` it returns the `PC` currently held by the `KSP`)

2141:   Level: developer

2143:   Note:
2144:   This routine is almost never used since `KSP` creates its own `PC` when needed.
2145:   Use `KSPGetPC()` to retrieve the preconditioner context instead of creating a new one.

2147: .seealso: [](ch_ksp), `KSPGetPC()`, `KSP`
2148: @*/
2149: PetscErrorCode KSPSetPC(KSP ksp, PC pc)
2150: {
2151:   PetscFunctionBegin;
2153:   if (pc) {
2155:     PetscCheckSameComm(ksp, 1, pc, 2);
2156:   }
2157:   if (ksp->pc != pc && ksp->setupstage) ksp->setupstage = KSP_SETUP_NEWMATRIX;
2158:   PetscCall(PetscObjectReference((PetscObject)pc));
2159:   PetscCall(PCDestroy(&ksp->pc));
2160:   ksp->pc = pc;
2161:   PetscFunctionReturn(PETSC_SUCCESS);
2162: }

2164: PETSC_INTERN PetscErrorCode PCCreate_MPI(PC);

2166: // PetscClangLinter pragma disable: -fdoc-internal-linkage
2167: /*@C
2168:    KSPCheckPCMPI - Checks if `-mpi_linear_solver_server` is active and the `PC` should be changed to `PCMPI`

2170:    Collective, No Fortran Support

2172:    Input Parameter:
2173: .  ksp - iterative solver obtained from `KSPCreate()`

2175:    Level: developer

2177: .seealso: [](ch_ksp), `KSPSetPC()`, `KSP`, `PCMPIServerBegin()`, `PCMPIServerEnd()`
2178: @*/
2179: PETSC_INTERN PetscErrorCode KSPCheckPCMPI(KSP ksp)
2180: {
2181:   PetscBool isPCMPI;

2183:   PetscFunctionBegin;
2185:   PetscCall(PetscObjectTypeCompare((PetscObject)ksp->pc, PCMPI, &isPCMPI));
2186:   if (PCMPIServerActive && ksp->nestlevel == 0 && !isPCMPI) {
2187:     const char *prefix;
2188:     char       *found = NULL;

2190:     PetscCall(KSPGetOptionsPrefix(ksp, &prefix));
2191:     if (prefix) PetscCall(PetscStrstr(prefix, "mpi_linear_solver_server_", &found));
2192:     if (!found) PetscCall(KSPAppendOptionsPrefix(ksp, "mpi_linear_solver_server_"));
2193:     PetscCall(PetscInfo(NULL, "In MPI Linear Solver Server and detected (root) PC that must be changed to PCMPI\n"));
2194:     PetscCall(PCSetType(ksp->pc, PCMPI));
2195:   }
2196:   PetscFunctionReturn(PETSC_SUCCESS);
2197: }

2199: /*@
2200:   KSPGetPC - Returns a pointer to the preconditioner context with the `KSP`

2202:   Not Collective

2204:   Input Parameter:
2205: . ksp - iterative solver obtained from `KSPCreate()`

2207:   Output Parameter:
2208: . pc - preconditioner context

2210:   Level: beginner

2212:   Note:
2213:   The `PC` is created if it does not already exist.

2215:   Developer Note:
2216:   Calls `KSPCheckPCMPI()` to check if the `KSP` is effected by `-mpi_linear_solver_server`

2218: .seealso: [](ch_ksp), `KSPSetPC()`, `KSP`, `PC`
2219: @*/
2220: PetscErrorCode KSPGetPC(KSP ksp, PC *pc)
2221: {
2222:   PetscFunctionBegin;
2224:   PetscAssertPointer(pc, 2);
2225:   if (!ksp->pc) {
2226:     PetscCall(PCCreate(PetscObjectComm((PetscObject)ksp), &ksp->pc));
2227:     PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp->pc, (PetscObject)ksp, 0));
2228:     PetscCall(PetscObjectSetOptions((PetscObject)ksp->pc, ((PetscObject)ksp)->options));
2229:     PetscCall(PCSetKSPNestLevel(ksp->pc, ksp->nestlevel));
2230:   }
2231:   PetscCall(KSPCheckPCMPI(ksp));
2232:   *pc = ksp->pc;
2233:   PetscFunctionReturn(PETSC_SUCCESS);
2234: }

2236: /*@
2237:   KSPMonitor - runs the user provided monitor routines, if they exist

2239:   Collective

2241:   Input Parameters:
2242: + ksp   - iterative solver obtained from `KSPCreate()`
2243: . it    - iteration number
2244: - rnorm - relative norm of the residual

2246:   Level: developer

2248:   Notes:
2249:   This routine is called by the `KSP` implementations.
2250:   It does not typically need to be called by the user.

2252:   For Krylov methods that do not keep a running value of the current solution (such as `KSPGMRES`) this
2253:   cannot be called after the `KSPConvergedReason` has been set but before the final solution has been computed.

2255: .seealso: [](ch_ksp), `KSPMonitorSet()`
2256: @*/
2257: PetscErrorCode KSPMonitor(KSP ksp, PetscInt it, PetscReal rnorm)
2258: {
2259:   PetscInt i, n = ksp->numbermonitors;

2261:   PetscFunctionBegin;
2262:   for (i = 0; i < n; i++) PetscCall((*ksp->monitor[i])(ksp, it, rnorm, ksp->monitorcontext[i]));
2263:   PetscFunctionReturn(PETSC_SUCCESS);
2264: }

2266: /*@C
2267:   KSPMonitorSet - Sets an ADDITIONAL function to be called at every iteration to monitor, i.e. display in some way, perhaps by printing in the terminal,
2268:   the residual norm computed in a `KSPSolve()`

2270:   Logically Collective

2272:   Input Parameters:
2273: + ksp            - iterative solver obtained from `KSPCreate()`
2274: . monitor        - pointer to function (if this is `NULL`, it turns off monitoring
2275: . ctx            - [optional] context for private data for the monitor routine (use `NULL` if no context is needed)
2276: - monitordestroy - [optional] routine that frees monitor context (may be `NULL`), see `PetscCtxDestroyFn` for the calling sequence

2278:   Calling sequence of `monitor`:
2279: + ksp   - iterative solver obtained from `KSPCreate()`
2280: . it    - iteration number
2281: . rnorm - (estimated) 2-norm of (preconditioned) residual
2282: - ctx   - optional monitoring context, as set by `KSPMonitorSet()`

2284:   Options Database Keys:
2285: + -ksp_monitor                             - sets `KSPMonitorResidual()`
2286: . -ksp_monitor draw                        - sets `KSPMonitorResidualDraw()` and plots residual
2287: . -ksp_monitor draw::draw_lg               - sets `KSPMonitorResidualDrawLG()` and plots residual
2288: . -ksp_monitor_pause_final                 - Pauses any graphics when the solve finishes (only works for internal monitors)
2289: . -ksp_monitor_true_residual               - sets `KSPMonitorTrueResidual()`
2290: . -ksp_monitor_true_residual draw::draw_lg - sets `KSPMonitorTrueResidualDrawLG()` and plots residual
2291: . -ksp_monitor_max                         - sets `KSPMonitorTrueResidualMax()`
2292: . -ksp_monitor_singular_value              - sets `KSPMonitorSingularValue()`
2293: - -ksp_monitor_cancel                      - cancels all monitors that have been hardwired into a code by calls to `KSPMonitorSet()`, but
2294:                                              does not cancel those set via the options database.

2296:   Level: beginner

2298:   Notes:
2299:   The options database option `-ksp_monitor` and related options are the easiest way to turn on `KSP` iteration monitoring

2301:   The default is to do no monitoring.  To print the residual, or preconditioned
2302:   residual if `KSPSetNormType`(ksp,`KSP_NORM_PRECONDITIONED`) was called, use
2303:   `KSPMonitorResidual()` as the monitoring routine, with a `PETSCVIEWERASCII` as the
2304:   context.

2306:   Several different monitoring routines may be set by calling
2307:   `KSPMonitorSet()` multiple times; all will be called in the
2308:   order in which they were set.

2310:   Fortran Note:
2311:   Only a single monitor function can be set for each `KSP` object

2313: .seealso: [](ch_ksp), `KSPMonitorResidual()`, `KSPMonitorCancel()`, `KSP`, `PetscCtxDestroyFn`
2314: @*/
2315: PetscErrorCode KSPMonitorSet(KSP ksp, PetscErrorCode (*monitor)(KSP ksp, PetscInt it, PetscReal rnorm, void *ctx), void *ctx, PetscCtxDestroyFn *monitordestroy)
2316: {
2317:   PetscInt  i;
2318:   PetscBool identical;

2320:   PetscFunctionBegin;
2322:   for (i = 0; i < ksp->numbermonitors; i++) {
2323:     PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))monitor, ctx, monitordestroy, (PetscErrorCode (*)(void))ksp->monitor[i], ksp->monitorcontext[i], ksp->monitordestroy[i], &identical));
2324:     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
2325:   }
2326:   PetscCheck(ksp->numbermonitors < MAXKSPMONITORS, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Too many KSP monitors set");
2327:   ksp->monitor[ksp->numbermonitors]          = monitor;
2328:   ksp->monitordestroy[ksp->numbermonitors]   = monitordestroy;
2329:   ksp->monitorcontext[ksp->numbermonitors++] = ctx;
2330:   PetscFunctionReturn(PETSC_SUCCESS);
2331: }

2333: /*@
2334:   KSPMonitorCancel - Clears all monitors for a `KSP` object.

2336:   Logically Collective

2338:   Input Parameter:
2339: . ksp - iterative solver obtained from `KSPCreate()`

2341:   Options Database Key:
2342: . -ksp_monitor_cancel - Cancels all monitors that have been hardwired into a code by calls to `KSPMonitorSet()`, but does not cancel those set via the options database.

2344:   Level: intermediate

2346: .seealso: [](ch_ksp), `KSPMonitorResidual()`, `KSPMonitorSet()`, `KSP`
2347: @*/
2348: PetscErrorCode KSPMonitorCancel(KSP ksp)
2349: {
2350:   PetscInt i;

2352:   PetscFunctionBegin;
2354:   for (i = 0; i < ksp->numbermonitors; i++) {
2355:     if (ksp->monitordestroy[i]) PetscCall((*ksp->monitordestroy[i])(&ksp->monitorcontext[i]));
2356:   }
2357:   ksp->numbermonitors = 0;
2358:   PetscFunctionReturn(PETSC_SUCCESS);
2359: }

2361: /*@C
2362:   KSPGetMonitorContext - Gets the monitoring context, as set by `KSPMonitorSet()` for the FIRST monitor only.

2364:   Not Collective

2366:   Input Parameter:
2367: . ksp - iterative solver obtained from `KSPCreate()`

2369:   Output Parameter:
2370: . ctx - monitoring context

2372:   Level: intermediate

2374: .seealso: [](ch_ksp), `KSPMonitorResidual()`, `KSP`
2375: @*/
2376: PetscErrorCode KSPGetMonitorContext(KSP ksp, void *ctx)
2377: {
2378:   PetscFunctionBegin;
2380:   *(void **)ctx = ksp->monitorcontext[0];
2381:   PetscFunctionReturn(PETSC_SUCCESS);
2382: }

2384: /*@
2385:   KSPSetResidualHistory - Sets the array used to hold the residual history.
2386:   If set, this array will contain the residual norms computed at each
2387:   iteration of the solver.

2389:   Not Collective

2391:   Input Parameters:
2392: + ksp   - iterative solver obtained from `KSPCreate()`
2393: . a     - array to hold history
2394: . na    - size of `a`
2395: - reset - `PETSC_TRUE` indicates the history counter is reset to zero
2396:            for each new linear solve

2398:   Level: advanced

2400:   Notes:
2401:   If provided, `a` is NOT freed by PETSc so the user needs to keep track of it and destroy once the `KSP` object is destroyed.
2402:   If 'a' is `NULL` then space is allocated for the history. If 'na' `PETSC_DECIDE` or (deprecated) `PETSC_DEFAULT` then a
2403:   default array of length 10,000 is allocated.

2405:   If the array is not long enough then once the iterations is longer than the array length `KSPSolve()` stops recording the history

2407: .seealso: [](ch_ksp), `KSPGetResidualHistory()`, `KSP`
2408: @*/
2409: PetscErrorCode KSPSetResidualHistory(KSP ksp, PetscReal a[], PetscCount na, PetscBool reset)
2410: {
2411:   PetscFunctionBegin;

2414:   PetscCall(PetscFree(ksp->res_hist_alloc));
2415:   if (na != PETSC_DECIDE && na != PETSC_DEFAULT && a) {
2416:     ksp->res_hist     = a;
2417:     ksp->res_hist_max = na;
2418:   } else {
2419:     if (na != PETSC_DECIDE && na != PETSC_DEFAULT) ksp->res_hist_max = (size_t)na;
2420:     else ksp->res_hist_max = 10000; /* like default ksp->max_it */
2421:     PetscCall(PetscCalloc1(ksp->res_hist_max, &ksp->res_hist_alloc));

2423:     ksp->res_hist = ksp->res_hist_alloc;
2424:   }
2425:   ksp->res_hist_len   = 0;
2426:   ksp->res_hist_reset = reset;
2427:   PetscFunctionReturn(PETSC_SUCCESS);
2428: }

2430: /*@C
2431:   KSPGetResidualHistory - Gets the array used to hold the residual history and the number of residuals it contains.

2433:   Not Collective

2435:   Input Parameter:
2436: . ksp - iterative solver obtained from `KSPCreate()`

2438:   Output Parameters:
2439: + a  - pointer to array to hold history (or `NULL`)
2440: - na - number of used entries in a (or `NULL`). Note this has different meanings depending on the `reset` argument to `KSPSetResidualHistory()`

2442:   Level: advanced

2444:   Note:
2445:   This array is borrowed and should not be freed by the caller.

2447:   Can only be called after a `KSPSetResidualHistory()` otherwise `a` and `na` are set to `NULL` and zero

2449:   When `reset` was `PETSC_TRUE` since a residual is computed before the first iteration, the value of `na` is generally one more than the value
2450:   returned with `KSPGetIterationNumber()`.

2452:   Some Krylov methods may not compute the final residual norm when convergence is declared because the maximum number of iterations allowed has been reached.
2453:   In this situation, when `reset` was `PETSC_TRUE`, `na` will then equal the number of iterations reported with `KSPGetIterationNumber()`

2455:   Some Krylov methods (such as `KSPSTCG`), under certain circumstances, do not compute the final residual norm. In this situation, when `reset` was `PETSC_TRUE`,
2456:   `na` will then equal the number of iterations reported with `KSPGetIterationNumber()`

2458:   `KSPBCGSL` does not record the residual norms for the "subiterations" hence the results from `KSPGetResidualHistory()` and `KSPGetIterationNumber()` will be different

2460:   Fortran Note:
2461:   Call `KSPRestoreResidualHistory()` when access to the history is no longer needed.

2463: .seealso: [](ch_ksp), `KSPSetResidualHistory()`, `KSP`, `KSPGetIterationNumber()`, `KSPSTCG`, `KSPBCGSL`
2464: @*/
2465: PetscErrorCode KSPGetResidualHistory(KSP ksp, const PetscReal *a[], PetscInt *na)
2466: {
2467:   PetscFunctionBegin;
2469:   if (a) *a = ksp->res_hist;
2470:   if (na) PetscCall(PetscIntCast(ksp->res_hist_len, na));
2471:   PetscFunctionReturn(PETSC_SUCCESS);
2472: }

2474: /*@
2475:   KSPSetErrorHistory - Sets the array used to hold the error history. If set, this array will contain the error norms computed at each iteration of the solver.

2477:   Not Collective

2479:   Input Parameters:
2480: + ksp   - iterative solver obtained from `KSPCreate()`
2481: . a     - array to hold history
2482: . na    - size of `a`
2483: - reset - `PETSC_TRUE` indicates the history counter is reset to zero for each new linear solve

2485:   Level: advanced

2487:   Notes:
2488:   If provided, `a` is NOT freed by PETSc so the user needs to keep track of it and destroy once the `KSP` object is destroyed.
2489:   If 'a' is `NULL` then space is allocated for the history. If 'na' is `PETSC_DECIDE` or (deprecated) `PETSC_DEFAULT` then a default array of length 1,0000 is allocated.

2491:   If the array is not long enough then once the iterations is longer than the array length `KSPSolve()` stops recording the history

2493: .seealso: [](ch_ksp), `KSPGetErrorHistory()`, `KSPSetResidualHistory()`, `KSP`
2494: @*/
2495: PetscErrorCode KSPSetErrorHistory(KSP ksp, PetscReal a[], PetscCount na, PetscBool reset)
2496: {
2497:   PetscFunctionBegin;

2500:   PetscCall(PetscFree(ksp->err_hist_alloc));
2501:   if (na != PETSC_DECIDE && na != PETSC_DEFAULT && a) {
2502:     ksp->err_hist     = a;
2503:     ksp->err_hist_max = na;
2504:   } else {
2505:     if (na != PETSC_DECIDE && na != PETSC_DEFAULT) ksp->err_hist_max = (size_t)na;
2506:     else ksp->err_hist_max = 10000; /* like default ksp->max_it */
2507:     PetscCall(PetscCalloc1(ksp->err_hist_max, &ksp->err_hist_alloc));
2508:     ksp->err_hist = ksp->err_hist_alloc;
2509:   }
2510:   ksp->err_hist_len   = 0;
2511:   ksp->err_hist_reset = reset;
2512:   PetscFunctionReturn(PETSC_SUCCESS);
2513: }

2515: /*@C
2516:   KSPGetErrorHistory - Gets the array used to hold the error history and the number of residuals it contains.

2518:   Not Collective

2520:   Input Parameter:
2521: . ksp - iterative solver obtained from `KSPCreate()`

2523:   Output Parameters:
2524: + a  - pointer to array to hold history (or `NULL`)
2525: - na - number of used entries in a (or `NULL`)

2527:   Level: advanced

2529:   Note:
2530:   This array is borrowed and should not be freed by the caller.
2531:   Can only be called after a `KSPSetErrorHistory()` otherwise `a` and `na` are set to `NULL` and zero

2533:   Fortran Note:
2534: .vb
2535:   PetscReal, pointer :: a(:)
2536: .ve

2538: .seealso: [](ch_ksp), `KSPSetErrorHistory()`, `KSPGetResidualHistory()`, `KSP`
2539: @*/
2540: PetscErrorCode KSPGetErrorHistory(KSP ksp, const PetscReal *a[], PetscInt *na)
2541: {
2542:   PetscFunctionBegin;
2544:   if (a) *a = ksp->err_hist;
2545:   if (na) PetscCall(PetscIntCast(ksp->err_hist_len, na));
2546:   PetscFunctionReturn(PETSC_SUCCESS);
2547: }

2549: /*@
2550:   KSPComputeConvergenceRate - Compute the convergence rate for the iteration <https:/en.wikipedia.org/wiki/Coefficient_of_determination>

2552:   Not Collective

2554:   Input Parameter:
2555: . ksp - The `KSP`

2557:   Output Parameters:
2558: + cr   - The residual contraction rate
2559: . rRsq - The coefficient of determination, $R^2$, indicating the linearity of the data
2560: . ce   - The error contraction rate
2561: - eRsq - The coefficient of determination, $R^2$, indicating the linearity of the data

2563:   Level: advanced

2565:   Note:
2566:   Suppose that the residual is reduced linearly, $r_k = c^k r_0$, which means $log r_k = log r_0 + k log c$. After linear regression,
2567:   the slope is $\log c$. The coefficient of determination is given by $1 - \frac{\sum_i (y_i - f(x_i))^2}{\sum_i (y_i - \bar y)}$,

2569: .seealso: [](ch_ksp), `KSP`, `KSPConvergedRateView()`
2570: @*/
2571: PetscErrorCode KSPComputeConvergenceRate(KSP ksp, PetscReal *cr, PetscReal *rRsq, PetscReal *ce, PetscReal *eRsq)
2572: {
2573:   PetscReal const *hist;
2574:   PetscReal       *x, *y, slope, intercept, mean = 0.0, var = 0.0, res = 0.0;
2575:   PetscInt         n, k;

2577:   PetscFunctionBegin;
2578:   if (cr || rRsq) {
2579:     PetscCall(KSPGetResidualHistory(ksp, &hist, &n));
2580:     if (!n) {
2581:       if (cr) *cr = 0.0;
2582:       if (rRsq) *rRsq = -1.0;
2583:     } else {
2584:       PetscCall(PetscMalloc2(n, &x, n, &y));
2585:       for (k = 0; k < n; ++k) {
2586:         x[k] = k;
2587:         y[k] = PetscLogReal(hist[k]);
2588:         mean += y[k];
2589:       }
2590:       mean /= n;
2591:       PetscCall(PetscLinearRegression(n, x, y, &slope, &intercept));
2592:       for (k = 0; k < n; ++k) {
2593:         res += PetscSqr(y[k] - (slope * x[k] + intercept));
2594:         var += PetscSqr(y[k] - mean);
2595:       }
2596:       PetscCall(PetscFree2(x, y));
2597:       if (cr) *cr = PetscExpReal(slope);
2598:       if (rRsq) *rRsq = var < PETSC_MACHINE_EPSILON ? 0.0 : 1.0 - (res / var);
2599:     }
2600:   }
2601:   if (ce || eRsq) {
2602:     PetscCall(KSPGetErrorHistory(ksp, &hist, &n));
2603:     if (!n) {
2604:       if (ce) *ce = 0.0;
2605:       if (eRsq) *eRsq = -1.0;
2606:     } else {
2607:       PetscCall(PetscMalloc2(n, &x, n, &y));
2608:       for (k = 0; k < n; ++k) {
2609:         x[k] = k;
2610:         y[k] = PetscLogReal(hist[k]);
2611:         mean += y[k];
2612:       }
2613:       mean /= n;
2614:       PetscCall(PetscLinearRegression(n, x, y, &slope, &intercept));
2615:       for (k = 0; k < n; ++k) {
2616:         res += PetscSqr(y[k] - (slope * x[k] + intercept));
2617:         var += PetscSqr(y[k] - mean);
2618:       }
2619:       PetscCall(PetscFree2(x, y));
2620:       if (ce) *ce = PetscExpReal(slope);
2621:       if (eRsq) *eRsq = var < PETSC_MACHINE_EPSILON ? 0.0 : 1.0 - (res / var);
2622:     }
2623:   }
2624:   PetscFunctionReturn(PETSC_SUCCESS);
2625: }

2627: /*@C
2628:   KSPSetConvergenceTest - Sets the function to be used to determine convergence of `KSPSolve()`

2630:   Logically Collective

2632:   Input Parameters:
2633: + ksp      - iterative solver obtained from `KSPCreate()`
2634: . converge - pointer to the function
2635: . ctx      - context for private data for the convergence routine (may be `NULL`)
2636: - destroy  - a routine for destroying the context (may be `NULL`)

2638:   Calling sequence of `converge`:
2639: + ksp    - iterative solver obtained from `KSPCreate()`
2640: . it     - iteration number
2641: . rnorm  - (estimated) 2-norm of (preconditioned) residual
2642: . reason - the reason why it has converged or diverged
2643: - ctx    - optional convergence context, as set by `KSPSetConvergenceTest()`

2645:   Calling sequence of `destroy`:
2646: . ctx - the context

2648:   Level: advanced

2650:   Notes:
2651:   Must be called after the `KSP` type has been set so put this after
2652:   a call to `KSPSetType()`, or `KSPSetFromOptions()`.

2654:   The default convergence test, `KSPConvergedDefault()`, aborts if the
2655:   residual grows to more than 10000 times the initial residual.

2657:   The default is a combination of relative and absolute tolerances.
2658:   The residual value that is tested may be an approximation; routines
2659:   that need exact values should compute them.

2661:   In the default PETSc convergence test, the precise values of reason
2662:   are macros such as `KSP_CONVERGED_RTOL`, which are defined in petscksp.h.

2664: .seealso: [](ch_ksp), `KSP`, `KSPConvergedDefault()`, `KSPGetConvergenceContext()`, `KSPSetTolerances()`, `KSPGetConvergenceTest()`, `KSPGetAndClearConvergenceTest()`
2665: @*/
2666: PetscErrorCode KSPSetConvergenceTest(KSP ksp, PetscErrorCode (*converge)(KSP ksp, PetscInt it, PetscReal rnorm, KSPConvergedReason *reason, void *ctx), void *ctx, PetscErrorCode (*destroy)(void *ctx))
2667: {
2668:   PetscFunctionBegin;
2670:   if (ksp->convergeddestroy) PetscCall((*ksp->convergeddestroy)(ksp->cnvP));
2671:   ksp->converged        = converge;
2672:   ksp->convergeddestroy = destroy;
2673:   ksp->cnvP             = ctx;
2674:   PetscFunctionReturn(PETSC_SUCCESS);
2675: }

2677: /*@C
2678:   KSPGetConvergenceTest - Gets the function to be used to determine convergence.

2680:   Logically Collective

2682:   Input Parameter:
2683: . ksp - iterative solver obtained from `KSPCreate()`

2685:   Output Parameters:
2686: + converge - pointer to convergence test function
2687: . ctx      - context for private data for the convergence routine (may be `NULL`)
2688: - destroy  - a routine for destroying the context (may be `NULL`)

2690:   Calling sequence of `converge`:
2691: + ksp    - iterative solver obtained from `KSPCreate()`
2692: . it     - iteration number
2693: . rnorm  - (estimated) 2-norm of (preconditioned) residual
2694: . reason - the reason why it has converged or diverged
2695: - ctx    - optional convergence context, as set by `KSPSetConvergenceTest()`

2697:   Calling sequence of `destroy`:
2698: . ctx - the convergence test context

2700:   Level: advanced

2702: .seealso: [](ch_ksp), `KSP`, `KSPConvergedDefault()`, `KSPGetConvergenceContext()`, `KSPSetTolerances()`, `KSPSetConvergenceTest()`, `KSPGetAndClearConvergenceTest()`
2703: @*/
2704: PetscErrorCode KSPGetConvergenceTest(KSP ksp, PetscErrorCode (**converge)(KSP ksp, PetscInt it, PetscReal rnorm, KSPConvergedReason *reason, void *ctx), void **ctx, PetscErrorCode (**destroy)(void *ctx))
2705: {
2706:   PetscFunctionBegin;
2708:   if (converge) *converge = ksp->converged;
2709:   if (destroy) *destroy = ksp->convergeddestroy;
2710:   if (ctx) *ctx = ksp->cnvP;
2711:   PetscFunctionReturn(PETSC_SUCCESS);
2712: }

2714: /*@C
2715:   KSPGetAndClearConvergenceTest - Gets the function to be used to determine convergence. Removes the current test without calling destroy on the test context

2717:   Logically Collective

2719:   Input Parameter:
2720: . ksp - iterative solver obtained from `KSPCreate()`

2722:   Output Parameters:
2723: + converge - pointer to convergence test function
2724: . ctx      - context for private data for the convergence routine
2725: - destroy  - a routine for destroying the context

2727:   Calling sequence of `converge`:
2728: + ksp    - iterative solver obtained from `KSPCreate()`
2729: . it     - iteration number
2730: . rnorm  - (estimated) 2-norm of (preconditioned) residual
2731: . reason - the reason why it has converged or diverged
2732: - ctx    - optional convergence context, as set by `KSPSetConvergenceTest()`

2734:   Calling sequence of `destroy`:
2735: . ctx - the convergence test context

2737:   Level: advanced

2739:   Note:
2740:   This is intended to be used to allow transferring the convergence test (and its context) to another testing object (for example another `KSP`)
2741:   and then calling `KSPSetConvergenceTest()` on this original `KSP`. If you just called `KSPGetConvergenceTest()` followed
2742:   by `KSPSetConvergenceTest()` the original context information
2743:   would be destroyed and hence the transferred context would be invalid and trigger a crash on use

2745: .seealso: [](ch_ksp), `KSP`, `KSPConvergedDefault()`, `KSPGetConvergenceContext()`, `KSPSetTolerances()`, `KSPSetConvergenceTest()`, `KSPGetConvergenceTest()`
2746: @*/
2747: PetscErrorCode KSPGetAndClearConvergenceTest(KSP ksp, PetscErrorCode (**converge)(KSP ksp, PetscInt it, PetscReal rnorm, KSPConvergedReason *reason, void *ctx), void **ctx, PetscErrorCode (**destroy)(void *ctx))
2748: {
2749:   PetscFunctionBegin;
2751:   *converge             = ksp->converged;
2752:   *destroy              = ksp->convergeddestroy;
2753:   *ctx                  = ksp->cnvP;
2754:   ksp->converged        = NULL;
2755:   ksp->cnvP             = NULL;
2756:   ksp->convergeddestroy = NULL;
2757:   PetscFunctionReturn(PETSC_SUCCESS);
2758: }

2760: /*@C
2761:   KSPGetConvergenceContext - Gets the convergence context set with `KSPSetConvergenceTest()`.

2763:   Not Collective

2765:   Input Parameter:
2766: . ksp - iterative solver obtained from `KSPCreate()`

2768:   Output Parameter:
2769: . ctx - monitoring context

2771:   Level: advanced

2773: .seealso: [](ch_ksp), `KSP`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSPGetConvergenceTest()`
2774: @*/
2775: PetscErrorCode KSPGetConvergenceContext(KSP ksp, void *ctx)
2776: {
2777:   PetscFunctionBegin;
2779:   *(void **)ctx = ksp->cnvP;
2780:   PetscFunctionReturn(PETSC_SUCCESS);
2781: }

2783: /*@
2784:   KSPBuildSolution - Builds the approximate solution in a vector provided.

2786:   Collective

2788:   Input Parameter:
2789: . ksp - iterative solver obtained from `KSPCreate()`

2791:   Output Parameter:
2792:    Provide exactly one of
2793: + v - location to stash solution, optional, otherwise pass `NULL`
2794: - V - the solution is returned in this location. This vector is created internally. This vector should NOT be destroyed by the user with `VecDestroy()`.

2796:   Level: developer

2798:   Notes:
2799:   This routine can be used in one of two ways
2800: .vb
2801:       KSPBuildSolution(ksp,NULL,&V);
2802:    or
2803:       KSPBuildSolution(ksp,v,NULL); or KSPBuildSolution(ksp,v,&v);
2804: .ve
2805:   In the first case an internal vector is allocated to store the solution
2806:   (the user cannot destroy this vector). In the second case the solution
2807:   is generated in the vector that the user provides. Note that for certain
2808:   methods, such as `KSPCG`, the second case requires a copy of the solution,
2809:   while in the first case the call is essentially free since it simply
2810:   returns the vector where the solution already is stored. For some methods
2811:   like `KSPGMRES` during the solve this is a reasonably expensive operation and should only be
2812:   used if truly needed.

2814: .seealso: [](ch_ksp), `KSPGetSolution()`, `KSPBuildResidual()`, `KSP`
2815: @*/
2816: PetscErrorCode KSPBuildSolution(KSP ksp, Vec v, Vec *V)
2817: {
2818:   PetscFunctionBegin;
2820:   PetscCheck(V || v, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_WRONG, "Must provide either v or V");
2821:   if (!V) V = &v;
2822:   if (ksp->reason != KSP_CONVERGED_ITERATING) {
2823:     if (!v) PetscCall(KSPGetSolution(ksp, V));
2824:     else PetscCall(VecCopy(ksp->vec_sol, v));
2825:   } else {
2826:     PetscUseTypeMethod(ksp, buildsolution, v, V);
2827:   }
2828:   PetscFunctionReturn(PETSC_SUCCESS);
2829: }

2831: /*@
2832:   KSPBuildResidual - Builds the residual in a vector provided.

2834:   Collective

2836:   Input Parameter:
2837: . ksp - iterative solver obtained from `KSPCreate()`

2839:   Output Parameters:
2840: + t - work vector.  If not provided then one is generated.
2841: . v - optional location to stash residual.  If `v` is not provided, then a location is generated.
2842: - V - the residual

2844:   Level: advanced

2846:   Note:
2847:   Regardless of whether or not `v` is provided, the residual is
2848:   returned in `V`.

2850: .seealso: [](ch_ksp), `KSP`, `KSPBuildSolution()`
2851: @*/
2852: PetscErrorCode KSPBuildResidual(KSP ksp, Vec t, Vec v, Vec *V)
2853: {
2854:   PetscBool flag = PETSC_FALSE;
2855:   Vec       w = v, tt = t;

2857:   PetscFunctionBegin;
2859:   if (!w) PetscCall(VecDuplicate(ksp->vec_rhs, &w));
2860:   if (!tt) {
2861:     PetscCall(VecDuplicate(ksp->vec_sol, &tt));
2862:     flag = PETSC_TRUE;
2863:   }
2864:   PetscUseTypeMethod(ksp, buildresidual, tt, w, V);
2865:   if (flag) PetscCall(VecDestroy(&tt));
2866:   PetscFunctionReturn(PETSC_SUCCESS);
2867: }

2869: /*@
2870:   KSPSetDiagonalScale - Tells `KSP` to symmetrically diagonally scale the system
2871:   before solving. This actually CHANGES the matrix (and right-hand side).

2873:   Logically Collective

2875:   Input Parameters:
2876: + ksp   - the `KSP` context
2877: - scale - `PETSC_TRUE` or `PETSC_FALSE`

2879:   Options Database Keys:
2880: + -ksp_diagonal_scale     - perform a diagonal scaling before the solve
2881: - -ksp_diagonal_scale_fix - scale the matrix back AFTER the solve

2883:   Level: advanced

2885:   Notes:
2886:   Scales the matrix by  $D^{-1/2}  A  D^{-1/2}  [D^{1/2} x ] = D^{-1/2} b $
2887:   where $D_{ii}$ is $1/abs(A_{ii}) $ unless $A_{ii}$ is zero and then it is 1.

2889:   BE CAREFUL with this routine: it actually scales the matrix and right
2890:   hand side that define the system. After the system is solved the matrix
2891:   and right-hand side remain scaled unless you use `KSPSetDiagonalScaleFix()`

2893:   This should NOT be used within the `SNES` solves if you are using a line
2894:   search.

2896:   If you use this with the `PCType` `PCEISENSTAT` preconditioner than you can
2897:   use the `PCEisenstatSetNoDiagonalScaling()` option, or `-pc_eisenstat_no_diagonal_scaling`
2898:   to save some unneeded, redundant flops.

2900: .seealso: [](ch_ksp), `KSPGetDiagonalScale()`, `KSPSetDiagonalScaleFix()`, `KSP`
2901: @*/
2902: PetscErrorCode KSPSetDiagonalScale(KSP ksp, PetscBool scale)
2903: {
2904:   PetscFunctionBegin;
2907:   ksp->dscale = scale;
2908:   PetscFunctionReturn(PETSC_SUCCESS);
2909: }

2911: /*@
2912:   KSPGetDiagonalScale - Checks if `KSP` solver scales the matrix and right-hand side, that is if `KSPSetDiagonalScale()` has been called

2914:   Not Collective

2916:   Input Parameter:
2917: . ksp - the `KSP` context

2919:   Output Parameter:
2920: . scale - `PETSC_TRUE` or `PETSC_FALSE`

2922:   Level: intermediate

2924: .seealso: [](ch_ksp), `KSP`, `KSPSetDiagonalScale()`, `KSPSetDiagonalScaleFix()`
2925: @*/
2926: PetscErrorCode KSPGetDiagonalScale(KSP ksp, PetscBool *scale)
2927: {
2928:   PetscFunctionBegin;
2930:   PetscAssertPointer(scale, 2);
2931:   *scale = ksp->dscale;
2932:   PetscFunctionReturn(PETSC_SUCCESS);
2933: }

2935: /*@
2936:   KSPSetDiagonalScaleFix - Tells `KSP` to diagonally scale the system back after solving.

2938:   Logically Collective

2940:   Input Parameters:
2941: + ksp - the `KSP` context
2942: - fix - `PETSC_TRUE` to scale back after the system solve, `PETSC_FALSE` to not
2943:          rescale (default)

2945:   Level: intermediate

2947:   Notes:
2948:   Must be called after `KSPSetDiagonalScale()`

2950:   Using this will slow things down, because it rescales the matrix before and
2951:   after each linear solve. This is intended mainly for testing to allow one
2952:   to easily get back the original system to make sure the solution computed is
2953:   accurate enough.

2955: .seealso: [](ch_ksp), `KSPGetDiagonalScale()`, `KSPSetDiagonalScale()`, `KSPGetDiagonalScaleFix()`, `KSP`
2956: @*/
2957: PetscErrorCode KSPSetDiagonalScaleFix(KSP ksp, PetscBool fix)
2958: {
2959:   PetscFunctionBegin;
2962:   ksp->dscalefix = fix;
2963:   PetscFunctionReturn(PETSC_SUCCESS);
2964: }

2966: /*@
2967:   KSPGetDiagonalScaleFix - Determines if `KSP` diagonally scales the system back after solving. That is `KSPSetDiagonalScaleFix()` has been called

2969:   Not Collective

2971:   Input Parameter:
2972: . ksp - the `KSP` context

2974:   Output Parameter:
2975: . fix - `PETSC_TRUE` to scale back after the system solve, `PETSC_FALSE` to not
2976:          rescale (default)

2978:   Level: intermediate

2980: .seealso: [](ch_ksp), `KSPGetDiagonalScale()`, `KSPSetDiagonalScale()`, `KSPSetDiagonalScaleFix()`, `KSP`
2981: @*/
2982: PetscErrorCode KSPGetDiagonalScaleFix(KSP ksp, PetscBool *fix)
2983: {
2984:   PetscFunctionBegin;
2986:   PetscAssertPointer(fix, 2);
2987:   *fix = ksp->dscalefix;
2988:   PetscFunctionReturn(PETSC_SUCCESS);
2989: }

2991: /*@C
2992:   KSPSetComputeOperators - set routine to compute the linear operators

2994:   Logically Collective

2996:   Input Parameters:
2997: + ksp  - the `KSP` context
2998: . func - function to compute the operators, see `KSPComputeOperatorsFn` for the calling sequence
2999: - ctx  - optional context

3001:   Level: beginner

3003:   Notes:
3004:   `func()` will be called automatically at the very next call to `KSPSolve()`. It will NOT be called at future `KSPSolve()` calls
3005:   unless either `KSPSetComputeOperators()` or `KSPSetOperators()` is called before that `KSPSolve()` is called. This allows the same system to be solved several times
3006:   with different right-hand side functions but is a confusing API since one might expect it to be called for each `KSPSolve()`

3008:   To reuse the same preconditioner for the next `KSPSolve()` and not compute a new one based on the most recently computed matrix call `KSPSetReusePreconditioner()`

3010:   Developer Note:
3011:   Perhaps this routine and `KSPSetComputeRHS()` could be combined into a new API that makes clear when new matrices are computing without requiring call this
3012:   routine to indicate when the new matrix should be computed.

3014: .seealso: [](ch_ksp), `KSP`, `KSPSetOperators()`, `KSPSetComputeRHS()`, `DMKSPSetComputeOperators()`, `KSPSetComputeInitialGuess()`, `KSPComputeOperatorsFn`
3015: @*/
3016: PetscErrorCode KSPSetComputeOperators(KSP ksp, KSPComputeOperatorsFn *func, void *ctx)
3017: {
3018:   DM dm;

3020:   PetscFunctionBegin;
3022:   PetscCall(KSPGetDM(ksp, &dm));
3023:   PetscCall(DMKSPSetComputeOperators(dm, func, ctx));
3024:   if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX;
3025:   PetscFunctionReturn(PETSC_SUCCESS);
3026: }

3028: /*@C
3029:   KSPSetComputeRHS - set routine to compute the right-hand side of the linear system

3031:   Logically Collective

3033:   Input Parameters:
3034: + ksp  - the `KSP` context
3035: . func - function to compute the right-hand side, see `KSPComputeRHSFn` for the calling sequence
3036: - ctx  - optional context

3038:   Level: beginner

3040:   Note:
3041:   The routine you provide will be called EACH you call `KSPSolve()` to prepare the new right-hand side for that solve

3043: .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `DMKSPSetComputeRHS()`, `KSPSetComputeOperators()`, `KSPSetOperators()`, `KSPComputeRHSFn`
3044: @*/
3045: PetscErrorCode KSPSetComputeRHS(KSP ksp, KSPComputeRHSFn *func, void *ctx)
3046: {
3047:   DM dm;

3049:   PetscFunctionBegin;
3051:   PetscCall(KSPGetDM(ksp, &dm));
3052:   PetscCall(DMKSPSetComputeRHS(dm, func, ctx));
3053:   PetscFunctionReturn(PETSC_SUCCESS);
3054: }

3056: /*@C
3057:   KSPSetComputeInitialGuess - set routine to compute the initial guess of the linear system

3059:   Logically Collective

3061:   Input Parameters:
3062: + ksp  - the `KSP` context
3063: . func - function to compute the initial guess, see `KSPComputeInitialGuessFn` for calling sequence
3064: - ctx  - optional context

3066:   Level: beginner

3068:   Note:
3069:   This should only be used in conjunction with `KSPSetComputeRHS()` and `KSPSetComputeOperators()`, otherwise
3070:   call `KSPSetInitialGuessNonzero()` and set the initial guess values in the solution vector passed to `KSPSolve()` before calling the solver

3072: .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `KSPSetComputeRHS()`, `KSPSetComputeOperators()`, `DMKSPSetComputeInitialGuess()`, `KSPSetInitialGuessNonzero()`,
3073:           `KSPComputeInitialGuessFn`
3074: @*/
3075: PetscErrorCode KSPSetComputeInitialGuess(KSP ksp, KSPComputeInitialGuessFn *func, void *ctx)
3076: {
3077:   DM dm;

3079:   PetscFunctionBegin;
3081:   PetscCall(KSPGetDM(ksp, &dm));
3082:   PetscCall(DMKSPSetComputeInitialGuess(dm, func, ctx));
3083:   PetscFunctionReturn(PETSC_SUCCESS);
3084: }

3086: /*@
3087:   KSPSetUseExplicitTranspose - Determines the explicit transpose of the operator is formed in `KSPSolveTranspose()`. In some configurations (like GPUs) it may
3088:   be explicitly formed since the solve is much more efficient.

3090:   Logically Collective

3092:   Input Parameter:
3093: . ksp - the `KSP` context

3095:   Output Parameter:
3096: . flg - `PETSC_TRUE` to transpose the system in `KSPSolveTranspose()`, `PETSC_FALSE` to not transpose (default)

3098:   Level: advanced

3100: .seealso: [](ch_ksp), `KSPSolveTranspose()`, `KSP`
3101: @*/
3102: PetscErrorCode KSPSetUseExplicitTranspose(KSP ksp, PetscBool flg)
3103: {
3104:   PetscFunctionBegin;
3107:   ksp->transpose.use_explicittranspose = flg;
3108:   PetscFunctionReturn(PETSC_SUCCESS);
3109: }