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: . ctx               - [optional] context for private data for the `KSPConvergedReason` view routine (use `NULL` if context is not needed)
513: - reasonviewdestroy - [optional] routine that frees `ctx` (may be `NULL`), see `PetscCtxDestroyFn` for the calling sequence

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

520:   Level: intermediate

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

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

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

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

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

553:   Collective

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

558:   Level: intermediate

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

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

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

579:   Collective

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

584:   Level: intermediate

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

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

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

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

609:   Collective

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

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

618:   Level: intermediate

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

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

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

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

672: #include <petscdraw.h>

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

923:   ksp->guess_zero = guess_zero;

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

928:   PetscCall(KSPConvergedReasonViewFromOptions(ksp));

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

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

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

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

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

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

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

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

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

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

1011:   Collective

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

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

1036:   Level: beginner

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1105: static PetscErrorCode KSPUseExplicitTranspose_Private(KSP ksp)
1106: {
1107:   Mat J, Jpre;

1109:   PetscFunctionBegin;
1110:   PetscCall(KSPGetOperators(ksp, &J, &Jpre));
1111:   if (!ksp->transpose.reuse_transpose) {
1112:     PetscCall(MatTranspose(J, MAT_INITIAL_MATRIX, &ksp->transpose.AT));
1113:     if (J != Jpre) PetscCall(MatTranspose(Jpre, MAT_INITIAL_MATRIX, &ksp->transpose.BT));
1114:     ksp->transpose.reuse_transpose = PETSC_TRUE;
1115:   } else {
1116:     PetscCall(MatTranspose(J, MAT_REUSE_MATRIX, &ksp->transpose.AT));
1117:     if (J != Jpre) PetscCall(MatTranspose(Jpre, MAT_REUSE_MATRIX, &ksp->transpose.BT));
1118:   }
1119:   if (J == Jpre && ksp->transpose.BT != ksp->transpose.AT) {
1120:     PetscCall(PetscObjectReference((PetscObject)ksp->transpose.AT));
1121:     ksp->transpose.BT = ksp->transpose.AT;
1122:   }
1123:   PetscCall(KSPSetOperators(ksp, ksp->transpose.AT, ksp->transpose.BT));
1124:   PetscFunctionReturn(PETSC_SUCCESS);
1125: }

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

1130:   Collective

1132:   Input Parameters:
1133: + ksp - iterative solver obtained from `KSPCreate()`
1134: . b   - right-hand side vector
1135: - x   - solution vector

1137:   Level: developer

1139:   Note:
1140:   For complex numbers, this solve the non-Hermitian transpose system.

1142:   Developer Note:
1143:   We need to implement a `KSPSolveHermitianTranspose()`

