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 general, be less than this.

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

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

 89:   Level: advanced

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

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

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

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

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

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

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

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

133:   Not Collective

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

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

146:   Level: advanced

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

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

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

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

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

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

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

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

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

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

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

192:   Collective

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

197:   Level: advanced

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

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

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

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

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

236:   Collective

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

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

245:   Level: intermediate

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

250:   Reusing the preconditioner reduces the time needed to form new preconditioners but may (significantly) increase the number
251:   of iterations needed for future solves depending on how much the matrix entries have changed.

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

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

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

270:   Collective

272:   Input Parameter:
273: . ksp - iterative solver obtained from `KSPCreate()`

275:   Output Parameter:
276: . flag - the boolean flag indicating if the current preconditioner should be reused

278:   Level: intermediate

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

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

296:   Collective

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

302:   Level: developer

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

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

318:   Collective

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

323:   Level: developer

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

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

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

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

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

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

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

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

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

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

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

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

408:     if (!ksp->diagonal) { /* allocate vector to hold diagonal */
409:       PetscCall(MatCreateVecs(pmat, &ksp->diagonal, NULL));
410:     }
411:     PetscCall(MatGetDiagonal(pmat, ksp->diagonal));
412:     PetscCall(VecGetLocalSize(ksp->diagonal, &n));
413:     PetscCall(VecGetArray(ksp->diagonal, &xx));
414:     for (i = 0; i < n; i++) {
415:       if (xx[i] != 0.0) xx[i] = 1.0 / PetscSqrtReal(PetscAbsScalar(xx[i]));
416:       else {
417:         xx[i]    = 1.0;
418:         zeroflag = PETSC_TRUE;
419:       }
420:     }
421:     PetscCall(VecRestoreArray(ksp->diagonal, &xx));
422:     if (zeroflag) PetscCall(PetscInfo(ksp, "Zero detected in diagonal of matrix, using 1 at those locations\n"));
423:     PetscCall(MatDiagonalScale(pmat, ksp->diagonal, ksp->diagonal));
424:     if (mat != pmat) PetscCall(MatDiagonalScale(mat, ksp->diagonal, ksp->diagonal));
425:     ksp->dscalefix2 = PETSC_FALSE;
426:   }
427:   PetscCall(PetscLogEventEnd(KSP_SetUp, ksp, ksp->vec_rhs, ksp->vec_sol, 0));
428:   PetscCall(PCSetErrorIfFailure(ksp->pc, ksp->errorifnotconverged));
429:   PetscCall(PCSetUp(ksp->pc));
430:   PetscCall(PCGetFailedReason(ksp->pc, &pcreason));
431:   /* 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 */
432:   if (pcreason) ksp->reason = KSP_DIVERGED_PC_FAILED;

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

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

448:   Collective

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

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

458:   Level: beginner

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

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

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

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

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

507:   Logically Collective

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

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

521:   Level: intermediate

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

528:   Developer Note:
529:   Should be named KSPConvergedReasonViewAdd().

531: .seealso: [](ch_ksp), `KSPConvergedReasonView()`, `KSPConvergedReasonViewFn`, `KSPConvergedReasonViewCancel()`, `PetscCtxDestroyFn`
532: @*/
533: PetscErrorCode KSPConvergedReasonViewSet(KSP ksp, KSPConvergedReasonViewFn *f, void *vctx, PetscCtxDestroyFn *reasonviewdestroy)
534: {
535:   PetscFunctionBegin;
537:   for (PetscInt i = 0; i < ksp->numberreasonviews; i++) {
538:     PetscBool identical;

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

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

554:   Collective

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

559:   Level: intermediate

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

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

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

580:   Collective

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

585:   Level: intermediate

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

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

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

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

610:   Collective

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

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

619:   Level: intermediate

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

624:   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,
625:   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)}$,

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

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

673: #include <petscdraw.h>

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

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

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

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

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

735:   PetscFunctionBegin;
736:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
737:   PetscCall(KSPGetIterationNumber(ksp, &nits));
738:   if (!nits) {
739:     PetscCall(PetscViewerASCIIPrintf(viewer, "Zero iterations in solver, cannot approximate any singular values\n"));
740:     PetscFunctionReturn(PETSC_SUCCESS);
741:   }
742:   PetscCall(KSPComputeExtremeSingularValues(ksp, &smax, &smin));
743:   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)));
744:   PetscFunctionReturn(PETSC_SUCCESS);
745: }

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

751:   PetscFunctionBegin;
752:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
753:   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");
754:   if (isascii) {
755:     Mat       A;
756:     Vec       t;
757:     PetscReal norm;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

922:   ksp->guess_zero = guess_zero;

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

927:   PetscCall(KSPConvergedReasonViewFromOptions(ksp));

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

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

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

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

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

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

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

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

999:     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]);
1000:     PetscCall(PCGetFailedReason(ksp->pc, &reason));
1001:     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]);
1002:   }
1003:   level--;
1004:   PetscFunctionReturn(PETSC_SUCCESS);
1005: }

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

1010:   Collective

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

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

1033:   Level: beginner

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

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

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

1042:   `KSPSolve()` will normally return without generating an error regardless of whether the linear system was solved or if constructing the preconditioner failed.
1043:   Call `KSPGetConvergedReason()` to determine if the solver converged or failed and why. The option -ksp_error_if_not_converged or function `KSPSetErrorIfNotConverged()`
1044:   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
1045:   it may be fine that inner solvers in the preconditioner do not converge during the solution process.

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

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

1052:   $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()`).

1054:   `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
1055:   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
1056:   direction thus the solution which is a linear combination of the search directions has no component in the nullspace(A).

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

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

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

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

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

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

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

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

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

1105:   Collective

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

1112:   Level: developer

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

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

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

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

1159:   PetscFunctionBegin;
1160:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &flg));
1161:   if (flg) {
1162:     PetscCall(PCGetOperators(ksp->pc, &A, NULL));
1163:     if (!ksp->transpose_solve) PetscCall(MatMatMult(A, X, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &R));
1164:     else PetscCall(MatTransposeMatMult(A, X, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &R));
1165:     PetscCall(MatAYPX(R, -1.0, B, SAME_NONZERO_PATTERN));
1166:     PetscCall(MatGetSize(R, NULL, &N));
1167:     PetscCall(PetscMalloc1(N, &norms));
1168:     PetscCall(MatGetColumnNorms(R, NORM_2, norms));
1169:     PetscCall(MatDestroy(&R));
1170:     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]));
1171:     PetscCall(PetscFree(norms));
1172:   }
1173:   PetscFunctionReturn(PETSC_SUCCESS);
1174: }

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

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

1225:       PetscCall(KSPConvergedReasonViewFromOptions(ksp));

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

1239:         PetscCall(KSPConvergedReasonViewFromOptions(ksp));

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

