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(PetscOptionsPushCreateViewerOff(PETSC_FALSE));
764:     PetscCall(VecViewFromOptions(t, (PetscObject)ksp, "-ksp_view_final_residual_vec"));
765:     PetscCall(PetscOptionsPopCreateViewerOff());
766:     PetscCall(VecNorm(t, NORM_2, &norm));
767:     PetscCall(VecDestroy(&t));
768:     PetscCall(PetscViewerASCIIPrintf(viewer, "KSP final norm of residual %g\n", (double)norm));
769:   }
770:   PetscFunctionReturn(PETSC_SUCCESS);
771: }

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

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

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

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

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

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

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

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

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

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

844:   if (ksp->guess) {
845:     PetscObjectState ostate, state;

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

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

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

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

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

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

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

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

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

924:   ksp->guess_zero = guess_zero;

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

929:   PetscCall(KSPConvergedReasonViewFromOptions(ksp));

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

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

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

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

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

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

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

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

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

1009: /*@
1010:   KSPSolve - Solves a linear system associated with `KSP` object

1012:   Collective

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

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

1037:   Level: beginner

1039:   Notes:
1040:   See `KSPSetFromOptions()` for additional options database keys that affect `KSPSolve()`

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

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

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

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

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

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

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

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

1066:   Developer Notes:
1067:   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
1068:   $ 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
1069:   such as diagonal scaling we cannot apply the inverse of the preconditioner to a vector and thus cannot compute $nullspace(AB)$.

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

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

1076:   Understanding Convergence\:
1077:   The manual pages `KSPMonitorSet()`, `KSPComputeEigenvalues()`, and
1078:   `KSPComputeEigenvaluesExplicitly()` provide information on additional
1079:   options to monitor convergence and print eigenvalue information.

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

1089:   PetscFunctionBegin;
1093:   ksp->transpose_solve = PETSC_FALSE;
1094:   PetscCall(KSPSolve_Private(ksp, b, x));
1095:   PetscCall(PetscObjectTypeCompare((PetscObject)ksp->pc, PCMPI, &isPCMPI));
1096:   if (PCMPIServerActive && isPCMPI) {
1097:     KSP subksp;

1099:     PetscCall(PCMPIGetKSP(ksp->pc, &subksp));
1100:     ksp->its    = subksp->its;
1101:     ksp->reason = subksp->reason;
1102:   }
1103:   PetscFunctionReturn(PETSC_SUCCESS);
1104: }

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

1109:   Collective

1111:   Input Parameters:
1112: + ksp - iterative solver obtained from `KSPCreate()`
1113: . b   - right-hand side vector
1114: - x   - solution vector

1116:   Level: developer

1118:   Note:
1119:   For complex numbers this solve the non-Hermitian transpose system.

1121:   Developer Note:
1122:   We need to implement a `KSPSolveHermitianTranspose()`

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

1156: static PetscErrorCode KSPViewFinalMatResidual_Internal(KSP ksp, Mat B, Mat X, PetscViewer viewer, PetscViewerFormat format, PetscInt shift)
1157: {
1158:   Mat        A, R;
1159:   PetscReal *norms;
1160:   PetscInt   i, N;
1161:   PetscBool  flg;

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

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

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

1229:       PetscCall(KSPConvergedReasonViewFromOptions(ksp));

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

1243:         PetscCall(KSPConvergedReasonViewFromOptions(ksp));

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

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

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

1284:   Input Parameters:
1285: + ksp - iterative solver
1286: - B   - block of right-hand sides

1288:   Output Parameter:
1289: . X - block of solutions

1291:   Level: intermediate

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

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

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

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

1311:   Input Parameters:
1312: + ksp - iterative solver
1313: - B   - block of right-hand sides

1315:   Output Parameter:
1316: . X - block of solutions

1318:   Level: intermediate

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

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

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

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

1339:   Logically Collective

1341:   Input Parameters:
1342: + ksp - the `KSP` iterative solver
1343: - bs  - batch size

1345:   Level: advanced

1347:   Note:
1348:   Using a larger block size can improve the efficiency of the solver.

1350: .seealso: [](ch_ksp), `KSPMatSolve()`, `KSPGetMatSolveBatchSize()`, `-mat_mumps_icntl_27`, `-matmatmult_Bbn`
1351: @*/
1352: PetscErrorCode KSPSetMatSolveBatchSize(KSP ksp, PetscInt bs)
1353: {
1354:   PetscFunctionBegin;
1357:   ksp->nmax = bs;
1358:   PetscFunctionReturn(PETSC_SUCCESS);
1359: }

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

1364:   Input Parameter:
1365: . ksp - iterative solver context

1367:   Output Parameter:
1368: . bs - batch size

1370:   Level: advanced

1372: .seealso: [](ch_ksp), `KSPMatSolve()`, `KSPSetMatSolveBatchSize()`, `-mat_mumps_icntl_27`, `-matmatmult_Bbn`
1373: @*/
1374: PetscErrorCode KSPGetMatSolveBatchSize(KSP ksp, PetscInt *bs)
1375: {
1376:   PetscFunctionBegin;
1378:   PetscAssertPointer(bs, 2);
1379:   *bs = ksp->nmax;
1380:   PetscFunctionReturn(PETSC_SUCCESS);
1381: }

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

1386:   Collective

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

1391:   Level: beginner

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

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

1433:   Collective

1435:   Input Parameter:
1436: . ksp - iterative solver obtained from `KSPCreate()`

1438:   Level: intermediate

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

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

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

1464:   ksp->setupstage = KSP_SETUP_NEW;
1465:   ksp->nmax       = PETSC_DECIDE;
1466:   PetscFunctionReturn(PETSC_SUCCESS);
1467: }