1145: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPDestroy()`, `KSPSetTolerances()`, `KSPConvergedDefault()`,
1146:           `KSPSolve()`, `KSP`, `KSPSetOperators()`
1147: @*/
1148: PetscErrorCode KSPSolveTranspose(KSP ksp, Vec b, Vec x)
1149: {
1150:   PetscFunctionBegin;
1154:   if (ksp->transpose.use_explicittranspose) {
1155:     Mat J, Jpre;
1156:     PetscCall(KSPGetOperators(ksp, &J, &Jpre));
1157:     if (!ksp->transpose.reuse_transpose) {
1158:       PetscCall(MatTranspose(J, MAT_INITIAL_MATRIX, &ksp->transpose.AT));
1159:       if (J != Jpre) PetscCall(MatTranspose(Jpre, MAT_INITIAL_MATRIX, &ksp->transpose.BT));
1160:       ksp->transpose.reuse_transpose = PETSC_TRUE;
1161:     } else {
1162:       PetscCall(MatTranspose(J, MAT_REUSE_MATRIX, &ksp->transpose.AT));
1163:       if (J != Jpre) PetscCall(MatTranspose(Jpre, MAT_REUSE_MATRIX, &ksp->transpose.BT));
1164:     }
1165:     if (J == Jpre && ksp->transpose.BT != ksp->transpose.AT) {
1166:       PetscCall(PetscObjectReference((PetscObject)ksp->transpose.AT));
1167:       ksp->transpose.BT = ksp->transpose.AT;
1168:     }
1169:     PetscCall(KSPSetOperators(ksp, ksp->transpose.AT, ksp->transpose.BT));
1170:   } else {
1171:     ksp->transpose_solve = PETSC_TRUE;
1172:   }
1173:   PetscCall(KSPSolve_Private(ksp, b, x));
1174:   PetscFunctionReturn(PETSC_SUCCESS);
1175: }

1177: static PetscErrorCode KSPViewFinalMatResidual_Internal(KSP ksp, Mat B, Mat X, PetscViewer viewer, PetscViewerFormat format, PetscInt shift)
1178: {
1179:   Mat        A, R;
1180:   PetscReal *norms;
1181:   PetscInt   i, N;
1182:   PetscBool  flg;

1184:   PetscFunctionBegin;
1185:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &flg));
1186:   if (flg) {
1187:     PetscCall(PCGetOperators(ksp->pc, &A, NULL));
1188:     if (!ksp->transpose_solve) PetscCall(MatMatMult(A, X, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &R));
1189:     else PetscCall(MatTransposeMatMult(A, X, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &R));
1190:     PetscCall(MatAYPX(R, -1.0, B, SAME_NONZERO_PATTERN));
1191:     PetscCall(MatGetSize(R, NULL, &N));
1192:     PetscCall(PetscMalloc1(N, &norms));
1193:     PetscCall(MatGetColumnNorms(R, NORM_2, norms));
1194:     PetscCall(MatDestroy(&R));
1195:     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]));
1196:     PetscCall(PetscFree(norms));
1197:   }
1198:   PetscFunctionReturn(PETSC_SUCCESS);
1199: }

1201: static PetscErrorCode KSPMatSolve_Private(KSP ksp, Mat B, Mat X)
1202: {
1203:   Mat       A, P, vB, vX;
1204:   Vec       cb, cx;
1205:   PetscInt  n1, N1, n2, N2, Bbn = PETSC_DECIDE;
1206:   PetscBool match;

1208:   PetscFunctionBegin;
1212:   PetscCheckSameComm(ksp, 1, B, 2);
1213:   PetscCheckSameComm(ksp, 1, X, 3);
1214:   PetscCheckSameType(B, 2, X, 3);
1215:   PetscCheck(B->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
1216:   MatCheckPreallocated(X, 3);
1217:   if (!X->assembled) {
1218:     PetscCall(MatSetOption(X, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
1219:     PetscCall(MatAssemblyBegin(X, MAT_FINAL_ASSEMBLY));
1220:     PetscCall(MatAssemblyEnd(X, MAT_FINAL_ASSEMBLY));
1221:   }
1222:   PetscCheck(B != X, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_IDN, "B and X must be different matrices");
1223:   PetscCall(KSPGetOperators(ksp, &A, &P));
1224:   PetscCall(MatGetLocalSize(B, NULL, &n2));
1225:   PetscCall(MatGetLocalSize(X, NULL, &n1));
1226:   PetscCall(MatGetSize(B, NULL, &N2));
1227:   PetscCall(MatGetSize(X, NULL, &N1));
1228:   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);
1229:   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)B, &match, MATSEQDENSE, MATMPIDENSE, ""));
1230:   PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of right-hand sides not stored in a dense Mat");
1231:   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)X, &match, MATSEQDENSE, MATMPIDENSE, ""));
1232:   PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of solutions not stored in a dense Mat");
1233:   PetscCall(KSPSetUp(ksp));
1234:   PetscCall(KSPSetUpOnBlocks(ksp));
1235:   if (ksp->ops->matsolve) {
1236:     level++;
1237:     if (ksp->guess_zero) PetscCall(MatZeroEntries(X));
1238:     PetscCall(PetscLogEventBegin(!ksp->transpose_solve ? KSP_MatSolve : KSP_MatSolveTranspose, ksp, B, X, 0));
1239:     PetscCall(KSPGetMatSolveBatchSize(ksp, &Bbn));
1240:     /* by default, do a single solve with all columns */
1241:     if (Bbn == PETSC_DECIDE) Bbn = N2;
1242:     else PetscCheck(Bbn >= 1, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "KSPMatSolve() batch size %" PetscInt_FMT " must be positive", Bbn);
1243:     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));
1244:     /* if -ksp_matsolve_batch_size is greater than the actual number of columns, do a single solve with all columns */
1245:     if (Bbn >= N2) {
1246:       PetscUseTypeMethod(ksp, matsolve, B, X);
1247:       if (ksp->viewFinalRes) PetscCall(KSPViewFinalMatResidual_Internal(ksp, B, X, ksp->viewerFinalRes, ksp->formatFinalRes, 0));

1249:       PetscCall(KSPConvergedReasonViewFromOptions(ksp));

1251:       if (ksp->viewRate) {
1252:         PetscCall(PetscViewerPushFormat(ksp->viewerRate, PETSC_VIEWER_DEFAULT));
1253:         PetscCall(KSPConvergedRateView(ksp, ksp->viewerRate));
1254:         PetscCall(PetscViewerPopFormat(ksp->viewerRate));
1255:       }
1256:     } else {
1257:       for (n2 = 0; n2 < N2; n2 += Bbn) {
1258:         PetscCall(MatDenseGetSubMatrix(B, PETSC_DECIDE, PETSC_DECIDE, n2, PetscMin(n2 + Bbn, N2), &vB));
1259:         PetscCall(MatDenseGetSubMatrix(X, PETSC_DECIDE, PETSC_DECIDE, n2, PetscMin(n2 + Bbn, N2), &vX));
1260:         PetscUseTypeMethod(ksp, matsolve, vB, vX);
1261:         if (ksp->viewFinalRes) PetscCall(KSPViewFinalMatResidual_Internal(ksp, vB, vX, ksp->viewerFinalRes, ksp->formatFinalRes, n2));

1263:         PetscCall(KSPConvergedReasonViewFromOptions(ksp));

1265:         if (ksp->viewRate) {
1266:           PetscCall(PetscViewerPushFormat(ksp->viewerRate, PETSC_VIEWER_DEFAULT));
1267:           PetscCall(KSPConvergedRateView(ksp, ksp->viewerRate));
1268:           PetscCall(PetscViewerPopFormat(ksp->viewerRate));
1269:         }
1270:         PetscCall(MatDenseRestoreSubMatrix(B, &vB));
1271:         PetscCall(MatDenseRestoreSubMatrix(X, &vX));
1272:       }
1273:     }
1274:     if (ksp->viewMat) PetscCall(ObjectView((PetscObject)A, ksp->viewerMat, ksp->formatMat));
1275:     if (ksp->viewPMat) PetscCall(ObjectView((PetscObject)P, ksp->viewerPMat, ksp->formatPMat));
1276:     if (ksp->viewRhs) PetscCall(ObjectView((PetscObject)B, ksp->viewerRhs, ksp->formatRhs));
1277:     if (ksp->viewSol) PetscCall(ObjectView((PetscObject)X, ksp->viewerSol, ksp->formatSol));
1278:     if (ksp->view) PetscCall(KSPView(ksp, ksp->viewer));
1279:     PetscCall(PetscLogEventEnd(!ksp->transpose_solve ? KSP_MatSolve : KSP_MatSolveTranspose, ksp, B, X, 0));
1280:     if (ksp->errorifnotconverged && ksp->reason < 0 && (level == 1 || ksp->reason != KSP_DIVERGED_ITS)) {
1281:       PCFailedReason reason;

1283:       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]);
1284:       PetscCall(PCGetFailedReason(ksp->pc, &reason));
1285:       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]);
1286:     }
1287:     level--;
1288:   } else {
1289:     PetscCall(PetscInfo(ksp, "KSP type %s solving column by column\n", ((PetscObject)ksp)->type_name));
1290:     for (n2 = 0; n2 < N2; ++n2) {
1291:       PetscCall(MatDenseGetColumnVecRead(B, n2, &cb));
1292:       PetscCall(MatDenseGetColumnVecWrite(X, n2, &cx));
1293:       PetscCall(KSPSolve_Private(ksp, cb, cx));
1294:       PetscCall(MatDenseRestoreColumnVecWrite(X, n2, &cx));
1295:       PetscCall(MatDenseRestoreColumnVecRead(B, n2, &cb));
1296:     }
1297:   }
1298:   PetscFunctionReturn(PETSC_SUCCESS);
1299: }

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

1304:   Input Parameters:
1305: + ksp - iterative solver
1306: - B   - block of right-hand sides

1308:   Output Parameter:
1309: . X - block of solutions

1311:   Level: intermediate

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

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

1318: .seealso: [](ch_ksp), `KSPSolve()`, `MatMatSolve()`, `KSPMatSolveTranspose()`, `MATDENSE`, `KSPHPDDM`, `PCBJACOBI`, `PCASM`, `KSPSetMatSolveBatchSize()`
1319: @*/
1320: PetscErrorCode KSPMatSolve(KSP ksp, Mat B, Mat X)
1321: {
1322:   PetscFunctionBegin;
1323:   ksp->transpose_solve = PETSC_FALSE;
1324:   PetscCall(KSPMatSolve_Private(ksp, B, X));
1325:   PetscFunctionReturn(PETSC_SUCCESS);
1326: }

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

1331:   Input Parameters:
1332: + ksp - iterative solver
1333: - B   - block of right-hand sides

1335:   Output Parameter:
1336: . X - block of solutions

1338:   Level: intermediate

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

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

1346: .seealso: [](ch_ksp), `KSPSolveTranspose()`, `MatMatTransposeSolve()`, `KSPMatSolve()`, `MATDENSE`, `KSPHPDDM`, `PCBJACOBI`, `PCASM`
1347: @*/
1348: PetscErrorCode KSPMatSolveTranspose(KSP ksp, Mat B, Mat X)
1349: {
1350:   PetscFunctionBegin;
1351:   if (ksp->transpose.use_explicittranspose) PetscCall(KSPUseExplicitTranspose_Private(ksp));
1352:   else ksp->transpose_solve = PETSC_TRUE;
1353:   PetscCall(KSPMatSolve_Private(ksp, B, X));
1354:   PetscFunctionReturn(PETSC_SUCCESS);
1355: }

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

1360:   Logically Collective

1362:   Input Parameters:
1363: + ksp - the `KSP` iterative solver
1364: - bs  - batch size

1366:   Level: advanced

1368:   Note:
1369:   Using a larger block size can improve the efficiency of the solver.

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

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

1385:   Input Parameter:
1386: . ksp - iterative solver context

1388:   Output Parameter:
1389: . bs - batch size

1391:   Level: advanced

1393: .seealso: [](ch_ksp), `KSPMatSolve()`, `KSPSetMatSolveBatchSize()`, `-mat_mumps_icntl_27`, `-matproduct_batch_size`
1394: @*/
1395: PetscErrorCode KSPGetMatSolveBatchSize(KSP ksp, PetscInt *bs)
1396: {
1397:   PetscFunctionBegin;
1399:   PetscAssertPointer(bs, 2);
1400:   *bs = ksp->nmax;
1401:   PetscFunctionReturn(PETSC_SUCCESS);
1402: }

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

1407:   Collective

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

1412:   Level: beginner

1414: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPSolve()`, `KSPSetFromOptions()`, `KSP`
1415: @*/
1416: PetscErrorCode KSPResetViewers(KSP ksp)
1417: {
1418:   PetscFunctionBegin;
1420:   if (!ksp) PetscFunctionReturn(PETSC_SUCCESS);
1421:   PetscCall(PetscViewerDestroy(&ksp->viewer));
1422:   PetscCall(PetscViewerDestroy(&ksp->viewerPre));
1423:   PetscCall(PetscViewerDestroy(&ksp->viewerRate));
1424:   PetscCall(PetscViewerDestroy(&ksp->viewerMat));
1425:   PetscCall(PetscViewerDestroy(&ksp->viewerPMat));
1426:   PetscCall(PetscViewerDestroy(&ksp->viewerRhs));
1427:   PetscCall(PetscViewerDestroy(&ksp->viewerSol));
1428:   PetscCall(PetscViewerDestroy(&ksp->viewerMatExp));
1429:   PetscCall(PetscViewerDestroy(&ksp->viewerEV));
1430:   PetscCall(PetscViewerDestroy(&ksp->viewerSV));
1431:   PetscCall(PetscViewerDestroy(&ksp->viewerEVExp));
1432:   PetscCall(PetscViewerDestroy(&ksp->viewerFinalRes));
1433:   PetscCall(PetscViewerDestroy(&ksp->viewerPOpExp));
1434:   PetscCall(PetscViewerDestroy(&ksp->viewerDScale));
1435:   ksp->view         = PETSC_FALSE;
1436:   ksp->viewPre      = PETSC_FALSE;
1437:   ksp->viewMat      = PETSC_FALSE;
1438:   ksp->viewPMat     = PETSC_FALSE;
1439:   ksp->viewRhs      = PETSC_FALSE;
1440:   ksp->viewSol      = PETSC_FALSE;
1441:   ksp->viewMatExp   = PETSC_FALSE;
1442:   ksp->viewEV       = PETSC_FALSE;
1443:   ksp->viewSV       = PETSC_FALSE;
1444:   ksp->viewEVExp    = PETSC_FALSE;
1445:   ksp->viewFinalRes = PETSC_FALSE;
1446:   ksp->viewPOpExp   = PETSC_FALSE;
1447:   ksp->viewDScale   = PETSC_FALSE;
1448:   PetscFunctionReturn(PETSC_SUCCESS);
1449: }

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