1259:       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]);
1260:       PetscCall(PCGetFailedReason(ksp->pc, &reason));
1261:       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]);
1262:     }
1263:     level--;
1264:   } else {
1265:     PetscCall(PetscInfo(ksp, "KSP type %s solving column by column\n", ((PetscObject)ksp)->type_name));
1266:     for (n2 = 0; n2 < N2; ++n2) {
1267:       PetscCall(MatDenseGetColumnVecRead(B, n2, &cb));
1268:       PetscCall(MatDenseGetColumnVecWrite(X, n2, &cx));
1269:       PetscCall(KSPSolve_Private(ksp, cb, cx));
1270:       PetscCall(MatDenseRestoreColumnVecWrite(X, n2, &cx));
1271:       PetscCall(MatDenseRestoreColumnVecRead(B, n2, &cb));
1272:     }
1273:   }
1274:   PetscFunctionReturn(PETSC_SUCCESS);
1275: }

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

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

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

1287:   Level: intermediate

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

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

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

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

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

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

1314:   Level: intermediate

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

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

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

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

1335:   Logically Collective

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

1341:   Level: advanced

1343:   Note:
1344:   Using a larger block size can improve the efficiency of the solver.

1346: .seealso: [](ch_ksp), `KSPMatSolve()`, `KSPGetMatSolveBatchSize()`, `-mat_mumps_icntl_27`, `-matmatmult_Bbn`
1347: @*/
1348: PetscErrorCode KSPSetMatSolveBatchSize(KSP ksp, PetscInt bs)
1349: {
1350:   PetscFunctionBegin;
1353:   ksp->nmax = bs;
1354:   PetscFunctionReturn(PETSC_SUCCESS);
1355: }

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

1360:   Input Parameter:
1361: . ksp - iterative solver context

1363:   Output Parameter:
1364: . bs - batch size

1366:   Level: advanced

1368: .seealso: [](ch_ksp), `KSPMatSolve()`, `KSPSetMatSolveBatchSize()`, `-mat_mumps_icntl_27`, `-matmatmult_Bbn`
1369: @*/
1370: PetscErrorCode KSPGetMatSolveBatchSize(KSP ksp, PetscInt *bs)
1371: {
1372:   PetscFunctionBegin;
1374:   PetscAssertPointer(bs, 2);
1375:   *bs = ksp->nmax;
1376:   PetscFunctionReturn(PETSC_SUCCESS);
1377: }

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

1382:   Collective

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

1387:   Level: beginner

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

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

1429:   Collective

1431:   Input Parameter:
1432: . ksp - iterative solver obtained from `KSPCreate()`

1434:   Level: intermediate

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

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

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

1460:   ksp->setupstage = KSP_SETUP_NEW;
1461:   ksp->nmax       = PETSC_DECIDE;
1462:   PetscFunctionReturn(PETSC_SUCCESS);
1463: }