1469: /*@
1470:   KSPDestroy - Destroys a `KSP` context.

1472:   Collective

1474:   Input Parameter:
1475: . ksp - iterative solver obtained from `KSPCreate()`

1477:   Level: beginner

1479: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPSolve()`, `KSP`
1480: @*/
1481: PetscErrorCode KSPDestroy(KSP *ksp)
1482: {
1483:   PC pc;

1485:   PetscFunctionBegin;
1486:   if (!*ksp) PetscFunctionReturn(PETSC_SUCCESS);
1488:   if (--((PetscObject)*ksp)->refct > 0) {
1489:     *ksp = NULL;
1490:     PetscFunctionReturn(PETSC_SUCCESS);
1491:   }

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

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

1507:   if ((*ksp)->transpose.use_explicittranspose) {
1508:     PetscCall(MatDestroy(&(*ksp)->transpose.AT));
1509:     PetscCall(MatDestroy(&(*ksp)->transpose.BT));
1510:     (*ksp)->transpose.reuse_transpose = PETSC_FALSE;
1511:   }

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

1525: /*@
1526:   KSPSetPCSide - Sets the preconditioning side.

1528:   Logically Collective

1530:   Input Parameter:
1531: . ksp - iterative solver obtained from `KSPCreate()`

1533:   Output Parameter:
1534: . side - the preconditioning side, where side is one of
1535: .vb
1536:   PC_LEFT      - left preconditioning (default)
1537:   PC_RIGHT     - right preconditioning
1538:   PC_SYMMETRIC - symmetric preconditioning
1539: .ve

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

1544:   Level: intermediate

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

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

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

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

1558: .seealso: [](ch_ksp), `KSPGetPCSide()`, `KSPSetNormType()`, `KSPGetNormType()`, `KSP`, `KSPSetPreSolve()`, `KSPSetPostSolve()`
1559: @*/
1560: PetscErrorCode KSPSetPCSide(KSP ksp, PCSide side)
1561: {
1562:   PetscFunctionBegin;
1565:   ksp->pc_side = ksp->pc_side_set = side;
1566:   PetscFunctionReturn(PETSC_SUCCESS);
1567: }

1569: /*@
1570:   KSPGetPCSide - Gets the preconditioning side.

1572:   Not Collective

1574:   Input Parameter:
1575: . ksp - iterative solver obtained from `KSPCreate()`

1577:   Output Parameter:
1578: . side - the preconditioning side, where side is one of
1579: .vb
1580:   PC_LEFT      - left preconditioning (default)
1581:   PC_RIGHT     - right preconditioning
1582:   PC_SYMMETRIC - symmetric preconditioning
1583: .ve

1585:   Level: intermediate

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

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

1603:   Not Collective

1605:   Input Parameter:
1606: . ksp - the Krylov subspace context

1608:   Output Parameters:
1609: + rtol   - the relative convergence tolerance
1610: . abstol - the absolute convergence tolerance
1611: . dtol   - the divergence tolerance
1612: - maxits - maximum number of iterations

1614:   Level: intermediate

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

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

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

1636:   Logically Collective

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

1645:   Options Database Keys:
1646: + -ksp_atol <abstol>   - Sets `abstol`
1647: . -ksp_rtol <rtol>     - Sets `rtol`
1648: . -ksp_divtol <dtol>   - Sets `dtol`
1649: - -ksp_max_it <maxits> - Sets `maxits`

1651:   Level: intermediate

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

1657:   All parameters must be non-negative.

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

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

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

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

1668:   Fortran Note:
1669:   Use `PETSC_CURRENT_INTEGER`, `PETSC_CURRENT_REAL`, `PETSC_DETERMINE_INTEGER`, or `PETSC_DETERMINE_REAL`

1671: .seealso: [](ch_ksp), `KSPGetTolerances()`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSP`, `KSPSetMinimumIterations()`
1672: @*/
1673: PetscErrorCode KSPSetTolerances(KSP ksp, PetscReal rtol, PetscReal abstol, PetscReal dtol, PetscInt maxits)
1674: {
1675:   PetscFunctionBegin;

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

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

1716:   Logically Collective

1718:   Input Parameters:
1719: + ksp   - the Krylov subspace context
1720: - minit - minimum number of iterations to use

1722:   Options Database Key:
1723: . -ksp_min_it <minits> - Sets `minit`

1725:   Level: intermediate

1727:   Notes:
1728:   Use `KSPSetTolerances()` to set a variety of other tolerances

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

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

1737: .seealso: [](ch_ksp), `KSPGetTolerances()`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSP`, `KSPSetTolerances()`, `KSPGetMinimumIterations()`
1738: @*/
1739: PetscErrorCode KSPSetMinimumIterations(KSP ksp, PetscInt minit)
1740: {
1741:   PetscFunctionBegin;

1745:   PetscCheck(minit >= 0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Minimum number of iterations %" PetscInt_FMT " must be non-negative", minit);
1746:   ksp->min_it = minit;
1747:   PetscFunctionReturn(PETSC_SUCCESS);
1748: }

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

1753:   Not Collective

1755:   Input Parameter:
1756: . ksp - the Krylov subspace context

1758:   Output Parameter:
1759: . minit - minimum number of iterations to use

1761:   Level: intermediate

1763: .seealso: [](ch_ksp), `KSPGetTolerances()`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSP`, `KSPSetTolerances()`, `KSPSetMinimumIterations()`
1764: @*/
1765: PetscErrorCode KSPGetMinimumIterations(KSP ksp, PetscInt *minit)
1766: {
1767:   PetscFunctionBegin;
1769:   PetscAssertPointer(minit, 2);

1771:   *minit = ksp->min_it;
1772:   PetscFunctionReturn(PETSC_SUCCESS);
1773: }

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

1780:   Logically Collective

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

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

1789:   Level: beginner

1791: .seealso: [](ch_ksp), `KSPGetInitialGuessNonzero()`, `KSPGuessSetType()`, `KSPGuessType`, `KSP`
1792: @*/
1793: PetscErrorCode KSPSetInitialGuessNonzero(KSP ksp, PetscBool flg)
1794: {
1795:   PetscFunctionBegin;
1798:   ksp->guess_zero = (PetscBool)!flg;
1799:   PetscFunctionReturn(PETSC_SUCCESS);
1800: }

1802: /*@
1803:   KSPGetInitialGuessNonzero - Determines whether the `KSP` solver is using
1804:   a zero initial guess.

1806:   Not Collective

1808:   Input Parameter:
1809: . ksp - iterative solver obtained from `KSPCreate()`

1811:   Output Parameter:
1812: . flag - `PETSC_TRUE` if guess is nonzero, else `PETSC_FALSE`

1814:   Level: intermediate

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

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

1831:   Logically Collective

1833:   Input Parameters:
1834: + ksp - iterative solver obtained from `KSPCreate()`
1835: - flg - `PETSC_TRUE` indicates you want the error generated

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

1840:   Level: intermediate

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

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

1849: .seealso: [](ch_ksp), `KSPGetErrorIfNotConverged()`, `KSP`
1850: @*/
1851: PetscErrorCode KSPSetErrorIfNotConverged(KSP ksp, PetscBool flg)
1852: {
1853:   PetscFunctionBegin;
1856:   ksp->errorifnotconverged = flg;
1857:   PetscFunctionReturn(PETSC_SUCCESS);
1858: }

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

1863:   Not Collective

1865:   Input Parameter:
1866: . ksp - iterative solver obtained from KSPCreate()

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

1871:   Level: intermediate

1873: .seealso: [](ch_ksp), `KSPSetErrorIfNotConverged()`, `KSP`
1874: @*/
1875: PetscErrorCode KSPGetErrorIfNotConverged(KSP ksp, PetscBool *flag)
1876: {
1877:   PetscFunctionBegin;
1879:   PetscAssertPointer(flag, 2);
1880:   *flag = ksp->errorifnotconverged;
1881:   PetscFunctionReturn(PETSC_SUCCESS);
1882: }

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

1887:   Logically Collective

1889:   Input Parameters:
1890: + ksp - iterative solver obtained from `KSPCreate()`
1891: - flg - `PETSC_TRUE` or `PETSC_FALSE`

1893:   Level: advanced

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

1899: .seealso: [](ch_ksp), `KSPGetInitialGuessKnoll()`, `KSPGuess`, `KSPSetInitialGuessNonzero()`, `KSPGetInitialGuessNonzero()`, `KSP`
1900: @*/
1901: PetscErrorCode KSPSetInitialGuessKnoll(KSP ksp, PetscBool flg)
1902: {
1903:   PetscFunctionBegin;
1906:   ksp->guess_knoll = flg;
1907:   PetscFunctionReturn(PETSC_SUCCESS);
1908: }

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

1914:   Not Collective

1916:   Input Parameter:
1917: . ksp - iterative solver obtained from `KSPCreate()`

1919:   Output Parameter:
1920: . flag - `PETSC_TRUE` if using Knoll trick, else `PETSC_FALSE`

1922:   Level: advanced

1924: .seealso: [](ch_ksp), `KSPSetInitialGuessKnoll()`, `KSPSetInitialGuessNonzero()`, `KSPGetInitialGuessNonzero()`, `KSP`
1925: @*/
1926: PetscErrorCode KSPGetInitialGuessKnoll(KSP ksp, PetscBool *flag)
1927: {
1928:   PetscFunctionBegin;
1930:   PetscAssertPointer(flag, 2);
1931:   *flag = ksp->guess_knoll;
1932:   PetscFunctionReturn(PETSC_SUCCESS);
1933: }

1935: /*@
1936:   KSPGetComputeSingularValues - Gets the flag indicating whether the extreme singular
1937:   values will be calculated via a Lanczos or Arnoldi process as the linear
1938:   system is solved.

1940:   Not Collective

1942:   Input Parameter:
1943: . ksp - iterative solver obtained from `KSPCreate()`

1945:   Output Parameter:
1946: . flg - `PETSC_TRUE` or `PETSC_FALSE`

1948:   Options Database Key:
1949: . -ksp_monitor_singular_value - Activates `KSPSetComputeSingularValues()`

1951:   Level: advanced

1953:   Notes:
1954:   This option is not valid for `KSPType`.

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

1960: .seealso: [](ch_ksp), `KSPComputeExtremeSingularValues()`, `KSPMonitorSingularValue()`, `KSP`
1961: @*/
1962: PetscErrorCode KSPGetComputeSingularValues(KSP ksp, PetscBool *flg)
1963: {
1964:   PetscFunctionBegin;
1966:   PetscAssertPointer(flg, 2);
1967:   *flg = ksp->calc_sings;
1968:   PetscFunctionReturn(PETSC_SUCCESS);
1969: }

1971: /*@
1972:   KSPSetComputeSingularValues - Sets a flag so that the extreme singular
1973:   values will be calculated via a Lanczos or Arnoldi process as the linear
1974:   system is solved.

1976:   Logically Collective

1978:   Input Parameters:
1979: + ksp - iterative solver obtained from `KSPCreate()`
1980: - flg - `PETSC_TRUE` or `PETSC_FALSE`

1982:   Options Database Key:
1983: . -ksp_monitor_singular_value - Activates `KSPSetComputeSingularValues()`

1985:   Level: advanced

1987:   Notes:
1988:   This option is not valid for all iterative methods.

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

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

1996: .seealso: [](ch_ksp), `KSPComputeExtremeSingularValues()`, `KSPMonitorSingularValue()`, `KSP`, `KSPSetComputeRitz()`
1997: @*/
1998: PetscErrorCode KSPSetComputeSingularValues(KSP ksp, PetscBool flg)
1999: {
2000:   PetscFunctionBegin;
2003:   ksp->calc_sings = flg;
2004:   PetscFunctionReturn(PETSC_SUCCESS);
2005: }

2007: /*@
2008:   KSPGetComputeEigenvalues - Gets the flag indicating that the extreme eigenvalues
2009:   values will be calculated via a Lanczos or Arnoldi process as the linear
2010:   system is solved.

2012:   Not Collective

2014:   Input Parameter:
2015: . ksp - iterative solver obtained from `KSPCreate()`

2017:   Output Parameter:
2018: . flg - `PETSC_TRUE` or `PETSC_FALSE`

2020:   Level: advanced

2022:   Note:
2023:   Currently this option is not valid for all iterative methods.

2025: .seealso: [](ch_ksp), `KSPComputeEigenvalues()`, `KSPComputeEigenvaluesExplicitly()`, `KSP`, `KSPSetComputeRitz()`
2026: @*/
2027: PetscErrorCode KSPGetComputeEigenvalues(KSP ksp, PetscBool *flg)
2028: {
2029:   PetscFunctionBegin;
2031:   PetscAssertPointer(flg, 2);
2032:   *flg = ksp->calc_sings;
2033:   PetscFunctionReturn(PETSC_SUCCESS);
2034: }

2036: /*@
2037:   KSPSetComputeEigenvalues - Sets a flag so that the extreme eigenvalues
2038:   values will be calculated via a Lanczos or Arnoldi process as the linear
2039:   system is solved.

2041:   Logically Collective

2043:   Input Parameters:
2044: + ksp - iterative solver obtained from `KSPCreate()`
2045: - flg - `PETSC_TRUE` or `PETSC_FALSE`

2047:   Level: advanced

2049:   Note:
2050:   Currently this option is not valid for all iterative methods.

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

2054: .seealso: [](ch_ksp), `KSPComputeEigenvalues()`, `KSPComputeEigenvaluesExplicitly()`, `KSP`, `KSPSetComputeRitz()`
2055: @*/
2056: PetscErrorCode KSPSetComputeEigenvalues(KSP ksp, PetscBool flg)
2057: {
2058:   PetscFunctionBegin;
2061:   ksp->calc_sings = flg;
2062:   PetscFunctionReturn(PETSC_SUCCESS);
2063: }

2065: /*@
2066:   KSPSetComputeRitz - Sets a flag so that the Ritz or harmonic Ritz pairs
2067:   will be calculated via a Lanczos or Arnoldi process as the linear
2068:   system is solved.

2070:   Logically Collective

2072:   Input Parameters:
2073: + ksp - iterative solver obtained from `KSPCreate()`
2074: - flg - `PETSC_TRUE` or `PETSC_FALSE`

2076:   Level: advanced

2078:   Note:
2079:   Currently this option is only valid for the `KSPGMRES` method.

2081: .seealso: [](ch_ksp), `KSPComputeRitz()`, `KSP`, `KSPComputeEigenvalues()`, `KSPComputeExtremeSingularValues()`
2082: @*/
2083: PetscErrorCode KSPSetComputeRitz(KSP ksp, PetscBool flg)
2084: {
2085:   PetscFunctionBegin;
2088:   ksp->calc_ritz = flg;
2089:   PetscFunctionReturn(PETSC_SUCCESS);
2090: }

2092: /*@
2093:   KSPGetRhs - Gets the right-hand-side vector for the linear system to
2094:   be solved.

2096:   Not Collective

2098:   Input Parameter:
2099: . ksp - iterative solver obtained from `KSPCreate()`

2101:   Output Parameter:
2102: . r - right-hand-side vector

2104:   Level: developer

2106: .seealso: [](ch_ksp), `KSPGetSolution()`, `KSPSolve()`, `KSP`
2107: @*/
2108: PetscErrorCode KSPGetRhs(KSP ksp, Vec *r)
2109: {
2110:   PetscFunctionBegin;
2112:   PetscAssertPointer(r, 2);
2113:   *r = ksp->vec_rhs;
2114:   PetscFunctionReturn(PETSC_SUCCESS);
2115: }

2117: /*@
2118:   KSPGetSolution - Gets the location of the solution for the
2119:   linear system to be solved.

2121:   Not Collective

2123:   Input Parameter:
2124: . ksp - iterative solver obtained from `KSPCreate()`

2126:   Output Parameter:
2127: . v - solution vector

2129:   Level: developer

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

2135: .seealso: [](ch_ksp), `KSPGetRhs()`, `KSPBuildSolution()`, `KSPSolve()`, `KSP`
2136: @*/
2137: PetscErrorCode KSPGetSolution(KSP ksp, Vec *v)
2138: {
2139:   PetscFunctionBegin;
2141:   PetscAssertPointer(v, 2);
2142:   *v = ksp->vec_sol;
2143:   PetscFunctionReturn(PETSC_SUCCESS);
2144: }

2146: /*@
2147:   KSPSetPC - Sets the preconditioner to be used to calculate the
2148:   application of the preconditioner on a vector into a `KSP`.

2150:   Collective

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

2156:   Level: developer

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

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

2179: PETSC_INTERN PetscErrorCode PCCreate_MPI(PC);

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

2185:    Collective, No Fortran Support

2187:    Input Parameter:
2188: .  ksp - iterative solver obtained from `KSPCreate()`

2190:    Level: developer

2192: .seealso: [](ch_ksp), `KSPSetPC()`, `KSP`, `PCMPIServerBegin()`, `PCMPIServerEnd()`
2193: @*/
2194: PETSC_INTERN PetscErrorCode KSPCheckPCMPI(KSP ksp)
2195: {
2196:   PetscBool isPCMPI;

2198:   PetscFunctionBegin;
2200:   PetscCall(PetscObjectTypeCompare((PetscObject)ksp->pc, PCMPI, &isPCMPI));
2201:   if (PCMPIServerActive && ksp->nestlevel == 0 && !isPCMPI) {
2202:     const char *prefix;
2203:     char       *found = NULL;

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

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

2217:   Not Collective

2219:   Input Parameter:
2220: . ksp - iterative solver obtained from `KSPCreate()`

2222:   Output Parameter:
2223: . pc - preconditioner context

2225:   Level: beginner

2227:   Note:
2228:   The `PC` is created if it does not already exist.

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

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

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

2254:   Collective

2256:   Input Parameters:
2257: + ksp   - iterative solver obtained from `KSPCreate()`
2258: . it    - iteration number
2259: - rnorm - relative norm of the residual

2261:   Level: developer

2263:   Notes:
2264:   This routine is called by the `KSP` implementations.
2265:   It does not typically need to be called by the user.

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

2270: .seealso: [](ch_ksp), `KSPMonitorSet()`
2271: @*/
2272: PetscErrorCode KSPMonitor(KSP ksp, PetscInt it, PetscReal rnorm)
2273: {
2274:   PetscInt i, n = ksp->numbermonitors;

2276:   PetscFunctionBegin;
2277:   for (i = 0; i < n; i++) PetscCall((*ksp->monitor[i])(ksp, it, rnorm, ksp->monitorcontext[i]));
2278:   PetscFunctionReturn(PETSC_SUCCESS);
2279: }

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

2285:   Logically Collective

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

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

2306:   Level: beginner

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

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

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

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

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

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

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

2344: /*@
2345:   KSPMonitorCancel - Clears all monitors for a `KSP` object.

2347:   Logically Collective

2349:   Input Parameter:
2350: . ksp - iterative solver obtained from `KSPCreate()`

2352:   Options Database Key:
2353: . -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.

2355:   Level: intermediate

2357: .seealso: [](ch_ksp), `KSPMonitorResidual()`, `KSPMonitorSet()`, `KSP`
2358: @*/
2359: PetscErrorCode KSPMonitorCancel(KSP ksp)
2360: {
2361:   PetscInt i;

2363:   PetscFunctionBegin;
2365:   for (i = 0; i < ksp->numbermonitors; i++) {
2366:     if (ksp->monitordestroy[i]) PetscCall((*ksp->monitordestroy[i])(&ksp->monitorcontext[i]));
2367:   }
2368:   ksp->numbermonitors = 0;
2369:   PetscFunctionReturn(PETSC_SUCCESS);
2370: }

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

2375:   Not Collective

2377:   Input Parameter:
2378: . ksp - iterative solver obtained from `KSPCreate()`

2380:   Output Parameter:
2381: . ctx - monitoring context

2383:   Level: intermediate

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

2395: /*@
2396:   KSPSetResidualHistory - Sets the array used to hold the residual history.
2397:   If set, this array will contain the residual norms computed at each
2398:   iteration of the solver.

2400:   Not Collective

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

2409:   Level: advanced

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

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

2418: .seealso: [](ch_ksp), `KSPGetResidualHistory()`, `KSP`
2419: @*/
2420: PetscErrorCode KSPSetResidualHistory(KSP ksp, PetscReal a[], PetscCount na, PetscBool reset)
2421: {
2422:   PetscFunctionBegin;

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

2434:     ksp->res_hist = ksp->res_hist_alloc;
2435:   }
2436:   ksp->res_hist_len   = 0;
2437:   ksp->res_hist_reset = reset;
2438:   PetscFunctionReturn(PETSC_SUCCESS);
2439: }

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

2444:   Not Collective

2446:   Input Parameter:
2447: . ksp - iterative solver obtained from `KSPCreate()`

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

2453:   Level: advanced

2455:   Note:
2456:   This array is borrowed and should not be freed by the caller.

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

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

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

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

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

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

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

2485: /*@
2486:   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.

2488:   Not Collective

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

2496:   Level: advanced

2498:   Notes:
2499:   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.
2500:   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.

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

2504: .seealso: [](ch_ksp), `KSPGetErrorHistory()`, `KSPSetResidualHistory()`, `KSP`
2505: @*/
2506: PetscErrorCode KSPSetErrorHistory(KSP ksp, PetscReal a[], PetscCount na, PetscBool reset)
2507: {
2508:   PetscFunctionBegin;

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

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

2529:   Not Collective

2531:   Input Parameter:
2532: . ksp - iterative solver obtained from `KSPCreate()`

2534:   Output Parameters:
2535: + a  - pointer to array to hold history (or `NULL`)
2536: - na - number of used entries in a (or `NULL`)

2538:   Level: advanced

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

2544:   Fortran Note:
2545: .vb
2546:   PetscReal, pointer :: a(:)
2547: .ve

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

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

2563:   Not Collective

2565:   Input Parameter:
2566: . ksp - The `KSP`

2568:   Output Parameters:
2569: + cr   - The residual contraction rate
2570: . rRsq - The coefficient of determination, $R^2$, indicating the linearity of the data
2571: . ce   - The error contraction rate
2572: - eRsq - The coefficient of determination, $R^2$, indicating the linearity of the data

2574:   Level: advanced

2576:   Note:
2577:   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,
2578:   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)}$,

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

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

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

2641:   Logically Collective

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

2649:   Level: advanced

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

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

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

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

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

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

2681:   Logically Collective

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

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

2691:   Level: advanced

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

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

2708:   Logically Collective

2710:   Input Parameter:
2711: . ksp - iterative solver obtained from `KSPCreate()`

2713:   Output Parameters:
2714: + converge - pointer to convergence test function, see `KSPConvergenceTestFn`
2715: . ctx      - context for private data for the convergence routine
2716: - destroy  - a routine for destroying the context

2718:   Level: advanced

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

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

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

2744:   Not Collective

2746:   Input Parameter:
2747: . ksp - iterative solver obtained from `KSPCreate()`

2749:   Output Parameter:
2750: . ctx - monitoring context

2752:   Level: advanced

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

2764: /*@
2765:   KSPBuildSolution - Builds the approximate solution in a vector provided.

2767:   Collective

2769:   Input Parameter:
2770: . ksp - iterative solver obtained from `KSPCreate()`

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

2777:   Level: developer

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

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

2812: /*@
2813:   KSPBuildResidual - Builds the residual in a vector provided.

2815:   Collective

2817:   Input Parameter:
2818: . ksp - iterative solver obtained from `KSPCreate()`

2820:   Output Parameters:
2821: + t - work vector.  If not provided then one is generated.
2822: . v - optional location to stash residual.  If `v` is not provided, then a location is generated.
2823: - V - the residual

2825:   Level: advanced

2827:   Note:
2828:   Regardless of whether or not `v` is provided, the residual is
2829:   returned in `V`.

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

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

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

2854:   Logically Collective

2856:   Input Parameters:
2857: + ksp   - the `KSP` context
2858: - scale - `PETSC_TRUE` or `PETSC_FALSE`

2860:   Options Database Keys:
2861: + -ksp_diagonal_scale     - perform a diagonal scaling before the solve
2862: - -ksp_diagonal_scale_fix - scale the matrix back AFTER the solve

2864:   Level: advanced

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

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

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

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

2881: .seealso: [](ch_ksp), `KSPGetDiagonalScale()`, `KSPSetDiagonalScaleFix()`, `KSP`
2882: @*/
2883: PetscErrorCode KSPSetDiagonalScale(KSP ksp, PetscBool scale)
2884: {
2885:   PetscFunctionBegin;
2888:   ksp->dscale = scale;
2889:   PetscFunctionReturn(PETSC_SUCCESS);
2890: }

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

2895:   Not Collective

2897:   Input Parameter:
2898: . ksp - the `KSP` context

2900:   Output Parameter:
2901: . scale - `PETSC_TRUE` or `PETSC_FALSE`

2903:   Level: intermediate

2905: .seealso: [](ch_ksp), `KSP`, `KSPSetDiagonalScale()`, `KSPSetDiagonalScaleFix()`
2906: @*/
2907: PetscErrorCode KSPGetDiagonalScale(KSP ksp, PetscBool *scale)
2908: {
2909:   PetscFunctionBegin;
2911:   PetscAssertPointer(scale, 2);
2912:   *scale = ksp->dscale;
2913:   PetscFunctionReturn(PETSC_SUCCESS);
2914: }

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

2919:   Logically Collective

2921:   Input Parameters:
2922: + ksp - the `KSP` context
2923: - fix - `PETSC_TRUE` to scale back after the system solve, `PETSC_FALSE` to not
2924:          rescale (default)

2926:   Level: intermediate

2928:   Notes:
2929:   Must be called after `KSPSetDiagonalScale()`

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

2936: .seealso: [](ch_ksp), `KSPGetDiagonalScale()`, `KSPSetDiagonalScale()`, `KSPGetDiagonalScaleFix()`, `KSP`
2937: @*/
2938: PetscErrorCode KSPSetDiagonalScaleFix(KSP ksp, PetscBool fix)
2939: {
2940:   PetscFunctionBegin;
2943:   ksp->dscalefix = fix;
2944:   PetscFunctionReturn(PETSC_SUCCESS);
2945: }

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

2950:   Not Collective

2952:   Input Parameter:
2953: . ksp - the `KSP` context

2955:   Output Parameter:
2956: . fix - `PETSC_TRUE` to scale back after the system solve, `PETSC_FALSE` to not
2957:          rescale (default)

2959:   Level: intermediate

2961: .seealso: [](ch_ksp), `KSPGetDiagonalScale()`, `KSPSetDiagonalScale()`, `KSPSetDiagonalScaleFix()`, `KSP`
2962: @*/
2963: PetscErrorCode KSPGetDiagonalScaleFix(KSP ksp, PetscBool *fix)
2964: {
2965:   PetscFunctionBegin;
2967:   PetscAssertPointer(fix, 2);
2968:   *fix = ksp->dscalefix;
2969:   PetscFunctionReturn(PETSC_SUCCESS);
2970: }

2972: /*@C
2973:   KSPSetComputeOperators - set routine to compute the linear operators

2975:   Logically Collective

2977:   Input Parameters:
2978: + ksp  - the `KSP` context
2979: . func - function to compute the operators, see `KSPComputeOperatorsFn` for the calling sequence
2980: - ctx  - optional context

2982:   Level: beginner

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

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

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

2995: .seealso: [](ch_ksp), `KSP`, `KSPSetOperators()`, `KSPSetComputeRHS()`, `DMKSPSetComputeOperators()`, `KSPSetComputeInitialGuess()`, `KSPComputeOperatorsFn`
2996: @*/
2997: PetscErrorCode KSPSetComputeOperators(KSP ksp, KSPComputeOperatorsFn *func, void *ctx)
2998: {
2999:   DM dm;

3001:   PetscFunctionBegin;
3003:   PetscCall(KSPGetDM(ksp, &dm));
3004:   PetscCall(DMKSPSetComputeOperators(dm, func, ctx));
3005:   if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX;
3006:   PetscFunctionReturn(PETSC_SUCCESS);
3007: }

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

3012:   Logically Collective

3014:   Input Parameters:
3015: + ksp  - the `KSP` context
3016: . func - function to compute the right-hand side, see `KSPComputeRHSFn` for the calling sequence
3017: - ctx  - optional context

3019:   Level: beginner

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

3024: .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `DMKSPSetComputeRHS()`, `KSPSetComputeOperators()`, `KSPSetOperators()`, `KSPComputeRHSFn`
3025: @*/
3026: PetscErrorCode KSPSetComputeRHS(KSP ksp, KSPComputeRHSFn *func, void *ctx)
3027: {
3028:   DM dm;

3030:   PetscFunctionBegin;
3032:   PetscCall(KSPGetDM(ksp, &dm));
3033:   PetscCall(DMKSPSetComputeRHS(dm, func, ctx));
3034:   PetscFunctionReturn(PETSC_SUCCESS);
3035: }

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

3040:   Logically Collective

3042:   Input Parameters:
3043: + ksp  - the `KSP` context
3044: . func - function to compute the initial guess, see `KSPComputeInitialGuessFn` for calling sequence
3045: - ctx  - optional context

3047:   Level: beginner

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

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

3060:   PetscFunctionBegin;
3062:   PetscCall(KSPGetDM(ksp, &dm));
3063:   PetscCall(DMKSPSetComputeInitialGuess(dm, func, ctx));
3064:   PetscFunctionReturn(PETSC_SUCCESS);
3065: }

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

3071:   Logically Collective

3073:   Input Parameter:
3074: . ksp - the `KSP` context

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

3079:   Level: advanced

3081: .seealso: [](ch_ksp), `KSPSolveTranspose()`, `KSP`
3082: @*/
3083: PetscErrorCode KSPSetUseExplicitTranspose(KSP ksp, PetscBool flg)
3084: {
3085:   PetscFunctionBegin;
3088:   ksp->transpose.use_explicittranspose = flg;
3089:   PetscFunctionReturn(PETSC_SUCCESS);
3090: }