1454:   Collective

1456:   Input Parameter:
1457: . ksp - iterative solver obtained from `KSPCreate()`

1459:   Level: intermediate

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

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

1466: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPSolve()`, `KSP`
1467: @*/
1468: PetscErrorCode KSPReset(KSP ksp)
1469: {
1470:   PetscFunctionBegin;
1472:   if (!ksp) PetscFunctionReturn(PETSC_SUCCESS);
1473:   PetscTryTypeMethod(ksp, reset);
1474:   if (ksp->pc) PetscCall(PCReset(ksp->pc));
1475:   if (ksp->guess) {
1476:     KSPGuess guess = ksp->guess;
1477:     PetscTryTypeMethod(guess, reset);
1478:   }
1479:   PetscCall(VecDestroyVecs(ksp->nwork, &ksp->work));
1480:   PetscCall(VecDestroy(&ksp->vec_rhs));
1481:   PetscCall(VecDestroy(&ksp->vec_sol));
1482:   PetscCall(VecDestroy(&ksp->diagonal));
1483:   PetscCall(VecDestroy(&ksp->truediagonal));

1485:   ksp->setupstage = KSP_SETUP_NEW;
1486:   ksp->nmax       = PETSC_DECIDE;
1487:   PetscFunctionReturn(PETSC_SUCCESS);
1488: }

1490: /*@
1491:   KSPDestroy - Destroys a `KSP` context.

1493:   Collective

1495:   Input Parameter:
1496: . ksp - iterative solver obtained from `KSPCreate()`

1498:   Level: beginner

1500: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPSolve()`, `KSP`
1501: @*/
1502: PetscErrorCode KSPDestroy(KSP *ksp)
1503: {
1504:   PC pc;

1506:   PetscFunctionBegin;
1507:   if (!*ksp) PetscFunctionReturn(PETSC_SUCCESS);
1509:   if (--((PetscObject)*ksp)->refct > 0) {
1510:     *ksp = NULL;
1511:     PetscFunctionReturn(PETSC_SUCCESS);
1512:   }

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

1516:   /*
1517:    Avoid a cascading call to PCReset(ksp->pc) from the following call:
1518:    PCReset() shouldn't be called from KSPDestroy() as it is unprotected by pc's
1519:    refcount (and may be shared, e.g., by other ksps).
1520:    */
1521:   pc         = (*ksp)->pc;
1522:   (*ksp)->pc = NULL;
1523:   PetscCall(KSPReset(*ksp));
1524:   PetscCall(KSPResetViewers(*ksp));
1525:   (*ksp)->pc = pc;
1526:   PetscTryTypeMethod(*ksp, destroy);

1528:   if ((*ksp)->transpose.use_explicittranspose) {
1529:     PetscCall(MatDestroy(&(*ksp)->transpose.AT));
1530:     PetscCall(MatDestroy(&(*ksp)->transpose.BT));
1531:     (*ksp)->transpose.reuse_transpose = PETSC_FALSE;
1532:   }

1534:   PetscCall(KSPGuessDestroy(&(*ksp)->guess));
1535:   PetscCall(DMDestroy(&(*ksp)->dm));
1536:   PetscCall(PCDestroy(&(*ksp)->pc));
1537:   PetscCall(PetscFree((*ksp)->res_hist_alloc));
1538:   PetscCall(PetscFree((*ksp)->err_hist_alloc));
1539:   if ((*ksp)->convergeddestroy) PetscCall((*(*ksp)->convergeddestroy)(&(*ksp)->cnvP));
1540:   PetscCall(KSPMonitorCancel(*ksp));
1541:   PetscCall(KSPConvergedReasonViewCancel(*ksp));
1542:   PetscCall(PetscHeaderDestroy(ksp));
1543:   PetscFunctionReturn(PETSC_SUCCESS);
1544: }

1546: /*@
1547:   KSPSetPCSide - Sets the preconditioning side.

1549:   Logically Collective

1551:   Input Parameter:
1552: . ksp - iterative solver obtained from `KSPCreate()`

1554:   Output Parameter:
1555: . side - the preconditioning side, where side is one of
1556: .vb
1557:   PC_LEFT      - left preconditioning (default)
1558:   PC_RIGHT     - right preconditioning
1559:   PC_SYMMETRIC - symmetric preconditioning
1560: .ve

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

1565:   Level: intermediate

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

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

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

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

1579: .seealso: [](ch_ksp), `KSPGetPCSide()`, `KSPSetNormType()`, `KSPGetNormType()`, `KSP`, `KSPSetPreSolve()`, `KSPSetPostSolve()`
1580: @*/
1581: PetscErrorCode KSPSetPCSide(KSP ksp, PCSide side)
1582: {
1583:   PetscFunctionBegin;
1586:   ksp->pc_side = ksp->pc_side_set = side;
1587:   PetscFunctionReturn(PETSC_SUCCESS);
1588: }

1590: /*@
1591:   KSPGetPCSide - Gets the preconditioning side.

1593:   Not Collective

1595:   Input Parameter:
1596: . ksp - iterative solver obtained from `KSPCreate()`

1598:   Output Parameter:
1599: . side - the preconditioning side, where side is one of
1600: .vb
1601:   PC_LEFT      - left preconditioning (default)
1602:   PC_RIGHT     - right preconditioning
1603:   PC_SYMMETRIC - symmetric preconditioning
1604: .ve

1606:   Level: intermediate

1608: .seealso: [](ch_ksp), `KSPSetPCSide()`, `KSP`
1609: @*/
1610: PetscErrorCode KSPGetPCSide(KSP ksp, PCSide *side)
1611: {
1612:   PetscFunctionBegin;
1614:   PetscAssertPointer(side, 2);
1615:   PetscCall(KSPSetUpNorms_Private(ksp, PETSC_TRUE, &ksp->normtype, &ksp->pc_side));
1616:   *side = ksp->pc_side;
1617:   PetscFunctionReturn(PETSC_SUCCESS);
1618: }

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

1624:   Not Collective

1626:   Input Parameter:
1627: . ksp - the Krylov subspace context

1629:   Output Parameters:
1630: + rtol   - the relative convergence tolerance
1631: . abstol - the absolute convergence tolerance
1632: . dtol   - the divergence tolerance
1633: - maxits - maximum number of iterations

1635:   Level: intermediate

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

1640: .seealso: [](ch_ksp), `KSPSetTolerances()`, `KSP`, `KSPSetMinimumIterations()`, `KSPGetMinimumIterations()`
1641: @*/
1642: PetscErrorCode KSPGetTolerances(KSP ksp, PeOp PetscReal *rtol, PeOp PetscReal *abstol, PeOp PetscReal *dtol, PeOp PetscInt *maxits)
1643: {
1644:   PetscFunctionBegin;
1646:   if (abstol) *abstol = ksp->abstol;
1647:   if (rtol) *rtol = ksp->rtol;
1648:   if (dtol) *dtol = ksp->divtol;
1649:   if (maxits) *maxits = ksp->max_it;
1650:   PetscFunctionReturn(PETSC_SUCCESS);
1651: }

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

1657:   Logically Collective

1659:   Input Parameters:
1660: + ksp    - the Krylov subspace context
1661: . rtol   - the relative convergence tolerance, relative decrease in the (possibly preconditioned) residual norm
1662: . abstol - the absolute convergence tolerance   absolute size of the (possibly preconditioned) residual norm
1663: . dtol   - the divergence tolerance,   amount (possibly preconditioned) residual norm can increase before `KSPConvergedDefault()` concludes that the method is diverging
1664: - maxits - maximum number of iterations to use

1666:   Options Database Keys:
1667: + -ksp_atol <abstol>   - Sets `abstol`
1668: . -ksp_rtol <rtol>     - Sets `rtol`
1669: . -ksp_divtol <dtol>   - Sets `dtol`
1670: - -ksp_max_it <maxits> - Sets `maxits`

1672:   Level: intermediate

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

1678:   All parameters must be non-negative.

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

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

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

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

1689:   Fortran Note:
1690:   Use `PETSC_CURRENT_INTEGER`, `PETSC_CURRENT_REAL`, `PETSC_DETERMINE_INTEGER`, or `PETSC_DETERMINE_REAL`