1465: /*@
1466:   KSPDestroy - Destroys a `KSP` context.

1468:   Collective

1470:   Input Parameter:
1471: . ksp - iterative solver obtained from `KSPCreate()`

1473:   Level: beginner

1475: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPSolve()`, `KSP`
1476: @*/
1477: PetscErrorCode KSPDestroy(KSP *ksp)
1478: {
1479:   PC pc;

1481:   PetscFunctionBegin;
1482:   if (!*ksp) PetscFunctionReturn(PETSC_SUCCESS);
1484:   if (--((PetscObject)*ksp)->refct > 0) {
1485:     *ksp = NULL;
1486:     PetscFunctionReturn(PETSC_SUCCESS);
1487:   }

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

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

1503:   if ((*ksp)->transpose.use_explicittranspose) {
1504:     PetscCall(MatDestroy(&(*ksp)->transpose.AT));
1505:     PetscCall(MatDestroy(&(*ksp)->transpose.BT));
1506:     (*ksp)->transpose.reuse_transpose = PETSC_FALSE;
1507:   }

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

1521: /*@
1522:   KSPSetPCSide - Sets the preconditioning side.

1524:   Logically Collective

1526:   Input Parameter:
1527: . ksp - iterative solver obtained from `KSPCreate()`

1529:   Output Parameter:
1530: . side - the preconditioning side, where side is one of
1531: .vb
1532:   PC_LEFT      - left preconditioning (default)
1533:   PC_RIGHT     - right preconditioning
1534:   PC_SYMMETRIC - symmetric preconditioning
1535: .ve

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

1540:   Level: intermediate

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

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

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

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

1554: .seealso: [](ch_ksp), `KSPGetPCSide()`, `KSPSetNormType()`, `KSPGetNormType()`, `KSP`, `KSPSetPreSolve()`, `KSPSetPostSolve()`
1555: @*/
1556: PetscErrorCode KSPSetPCSide(KSP ksp, PCSide side)
1557: {
1558:   PetscFunctionBegin;
1561:   ksp->pc_side = ksp->pc_side_set = side;
1562:   PetscFunctionReturn(PETSC_SUCCESS);
1563: }

1565: /*@
1566:   KSPGetPCSide - Gets the preconditioning side.

1568:   Not Collective

1570:   Input Parameter:
1571: . ksp - iterative solver obtained from `KSPCreate()`

1573:   Output Parameter:
1574: . side - the preconditioning side, where side is one of
1575: .vb
1576:   PC_LEFT      - left preconditioning (default)
1577:   PC_RIGHT     - right preconditioning
1578:   PC_SYMMETRIC - symmetric preconditioning
1579: .ve

1581:   Level: intermediate

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

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

1599:   Not Collective

1601:   Input Parameter:
1602: . ksp - the Krylov subspace context

1604:   Output Parameters:
1605: + rtol   - the relative convergence tolerance
1606: . abstol - the absolute convergence tolerance
1607: . dtol   - the divergence tolerance
1608: - maxits - maximum number of iterations

1610:   Level: intermediate

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

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

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

1632:   Logically Collective

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

1641:   Options Database Keys:
1642: + -ksp_atol <abstol>   - Sets `abstol`
1643: . -ksp_rtol <rtol>     - Sets `rtol`
1644: . -ksp_divtol <dtol>   - Sets `dtol`
1645: - -ksp_max_it <maxits> - Sets `maxits`

1647:   Level: intermediate

1649:   Notes:
1650:   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.
1651:   The norm used depends on the `KSPNormType` that has been set with `KSPSetNormType()`, the default depends on the `KSPType` used.

1653:   All parameters must be non-negative.

1655:   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).

1657:   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.

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

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

1664:   Fortran Note:
1665:   Use `PETSC_CURRENT_INTEGER`, `PETSC_CURRENT_REAL`, `PETSC_DETERMINE_INTEGER`, or `PETSC_DETERMINE_REAL`

1667: .seealso: [](ch_ksp), `KSPGetTolerances()`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSP`, `KSPSetMinimumIterations()`
1668: @*/
1669: PetscErrorCode KSPSetTolerances(KSP ksp, PetscReal rtol, PetscReal abstol, PetscReal dtol, PetscInt maxits)
1670: {
1671:   PetscFunctionBegin;

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

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

1712:   Logically Collective

1714:   Input Parameters:
1715: + ksp   - the Krylov subspace context
1716: - minit - minimum number of iterations to use

1718:   Options Database Key:
1719: . -ksp_min_it <minits> - Sets `minit`

1721:   Level: intermediate

1723:   Notes:
1724:   Use `KSPSetTolerances()` to set a variety of other tolerances

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

1729:   If the initial residual norm is small enough solvers may return immediately without computing any improvement to the solution. Using this routine
1730:   prevents that which usually ensures the solution is changed (often minimally) from the previous solution. This option may be used with ODE integrators
1731:   to ensure the integrator does not fall into a false steady-state solution of the ODE.

1733: .seealso: [](ch_ksp), `KSPGetTolerances()`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSP`, `KSPSetTolerances()`, `KSPGetMinimumIterations()`
1734: @*/
1735: PetscErrorCode KSPSetMinimumIterations(KSP ksp, PetscInt minit)
1736: {
1737:   PetscFunctionBegin;

1741:   PetscCheck(minit >= 0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Minimum number of iterations %" PetscInt_FMT " must be non-negative", minit);
1742:   ksp->min_it = minit;
1743:   PetscFunctionReturn(PETSC_SUCCESS);
1744: }

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

1749:   Not Collective

1751:   Input Parameter:
1752: . ksp - the Krylov subspace context

1754:   Output Parameter:
1755: . minit - minimum number of iterations to use

1757:   Level: intermediate

1759: .seealso: [](ch_ksp), `KSPGetTolerances()`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSP`, `KSPSetTolerances()`, `KSPSetMinimumIterations()`
1760: @*/
1761: PetscErrorCode KSPGetMinimumIterations(KSP ksp, PetscInt *minit)
1762: {
1763:   PetscFunctionBegin;
1765:   PetscAssertPointer(minit, 2);

1767:   *minit = ksp->min_it;
1768:   PetscFunctionReturn(PETSC_SUCCESS);
1769: }

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

1776:   Logically Collective

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

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

1785:   Level: beginner

1787: .seealso: [](ch_ksp), `KSPGetInitialGuessNonzero()`, `KSPGuessSetType()`, `KSPGuessType`, `KSP`
1788: @*/
1789: PetscErrorCode KSPSetInitialGuessNonzero(KSP ksp, PetscBool flg)
1790: {
1791:   PetscFunctionBegin;
1794:   ksp->guess_zero = (PetscBool)!flg;
1795:   PetscFunctionReturn(PETSC_SUCCESS);
1796: }

1798: /*@
1799:   KSPGetInitialGuessNonzero - Determines whether the `KSP` solver is using
1800:   a zero initial guess.

1802:   Not Collective

1804:   Input Parameter:
1805: . ksp - iterative solver obtained from `KSPCreate()`

1807:   Output Parameter:
1808: . flag - `PETSC_TRUE` if guess is nonzero, else `PETSC_FALSE`

1810:   Level: intermediate

1812: .seealso: [](ch_ksp), `KSPSetInitialGuessNonzero()`, `KSP`
1813: @*/
1814: PetscErrorCode KSPGetInitialGuessNonzero(KSP ksp, PetscBool *flag)
1815: {
1816:   PetscFunctionBegin;
1818:   PetscAssertPointer(flag, 2);
1819:   if (ksp->guess_zero) *flag = PETSC_FALSE;
1820:   else *flag = PETSC_TRUE;
1821:   PetscFunctionReturn(PETSC_SUCCESS);
1822: }

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

1827:   Logically Collective

1829:   Input Parameters:
1830: + ksp - iterative solver obtained from `KSPCreate()`
1831: - flg - `PETSC_TRUE` indicates you want the error generated

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

1836:   Level: intermediate

1838:   Notes:
1839:   Normally PETSc continues if a linear solver fails to converge, you can call `KSPGetConvergedReason()` after a `KSPSolve()`
1840:   to determine if it has converged. This functionality is mostly helpful while running in a debugger (`-start_in_debugger`) to determine exactly where
1841:   the failure occurs and why.

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

1845: .seealso: [](ch_ksp), `KSPGetErrorIfNotConverged()`, `KSP`
1846: @*/
1847: PetscErrorCode KSPSetErrorIfNotConverged(KSP ksp, PetscBool flg)
1848: {
1849:   PetscFunctionBegin;
1852:   ksp->errorifnotconverged = flg;
1853:   PetscFunctionReturn(PETSC_SUCCESS);
1854: }

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

1859:   Not Collective

1861:   Input Parameter:
1862: . ksp - iterative solver obtained from KSPCreate()

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

1867:   Level: intermediate

1869: .seealso: [](ch_ksp), `KSPSetErrorIfNotConverged()`, `KSP`
1870: @*/
1871: PetscErrorCode KSPGetErrorIfNotConverged(KSP ksp, PetscBool *flag)
1872: {
1873:   PetscFunctionBegin;
1875:   PetscAssertPointer(flag, 2);
1876:   *flag = ksp->errorifnotconverged;
1877:   PetscFunctionReturn(PETSC_SUCCESS);
1878: }

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

1883:   Logically Collective

1885:   Input Parameters:
1886: + ksp - iterative solver obtained from `KSPCreate()`
1887: - flg - `PETSC_TRUE` or `PETSC_FALSE`

1889:   Level: advanced

1891:   Developer Note:
1892:   The Knoll trick is not currently implemented using the `KSPGuess` class which provides a variety of ways of computing
1893:   an initial guess based on previous solves.

1895: .seealso: [](ch_ksp), `KSPGetInitialGuessKnoll()`, `KSPGuess`, `KSPSetInitialGuessNonzero()`, `KSPGetInitialGuessNonzero()`, `KSP`
1896: @*/
1897: PetscErrorCode KSPSetInitialGuessKnoll(KSP ksp, PetscBool flg)
1898: {
1899:   PetscFunctionBegin;
1902:   ksp->guess_knoll = flg;
1903:   PetscFunctionReturn(PETSC_SUCCESS);
1904: }

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

1910:   Not Collective

1912:   Input Parameter:
1913: . ksp - iterative solver obtained from `KSPCreate()`

1915:   Output Parameter:
1916: . flag - `PETSC_TRUE` if using Knoll trick, else `PETSC_FALSE`

1918:   Level: advanced

1920: .seealso: [](ch_ksp), `KSPSetInitialGuessKnoll()`, `KSPSetInitialGuessNonzero()`, `KSPGetInitialGuessNonzero()`, `KSP`
1921: @*/
1922: PetscErrorCode KSPGetInitialGuessKnoll(KSP ksp, PetscBool *flag)
1923: {
1924:   PetscFunctionBegin;
1926:   PetscAssertPointer(flag, 2);
1927:   *flag = ksp->guess_knoll;
1928:   PetscFunctionReturn(PETSC_SUCCESS);
1929: }

1931: /*@
1932:   KSPGetComputeSingularValues - Gets the flag indicating whether the extreme singular
1933:   values will be calculated via a Lanczos or Arnoldi process as the linear
1934:   system is solved.

1936:   Not Collective

1938:   Input Parameter:
1939: . ksp - iterative solver obtained from `KSPCreate()`

1941:   Output Parameter:
1942: . flg - `PETSC_TRUE` or `PETSC_FALSE`

1944:   Options Database Key:
1945: . -ksp_monitor_singular_value - Activates `KSPSetComputeSingularValues()`

1947:   Level: advanced

1949:   Notes:
1950:   This option is not valid for `KSPType`.

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

1956: .seealso: [](ch_ksp), `KSPComputeExtremeSingularValues()`, `KSPMonitorSingularValue()`, `KSP`
1957: @*/
1958: PetscErrorCode KSPGetComputeSingularValues(KSP ksp, PetscBool *flg)
1959: {
1960:   PetscFunctionBegin;
1962:   PetscAssertPointer(flg, 2);
1963:   *flg = ksp->calc_sings;
1964:   PetscFunctionReturn(PETSC_SUCCESS);
1965: }

1967: /*@
1968:   KSPSetComputeSingularValues - Sets a flag so that the extreme singular
1969:   values will be calculated via a Lanczos or Arnoldi process as the linear
1970:   system is solved.

1972:   Logically Collective

1974:   Input Parameters:
1975: + ksp - iterative solver obtained from `KSPCreate()`
1976: - flg - `PETSC_TRUE` or `PETSC_FALSE`

1978:   Options Database Key:
1979: . -ksp_monitor_singular_value - Activates `KSPSetComputeSingularValues()`

1981:   Level: advanced

1983:   Notes:
1984:   This option is not valid for all iterative methods.

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

1990:   Consider using the excellent package SLEPc for accurate efficient computations of singular or eigenvalues.

1992: .seealso: [](ch_ksp), `KSPComputeExtremeSingularValues()`, `KSPMonitorSingularValue()`, `KSP`, `KSPSetComputeRitz()`
1993: @*/
1994: PetscErrorCode KSPSetComputeSingularValues(KSP ksp, PetscBool flg)
1995: {
1996:   PetscFunctionBegin;
1999:   ksp->calc_sings = flg;
2000:   PetscFunctionReturn(PETSC_SUCCESS);
2001: }

2003: /*@
2004:   KSPGetComputeEigenvalues - Gets the flag indicating that the extreme eigenvalues
2005:   values will be calculated via a Lanczos or Arnoldi process as the linear
2006:   system is solved.

2008:   Not Collective

2010:   Input Parameter:
2011: . ksp - iterative solver obtained from `KSPCreate()`

2013:   Output Parameter:
2014: . flg - `PETSC_TRUE` or `PETSC_FALSE`

2016:   Level: advanced

2018:   Note:
2019:   Currently this option is not valid for all iterative methods.

2021: .seealso: [](ch_ksp), `KSPComputeEigenvalues()`, `KSPComputeEigenvaluesExplicitly()`, `KSP`, `KSPSetComputeRitz()`
2022: @*/
2023: PetscErrorCode KSPGetComputeEigenvalues(KSP ksp, PetscBool *flg)
2024: {
2025:   PetscFunctionBegin;
2027:   PetscAssertPointer(flg, 2);
2028:   *flg = ksp->calc_sings;
2029:   PetscFunctionReturn(PETSC_SUCCESS);
2030: }

2032: /*@
2033:   KSPSetComputeEigenvalues - Sets a flag so that the extreme eigenvalues
2034:   values will be calculated via a Lanczos or Arnoldi process as the linear
2035:   system is solved.

2037:   Logically Collective

2039:   Input Parameters:
2040: + ksp - iterative solver obtained from `KSPCreate()`
2041: - flg - `PETSC_TRUE` or `PETSC_FALSE`

2043:   Level: advanced

2045:   Note:
2046:   Currently this option is not valid for all iterative methods.

2048:   Consider using the excellent package SLEPc for accurate efficient computations of singular or eigenvalues.

2050: .seealso: [](ch_ksp), `KSPComputeEigenvalues()`, `KSPComputeEigenvaluesExplicitly()`, `KSP`, `KSPSetComputeRitz()`
2051: @*/
2052: PetscErrorCode KSPSetComputeEigenvalues(KSP ksp, PetscBool flg)
2053: {
2054:   PetscFunctionBegin;
2057:   ksp->calc_sings = flg;
2058:   PetscFunctionReturn(PETSC_SUCCESS);
2059: }

2061: /*@
2062:   KSPSetComputeRitz - Sets a flag so that the Ritz or harmonic Ritz pairs
2063:   will be calculated via a Lanczos or Arnoldi process as the linear
2064:   system is solved.

2066:   Logically Collective

2068:   Input Parameters:
2069: + ksp - iterative solver obtained from `KSPCreate()`
2070: - flg - `PETSC_TRUE` or `PETSC_FALSE`

2072:   Level: advanced

2074:   Note:
2075:   Currently this option is only valid for the `KSPGMRES` method.

2077: .seealso: [](ch_ksp), `KSPComputeRitz()`, `KSP`, `KSPComputeEigenvalues()`, `KSPComputeExtremeSingularValues()`
2078: @*/
2079: PetscErrorCode KSPSetComputeRitz(KSP ksp, PetscBool flg)
2080: {
2081:   PetscFunctionBegin;
2084:   ksp->calc_ritz = flg;
2085:   PetscFunctionReturn(PETSC_SUCCESS);
2086: }

2088: /*@
2089:   KSPGetRhs - Gets the right-hand-side vector for the linear system to
2090:   be solved.

2092:   Not Collective

2094:   Input Parameter:
2095: . ksp - iterative solver obtained from `KSPCreate()`

2097:   Output Parameter:
2098: . r - right-hand-side vector

2100:   Level: developer

2102: .seealso: [](ch_ksp), `KSPGetSolution()`, `KSPSolve()`, `KSP`
2103: @*/
2104: PetscErrorCode KSPGetRhs(KSP ksp, Vec *r)
2105: {
2106:   PetscFunctionBegin;
2108:   PetscAssertPointer(r, 2);
2109:   *r = ksp->vec_rhs;
2110:   PetscFunctionReturn(PETSC_SUCCESS);
2111: }

2113: /*@
2114:   KSPGetSolution - Gets the location of the solution for the
2115:   linear system to be solved.

2117:   Not Collective

2119:   Input Parameter:
2120: . ksp - iterative solver obtained from `KSPCreate()`

2122:   Output Parameter:
2123: . v - solution vector

2125:   Level: developer

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

2131: .seealso: [](ch_ksp), `KSPGetRhs()`, `KSPBuildSolution()`, `KSPSolve()`, `KSP`
2132: @*/
2133: PetscErrorCode KSPGetSolution(KSP ksp, Vec *v)
2134: {
2135:   PetscFunctionBegin;
2137:   PetscAssertPointer(v, 2);
2138:   *v = ksp->vec_sol;
2139:   PetscFunctionReturn(PETSC_SUCCESS);
2140: }

2142: /*@
2143:   KSPSetPC - Sets the preconditioner to be used to calculate the
2144:   application of the preconditioner on a vector into a `KSP`.

2146:   Collective

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

2152:   Level: developer

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

2158: .seealso: [](ch_ksp), `KSPGetPC()`, `KSP`
2159: @*/
2160: PetscErrorCode KSPSetPC(KSP ksp, PC pc)
2161: {
2162:   PetscFunctionBegin;
2164:   if (pc) {
2166:     PetscCheckSameComm(ksp, 1, pc, 2);
2167:   }
2168:   if (ksp->pc != pc && ksp->setupstage) ksp->setupstage = KSP_SETUP_NEWMATRIX;
2169:   PetscCall(PetscObjectReference((PetscObject)pc));
2170:   PetscCall(PCDestroy(&ksp->pc));
2171:   ksp->pc = pc;
2172:   PetscFunctionReturn(PETSC_SUCCESS);
2173: }

2175: PETSC_INTERN PetscErrorCode PCCreate_MPI(PC);

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

2181:    Collective, No Fortran Support

2183:    Input Parameter:
2184: .  ksp - iterative solver obtained from `KSPCreate()`

2186:    Level: developer

2188: .seealso: [](ch_ksp), `KSPSetPC()`, `KSP`, `PCMPIServerBegin()`, `PCMPIServerEnd()`
2189: @*/
2190: PETSC_INTERN PetscErrorCode KSPCheckPCMPI(KSP ksp)
2191: {
2192:   PetscBool isPCMPI;

2194:   PetscFunctionBegin;
2196:   PetscCall(PetscObjectTypeCompare((PetscObject)ksp->pc, PCMPI, &isPCMPI));
2197:   if (PCMPIServerActive && ksp->nestlevel == 0 && !isPCMPI) {
2198:     const char *prefix;
2199:     char       *found = NULL;

2201:     PetscCall(KSPGetOptionsPrefix(ksp, &prefix));
2202:     if (prefix) PetscCall(PetscStrstr(prefix, "mpi_linear_solver_server_", &found));
2203:     if (!found) PetscCall(KSPAppendOptionsPrefix(ksp, "mpi_linear_solver_server_"));
2204:     PetscCall(PetscInfo(NULL, "In MPI Linear Solver Server and detected (root) PC that must be changed to PCMPI\n"));
2205:     PetscCall(PCSetType(ksp->pc, PCMPI));
2206:   }
2207:   PetscFunctionReturn(PETSC_SUCCESS);
2208: }

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

2213:   Not Collective

2215:   Input Parameter:
2216: . ksp - iterative solver obtained from `KSPCreate()`

2218:   Output Parameter:
2219: . pc - preconditioner context

2221:   Level: beginner

2223:   Note:
2224:   The `PC` is created if it does not already exist.

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

2229: .seealso: [](ch_ksp), `KSPSetPC()`, `KSP`, `PC`
2230: @*/
2231: PetscErrorCode KSPGetPC(KSP ksp, PC *pc)
2232: {
2233:   PetscFunctionBegin;
2235:   PetscAssertPointer(pc, 2);
2236:   if (!ksp->pc) {
2237:     PetscCall(PCCreate(PetscObjectComm((PetscObject)ksp), &ksp->pc));
2238:     PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp->pc, (PetscObject)ksp, 0));
2239:     PetscCall(PetscObjectSetOptions((PetscObject)ksp->pc, ((PetscObject)ksp)->options));
2240:     PetscCall(PCSetKSPNestLevel(ksp->pc, ksp->nestlevel));
2241:   }
2242:   PetscCall(KSPCheckPCMPI(ksp));
2243:   *pc = ksp->pc;
2244:   PetscFunctionReturn(PETSC_SUCCESS);
2245: }

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