1692: .seealso: [](ch_ksp), `KSPGetTolerances()`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSP`, `KSPSetMinimumIterations()`
1693: @*/
1694: PetscErrorCode KSPSetTolerances(KSP ksp, PetscReal rtol, PetscReal abstol, PetscReal dtol, PetscInt maxits)
1695: {
1696:   PetscFunctionBegin;

1703:   if (rtol == (PetscReal)PETSC_DETERMINE) {
1704:     ksp->rtol = ksp->default_rtol;
1705:   } else if (rtol != (PetscReal)PETSC_CURRENT) {
1706:     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);
1707:     ksp->rtol = rtol;
1708:   }
1709:   if (abstol == (PetscReal)PETSC_DETERMINE) {
1710:     ksp->abstol = ksp->default_abstol;
1711:   } else if (abstol != (PetscReal)PETSC_CURRENT) {
1712:     PetscCheck(abstol >= 0.0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %g must be non-negative", (double)abstol);
1713:     ksp->abstol = abstol;
1714:   }
1715:   if (dtol == (PetscReal)PETSC_DETERMINE) {
1716:     ksp->divtol = ksp->default_divtol;
1717:   } else if (dtol == (PetscReal)PETSC_UNLIMITED) {
1718:     ksp->divtol = PETSC_MAX_REAL;
1719:   } else if (dtol != (PetscReal)PETSC_CURRENT) {
1720:     PetscCheck(dtol >= 0.0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Divergence tolerance %g must be larger than 1.0", (double)dtol);
1721:     ksp->divtol = dtol;
1722:   }
1723:   if (maxits == PETSC_DETERMINE) {
1724:     ksp->max_it = ksp->default_max_it;
1725:   } else if (maxits == PETSC_UNLIMITED) {
1726:     ksp->max_it = PETSC_INT_MAX;
1727:   } else if (maxits != PETSC_CURRENT) {
1728:     PetscCheck(maxits >= 0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", maxits);
1729:     ksp->max_it = maxits;
1730:   }
1731:   PetscFunctionReturn(PETSC_SUCCESS);
1732: }

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

1737:   Logically Collective

1739:   Input Parameters:
1740: + ksp   - the Krylov subspace context
1741: - minit - minimum number of iterations to use

1743:   Options Database Key:
1744: . -ksp_min_it <minits> - Sets `minit`

1746:   Level: intermediate

1748:   Notes:
1749:   Use `KSPSetTolerances()` to set a variety of other tolerances

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

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

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

1766:   PetscCheck(minit >= 0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Minimum number of iterations %" PetscInt_FMT " must be non-negative", minit);
1767:   ksp->min_it = minit;
1768:   PetscFunctionReturn(PETSC_SUCCESS);
1769: }

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

1774:   Not Collective

1776:   Input Parameter:
1777: . ksp - the Krylov subspace context

1779:   Output Parameter:
1780: . minit - minimum number of iterations to use

1782:   Level: intermediate

1784: .seealso: [](ch_ksp), `KSPGetTolerances()`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSP`, `KSPSetTolerances()`, `KSPSetMinimumIterations()`
1785: @*/
1786: PetscErrorCode KSPGetMinimumIterations(KSP ksp, PetscInt *minit)
1787: {
1788:   PetscFunctionBegin;
1790:   PetscAssertPointer(minit, 2);

1792:   *minit = ksp->min_it;
1793:   PetscFunctionReturn(PETSC_SUCCESS);
1794: }

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

1801:   Logically Collective

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

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

1810:   Level: beginner

1812: .seealso: [](ch_ksp), `KSPGetInitialGuessNonzero()`, `KSPGuessSetType()`, `KSPGuessType`, `KSP`
1813: @*/
1814: PetscErrorCode KSPSetInitialGuessNonzero(KSP ksp, PetscBool flg)
1815: {
1816:   PetscFunctionBegin;
1819:   ksp->guess_zero = (PetscBool)!flg;
1820:   PetscFunctionReturn(PETSC_SUCCESS);
1821: }

1823: /*@
1824:   KSPGetInitialGuessNonzero - Determines whether the `KSP` solver is using
1825:   a zero initial guess.

1827:   Not Collective

1829:   Input Parameter:
1830: . ksp - iterative solver obtained from `KSPCreate()`

1832:   Output Parameter:
1833: . flag - `PETSC_TRUE` if guess is nonzero, else `PETSC_FALSE`

1835:   Level: intermediate

1837: .seealso: [](ch_ksp), `KSPSetInitialGuessNonzero()`, `KSP`
1838: @*/
1839: PetscErrorCode KSPGetInitialGuessNonzero(KSP ksp, PetscBool *flag)
1840: {
1841:   PetscFunctionBegin;
1843:   PetscAssertPointer(flag, 2);
1844:   if (ksp->guess_zero) *flag = PETSC_FALSE;
1845:   else *flag = PETSC_TRUE;
1846:   PetscFunctionReturn(PETSC_SUCCESS);
1847: }

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

1852:   Logically Collective

1854:   Input Parameters:
1855: + ksp - iterative solver obtained from `KSPCreate()`
1856: - flg - `PETSC_TRUE` indicates you want the error generated

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

1861:   Level: intermediate

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

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

1870: .seealso: [](ch_ksp), `KSPGetErrorIfNotConverged()`, `KSP`
1871: @*/
1872: PetscErrorCode KSPSetErrorIfNotConverged(KSP ksp, PetscBool flg)
1873: {
1874:   PetscFunctionBegin;
1877:   ksp->errorifnotconverged = flg;
1878:   PetscFunctionReturn(PETSC_SUCCESS);
1879: }

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

1884:   Not Collective

1886:   Input Parameter:
1887: . ksp - iterative solver obtained from KSPCreate()

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

1892:   Level: intermediate

1894: .seealso: [](ch_ksp), `KSPSetErrorIfNotConverged()`, `KSP`
1895: @*/
1896: PetscErrorCode KSPGetErrorIfNotConverged(KSP ksp, PetscBool *flag)
1897: {
1898:   PetscFunctionBegin;
1900:   PetscAssertPointer(flag, 2);
1901:   *flag = ksp->errorifnotconverged;
1902:   PetscFunctionReturn(PETSC_SUCCESS);
1903: }

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

1908:   Logically Collective

1910:   Input Parameters:
1911: + ksp - iterative solver obtained from `KSPCreate()`
1912: - flg - `PETSC_TRUE` or `PETSC_FALSE`

1914:   Level: advanced

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

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

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

1935:   Not Collective

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

1940:   Output Parameter:
1941: . flag - `PETSC_TRUE` if using Knoll trick, else `PETSC_FALSE`

1943:   Level: advanced

1945: .seealso: [](ch_ksp), `KSPSetInitialGuessKnoll()`, `KSPSetInitialGuessNonzero()`, `KSPGetInitialGuessNonzero()`, `KSP`
1946: @*/
1947: PetscErrorCode KSPGetInitialGuessKnoll(KSP ksp, PetscBool *flag)
1948: {
1949:   PetscFunctionBegin;
1951:   PetscAssertPointer(flag, 2);
1952:   *flag = ksp->guess_knoll;
1953:   PetscFunctionReturn(PETSC_SUCCESS);
1954: }

1956: /*@
1957:   KSPGetComputeSingularValues - Gets the flag indicating whether the extreme singular
1958:   values will be calculated via a Lanczos or Arnoldi process as the linear
1959:   system is solved.

1961:   Not Collective

1963:   Input Parameter:
1964: . ksp - iterative solver obtained from `KSPCreate()`

1966:   Output Parameter:
1967: . flg - `PETSC_TRUE` or `PETSC_FALSE`

1969:   Options Database Key:
1970: . -ksp_monitor_singular_value - Activates `KSPSetComputeSingularValues()`

1972:   Level: advanced

1974:   Notes:
1975:   This option is not valid for `KSPType`.

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

1981: .seealso: [](ch_ksp), `KSPComputeExtremeSingularValues()`, `KSPMonitorSingularValue()`, `KSP`
1982: @*/
1983: PetscErrorCode KSPGetComputeSingularValues(KSP ksp, PetscBool *flg)
1984: {
1985:   PetscFunctionBegin;
1987:   PetscAssertPointer(flg, 2);
1988:   *flg = ksp->calc_sings;
1989:   PetscFunctionReturn(PETSC_SUCCESS);
1990: }

1992: /*@
1993:   KSPSetComputeSingularValues - Sets a flag so that the extreme singular
1994:   values will be calculated via a Lanczos or Arnoldi process as the linear
1995:   system is solved.

1997:   Logically Collective

1999:   Input Parameters:
2000: + ksp - iterative solver obtained from `KSPCreate()`
2001: - flg - `PETSC_TRUE` or `PETSC_FALSE`

2003:   Options Database Key:
2004: . -ksp_monitor_singular_value - Activates `KSPSetComputeSingularValues()`

2006:   Level: advanced

2008:   Notes:
2009:   This option is not valid for all iterative methods.

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

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

2017: .seealso: [](ch_ksp), `KSPComputeExtremeSingularValues()`, `KSPMonitorSingularValue()`, `KSP`, `KSPSetComputeRitz()`
2018: @*/
2019: PetscErrorCode KSPSetComputeSingularValues(KSP ksp, PetscBool flg)
2020: {
2021:   PetscFunctionBegin;
2024:   ksp->calc_sings = flg;
2025:   PetscFunctionReturn(PETSC_SUCCESS);
2026: }

2028: /*@
2029:   KSPGetComputeEigenvalues - Gets the flag indicating that the extreme eigenvalues
2030:   values will be calculated via a Lanczos or Arnoldi process as the linear
2031:   system is solved.

2033:   Not Collective

2035:   Input Parameter:
2036: . ksp - iterative solver obtained from `KSPCreate()`

2038:   Output Parameter:
2039: . flg - `PETSC_TRUE` or `PETSC_FALSE`

2041:   Level: advanced

2043:   Note:
2044:   Currently this option is not valid for all iterative methods.

2046: .seealso: [](ch_ksp), `KSPComputeEigenvalues()`, `KSPComputeEigenvaluesExplicitly()`, `KSP`, `KSPSetComputeRitz()`
2047: @*/
2048: PetscErrorCode KSPGetComputeEigenvalues(KSP ksp, PetscBool *flg)
2049: {
2050:   PetscFunctionBegin;
2052:   PetscAssertPointer(flg, 2);
2053:   *flg = ksp->calc_sings;
2054:   PetscFunctionReturn(PETSC_SUCCESS);
2055: }

2057: /*@
2058:   KSPSetComputeEigenvalues - Sets a flag so that the extreme eigenvalues
2059:   values will be calculated via a Lanczos or Arnoldi process as the linear
2060:   system is solved.

2062:   Logically Collective

2064:   Input Parameters:
2065: + ksp - iterative solver obtained from `KSPCreate()`
2066: - flg - `PETSC_TRUE` or `PETSC_FALSE`

2068:   Level: advanced

2070:   Note:
2071:   Currently this option is not valid for all iterative methods.

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

2075: .seealso: [](ch_ksp), `KSPComputeEigenvalues()`, `KSPComputeEigenvaluesExplicitly()`, `KSP`, `KSPSetComputeRitz()`
2076: @*/
2077: PetscErrorCode KSPSetComputeEigenvalues(KSP ksp, PetscBool flg)
2078: {
2079:   PetscFunctionBegin;
2082:   ksp->calc_sings = flg;
2083:   PetscFunctionReturn(PETSC_SUCCESS);
2084: }

2086: /*@
2087:   KSPSetComputeRitz - Sets a flag so that the Ritz or harmonic Ritz pairs
2088:   will be calculated via a Lanczos or Arnoldi process as the linear
2089:   system is solved.

2091:   Logically Collective

2093:   Input Parameters:
2094: + ksp - iterative solver obtained from `KSPCreate()`
2095: - flg - `PETSC_TRUE` or `PETSC_FALSE`

2097:   Level: advanced

2099:   Note:
2100:   Currently this option is only valid for the `KSPGMRES` method.

2102: .seealso: [](ch_ksp), `KSPComputeRitz()`, `KSP`, `KSPComputeEigenvalues()`, `KSPComputeExtremeSingularValues()`
2103: @*/
2104: PetscErrorCode KSPSetComputeRitz(KSP ksp, PetscBool flg)
2105: {
2106:   PetscFunctionBegin;
2109:   ksp->calc_ritz = flg;
2110:   PetscFunctionReturn(PETSC_SUCCESS);
2111: }

2113: /*@
2114:   KSPGetRhs - Gets the right-hand-side vector for the linear system to
2115:   be solved.

2117:   Not Collective

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

2122:   Output Parameter:
2123: . r - right-hand-side vector

2125:   Level: developer

2127: .seealso: [](ch_ksp), `KSPGetSolution()`, `KSPSolve()`, `KSP`
2128: @*/
2129: PetscErrorCode KSPGetRhs(KSP ksp, Vec *r)
2130: {
2131:   PetscFunctionBegin;
2133:   PetscAssertPointer(r, 2);
2134:   *r = ksp->vec_rhs;
2135:   PetscFunctionReturn(PETSC_SUCCESS);
2136: }

2138: /*@
2139:   KSPGetSolution - Gets the location of the solution for the
2140:   linear system to be solved.

2142:   Not Collective

2144:   Input Parameter:
2145: . ksp - iterative solver obtained from `KSPCreate()`

2147:   Output Parameter:
2148: . v - solution vector

2150:   Level: developer

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

2156: .seealso: [](ch_ksp), `KSPGetRhs()`, `KSPBuildSolution()`, `KSPSolve()`, `KSP`
2157: @*/
2158: PetscErrorCode KSPGetSolution(KSP ksp, Vec *v)
2159: {
2160:   PetscFunctionBegin;
2162:   PetscAssertPointer(v, 2);
2163:   *v = ksp->vec_sol;
2164:   PetscFunctionReturn(PETSC_SUCCESS);
2165: }

2167: /*@
2168:   KSPSetPC - Sets the preconditioner to be used to calculate the
2169:   application of the preconditioner on a vector into a `KSP`.

2171:   Collective

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

2177:   Level: developer

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

2183: .seealso: [](ch_ksp), `KSPGetPC()`, `KSP`
2184: @*/
2185: PetscErrorCode KSPSetPC(KSP ksp, PC pc)
2186: {
2187:   PetscFunctionBegin;
2189:   if (pc) {
2191:     PetscCheckSameComm(ksp, 1, pc, 2);
2192:   }
2193:   if (ksp->pc != pc && ksp->setupstage) ksp->setupstage = KSP_SETUP_NEWMATRIX;
2194:   PetscCall(PetscObjectReference((PetscObject)pc));
2195:   PetscCall(PCDestroy(&ksp->pc));
2196:   ksp->pc = pc;
2197:   PetscFunctionReturn(PETSC_SUCCESS);
2198: }

2200: PETSC_INTERN PetscErrorCode PCCreate_MPI(PC);

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

2206:    Collective, No Fortran Support

2208:    Input Parameter:
2209: .  ksp - iterative solver obtained from `KSPCreate()`

2211:    Level: developer

2213: .seealso: [](ch_ksp), `KSPSetPC()`, `KSP`, `PCMPIServerBegin()`, `PCMPIServerEnd()`
2214: @*/
2215: PETSC_INTERN PetscErrorCode KSPCheckPCMPI(KSP ksp)
2216: {
2217:   PetscBool isPCMPI;

2219:   PetscFunctionBegin;
2221:   PetscCall(PetscObjectTypeCompare((PetscObject)ksp->pc, PCMPI, &isPCMPI));
2222:   if (PCMPIServerActive && ksp->nestlevel == 0 && !isPCMPI) {
2223:     const char *prefix;
2224:     char       *found = NULL;

2226:     PetscCall(KSPGetOptionsPrefix(ksp, &prefix));
2227:     if (prefix) PetscCall(PetscStrstr(prefix, "mpi_linear_solver_server_", &found));
2228:     if (!found) PetscCall(KSPAppendOptionsPrefix(ksp, "mpi_linear_solver_server_"));
2229:     PetscCall(PetscInfo(NULL, "In MPI Linear Solver Server and detected (root) PC that must be changed to PCMPI\n"));
2230:     PetscCall(PCSetType(ksp->pc, PCMPI));
2231:   }
2232:   PetscFunctionReturn(PETSC_SUCCESS);
2233: }

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

2238:   Not Collective

2240:   Input Parameter:
2241: . ksp - iterative solver obtained from `KSPCreate()`

2243:   Output Parameter:
2244: . pc - preconditioner context

2246:   Level: beginner

2248:   Note:
2249:   The `PC` is created if it does not already exist.

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

2254: .seealso: [](ch_ksp), `KSPSetPC()`, `KSP`, `PC`
2255: @*/
2256: PetscErrorCode KSPGetPC(KSP ksp, PC *pc)
2257: {
2258:   PetscFunctionBegin;
2260:   PetscAssertPointer(pc, 2);
2261:   if (!ksp->pc) {
2262:     PetscCall(PCCreate(PetscObjectComm((PetscObject)ksp), &ksp->pc));
2263:     PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp->pc, (PetscObject)ksp, 0));
2264:     PetscCall(PetscObjectSetOptions((PetscObject)ksp->pc, ((PetscObject)ksp)->options));
2265:     PetscCall(PCSetKSPNestLevel(ksp->pc, ksp->nestlevel));
2266:   }
2267:   PetscCall(KSPCheckPCMPI(ksp));
2268:   *pc = ksp->pc;
2269:   PetscFunctionReturn(PETSC_SUCCESS);
2270: }

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

2275:   Collective

2277:   Input Parameters:
2278: + ksp   - iterative solver obtained from `KSPCreate()`
2279: . it    - iteration number
2280: - rnorm - relative norm of the residual

2282:   Level: developer

2284:   Notes:
2285:   This routine is called by the `KSP` implementations.
2286:   It does not typically need to be called by the user.

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