2250:   Collective

2252:   Input Parameters:
2253: + ksp   - iterative solver obtained from `KSPCreate()`
2254: . it    - iteration number
2255: - rnorm - relative norm of the residual

2257:   Level: developer

2259:   Notes:
2260:   This routine is called by the `KSP` implementations.
2261:   It does not typically need to be called by the user.

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

2266: .seealso: [](ch_ksp), `KSPMonitorSet()`
2267: @*/
2268: PetscErrorCode KSPMonitor(KSP ksp, PetscInt it, PetscReal rnorm)
2269: {
2270:   PetscInt i, n = ksp->numbermonitors;

2272:   PetscFunctionBegin;
2273:   for (i = 0; i < n; i++) PetscCall((*ksp->monitor[i])(ksp, it, rnorm, ksp->monitorcontext[i]));
2274:   PetscFunctionReturn(PETSC_SUCCESS);
2275: }

2277: /*@C
2278:   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,
2279:   the residual norm computed in a `KSPSolve()`

2281:   Logically Collective

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

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

2302:   Level: beginner

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

2307:   `KSPMonitorRegister()` provides a way to associate an options database key with `KSP` monitor function.

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

2314:   Several different monitoring routines may be set by calling
2315:   `KSPMonitorSet()` multiple times; they will be called in the
2316:   order in which they were set.

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

2321: .seealso: [](ch_ksp), `KSPMonitorResidual()`, `KSPMonitorRegister()`, `KSPMonitorCancel()`, `KSP`, `PetscCtxDestroyFn`
2322: @*/
2323: PetscErrorCode KSPMonitorSet(KSP ksp, KSPMonitorFn *monitor, void *ctx, PetscCtxDestroyFn *monitordestroy)
2324: {
2325:   PetscFunctionBegin;
2327:   for (PetscInt i = 0; i < ksp->numbermonitors; i++) {
2328:     PetscBool identical;

2330:     PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))(PetscVoidFn *)monitor, ctx, monitordestroy, (PetscErrorCode (*)(void))(PetscVoidFn *)ksp->monitor[i], ksp->monitorcontext[i], ksp->monitordestroy[i], &identical));
2331:     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
2332:   }
2333:   PetscCheck(ksp->numbermonitors < MAXKSPMONITORS, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Too many KSP monitors set");
2334:   ksp->monitor[ksp->numbermonitors]          = monitor;
2335:   ksp->monitordestroy[ksp->numbermonitors]   = monitordestroy;
2336:   ksp->monitorcontext[ksp->numbermonitors++] = ctx;
2337:   PetscFunctionReturn(PETSC_SUCCESS);
2338: }