2291: .seealso: [](ch_ksp), `KSPMonitorSet()`
2292: @*/
2293: PetscErrorCode KSPMonitor(KSP ksp, PetscInt it, PetscReal rnorm)
2294: {
2295:   PetscInt i, n = ksp->numbermonitors;

2297:   PetscFunctionBegin;
2298:   for (i = 0; i < n; i++) PetscCall((*ksp->monitor[i])(ksp, it, rnorm, ksp->monitorcontext[i]));
2299:   PetscFunctionReturn(PETSC_SUCCESS);
2300: }

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

2306:   Logically Collective

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

2314:   Options Database Keys:
2315: + -ksp_monitor                             - sets `KSPMonitorResidual()`
2316: . -ksp_monitor hdf5:filename               - sets `KSPMonitorResidualView()` and saves residual
2317: . -ksp_monitor draw                        - sets `KSPMonitorResidualView()` and plots residual
2318: . -ksp_monitor draw::draw_lg               - sets `KSPMonitorResidualDrawLG()` and plots residual
2319: . -ksp_monitor_pause_final                 - Pauses any graphics when the solve finishes (only works for internal monitors)
2320: . -ksp_monitor_true_residual               - sets `KSPMonitorTrueResidual()`
2321: . -ksp_monitor_true_residual draw::draw_lg - sets `KSPMonitorTrueResidualDrawLG()` and plots residual
2322: . -ksp_monitor_max                         - sets `KSPMonitorTrueResidualMax()`
2323: . -ksp_monitor_singular_value              - sets `KSPMonitorSingularValue()`
2324: - -ksp_monitor_cancel                      - cancels all monitors that have been hardwired into a code by calls to `KSPMonitorSet()`, but
2325:                                              does not cancel those set via the options database.

2327:   Level: beginner

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

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

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

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

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

2346: .seealso: [](ch_ksp), `KSPMonitorResidual()`, `KSPMonitorRegister()`, `KSPMonitorCancel()`, `KSP`, `PetscCtxDestroyFn`
2347: @*/
2348: PetscErrorCode KSPMonitorSet(KSP ksp, KSPMonitorFn *monitor, PetscCtx ctx, PetscCtxDestroyFn *monitordestroy)
2349: {
2350:   PetscFunctionBegin;
2352:   for (PetscInt i = 0; i < ksp->numbermonitors; i++) {
2353:     PetscBool identical;

2355:     PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))(PetscVoidFn *)monitor, ctx, monitordestroy, (PetscErrorCode (*)(void))(PetscVoidFn *)ksp->monitor[i], ksp->monitorcontext[i], ksp->monitordestroy[i], &identical));
2356:     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
2357:   }
2358:   PetscCheck(ksp->numbermonitors < MAXKSPMONITORS, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Too many KSP monitors set");
2359:   ksp->monitor[ksp->numbermonitors]          = monitor;
2360:   ksp->monitordestroy[ksp->numbermonitors]   = monitordestroy;
2361:   ksp->monitorcontext[ksp->numbermonitors++] = ctx;
2362:   PetscFunctionReturn(PETSC_SUCCESS);
2363: }

2365: /*@
2366:   KSPMonitorCancel - Clears all monitors for a `KSP` object.

2368:   Logically Collective

2370:   Input Parameter:
2371: . ksp - iterative solver obtained from `KSPCreate()`

2373:   Options Database Key:
2374: . -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.

2376:   Level: intermediate

2378: .seealso: [](ch_ksp), `KSPMonitorResidual()`, `KSPMonitorSet()`, `KSP`
2379: @*/
2380: PetscErrorCode KSPMonitorCancel(KSP ksp)
2381: {
2382:   PetscInt i;

2384:   PetscFunctionBegin;
2386:   for (i = 0; i < ksp->numbermonitors; i++) {
2387:     if (ksp->monitordestroy[i]) PetscCall((*ksp->monitordestroy[i])(&ksp->monitorcontext[i]));
2388:   }
2389:   ksp->numbermonitors = 0;
2390:   PetscFunctionReturn(PETSC_SUCCESS);
2391: }

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

2396:   Not Collective

2398:   Input Parameter:
2399: . ksp - iterative solver obtained from `KSPCreate()`

2401:   Output Parameter:
2402: . ctx - monitoring context

2404:   Level: intermediate

2406:   Fortran Notes:
2407:   This only works when the context is a Fortran derived type or a `PetscObject`. Declare `ctx` with
2408: .vb
2409:   type(tUsertype), pointer :: ctx
2410: .ve

2412: .seealso: [](ch_ksp), `KSPMonitorResidual()`, `KSP`
2413: @*/
2414: PetscErrorCode KSPGetMonitorContext(KSP ksp, PetscCtxRt ctx)
2415: {
2416:   PetscFunctionBegin;
2418:   *(void **)ctx = ksp->monitorcontext[0];
2419:   PetscFunctionReturn(PETSC_SUCCESS);
2420: }

2422: /*@
2423:   KSPSetResidualHistory - Sets the array used to hold the residual history.
2424:   If set, this array will contain the residual norms computed at each
2425:   iteration of the solver.

2427:   Not Collective

2429:   Input Parameters:
2430: + ksp   - iterative solver obtained from `KSPCreate()`
2431: . a     - array to hold history
2432: . na    - size of `a`
2433: - reset - `PETSC_TRUE` indicates the history counter is reset to zero
2434:            for each new linear solve

2436:   Level: advanced

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

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

2445: .seealso: [](ch_ksp), `KSPGetResidualHistory()`, `KSP`
2446: @*/
2447: PetscErrorCode KSPSetResidualHistory(KSP ksp, PetscReal a[], PetscCount na, PetscBool reset)
2448: {
2449:   PetscFunctionBegin;

2452:   PetscCall(PetscFree(ksp->res_hist_alloc));
2453:   if (na != PETSC_DECIDE && na != PETSC_DEFAULT && a) {
2454:     ksp->res_hist     = a;
2455:     ksp->res_hist_max = na;
2456:   } else {
2457:     if (na != PETSC_DECIDE && na != PETSC_DEFAULT) ksp->res_hist_max = (size_t)na;
2458:     else ksp->res_hist_max = 10000; /* like default ksp->max_it */
2459:     PetscCall(PetscCalloc1(ksp->res_hist_max, &ksp->res_hist_alloc));

2461:     ksp->res_hist = ksp->res_hist_alloc;
2462:   }
2463:   ksp->res_hist_len   = 0;
2464:   ksp->res_hist_reset = reset;
2465:   PetscFunctionReturn(PETSC_SUCCESS);
2466: }

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

2471:   Not Collective

2473:   Input Parameter:
2474: . ksp - iterative solver obtained from `KSPCreate()`

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

2480:   Level: advanced

2482:   Note:
2483:   This array is borrowed and should not be freed by the caller.

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

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

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

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

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

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

2501: .seealso: [](ch_ksp), `KSPSetResidualHistory()`, `KSP`, `KSPGetIterationNumber()`, `KSPSTCG`, `KSPBCGSL`
2502: @*/
2503: PetscErrorCode KSPGetResidualHistory(KSP ksp, const PetscReal *a[], PetscInt *na)
2504: {
2505:   PetscFunctionBegin;
2507:   if (a) *a = ksp->res_hist;
2508:   if (na) PetscCall(PetscIntCast(ksp->res_hist_len, na));
2509:   PetscFunctionReturn(PETSC_SUCCESS);
2510: }

2512: /*@
2513:   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.

2515:   Not Collective

2517:   Input Parameters:
2518: + ksp   - iterative solver obtained from `KSPCreate()`
2519: . a     - array to hold history
2520: . na    - size of `a`
2521: - reset - `PETSC_TRUE` indicates the history counter is reset to zero for each new linear solve

2523:   Level: advanced

2525:   Notes:
2526:   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.
2527:   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.

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

2531: .seealso: [](ch_ksp), `KSPGetErrorHistory()`, `KSPSetResidualHistory()`, `KSP`
2532: @*/
2533: PetscErrorCode KSPSetErrorHistory(KSP ksp, PetscReal a[], PetscCount na, PetscBool reset)
2534: {
2535:   PetscFunctionBegin;

2538:   PetscCall(PetscFree(ksp->err_hist_alloc));
2539:   if (na != PETSC_DECIDE && na != PETSC_DEFAULT && a) {
2540:     ksp->err_hist     = a;
2541:     ksp->err_hist_max = na;
2542:   } else {
2543:     if (na != PETSC_DECIDE && na != PETSC_DEFAULT) ksp->err_hist_max = (size_t)na;
2544:     else ksp->err_hist_max = 10000; /* like default ksp->max_it */
2545:     PetscCall(PetscCalloc1(ksp->err_hist_max, &ksp->err_hist_alloc));
2546:     ksp->err_hist = ksp->err_hist_alloc;
2547:   }
2548:   ksp->err_hist_len   = 0;
2549:   ksp->err_hist_reset = reset;
2550:   PetscFunctionReturn(PETSC_SUCCESS);
2551: }

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

2556:   Not Collective

2558:   Input Parameter:
2559: . ksp - iterative solver obtained from `KSPCreate()`

2561:   Output Parameters:
2562: + a  - pointer to array to hold history (or `NULL`)
2563: - na - number of used entries in a (or `NULL`)

2565:   Level: advanced

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

2571:   Fortran Note:
2572: .vb
2573:   PetscReal, pointer :: a(:)
2574: .ve

2576: .seealso: [](ch_ksp), `KSPSetErrorHistory()`, `KSPGetResidualHistory()`, `KSP`
2577: @*/
2578: PetscErrorCode KSPGetErrorHistory(KSP ksp, const PetscReal *a[], PetscInt *na)
2579: {
2580:   PetscFunctionBegin;
2582:   if (a) *a = ksp->err_hist;
2583:   if (na) PetscCall(PetscIntCast(ksp->err_hist_len, na));
2584:   PetscFunctionReturn(PETSC_SUCCESS);
2585: }

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

2590:   Not Collective

2592:   Input Parameter:
2593: . ksp - The `KSP`

2595:   Output Parameters:
2596: + cr   - The residual contraction rate
2597: . rRsq - The coefficient of determination, $R^2$, indicating the linearity of the data
2598: . ce   - The error contraction rate
2599: - eRsq - The coefficient of determination, $R^2$, indicating the linearity of the data

2601:   Level: advanced

2603:   Note:
2604:   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,
2605:   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)}$,

2607: .seealso: [](ch_ksp), `KSP`, `KSPConvergedRateView()`
2608: @*/
2609: PetscErrorCode KSPComputeConvergenceRate(KSP ksp, PetscReal *cr, PetscReal *rRsq, PetscReal *ce, PetscReal *eRsq)
2610: {
2611:   PetscReal const *hist;
2612:   PetscReal       *x, *y, slope, intercept, mean = 0.0, var = 0.0, res = 0.0;
2613:   PetscInt         n, k;

2615:   PetscFunctionBegin;
2616:   if (cr || rRsq) {
2617:     PetscCall(KSPGetResidualHistory(ksp, &hist, &n));
2618:     if (!n) {
2619:       if (cr) *cr = 0.0;
2620:       if (rRsq) *rRsq = -1.0;
2621:     } else {
2622:       PetscCall(PetscMalloc2(n, &x, n, &y));
2623:       for (k = 0; k < n; ++k) {
2624:         x[k] = k;
2625:         y[k] = PetscLogReal(hist[k]);
2626:         mean += y[k];
2627:       }
2628:       mean /= n;
2629:       PetscCall(PetscLinearRegression(n, x, y, &slope, &intercept));
2630:       for (k = 0; k < n; ++k) {
2631:         res += PetscSqr(y[k] - (slope * x[k] + intercept));
2632:         var += PetscSqr(y[k] - mean);
2633:       }
2634:       PetscCall(PetscFree2(x, y));
2635:       if (cr) *cr = PetscExpReal(slope);
2636:       if (rRsq) *rRsq = var < PETSC_MACHINE_EPSILON ? 0.0 : 1.0 - (res / var);
2637:     }
2638:   }
2639:   if (ce || eRsq) {
2640:     PetscCall(KSPGetErrorHistory(ksp, &hist, &n));
2641:     if (!n) {
2642:       if (ce) *ce = 0.0;
2643:       if (eRsq) *eRsq = -1.0;
2644:     } else {
2645:       PetscCall(PetscMalloc2(n, &x, n, &y));
2646:       for (k = 0; k < n; ++k) {
2647:         x[k] = k;
2648:         y[k] = PetscLogReal(hist[k]);
2649:         mean += y[k];
2650:       }
2651:       mean /= n;
2652:       PetscCall(PetscLinearRegression(n, x, y, &slope, &intercept));
2653:       for (k = 0; k < n; ++k) {
2654:         res += PetscSqr(y[k] - (slope * x[k] + intercept));
2655:         var += PetscSqr(y[k] - mean);
2656:       }
2657:       PetscCall(PetscFree2(x, y));
2658:       if (ce) *ce = PetscExpReal(slope);
2659:       if (eRsq) *eRsq = var < PETSC_MACHINE_EPSILON ? 0.0 : 1.0 - (res / var);
2660:     }
2661:   }
2662:   PetscFunctionReturn(PETSC_SUCCESS);
2663: }

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

2668:   Logically Collective

2670:   Input Parameters:
2671: + ksp      - iterative solver obtained from `KSPCreate()`
2672: . converge - pointer to the function, see `KSPConvergenceTestFn`
2673: . ctx      - context for private data for the convergence routine (may be `NULL`)
2674: - destroy  - a routine for destroying the context (may be `NULL`)

2676:   Level: advanced

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

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

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

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

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

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

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 (may be `NULL`)
2716: - destroy  - a routine for destroying the context (may be `NULL`)

2718:   Level: advanced

2720: .seealso: [](ch_ksp), `KSP`, `KSPConvergedDefault()`, `KSPGetConvergenceContext()`, `KSPSetTolerances()`, `KSPSetConvergenceTest()`, `KSPGetAndClearConvergenceTest()`
2721: @*/
2722: PetscErrorCode KSPGetConvergenceTest(KSP ksp, KSPConvergenceTestFn **converge, PetscCtxRt ctx, PetscCtxDestroyFn **destroy)
2723: {
2724:   PetscFunctionBegin;
2726:   if (converge) *converge = ksp->converged;
2727:   if (destroy) *destroy = ksp->convergeddestroy;
2728:   if (ctx) *(void **)ctx = ksp->cnvP;
2729:   PetscFunctionReturn(PETSC_SUCCESS);
2730: }

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

2735:   Logically Collective

2737:   Input Parameter:
2738: . ksp - iterative solver obtained from `KSPCreate()`

2740:   Output Parameters:
2741: + converge - pointer to convergence test function, see `KSPConvergenceTestFn`
2742: . ctx      - context for private data for the convergence routine
2743: - destroy  - a routine for destroying the context

2745:   Level: advanced

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

2753: .seealso: [](ch_ksp), `KSP`, `KSPConvergedDefault()`, `KSPGetConvergenceContext()`, `KSPSetTolerances()`, `KSPSetConvergenceTest()`, `KSPGetConvergenceTest()`
2754: @*/
2755: PetscErrorCode KSPGetAndClearConvergenceTest(KSP ksp, KSPConvergenceTestFn **converge, PetscCtxRt ctx, PetscCtxDestroyFn **destroy)
2756: {
2757:   PetscFunctionBegin;
2759:   *converge             = ksp->converged;
2760:   *destroy              = ksp->convergeddestroy;
2761:   *(void **)ctx         = ksp->cnvP;
2762:   ksp->converged        = NULL;
2763:   ksp->cnvP             = NULL;
2764:   ksp->convergeddestroy = NULL;
2765:   PetscFunctionReturn(PETSC_SUCCESS);
2766: }

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

2771:   Not Collective

2773:   Input Parameter:
2774: . ksp - iterative solver obtained from `KSPCreate()`

2776:   Output Parameter:
2777: . ctx - monitoring context

2779:   Level: advanced

2781:   Fortran Note:
2782:   This only works when the context is a Fortran derived type or a `PetscObject`. Declare `ctx` with
2783: .vb
2784:   type(tUsertype), pointer :: ctx
2785: .ve