2340: /*@
2341:   KSPMonitorCancel - Clears all monitors for a `KSP` object.

2343:   Logically Collective

2345:   Input Parameter:
2346: . ksp - iterative solver obtained from `KSPCreate()`

2348:   Options Database Key:
2349: . -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.

2351:   Level: intermediate

2353: .seealso: [](ch_ksp), `KSPMonitorResidual()`, `KSPMonitorSet()`, `KSP`
2354: @*/
2355: PetscErrorCode KSPMonitorCancel(KSP ksp)
2356: {
2357:   PetscInt i;

2359:   PetscFunctionBegin;
2361:   for (i = 0; i < ksp->numbermonitors; i++) {
2362:     if (ksp->monitordestroy[i]) PetscCall((*ksp->monitordestroy[i])(&ksp->monitorcontext[i]));
2363:   }
2364:   ksp->numbermonitors = 0;
2365:   PetscFunctionReturn(PETSC_SUCCESS);
2366: }

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

2371:   Not Collective

2373:   Input Parameter:
2374: . ksp - iterative solver obtained from `KSPCreate()`

2376:   Output Parameter:
2377: . ctx - monitoring context

2379:   Level: intermediate

2381: .seealso: [](ch_ksp), `KSPMonitorResidual()`, `KSP`
2382: @*/
2383: PetscErrorCode KSPGetMonitorContext(KSP ksp, void *ctx)
2384: {
2385:   PetscFunctionBegin;
2387:   *(void **)ctx = ksp->monitorcontext[0];
2388:   PetscFunctionReturn(PETSC_SUCCESS);
2389: }