2787: .seealso: [](ch_ksp), `KSP`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSPGetConvergenceTest()`
2788: @*/
2789: PetscErrorCode KSPGetConvergenceContext(KSP ksp, PetscCtxRt ctx)
2790: {
2791:   PetscFunctionBegin;
2793:   *(void **)ctx = ksp->cnvP;
2794:   PetscFunctionReturn(PETSC_SUCCESS);
2795: }

2797: /*@
2798:   KSPBuildSolution - Builds the approximate solution in a vector provided.

2800:   Collective

2802:   Input Parameter:
2803: . ksp - iterative solver obtained from `KSPCreate()`

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

2810:   Level: developer

2812:   Notes:
2813:   This routine can be used in one of two ways
2814: .vb
2815:       KSPBuildSolution(ksp,NULL,&V);
2816:    or
2817:       KSPBuildSolution(ksp,v,NULL); or KSPBuildSolution(ksp,v,&v);
2818: .ve
2819:   In the first case an internal vector is allocated to store the solution
2820:   (the user cannot destroy this vector). In the second case the solution
2821:   is generated in the vector that the user provides. Note that for certain
2822:   methods, such as `KSPCG`, the second case requires a copy of the solution,
2823:   while in the first case the call is essentially free since it simply
2824:   returns the vector where the solution already is stored. For some methods
2825:   like `KSPGMRES` during the solve this is a reasonably expensive operation and should only be
2826:   used if truly needed.

2828: .seealso: [](ch_ksp), `KSPGetSolution()`, `KSPBuildResidual()`, `KSP`
2829: @*/
2830: PetscErrorCode KSPBuildSolution(KSP ksp, Vec v, Vec *V)
2831: {
2832:   PetscFunctionBegin;
2834:   PetscCheck(V || v, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_WRONG, "Must provide either v or V");
2835:   if (!V) V = &v;
2836:   if (ksp->reason != KSP_CONVERGED_ITERATING) {
2837:     if (!v) PetscCall(KSPGetSolution(ksp, V));
2838:     else PetscCall(VecCopy(ksp->vec_sol, v));
2839:   } else {
2840:     PetscUseTypeMethod(ksp, buildsolution, v, V);
2841:   }
2842:   PetscFunctionReturn(PETSC_SUCCESS);
2843: }

2845: /*@
2846:   KSPBuildResidual - Builds the residual in a vector provided.

2848:   Collective

2850:   Input Parameter:
2851: . ksp - iterative solver obtained from `KSPCreate()`

2853:   Output Parameters:
2854: + t - work vector.  If not provided then one is generated.
2855: . v - optional location to stash residual.  If `v` is not provided, then a location is generated.
2856: - V - the residual

2858:   Level: advanced

2860:   Note:
2861:   Regardless of whether or not `v` is provided, the residual is
2862:   returned in `V`.

2864: .seealso: [](ch_ksp), `KSP`, `KSPBuildSolution()`
2865: @*/
2866: PetscErrorCode KSPBuildResidual(KSP ksp, Vec t, Vec v, Vec *V)
2867: {
2868:   PetscBool flag = PETSC_FALSE;
2869:   Vec       w = v, tt = t;

2871:   PetscFunctionBegin;
2873:   if (!w) PetscCall(VecDuplicate(ksp->vec_rhs, &w));
2874:   if (!tt) {
2875:     PetscCall(VecDuplicate(ksp->vec_sol, &tt));
2876:     flag = PETSC_TRUE;
2877:   }
2878:   PetscUseTypeMethod(ksp, buildresidual, tt, w, V);
2879:   if (flag) PetscCall(VecDestroy(&tt));
2880:   PetscFunctionReturn(PETSC_SUCCESS);
2881: }

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

2887:   Logically Collective

2889:   Input Parameters:
2890: + ksp   - the `KSP` context
2891: - scale - `PETSC_TRUE` or `PETSC_FALSE`

2893:   Options Database Keys:
2894: + -ksp_diagonal_scale     - perform a diagonal scaling before the solve
2895: - -ksp_diagonal_scale_fix - scale the matrix back AFTER the solve

2897:   Level: advanced

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

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

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

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

2914: .seealso: [](ch_ksp), `KSPGetDiagonalScale()`, `KSPSetDiagonalScaleFix()`, `KSP`
2915: @*/
2916: PetscErrorCode KSPSetDiagonalScale(KSP ksp, PetscBool scale)
2917: {
2918:   PetscFunctionBegin;
2921:   ksp->dscale = scale;
2922:   PetscFunctionReturn(PETSC_SUCCESS);
2923: }

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

2928:   Not Collective

2930:   Input Parameter:
2931: . ksp - the `KSP` context

2933:   Output Parameter:
2934: . scale - `PETSC_TRUE` or `PETSC_FALSE`

2936:   Level: intermediate

2938: .seealso: [](ch_ksp), `KSP`, `KSPSetDiagonalScale()`, `KSPSetDiagonalScaleFix()`
2939: @*/
2940: PetscErrorCode KSPGetDiagonalScale(KSP ksp, PetscBool *scale)
2941: {
2942:   PetscFunctionBegin;
2944:   PetscAssertPointer(scale, 2);
2945:   *scale = ksp->dscale;
2946:   PetscFunctionReturn(PETSC_SUCCESS);
2947: }

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

2952:   Logically Collective

2954:   Input Parameters:
2955: + ksp - the `KSP` context
2956: - fix - `PETSC_TRUE` to scale back after the system solve, `PETSC_FALSE` to not
2957:          rescale (default)

2959:   Level: intermediate

2961:   Notes:
2962:   Must be called after `KSPSetDiagonalScale()`

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

2969: .seealso: [](ch_ksp), `KSPGetDiagonalScale()`, `KSPSetDiagonalScale()`, `KSPGetDiagonalScaleFix()`, `KSP`
2970: @*/
2971: PetscErrorCode KSPSetDiagonalScaleFix(KSP ksp, PetscBool fix)
2972: {
2973:   PetscFunctionBegin;
2976:   ksp->dscalefix = fix;
2977:   PetscFunctionReturn(PETSC_SUCCESS);
2978: }

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

2983:   Not Collective

2985:   Input Parameter:
2986: . ksp - the `KSP` context

2988:   Output Parameter:
2989: . fix - `PETSC_TRUE` to scale back after the system solve, `PETSC_FALSE` to not
2990:          rescale (default)

2992:   Level: intermediate

2994: .seealso: [](ch_ksp), `KSPGetDiagonalScale()`, `KSPSetDiagonalScale()`, `KSPSetDiagonalScaleFix()`, `KSP`
2995: @*/
2996: PetscErrorCode KSPGetDiagonalScaleFix(KSP ksp, PetscBool *fix)
2997: {
2998:   PetscFunctionBegin;
3000:   PetscAssertPointer(fix, 2);
3001:   *fix = ksp->dscalefix;
3002:   PetscFunctionReturn(PETSC_SUCCESS);
3003: }

3005: /*@C
3006:   KSPSetComputeOperators - set routine to compute the linear operators

3008:   Logically Collective

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

3015:   Level: beginner

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

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

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

3028: .seealso: [](ch_ksp), `KSP`, `KSPSetOperators()`, `KSPSetComputeRHS()`, `DMKSPSetComputeOperators()`, `KSPSetComputeInitialGuess()`, `KSPComputeOperatorsFn`
3029: @*/
3030: PetscErrorCode KSPSetComputeOperators(KSP ksp, KSPComputeOperatorsFn *func, PetscCtx ctx)
3031: {
3032:   DM dm;

3034:   PetscFunctionBegin;
3036:   PetscCall(KSPGetDM(ksp, &dm));
3037:   PetscCall(DMKSPSetComputeOperators(dm, func, ctx));
3038:   if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX;
3039:   PetscFunctionReturn(PETSC_SUCCESS);
3040: }

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

3045:   Logically Collective

3047:   Input Parameters:
3048: + ksp  - the `KSP` context
3049: . func - function to compute the right-hand side, see `KSPComputeRHSFn` for the calling sequence
3050: - ctx  - optional context

3052:   Level: beginner

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

3057: .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `DMKSPSetComputeRHS()`, `KSPSetComputeOperators()`, `KSPSetOperators()`, `KSPComputeRHSFn`
3058: @*/
3059: PetscErrorCode KSPSetComputeRHS(KSP ksp, KSPComputeRHSFn *func, PetscCtx ctx)
3060: {
3061:   DM dm;

3063:   PetscFunctionBegin;
3065:   PetscCall(KSPGetDM(ksp, &dm));
3066:   PetscCall(DMKSPSetComputeRHS(dm, func, ctx));
3067:   PetscFunctionReturn(PETSC_SUCCESS);
3068: }

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

3073:   Logically Collective

3075:   Input Parameters:
3076: + ksp  - the `KSP` context
3077: . func - function to compute the initial guess, see `KSPComputeInitialGuessFn` for calling sequence
3078: - ctx  - optional context

3080:   Level: beginner

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

3086: .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `KSPSetComputeRHS()`, `KSPSetComputeOperators()`, `DMKSPSetComputeInitialGuess()`, `KSPSetInitialGuessNonzero()`,
3087:           `KSPComputeInitialGuessFn`
3088: @*/
3089: PetscErrorCode KSPSetComputeInitialGuess(KSP ksp, KSPComputeInitialGuessFn *func, PetscCtx ctx)
3090: {
3091:   DM dm;

3093:   PetscFunctionBegin;
3095:   PetscCall(KSPGetDM(ksp, &dm));
3096:   PetscCall(DMKSPSetComputeInitialGuess(dm, func, ctx));
3097:   PetscFunctionReturn(PETSC_SUCCESS);
3098: }

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

3104:   Logically Collective

3106:   Input Parameter:
3107: . ksp - the `KSP` context

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

3112:   Level: advanced

3114: .seealso: [](ch_ksp), `KSPSolveTranspose()`, `KSP`
3115: @*/
3116: PetscErrorCode KSPSetUseExplicitTranspose(KSP ksp, PetscBool flg)
3117: {
3118:   PetscFunctionBegin;
3121:   ksp->transpose.use_explicittranspose = flg;
3122:   PetscFunctionReturn(PETSC_SUCCESS);
3123: }