2391: /*@
2392:   KSPSetResidualHistory - Sets the array used to hold the residual history.
2393:   If set, this array will contain the residual norms computed at each
2394:   iteration of the solver.

2396:   Not Collective

2398:   Input Parameters:
2399: + ksp   - iterative solver obtained from `KSPCreate()`
2400: . a     - array to hold history
2401: . na    - size of `a`
2402: - reset - `PETSC_TRUE` indicates the history counter is reset to zero
2403:            for each new linear solve

2405:   Level: advanced

2407:   Notes:
2408:   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.
2409:   If 'a' is `NULL` then space is allocated for the history. If 'na' `PETSC_DECIDE` or (deprecated) `PETSC_DEFAULT` then a
2410:   default array of length 10,000 is allocated.

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

2414: .seealso: [](ch_ksp), `KSPGetResidualHistory()`, `KSP`
2415: @*/
2416: PetscErrorCode KSPSetResidualHistory(KSP ksp, PetscReal a[], PetscCount na, PetscBool reset)
2417: {
2418:   PetscFunctionBegin;

2421:   PetscCall(PetscFree(ksp->res_hist_alloc));
2422:   if (na != PETSC_DECIDE && na != PETSC_DEFAULT && a) {
2423:     ksp->res_hist     = a;
2424:     ksp->res_hist_max = na;
2425:   } else {
2426:     if (na != PETSC_DECIDE && na != PETSC_DEFAULT) ksp->res_hist_max = (size_t)na;
2427:     else ksp->res_hist_max = 10000; /* like default ksp->max_it */
2428:     PetscCall(PetscCalloc1(ksp->res_hist_max, &ksp->res_hist_alloc));

2430:     ksp->res_hist = ksp->res_hist_alloc;
2431:   }
2432:   ksp->res_hist_len   = 0;
2433:   ksp->res_hist_reset = reset;
2434:   PetscFunctionReturn(PETSC_SUCCESS);
2435: }

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

2440:   Not Collective

2442:   Input Parameter:
2443: . ksp - iterative solver obtained from `KSPCreate()`

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

2449:   Level: advanced

2451:   Note:
2452:   This array is borrowed and should not be freed by the caller.

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

2456:   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
2457:   returned with `KSPGetIterationNumber()`.

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

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

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

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

2470: .seealso: [](ch_ksp), `KSPSetResidualHistory()`, `KSP`, `KSPGetIterationNumber()`, `KSPSTCG`, `KSPBCGSL`
2471: @*/
2472: PetscErrorCode KSPGetResidualHistory(KSP ksp, const PetscReal *a[], PetscInt *na)
2473: {
2474:   PetscFunctionBegin;
2476:   if (a) *a = ksp->res_hist;
2477:   if (na) PetscCall(PetscIntCast(ksp->res_hist_len, na));
2478:   PetscFunctionReturn(PETSC_SUCCESS);
2479: }

2481: /*@
2482:   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.

2484:   Not Collective

2486:   Input Parameters:
2487: + ksp   - iterative solver obtained from `KSPCreate()`
2488: . a     - array to hold history
2489: . na    - size of `a`
2490: - reset - `PETSC_TRUE` indicates the history counter is reset to zero for each new linear solve

2492:   Level: advanced

2494:   Notes:
2495:   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.
2496:   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.

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

2500: .seealso: [](ch_ksp), `KSPGetErrorHistory()`, `KSPSetResidualHistory()`, `KSP`
2501: @*/
2502: PetscErrorCode KSPSetErrorHistory(KSP ksp, PetscReal a[], PetscCount na, PetscBool reset)
2503: {
2504:   PetscFunctionBegin;

2507:   PetscCall(PetscFree(ksp->err_hist_alloc));
2508:   if (na != PETSC_DECIDE && na != PETSC_DEFAULT && a) {
2509:     ksp->err_hist     = a;
2510:     ksp->err_hist_max = na;
2511:   } else {
2512:     if (na != PETSC_DECIDE && na != PETSC_DEFAULT) ksp->err_hist_max = (size_t)na;
2513:     else ksp->err_hist_max = 10000; /* like default ksp->max_it */
2514:     PetscCall(PetscCalloc1(ksp->err_hist_max, &ksp->err_hist_alloc));
2515:     ksp->err_hist = ksp->err_hist_alloc;
2516:   }
2517:   ksp->err_hist_len   = 0;
2518:   ksp->err_hist_reset = reset;
2519:   PetscFunctionReturn(PETSC_SUCCESS);
2520: }

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

2525:   Not Collective

2527:   Input Parameter:
2528: . ksp - iterative solver obtained from `KSPCreate()`

2530:   Output Parameters:
2531: + a  - pointer to array to hold history (or `NULL`)
2532: - na - number of used entries in a (or `NULL`)

2534:   Level: advanced

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

2540:   Fortran Note:
2541: .vb
2542:   PetscReal, pointer :: a(:)
2543: .ve

2545: .seealso: [](ch_ksp), `KSPSetErrorHistory()`, `KSPGetResidualHistory()`, `KSP`
2546: @*/
2547: PetscErrorCode KSPGetErrorHistory(KSP ksp, const PetscReal *a[], PetscInt *na)
2548: {
2549:   PetscFunctionBegin;
2551:   if (a) *a = ksp->err_hist;
2552:   if (na) PetscCall(PetscIntCast(ksp->err_hist_len, na));
2553:   PetscFunctionReturn(PETSC_SUCCESS);
2554: }

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

2559:   Not Collective

2561:   Input Parameter:
2562: . ksp - The `KSP`

2564:   Output Parameters:
2565: + cr   - The residual contraction rate
2566: . rRsq - The coefficient of determination, $R^2$, indicating the linearity of the data
2567: . ce   - The error contraction rate
2568: - eRsq - The coefficient of determination, $R^2$, indicating the linearity of the data

2570:   Level: advanced

2572:   Note:
2573:   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,
2574:   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)}$,

2576: .seealso: [](ch_ksp), `KSP`, `KSPConvergedRateView()`
2577: @*/
2578: PetscErrorCode KSPComputeConvergenceRate(KSP ksp, PetscReal *cr, PetscReal *rRsq, PetscReal *ce, PetscReal *eRsq)
2579: {
2580:   PetscReal const *hist;
2581:   PetscReal       *x, *y, slope, intercept, mean = 0.0, var = 0.0, res = 0.0;
2582:   PetscInt         n, k;

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

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

2637:   Logically Collective

2639:   Input Parameters:
2640: + ksp      - iterative solver obtained from `KSPCreate()`
2641: . converge - pointer to the function, see `KSPConvergenceTestFn`
2642: . ctx      - context for private data for the convergence routine (may be `NULL`)
2643: - destroy  - a routine for destroying the context (may be `NULL`)

2645:   Level: advanced

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

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

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

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

2661: .seealso: [](ch_ksp), `KSP`, `KSPConvergenceTestFn`, `KSPConvergedDefault()`, `KSPGetConvergenceContext()`, `KSPSetTolerances()`, `KSPGetConvergenceTest()`, `KSPGetAndClearConvergenceTest()`
2662: @*/
2663: PetscErrorCode KSPSetConvergenceTest(KSP ksp, KSPConvergenceTestFn *converge, void *ctx, PetscCtxDestroyFn *destroy)
2664: {
2665:   PetscFunctionBegin;
2667:   if (ksp->convergeddestroy) PetscCall((*ksp->convergeddestroy)(&ksp->cnvP));
2668:   ksp->converged        = converge;
2669:   ksp->convergeddestroy = destroy;
2670:   ksp->cnvP             = ctx;
2671:   PetscFunctionReturn(PETSC_SUCCESS);
2672: }

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

2677:   Logically Collective

2679:   Input Parameter:
2680: . ksp - iterative solver obtained from `KSPCreate()`

2682:   Output Parameters:
2683: + converge - pointer to convergence test function, see `KSPConvergenceTestFn`
2684: . ctx      - context for private data for the convergence routine (may be `NULL`)
2685: - destroy  - a routine for destroying the context (may be `NULL`)

2687:   Level: advanced

2689: .seealso: [](ch_ksp), `KSP`, `KSPConvergedDefault()`, `KSPGetConvergenceContext()`, `KSPSetTolerances()`, `KSPSetConvergenceTest()`, `KSPGetAndClearConvergenceTest()`
2690: @*/
2691: PetscErrorCode KSPGetConvergenceTest(KSP ksp, KSPConvergenceTestFn **converge, void **ctx, PetscCtxDestroyFn **destroy)
2692: {
2693:   PetscFunctionBegin;
2695:   if (converge) *converge = ksp->converged;
2696:   if (destroy) *destroy = ksp->convergeddestroy;
2697:   if (ctx) *ctx = ksp->cnvP;
2698:   PetscFunctionReturn(PETSC_SUCCESS);
2699: }

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

2704:   Logically Collective

2706:   Input Parameter:
2707: . ksp - iterative solver obtained from `KSPCreate()`

2709:   Output Parameters:
2710: + converge - pointer to convergence test function, see `KSPConvergenceTestFn`
2711: . ctx      - context for private data for the convergence routine
2712: - destroy  - a routine for destroying the context

2714:   Level: advanced

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

2722: .seealso: [](ch_ksp), `KSP`, `KSPConvergedDefault()`, `KSPGetConvergenceContext()`, `KSPSetTolerances()`, `KSPSetConvergenceTest()`, `KSPGetConvergenceTest()`
2723: @*/
2724: PetscErrorCode KSPGetAndClearConvergenceTest(KSP ksp, KSPConvergenceTestFn **converge, void **ctx, PetscCtxDestroyFn **destroy)
2725: {
2726:   PetscFunctionBegin;
2728:   *converge             = ksp->converged;
2729:   *destroy              = ksp->convergeddestroy;
2730:   *ctx                  = ksp->cnvP;
2731:   ksp->converged        = NULL;
2732:   ksp->cnvP             = NULL;
2733:   ksp->convergeddestroy = NULL;
2734:   PetscFunctionReturn(PETSC_SUCCESS);
2735: }

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

2740:   Not Collective

2742:   Input Parameter:
2743: . ksp - iterative solver obtained from `KSPCreate()`

2745:   Output Parameter:
2746: . ctx - monitoring context

2748:   Level: advanced

2750: .seealso: [](ch_ksp), `KSP`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSPGetConvergenceTest()`
2751: @*/
2752: PetscErrorCode KSPGetConvergenceContext(KSP ksp, void *ctx)
2753: {
2754:   PetscFunctionBegin;
2756:   *(void **)ctx = ksp->cnvP;
2757:   PetscFunctionReturn(PETSC_SUCCESS);
2758: }

2760: /*@
2761:   KSPBuildSolution - Builds the approximate solution in a vector provided.

2763:   Collective

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

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

2773:   Level: developer

2775:   Notes:
2776:   This routine can be used in one of two ways
2777: .vb
2778:       KSPBuildSolution(ksp,NULL,&V);
2779:    or
2780:       KSPBuildSolution(ksp,v,NULL); or KSPBuildSolution(ksp,v,&v);
2781: .ve
2782:   In the first case an internal vector is allocated to store the solution
2783:   (the user cannot destroy this vector). In the second case the solution
2784:   is generated in the vector that the user provides. Note that for certain
2785:   methods, such as `KSPCG`, the second case requires a copy of the solution,
2786:   while in the first case the call is essentially free since it simply
2787:   returns the vector where the solution already is stored. For some methods
2788:   like `KSPGMRES` during the solve this is a reasonably expensive operation and should only be
2789:   used if truly needed.

2791: .seealso: [](ch_ksp), `KSPGetSolution()`, `KSPBuildResidual()`, `KSP`
2792: @*/
2793: PetscErrorCode KSPBuildSolution(KSP ksp, Vec v, Vec *V)
2794: {
2795:   PetscFunctionBegin;
2797:   PetscCheck(V || v, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_WRONG, "Must provide either v or V");
2798:   if (!V) V = &v;
2799:   if (ksp->reason != KSP_CONVERGED_ITERATING) {
2800:     if (!v) PetscCall(KSPGetSolution(ksp, V));
2801:     else PetscCall(VecCopy(ksp->vec_sol, v));
2802:   } else {
2803:     PetscUseTypeMethod(ksp, buildsolution, v, V);
2804:   }
2805:   PetscFunctionReturn(PETSC_SUCCESS);
2806: }

2808: /*@
2809:   KSPBuildResidual - Builds the residual in a vector provided.

2811:   Collective

2813:   Input Parameter:
2814: . ksp - iterative solver obtained from `KSPCreate()`

2816:   Output Parameters:
2817: + t - work vector.  If not provided then one is generated.
2818: . v - optional location to stash residual.  If `v` is not provided, then a location is generated.
2819: - V - the residual

2821:   Level: advanced

2823:   Note:
2824:   Regardless of whether or not `v` is provided, the residual is
2825:   returned in `V`.

2827: .seealso: [](ch_ksp), `KSP`, `KSPBuildSolution()`
2828: @*/
2829: PetscErrorCode KSPBuildResidual(KSP ksp, Vec t, Vec v, Vec *V)
2830: {
2831:   PetscBool flag = PETSC_FALSE;
2832:   Vec       w = v, tt = t;

2834:   PetscFunctionBegin;
2836:   if (!w) PetscCall(VecDuplicate(ksp->vec_rhs, &w));
2837:   if (!tt) {
2838:     PetscCall(VecDuplicate(ksp->vec_sol, &tt));
2839:     flag = PETSC_TRUE;
2840:   }
2841:   PetscUseTypeMethod(ksp, buildresidual, tt, w, V);
2842:   if (flag) PetscCall(VecDestroy(&tt));
2843:   PetscFunctionReturn(PETSC_SUCCESS);
2844: }

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

2850:   Logically Collective

2852:   Input Parameters:
2853: + ksp   - the `KSP` context
2854: - scale - `PETSC_TRUE` or `PETSC_FALSE`

2856:   Options Database Keys:
2857: + -ksp_diagonal_scale     - perform a diagonal scaling before the solve
2858: - -ksp_diagonal_scale_fix - scale the matrix back AFTER the solve

2860:   Level: advanced

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

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

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

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

2877: .seealso: [](ch_ksp), `KSPGetDiagonalScale()`, `KSPSetDiagonalScaleFix()`, `KSP`
2878: @*/
2879: PetscErrorCode KSPSetDiagonalScale(KSP ksp, PetscBool scale)
2880: {
2881:   PetscFunctionBegin;
2884:   ksp->dscale = scale;
2885:   PetscFunctionReturn(PETSC_SUCCESS);
2886: }

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

2891:   Not Collective

2893:   Input Parameter:
2894: . ksp - the `KSP` context

2896:   Output Parameter:
2897: . scale - `PETSC_TRUE` or `PETSC_FALSE`

2899:   Level: intermediate

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

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

2915:   Logically Collective

2917:   Input Parameters:
2918: + ksp - the `KSP` context
2919: - fix - `PETSC_TRUE` to scale back after the system solve, `PETSC_FALSE` to not
2920:          rescale (default)

2922:   Level: intermediate

2924:   Notes:
2925:   Must be called after `KSPSetDiagonalScale()`

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

2932: .seealso: [](ch_ksp), `KSPGetDiagonalScale()`, `KSPSetDiagonalScale()`, `KSPGetDiagonalScaleFix()`, `KSP`
2933: @*/
2934: PetscErrorCode KSPSetDiagonalScaleFix(KSP ksp, PetscBool fix)
2935: {
2936:   PetscFunctionBegin;
2939:   ksp->dscalefix = fix;
2940:   PetscFunctionReturn(PETSC_SUCCESS);
2941: }

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

2946:   Not Collective

2948:   Input Parameter:
2949: . ksp - the `KSP` context

2951:   Output Parameter:
2952: . fix - `PETSC_TRUE` to scale back after the system solve, `PETSC_FALSE` to not
2953:          rescale (default)

2955:   Level: intermediate

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

2968: /*@C
2969:   KSPSetComputeOperators - set routine to compute the linear operators

2971:   Logically Collective

2973:   Input Parameters:
2974: + ksp  - the `KSP` context
2975: . func - function to compute the operators, see `KSPComputeOperatorsFn` for the calling sequence
2976: - ctx  - optional context

2978:   Level: beginner

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

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

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

2991: .seealso: [](ch_ksp), `KSP`, `KSPSetOperators()`, `KSPSetComputeRHS()`, `DMKSPSetComputeOperators()`, `KSPSetComputeInitialGuess()`, `KSPComputeOperatorsFn`
2992: @*/
2993: PetscErrorCode KSPSetComputeOperators(KSP ksp, KSPComputeOperatorsFn *func, void *ctx)
2994: {
2995:   DM dm;

2997:   PetscFunctionBegin;
2999:   PetscCall(KSPGetDM(ksp, &dm));
3000:   PetscCall(DMKSPSetComputeOperators(dm, func, ctx));
3001:   if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX;
3002:   PetscFunctionReturn(PETSC_SUCCESS);
3003: }

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

3008:   Logically Collective

3010:   Input Parameters:
3011: + ksp  - the `KSP` context
3012: . func - function to compute the right-hand side, see `KSPComputeRHSFn` for the calling sequence
3013: - ctx  - optional context

3015:   Level: beginner

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

3020: .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `DMKSPSetComputeRHS()`, `KSPSetComputeOperators()`, `KSPSetOperators()`, `KSPComputeRHSFn`
3021: @*/
3022: PetscErrorCode KSPSetComputeRHS(KSP ksp, KSPComputeRHSFn *func, void *ctx)
3023: {
3024:   DM dm;

3026:   PetscFunctionBegin;
3028:   PetscCall(KSPGetDM(ksp, &dm));
3029:   PetscCall(DMKSPSetComputeRHS(dm, func, ctx));
3030:   PetscFunctionReturn(PETSC_SUCCESS);
3031: }

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

3036:   Logically Collective

3038:   Input Parameters:
3039: + ksp  - the `KSP` context
3040: . func - function to compute the initial guess, see `KSPComputeInitialGuessFn` for calling sequence
3041: - ctx  - optional context

3043:   Level: beginner

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

3049: .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `KSPSetComputeRHS()`, `KSPSetComputeOperators()`, `DMKSPSetComputeInitialGuess()`, `KSPSetInitialGuessNonzero()`,
3050:           `KSPComputeInitialGuessFn`
3051: @*/
3052: PetscErrorCode KSPSetComputeInitialGuess(KSP ksp, KSPComputeInitialGuessFn *func, void *ctx)
3053: {
3054:   DM dm;

3056:   PetscFunctionBegin;
3058:   PetscCall(KSPGetDM(ksp, &dm));
3059:   PetscCall(DMKSPSetComputeInitialGuess(dm, func, ctx));
3060:   PetscFunctionReturn(PETSC_SUCCESS);
3061: }

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

3067:   Logically Collective

3069:   Input Parameter:
3070: . ksp - the `KSP` context

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

3075:   Level: advanced

3077: .seealso: [](ch_ksp), `KSPSolveTranspose()`, `KSP`
3078: @*/
3079: PetscErrorCode KSPSetUseExplicitTranspose(KSP ksp, PetscBool flg)
3080: {
3081:   PetscFunctionBegin;
3084:   ksp->transpose.use_explicittranspose = flg;
3085:   PetscFunctionReturn(PETSC_SUCCESS);
3086: }