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 context 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_eigenvalues) 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 context obtained from `KSPCreate()`
 79: - n   - size of arrays `r` and `c`. The number of eigenvalues computed `neig` will, in
 80:        general, be less than this.

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

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

 90:   Level: advanced

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

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

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

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

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

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

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

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

134:   Not Collective

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

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

147:   Level: advanced

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

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

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

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

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

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

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

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

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

178: .seealso: [](ch_ksp), `KSPSetComputeRitz()`, `KSP`, `KSPGMRES`, `KSPComputeEigenvalues()`, `KSPSetComputeSingularValues()`, `KSPMonitorSingularValue()`
179: @*/
180: PetscErrorCode KSPComputeRitz(KSP ksp, PetscBool ritz, PetscBool small, PetscInt *nrit, Vec S[], PetscReal tetar[], PetscReal tetai[])
181: {
182:   PetscFunctionBegin;
184:   PetscCheck(ksp->calc_ritz, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_WRONGSTATE, "Ritz pairs not requested before KSPSetUp()");
185:   PetscTryTypeMethod(ksp, computeritz, ritz, small, nrit, S, tetar, tetai);
186:   PetscFunctionReturn(PETSC_SUCCESS);
187: }
188: /*@
189:   KSPSetUpOnBlocks - Sets up the preconditioner for each block in
190:   the block Jacobi, block Gauss-Seidel, and overlapping Schwarz
191:   methods.

193:   Collective

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

198:   Level: advanced

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

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

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

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

233: /*@
234:   KSPSetReusePreconditioner - reuse the current preconditioner, do not construct a new one even if the operator changes

236:   Collective

238:   Input Parameters:
239: + ksp  - iterative context obtained from `KSPCreate()`
240: - flag - `PETSC_TRUE` to reuse the current preconditioner

242:   Level: intermediate

244: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSolve()`, `KSPDestroy()`, `PCSetReusePreconditioner()`, `KSP`
245: @*/
246: PetscErrorCode KSPSetReusePreconditioner(KSP ksp, PetscBool flag)
247: {
248:   PC pc;

250:   PetscFunctionBegin;
252:   PetscCall(KSPGetPC(ksp, &pc));
253:   PetscCall(PCSetReusePreconditioner(pc, flag));
254:   PetscFunctionReturn(PETSC_SUCCESS);
255: }

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

260:   Collective

262:   Input Parameter:
263: . ksp - iterative context obtained from `KSPCreate()`

265:   Output Parameter:
266: . flag - the boolean flag

268:   Level: intermediate

270: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSolve()`, `KSPDestroy()`, `KSPSetReusePreconditioner()`, `KSP`
271: @*/
272: PetscErrorCode KSPGetReusePreconditioner(KSP ksp, PetscBool *flag)
273: {
274:   PetscFunctionBegin;
276:   PetscAssertPointer(flag, 2);
277:   *flag = PETSC_FALSE;
278:   if (ksp->pc) PetscCall(PCGetReusePreconditioner(ksp->pc, flag));
279:   PetscFunctionReturn(PETSC_SUCCESS);
280: }

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

285:   Collective

287:   Input Parameters:
288: + ksp  - iterative context obtained from `KSPCreate()`
289: - flag - `PETSC_TRUE` to skip calling the `PCSetFromOptions()`

291:   Level: intermediate

293: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSolve()`, `KSPDestroy()`, `PCSetReusePreconditioner()`, `KSP`
294: @*/
295: PetscErrorCode KSPSetSkipPCSetFromOptions(KSP ksp, PetscBool flag)
296: {
297:   PetscFunctionBegin;
299:   ksp->skippcsetfromoptions = flag;
300:   PetscFunctionReturn(PETSC_SUCCESS);
301: }

303: /*@
304:   KSPSetUp - Sets up the internal data structures for the
305:   later use of an iterative solver.

307:   Collective

309:   Input Parameter:
310: . ksp - iterative context obtained from `KSPCreate()`

312:   Level: developer

314: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSolve()`, `KSPDestroy()`, `KSP`
315: @*/
316: PetscErrorCode KSPSetUp(KSP ksp)
317: {
318:   Mat            A, B;
319:   Mat            mat, pmat;
320:   MatNullSpace   nullsp;
321:   PCFailedReason pcreason;
322:   PC             pc;
323:   PetscBool      pcmpi;

325:   PetscFunctionBegin;
327:   PetscCall(KSPGetPC(ksp, &pc));
328:   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCMPI, &pcmpi));
329:   if (pcmpi) {
330:     PetscBool ksppreonly;
331:     PetscCall(PetscObjectTypeCompare((PetscObject)ksp, KSPPREONLY, &ksppreonly));
332:     if (!ksppreonly) PetscCall(KSPSetType(ksp, KSPPREONLY));
333:   }
334:   level++;

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

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

342:   if (ksp->dmActive && !ksp->setupstage) {
343:     /* first time in so build matrix and vector data structures using DM */
344:     if (!ksp->vec_rhs) PetscCall(DMCreateGlobalVector(ksp->dm, &ksp->vec_rhs));
345:     if (!ksp->vec_sol) PetscCall(DMCreateGlobalVector(ksp->dm, &ksp->vec_sol));
346:     PetscCall(DMCreateMatrix(ksp->dm, &A));
347:     PetscCall(KSPSetOperators(ksp, A, A));
348:     PetscCall(PetscObjectDereference((PetscObject)A));
349:   }

351:   if (ksp->dmActive) {
352:     DMKSP kdm;
353:     PetscCall(DMGetDMKSP(ksp->dm, &kdm));

355:     if (kdm->ops->computeinitialguess && ksp->setupstage != KSP_SETUP_NEWRHS) {
356:       /* only computes initial guess the first time through */
357:       PetscCallBack("KSP callback initial guess", (*kdm->ops->computeinitialguess)(ksp, ksp->vec_sol, kdm->initialguessctx));
358:       PetscCall(KSPSetInitialGuessNonzero(ksp, PETSC_TRUE));
359:     }
360:     if (kdm->ops->computerhs) PetscCallBack("KSP callback rhs", (*kdm->ops->computerhs)(ksp, ksp->vec_rhs, kdm->rhsctx));

362:     if (ksp->setupstage != KSP_SETUP_NEWRHS) {
363:       PetscCheck(kdm->ops->computeoperators, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_WRONGSTATE, "You called KSPSetDM() but did not use DMKSPSetComputeOperators() or KSPSetDMActive(ksp,PETSC_FALSE);");
364:       PetscCall(KSPGetOperators(ksp, &A, &B));
365:       PetscCallBack("KSP callback operators", (*kdm->ops->computeoperators)(ksp, A, B, kdm->operatorsctx));
366:     }
367:   }

369:   if (ksp->setupstage == KSP_SETUP_NEWRHS) {
370:     level--;
371:     PetscFunctionReturn(PETSC_SUCCESS);
372:   }
373:   PetscCall(PetscLogEventBegin(KSP_SetUp, ksp, ksp->vec_rhs, ksp->vec_sol, 0));

375:   switch (ksp->setupstage) {
376:   case KSP_SETUP_NEW:
377:     PetscUseTypeMethod(ksp, setup);
378:     break;
379:   case KSP_SETUP_NEWMATRIX: /* This should be replaced with a more general mechanism */
380:     if (ksp->setupnewmatrix) PetscUseTypeMethod(ksp, setup);
381:     break;
382:   default:
383:     break;
384:   }

386:   if (!ksp->pc) PetscCall(KSPGetPC(ksp, &ksp->pc));
387:   PetscCall(PCGetOperators(ksp->pc, &mat, &pmat));
388:   /* scale the matrix if requested */
389:   if (ksp->dscale) {
390:     PetscScalar *xx;
391:     PetscInt     i, n;
392:     PetscBool    zeroflag = PETSC_FALSE;

394:     if (!ksp->diagonal) { /* allocate vector to hold diagonal */
395:       PetscCall(MatCreateVecs(pmat, &ksp->diagonal, NULL));
396:     }
397:     PetscCall(MatGetDiagonal(pmat, ksp->diagonal));
398:     PetscCall(VecGetLocalSize(ksp->diagonal, &n));
399:     PetscCall(VecGetArray(ksp->diagonal, &xx));
400:     for (i = 0; i < n; i++) {
401:       if (xx[i] != 0.0) xx[i] = 1.0 / PetscSqrtReal(PetscAbsScalar(xx[i]));
402:       else {
403:         xx[i]    = 1.0;
404:         zeroflag = PETSC_TRUE;
405:       }
406:     }
407:     PetscCall(VecRestoreArray(ksp->diagonal, &xx));
408:     if (zeroflag) PetscCall(PetscInfo(ksp, "Zero detected in diagonal of matrix, using 1 at those locations\n"));
409:     PetscCall(MatDiagonalScale(pmat, ksp->diagonal, ksp->diagonal));
410:     if (mat != pmat) PetscCall(MatDiagonalScale(mat, ksp->diagonal, ksp->diagonal));
411:     ksp->dscalefix2 = PETSC_FALSE;
412:   }
413:   PetscCall(PetscLogEventEnd(KSP_SetUp, ksp, ksp->vec_rhs, ksp->vec_sol, 0));
414:   PetscCall(PCSetErrorIfFailure(ksp->pc, ksp->errorifnotconverged));
415:   PetscCall(PCSetUp(ksp->pc));
416:   PetscCall(PCGetFailedReasonRank(ksp->pc, &pcreason));
417:   /* 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 */
418:   if (pcreason) ksp->reason = KSP_DIVERGED_PC_FAILED;

420:   PetscCall(MatGetNullSpace(mat, &nullsp));
421:   if (nullsp) {
422:     PetscBool test = PETSC_FALSE;
423:     PetscCall(PetscOptionsGetBool(((PetscObject)ksp)->options, ((PetscObject)ksp)->prefix, "-ksp_test_null_space", &test, NULL));
424:     if (test) PetscCall(MatNullSpaceTest(nullsp, mat, NULL));
425:   }
426:   ksp->setupstage = KSP_SETUP_NEWRHS;
427:   level--;
428:   PetscFunctionReturn(PETSC_SUCCESS);
429: }

431: /*@C
432:   KSPConvergedReasonView - Displays the reason a `KSP` solve converged or diverged to a viewer

434:   Collective

436:   Input Parameters:
437: + ksp    - iterative context obtained from `KSPCreate()`
438: - viewer - the viewer to display the reason

440:   Options Database Keys:
441: + -ksp_converged_reason          - print reason for converged or diverged, also prints number of iterations
442: - -ksp_converged_reason ::failed - only print reason and number of iterations when diverged

444:   Level: beginner

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

450: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPDestroy()`, `KSPSetTolerances()`, `KSPConvergedDefault()`,
451:           `KSPSolveTranspose()`, `KSPGetIterationNumber()`, `KSP`, `KSPGetConvergedReason()`, `PetscViewerPushFormat()`, `PetscViewerPopFormat()`
452: @*/
453: PetscErrorCode KSPConvergedReasonView(KSP ksp, PetscViewer viewer)
454: {
455:   PetscBool         isAscii;
456:   PetscViewerFormat format;

458:   PetscFunctionBegin;
459:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)ksp));
460:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isAscii));
461:   if (isAscii) {
462:     PetscCall(PetscViewerGetFormat(viewer, &format));
463:     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ksp)->tablevel));
464:     if (ksp->reason > 0 && format != PETSC_VIEWER_FAILED) {
465:       if (((PetscObject)ksp)->prefix) {
466:         PetscCall(PetscViewerASCIIPrintf(viewer, "Linear %s solve converged due to %s iterations %" PetscInt_FMT "\n", ((PetscObject)ksp)->prefix, KSPConvergedReasons[ksp->reason], ksp->its));
467:       } else {
468:         PetscCall(PetscViewerASCIIPrintf(viewer, "Linear solve converged due to %s iterations %" PetscInt_FMT "\n", KSPConvergedReasons[ksp->reason], ksp->its));
469:       }
470:     } else if (ksp->reason <= 0) {
471:       if (((PetscObject)ksp)->prefix) {
472:         PetscCall(PetscViewerASCIIPrintf(viewer, "Linear %s solve did not converge due to %s iterations %" PetscInt_FMT "\n", ((PetscObject)ksp)->prefix, KSPConvergedReasons[ksp->reason], ksp->its));
473:       } else {
474:         PetscCall(PetscViewerASCIIPrintf(viewer, "Linear solve did not converge due to %s iterations %" PetscInt_FMT "\n", KSPConvergedReasons[ksp->reason], ksp->its));
475:       }
476:       if (ksp->reason == KSP_DIVERGED_PC_FAILED) {
477:         PCFailedReason reason;
478:         PetscCall(PCGetFailedReason(ksp->pc, &reason));
479:         PetscCall(PetscViewerASCIIPrintf(viewer, "               PC failed due to %s \n", PCFailedReasons[reason]));
480:       }
481:     }
482:     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ksp)->tablevel));
483:   }
484:   PetscFunctionReturn(PETSC_SUCCESS);
485: }

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

491:   Logically Collective

493:   Input Parameters:
494: + ksp               - the `KSP` context
495: . f                 - the ksp converged reason view function
496: . vctx              - [optional] user-defined context for private data for the
497:           ksp converged reason view routine (use `NULL` if no context is desired)
498: - reasonviewdestroy - [optional] routine that frees reasonview context
499:           (may be `NULL`)

501:   Options Database Keys:
502: + -ksp_converged_reason             - sets a default `KSPConvergedReasonView()`
503: - -ksp_converged_reason_view_cancel - cancels all converged reason viewers that have
504:                             been hardwired into a code by
505:                             calls to `KSPConvergedReasonViewSet()`, but
506:                             does not cancel those set via
507:                             the options database.

509:   Level: intermediate

511:   Note:
512:   Several different converged reason view routines may be set by calling
513:   `KSPConvergedReasonViewSet()` multiple times; all will be called in the
514:   order in which they were set.

516: .seealso: [](ch_ksp), `KSPConvergedReasonView()`, `KSPConvergedReasonViewCancel()`
517: @*/
518: PetscErrorCode KSPConvergedReasonViewSet(KSP ksp, PetscErrorCode (*f)(KSP, void *), void *vctx, PetscErrorCode (*reasonviewdestroy)(void **))
519: {
520:   PetscInt  i;
521:   PetscBool identical;

523:   PetscFunctionBegin;
525:   for (i = 0; i < ksp->numberreasonviews; i++) {
526:     PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))f, vctx, reasonviewdestroy, (PetscErrorCode(*)(void))ksp->reasonview[i], ksp->reasonviewcontext[i], ksp->reasonviewdestroy[i], &identical));
527:     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
528:   }
529:   PetscCheck(ksp->numberreasonviews < MAXKSPREASONVIEWS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many KSP reasonview set");
530:   ksp->reasonview[ksp->numberreasonviews]          = f;
531:   ksp->reasonviewdestroy[ksp->numberreasonviews]   = reasonviewdestroy;
532:   ksp->reasonviewcontext[ksp->numberreasonviews++] = (void *)vctx;
533:   PetscFunctionReturn(PETSC_SUCCESS);
534: }

536: /*@
537:   KSPConvergedReasonViewCancel - Clears all the reasonview functions for a `KSP` object set with `KSPConvergedReasonViewSet()`.

539:   Collective

541:   Input Parameter:
542: . ksp - iterative context obtained from `KSPCreate()`

544:   Level: intermediate

546: .seealso: [](ch_ksp), `KSPCreate()`, `KSPDestroy()`, `KSPReset()`, `KSPConvergedReasonViewSet()`
547: @*/
548: PetscErrorCode KSPConvergedReasonViewCancel(KSP ksp)
549: {
550:   PetscInt i;

552:   PetscFunctionBegin;
554:   for (i = 0; i < ksp->numberreasonviews; i++) {
555:     if (ksp->reasonviewdestroy[i]) PetscCall((*ksp->reasonviewdestroy[i])(&ksp->reasonviewcontext[i]));
556:   }
557:   ksp->numberreasonviews = 0;
558:   PetscFunctionReturn(PETSC_SUCCESS);
559: }

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

564:   Collective

566:   Input Parameter:
567: . ksp - the `KSP` object

569:   Level: intermediate

571: .seealso: [](ch_ksp), `KSPConvergedReasonView()`, `KSPConvergedReasonViewSet()`
572: @*/
573: PetscErrorCode KSPConvergedReasonViewFromOptions(KSP ksp)
574: {
575:   PetscViewer       viewer;
576:   PetscBool         flg;
577:   PetscViewerFormat format;
578:   PetscInt          i;

580:   PetscFunctionBegin;

582:   /* Call all user-provided reason review routines */
583:   for (i = 0; i < ksp->numberreasonviews; i++) PetscCall((*ksp->reasonview[i])(ksp, ksp->reasonviewcontext[i]));

585:   /* Call the default PETSc routine */
586:   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp), ((PetscObject)ksp)->options, ((PetscObject)ksp)->prefix, "-ksp_converged_reason", &viewer, &format, &flg));
587:   if (flg) {
588:     PetscCall(PetscViewerPushFormat(viewer, format));
589:     PetscCall(KSPConvergedReasonView(ksp, viewer));
590:     PetscCall(PetscViewerPopFormat(viewer));
591:     PetscCall(PetscViewerDestroy(&viewer));
592:   }
593:   PetscFunctionReturn(PETSC_SUCCESS);
594: }

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

599:   Collective

601:   Input Parameters:
602: + ksp    - iterative context obtained from `KSPCreate()`
603: - viewer - the viewer to display the reason

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

608:   Level: intermediate

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

613:   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,
614:   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)}$,

616: .seealso: [](ch_ksp), `KSPConvergedReasonView()`, `KSPGetConvergedRate()`, `KSPSetTolerances()`, `KSPConvergedDefault()`
617: @*/
618: PetscErrorCode KSPConvergedRateView(KSP ksp, PetscViewer viewer)
619: {
620:   PetscViewerFormat format;
621:   PetscBool         isAscii;
622:   PetscReal         rrate, rRsq, erate = 0.0, eRsq = 0.0;
623:   PetscInt          its;
624:   const char       *prefix, *reason = KSPConvergedReasons[ksp->reason];

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

662: #include <petscdraw.h>

664: static PetscErrorCode KSPViewEigenvalues_Internal(KSP ksp, PetscBool isExplicit, PetscViewer viewer, PetscViewerFormat format)
665: {
666:   PetscReal  *r, *c;
667:   PetscInt    n, i, neig;
668:   PetscBool   isascii, isdraw;
669:   PetscMPIInt rank;

671:   PetscFunctionBegin;
672:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)ksp), &rank));
673:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
674:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
675:   if (isExplicit) {
676:     PetscCall(VecGetSize(ksp->vec_sol, &n));
677:     PetscCall(PetscMalloc2(n, &r, n, &c));
678:     PetscCall(KSPComputeEigenvaluesExplicitly(ksp, n, r, c));
679:     neig = n;
680:   } else {
681:     PetscInt nits;

683:     PetscCall(KSPGetIterationNumber(ksp, &nits));
684:     n = nits + 2;
685:     if (!nits) {
686:       PetscCall(PetscViewerASCIIPrintf(viewer, "Zero iterations in solver, cannot approximate any eigenvalues\n"));
687:       PetscFunctionReturn(PETSC_SUCCESS);
688:     }
689:     PetscCall(PetscMalloc2(n, &r, n, &c));
690:     PetscCall(KSPComputeEigenvalues(ksp, n, r, c, &neig));
691:   }
692:   if (isascii) {
693:     PetscCall(PetscViewerASCIIPrintf(viewer, "%s computed eigenvalues\n", isExplicit ? "Explicitly" : "Iteratively"));
694:     for (i = 0; i < neig; ++i) {
695:       if (c[i] >= 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, "%g + %gi\n", (double)r[i], (double)c[i]));
696:       else PetscCall(PetscViewerASCIIPrintf(viewer, "%g - %gi\n", (double)r[i], -(double)c[i]));
697:     }
698:   } else if (isdraw && rank == 0) {
699:     PetscDraw   draw;
700:     PetscDrawSP drawsp;

702:     if (format == PETSC_VIEWER_DRAW_CONTOUR) {
703:       PetscCall(KSPPlotEigenContours_Private(ksp, neig, r, c));
704:     } else {
705:       PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
706:       PetscCall(PetscDrawSPCreate(draw, 1, &drawsp));
707:       PetscCall(PetscDrawSPReset(drawsp));
708:       for (i = 0; i < neig; ++i) PetscCall(PetscDrawSPAddPoint(drawsp, r + i, c + i));
709:       PetscCall(PetscDrawSPDraw(drawsp, PETSC_TRUE));
710:       PetscCall(PetscDrawSPSave(drawsp));
711:       PetscCall(PetscDrawSPDestroy(&drawsp));
712:     }
713:   }
714:   PetscCall(PetscFree2(r, c));
715:   PetscFunctionReturn(PETSC_SUCCESS);
716: }

718: static PetscErrorCode KSPViewSingularvalues_Internal(KSP ksp, PetscViewer viewer, PetscViewerFormat format)
719: {
720:   PetscReal smax, smin;
721:   PetscInt  nits;
722:   PetscBool isascii;

724:   PetscFunctionBegin;
725:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
726:   PetscCall(KSPGetIterationNumber(ksp, &nits));
727:   if (!nits) {
728:     PetscCall(PetscViewerASCIIPrintf(viewer, "Zero iterations in solver, cannot approximate any singular values\n"));
729:     PetscFunctionReturn(PETSC_SUCCESS);
730:   }
731:   PetscCall(KSPComputeExtremeSingularValues(ksp, &smax, &smin));
732:   if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, "Iteratively computed extreme singular values: max %g min %g max/min %g\n", (double)smax, (double)smin, (double)(smax / smin)));
733:   PetscFunctionReturn(PETSC_SUCCESS);
734: }

736: static PetscErrorCode KSPViewFinalResidual_Internal(KSP ksp, PetscViewer viewer, PetscViewerFormat format)
737: {
738:   PetscBool isascii;

740:   PetscFunctionBegin;
741:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
742:   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");
743:   if (isascii) {
744:     Mat       A;
745:     Vec       t;
746:     PetscReal norm;

748:     PetscCall(PCGetOperators(ksp->pc, &A, NULL));
749:     PetscCall(VecDuplicate(ksp->vec_rhs, &t));
750:     PetscCall(KSP_MatMult(ksp, A, ksp->vec_sol, t));
751:     PetscCall(VecAYPX(t, -1.0, ksp->vec_rhs));
752:     PetscCall(VecViewFromOptions(t, (PetscObject)ksp, "-ksp_view_final_residual_vec"));
753:     PetscCall(VecNorm(t, NORM_2, &norm));
754:     PetscCall(VecDestroy(&t));
755:     PetscCall(PetscViewerASCIIPrintf(viewer, "KSP final norm of residual %g\n", (double)norm));
756:   }
757:   PetscFunctionReturn(PETSC_SUCCESS);
758: }

760: static PetscErrorCode KSPMonitorPauseFinal_Internal(KSP ksp)
761: {
762:   PetscInt i;

764:   PetscFunctionBegin;
765:   if (!ksp->pauseFinal) PetscFunctionReturn(PETSC_SUCCESS);
766:   for (i = 0; i < ksp->numbermonitors; ++i) {
767:     PetscViewerAndFormat *vf = (PetscViewerAndFormat *)ksp->monitorcontext[i];
768:     PetscDraw             draw;
769:     PetscReal             lpause;

771:     if (!vf) continue;
772:     if (vf->lg) {
773:       if (!PetscCheckPointer(vf->lg, PETSC_OBJECT)) continue;
774:       if (((PetscObject)vf->lg)->classid != PETSC_DRAWLG_CLASSID) continue;
775:       PetscCall(PetscDrawLGGetDraw(vf->lg, &draw));
776:       PetscCall(PetscDrawGetPause(draw, &lpause));
777:       PetscCall(PetscDrawSetPause(draw, -1.0));
778:       PetscCall(PetscDrawPause(draw));
779:       PetscCall(PetscDrawSetPause(draw, lpause));
780:     } else {
781:       PetscBool isdraw;

783:       if (!PetscCheckPointer(vf->viewer, PETSC_OBJECT)) continue;
784:       if (((PetscObject)vf->viewer)->classid != PETSC_VIEWER_CLASSID) continue;
785:       PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERDRAW, &isdraw));
786:       if (!isdraw) continue;
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:   }
794:   PetscFunctionReturn(PETSC_SUCCESS);
795: }

797: static PetscErrorCode KSPSolve_Private(KSP ksp, Vec b, Vec x)
798: {
799:   PetscBool    flg = PETSC_FALSE, inXisinB = PETSC_FALSE, guess_zero;
800:   Mat          mat, pmat;
801:   MPI_Comm     comm;
802:   MatNullSpace nullsp;
803:   Vec          btmp, vec_rhs = NULL;

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

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

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

828:   /* reset the residual history list if requested */
829:   if (ksp->res_hist_reset) ksp->res_hist_len = 0;
830:   if (ksp->err_hist_reset) ksp->err_hist_len = 0;

832:   /* KSPSetUp() scales the matrix if needed */
833:   PetscCall(KSPSetUp(ksp));
834:   PetscCall(KSPSetUpOnBlocks(ksp));

836:   if (ksp->guess) {
837:     PetscObjectState ostate, state;

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

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

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

862:       PetscCall(PCGetOperators(ksp->pc, &mat, &pmat));
863:       PetscCall(MatDiagonalScale(pmat, ksp->diagonal, ksp->diagonal));
864:       if (mat != pmat) PetscCall(MatDiagonalScale(mat, ksp->diagonal, ksp->diagonal));
865:     }

867:     /* scale initial guess */
868:     if (!ksp->guess_zero) {
869:       if (!ksp->truediagonal) {
870:         PetscCall(VecDuplicate(ksp->diagonal, &ksp->truediagonal));
871:         PetscCall(VecCopy(ksp->diagonal, ksp->truediagonal));
872:         PetscCall(VecReciprocal(ksp->truediagonal));
873:       }
874:       PetscCall(VecPointwiseMult(ksp->vec_sol, ksp->vec_sol, ksp->truediagonal));
875:     }
876:   }
877:   PetscCall(PCPreSolve(ksp->pc, ksp));

879:   if (ksp->guess_zero && !ksp->guess_not_read) PetscCall(VecSet(ksp->vec_sol, 0.0));
880:   if (ksp->guess_knoll) { /* The Knoll trick is independent on the KSPGuess specified */
881:     PetscCall(PCApply(ksp->pc, ksp->vec_rhs, ksp->vec_sol));
882:     PetscCall(KSP_RemoveNullSpace(ksp, ksp->vec_sol));
883:     ksp->guess_zero = PETSC_FALSE;
884:   }

886:   /* can we mark the initial guess as zero for this solve? */
887:   guess_zero = ksp->guess_zero;
888:   if (!ksp->guess_zero) {
889:     PetscReal norm;

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

910:   PetscCall(VecLockReadPop(ksp->vec_rhs));
911:   if (nullsp) {
912:     ksp->vec_rhs = vec_rhs;
913:     PetscCall(VecDestroy(&btmp));
914:   }

916:   ksp->guess_zero = guess_zero;

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

921:   PetscCall(KSPConvergedReasonViewFromOptions(ksp));

923:   if (ksp->viewRate) {
924:     PetscCall(PetscViewerPushFormat(ksp->viewerRate, ksp->formatRate));
925:     PetscCall(KSPConvergedRateView(ksp, ksp->viewerRate));
926:     PetscCall(PetscViewerPopFormat(ksp->viewerRate));
927:   }
928:   PetscCall(PCPostSolve(ksp->pc, ksp));

930:   /* diagonal scale solution if called for */
931:   if (ksp->dscale) {
932:     PetscCall(VecPointwiseMult(ksp->vec_sol, ksp->vec_sol, ksp->diagonal));
933:     /* unscale right hand side and matrix */
934:     if (ksp->dscalefix) {
935:       Mat mat, pmat;

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

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

964:     PetscCall(PCGetOperators(ksp->pc, &A, NULL));
965:     if (ksp->transpose_solve) {
966:       Mat AT;

968:       PetscCall(MatCreateTranspose(A, &AT));
969:       PetscCall(MatComputeOperator(AT, MATAIJ, &B));
970:       PetscCall(MatDestroy(&AT));
971:     } else {
972:       PetscCall(MatComputeOperator(A, MATAIJ, &B));
973:     }
974:     PetscCall(ObjectView((PetscObject)B, ksp->viewerMatExp, ksp->formatMatExp));
975:     PetscCall(MatDestroy(&B));
976:   }
977:   if (ksp->viewPOpExp) {
978:     Mat B;

980:     PetscCall(KSPComputeOperator(ksp, MATAIJ, &B));
981:     PetscCall(ObjectView((PetscObject)B, ksp->viewerPOpExp, ksp->formatPOpExp));
982:     PetscCall(MatDestroy(&B));
983:   }

985:   if (inXisinB) {
986:     PetscCall(VecCopy(x, b));
987:     PetscCall(VecDestroy(&x));
988:   }
989:   PetscCall(PetscObjectSAWsBlock((PetscObject)ksp));
990:   if (ksp->errorifnotconverged && ksp->reason < 0 && ((level == 1) || (ksp->reason != KSP_DIVERGED_ITS))) {
991:     PCFailedReason reason;

993:     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]);
994:     PetscCall(PCGetFailedReason(ksp->pc, &reason));
995:     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]);
996:   }
997:   level--;
998:   PetscFunctionReturn(PETSC_SUCCESS);
999: }

1001: /*@
1002:   KSPSolve - Solves linear system.

1004:   Collective

1006:   Input Parameters:
1007: + ksp - iterative context obtained from `KSPCreate()`
1008: . b   - the right hand side vector
1009: - x   - the solution (this may be the same vector as `b`, then `b` will be overwritten with answer)

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

1027:   Level: beginner

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

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

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

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

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

1044:   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()`

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

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

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

1059:   If using a direct method (e.g., via the `KSP` solver
1060:   `KSPPREONLY` and a preconditioner such as `PCLU` or `PCILU`,
1061:   then its=1.  See `KSPSetTolerances()` and `KSPConvergedDefault()`
1062:   for more details.

1064:   Understanding Convergence\:
1065:   The routines `KSPMonitorSet()`, `KSPComputeEigenvalues()`, and
1066:   `KSPComputeEigenvaluesExplicitly()` provide information on additional
1067:   options to monitor convergence and print eigenvalue information.

1069: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPDestroy()`, `KSPSetTolerances()`, `KSPConvergedDefault()`,
1070:           `KSPSolveTranspose()`, `KSPGetIterationNumber()`, `MatNullSpaceCreate()`, `MatSetNullSpace()`, `MatSetTransposeNullSpace()`, `KSP`,
1071:           `KSPConvergedReasonView()`, `KSPCheckSolve()`, `KSPSetErrorIfNotConverged()`
1072: @*/
1073: PetscErrorCode KSPSolve(KSP ksp, Vec b, Vec x)
1074: {
1075:   PetscFunctionBegin;
1079:   ksp->transpose_solve = PETSC_FALSE;
1080:   PetscCall(KSPSolve_Private(ksp, b, x));
1081:   PetscFunctionReturn(PETSC_SUCCESS);
1082: }

1084: /*@
1085:   KSPSolveTranspose - Solves a linear system with the transpose of the matrix, $ A^T x = b$.

1087:   Collective

1089:   Input Parameters:
1090: + ksp - iterative context obtained from `KSPCreate()`
1091: . b   - right hand side vector
1092: - x   - solution vector

1094:   Level: developer

1096:   Note:
1097:   For complex numbers this solve the non-Hermitian transpose system.

1099:   Developer Note:
1100:   We need to implement a `KSPSolveHermitianTranspose()`

1102: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPDestroy()`, `KSPSetTolerances()`, `KSPConvergedDefault()`,
1103:           `KSPSolve()`, `KSP`
1104: @*/
1105: PetscErrorCode KSPSolveTranspose(KSP ksp, Vec b, Vec x)
1106: {
1107:   PetscFunctionBegin;
1111:   if (ksp->transpose.use_explicittranspose) {
1112:     Mat J, Jpre;
1113:     PetscCall(KSPGetOperators(ksp, &J, &Jpre));
1114:     if (!ksp->transpose.reuse_transpose) {
1115:       PetscCall(MatTranspose(J, MAT_INITIAL_MATRIX, &ksp->transpose.AT));
1116:       if (J != Jpre) PetscCall(MatTranspose(Jpre, MAT_INITIAL_MATRIX, &ksp->transpose.BT));
1117:       ksp->transpose.reuse_transpose = PETSC_TRUE;
1118:     } else {
1119:       PetscCall(MatTranspose(J, MAT_REUSE_MATRIX, &ksp->transpose.AT));
1120:       if (J != Jpre) PetscCall(MatTranspose(Jpre, MAT_REUSE_MATRIX, &ksp->transpose.BT));
1121:     }
1122:     if (J == Jpre && ksp->transpose.BT != ksp->transpose.AT) {
1123:       PetscCall(PetscObjectReference((PetscObject)ksp->transpose.AT));
1124:       ksp->transpose.BT = ksp->transpose.AT;
1125:     }
1126:     PetscCall(KSPSetOperators(ksp, ksp->transpose.AT, ksp->transpose.BT));
1127:   } else {
1128:     ksp->transpose_solve = PETSC_TRUE;
1129:   }
1130:   PetscCall(KSPSolve_Private(ksp, b, x));
1131:   PetscFunctionReturn(PETSC_SUCCESS);
1132: }

1134: static PetscErrorCode KSPViewFinalMatResidual_Internal(KSP ksp, Mat B, Mat X, PetscViewer viewer, PetscViewerFormat format, PetscInt shift)
1135: {
1136:   Mat        A, R;
1137:   PetscReal *norms;
1138:   PetscInt   i, N;
1139:   PetscBool  flg;

1141:   PetscFunctionBegin;
1142:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &flg));
1143:   if (flg) {
1144:     PetscCall(PCGetOperators(ksp->pc, &A, NULL));
1145:     if (!ksp->transpose_solve) PetscCall(MatMatMult(A, X, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &R));
1146:     else PetscCall(MatTransposeMatMult(A, X, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &R));
1147:     PetscCall(MatAYPX(R, -1.0, B, SAME_NONZERO_PATTERN));
1148:     PetscCall(MatGetSize(R, NULL, &N));
1149:     PetscCall(PetscMalloc1(N, &norms));
1150:     PetscCall(MatGetColumnNorms(R, NORM_2, norms));
1151:     PetscCall(MatDestroy(&R));
1152:     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]));
1153:     PetscCall(PetscFree(norms));
1154:   }
1155:   PetscFunctionReturn(PETSC_SUCCESS);
1156: }

1158: static PetscErrorCode KSPMatSolve_Private(KSP ksp, Mat B, Mat X)
1159: {
1160:   Mat       A, P, vB, vX;
1161:   Vec       cb, cx;
1162:   PetscInt  n1, N1, n2, N2, Bbn = PETSC_DECIDE;
1163:   PetscBool match;

1165:   PetscFunctionBegin;
1169:   PetscCheckSameComm(ksp, 1, B, 2);
1170:   PetscCheckSameComm(ksp, 1, X, 3);
1171:   PetscCheckSameType(B, 2, X, 3);
1172:   PetscCheck(B->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
1173:   MatCheckPreallocated(X, 3);
1174:   if (!X->assembled) {
1175:     PetscCall(MatSetOption(X, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
1176:     PetscCall(MatAssemblyBegin(X, MAT_FINAL_ASSEMBLY));
1177:     PetscCall(MatAssemblyEnd(X, MAT_FINAL_ASSEMBLY));
1178:   }
1179:   PetscCheck(B != X, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_IDN, "B and X must be different matrices");
1180:   PetscCheck(!ksp->transpose_solve || !ksp->transpose.use_explicittranspose, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "KSPMatSolveTranspose() does not support -ksp_use_explicittranspose");
1181:   PetscCall(KSPGetOperators(ksp, &A, &P));
1182:   PetscCall(MatGetLocalSize(B, NULL, &n2));
1183:   PetscCall(MatGetLocalSize(X, NULL, &n1));
1184:   PetscCall(MatGetSize(B, NULL, &N2));
1185:   PetscCall(MatGetSize(X, NULL, &N1));
1186:   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);
1187:   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)B, &match, MATSEQDENSE, MATMPIDENSE, ""));
1188:   PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of right-hand sides not stored in a dense Mat");
1189:   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)X, &match, MATSEQDENSE, MATMPIDENSE, ""));
1190:   PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of solutions not stored in a dense Mat");
1191:   PetscCall(KSPSetUp(ksp));
1192:   PetscCall(KSPSetUpOnBlocks(ksp));
1193:   if (ksp->ops->matsolve) {
1194:     level++;
1195:     if (ksp->guess_zero) PetscCall(MatZeroEntries(X));
1196:     PetscCall(PetscLogEventBegin(!ksp->transpose_solve ? KSP_MatSolve : KSP_MatSolveTranspose, ksp, B, X, 0));
1197:     PetscCall(KSPGetMatSolveBatchSize(ksp, &Bbn));
1198:     /* by default, do a single solve with all columns */
1199:     if (Bbn == PETSC_DECIDE) Bbn = N2;
1200:     else PetscCheck(Bbn >= 1, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "KSPMatSolve() batch size %" PetscInt_FMT " must be positive", Bbn);
1201:     PetscCall(PetscInfo(ksp, "KSP type %s solving using batches of width at most %" PetscInt_FMT "\n", ((PetscObject)ksp)->type_name, Bbn));
1202:     /* if -ksp_matsolve_batch_size is greater than the actual number of columns, do a single solve with all columns */
1203:     if (Bbn >= N2) {
1204:       PetscUseTypeMethod(ksp, matsolve, B, X);
1205:       if (ksp->viewFinalRes) PetscCall(KSPViewFinalMatResidual_Internal(ksp, B, X, ksp->viewerFinalRes, ksp->formatFinalRes, 0));

1207:       PetscCall(KSPConvergedReasonViewFromOptions(ksp));

1209:       if (ksp->viewRate) {
1210:         PetscCall(PetscViewerPushFormat(ksp->viewerRate, PETSC_VIEWER_DEFAULT));
1211:         PetscCall(KSPConvergedRateView(ksp, ksp->viewerRate));
1212:         PetscCall(PetscViewerPopFormat(ksp->viewerRate));
1213:       }
1214:     } else {
1215:       for (n2 = 0; n2 < N2; n2 += Bbn) {
1216:         PetscCall(MatDenseGetSubMatrix(B, PETSC_DECIDE, PETSC_DECIDE, n2, PetscMin(n2 + Bbn, N2), &vB));
1217:         PetscCall(MatDenseGetSubMatrix(X, PETSC_DECIDE, PETSC_DECIDE, n2, PetscMin(n2 + Bbn, N2), &vX));
1218:         PetscUseTypeMethod(ksp, matsolve, vB, vX);
1219:         if (ksp->viewFinalRes) PetscCall(KSPViewFinalMatResidual_Internal(ksp, vB, vX, ksp->viewerFinalRes, ksp->formatFinalRes, n2));

1221:         PetscCall(KSPConvergedReasonViewFromOptions(ksp));

1223:         if (ksp->viewRate) {
1224:           PetscCall(PetscViewerPushFormat(ksp->viewerRate, PETSC_VIEWER_DEFAULT));
1225:           PetscCall(KSPConvergedRateView(ksp, ksp->viewerRate));
1226:           PetscCall(PetscViewerPopFormat(ksp->viewerRate));
1227:         }
1228:         PetscCall(MatDenseRestoreSubMatrix(B, &vB));
1229:         PetscCall(MatDenseRestoreSubMatrix(X, &vX));
1230:       }
1231:     }
1232:     if (ksp->viewMat) PetscCall(ObjectView((PetscObject)A, ksp->viewerMat, ksp->formatMat));
1233:     if (ksp->viewPMat) PetscCall(ObjectView((PetscObject)P, ksp->viewerPMat, ksp->formatPMat));
1234:     if (ksp->viewRhs) PetscCall(ObjectView((PetscObject)B, ksp->viewerRhs, ksp->formatRhs));
1235:     if (ksp->viewSol) PetscCall(ObjectView((PetscObject)X, ksp->viewerSol, ksp->formatSol));
1236:     if (ksp->view) PetscCall(KSPView(ksp, ksp->viewer));
1237:     PetscCall(PetscLogEventEnd(!ksp->transpose_solve ? KSP_MatSolve : KSP_MatSolveTranspose, ksp, B, X, 0));
1238:     if (ksp->errorifnotconverged && ksp->reason < 0 && (level == 1 || ksp->reason != KSP_DIVERGED_ITS)) {
1239:       PCFailedReason reason;

1241:       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]);
1242:       PetscCall(PCGetFailedReason(ksp->pc, &reason));
1243:       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]);
1244:     }
1245:     level--;
1246:   } else {
1247:     PetscCall(PetscInfo(ksp, "KSP type %s solving column by column\n", ((PetscObject)ksp)->type_name));
1248:     for (n2 = 0; n2 < N2; ++n2) {
1249:       PetscCall(MatDenseGetColumnVecRead(B, n2, &cb));
1250:       PetscCall(MatDenseGetColumnVecWrite(X, n2, &cx));
1251:       PetscCall(KSPSolve_Private(ksp, cb, cx));
1252:       PetscCall(MatDenseRestoreColumnVecWrite(X, n2, &cx));
1253:       PetscCall(MatDenseRestoreColumnVecRead(B, n2, &cb));
1254:     }
1255:   }
1256:   PetscFunctionReturn(PETSC_SUCCESS);
1257: }

1259: /*@
1260:   KSPMatSolve - Solves a linear system with multiple right-hand sides stored as a `MATDENSE`. Unlike `KSPSolve()`, `B` and `X` must be different matrices.

1262:   Input Parameters:
1263: + ksp - iterative context
1264: - B   - block of right-hand sides

1266:   Output Parameter:
1267: . X - block of solutions

1269:   Level: intermediate

1271:   Note:
1272:   This is a stripped-down version of `KSPSolve()`, which only handles `-ksp_view`, `-ksp_converged_reason`, `-ksp_converged_rate`, and `-ksp_view_final_residual`.

1274: .seealso: [](ch_ksp), `KSPSolve()`, `MatMatSolve()`, `KSPMatSolveTranspose()`, `MATDENSE`, `KSPHPDDM`, `PCBJACOBI`, `PCASM`
1275: @*/
1276: PetscErrorCode KSPMatSolve(KSP ksp, Mat B, Mat X)
1277: {
1278:   PetscFunctionBegin;
1279:   ksp->transpose_solve = PETSC_FALSE;
1280:   PetscCall(KSPMatSolve_Private(ksp, B, X));
1281:   PetscFunctionReturn(PETSC_SUCCESS);
1282: }

1284: /*@
1285:   KSPMatSolveTranspose - Solves a linear system with the transposed matrix with multiple right-hand sides stored as a `MATDENSE`. Unlike `KSPSolveTranspose()`,
1286:   `B` and `X` must be different matrices and the transposed matrix cannot be assembled explicitly for the user.

1288:   Input Parameters:
1289: + ksp - iterative context
1290: - B   - block of right-hand sides

1292:   Output Parameter:
1293: . X - block of solutions

1295:   Level: intermediate

1297:   Note:
1298:   This is a stripped-down version of `KSPSolveTranspose()`, which only handles `-ksp_view`, `-ksp_converged_reason`, `-ksp_converged_rate`, and `-ksp_view_final_residual`.

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

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

1313:   Logically Collective

1315:   Input Parameters:
1316: + ksp - iterative context
1317: - bs  - batch size

1319:   Level: advanced

1321: .seealso: [](ch_ksp), `KSPMatSolve()`, `KSPGetMatSolveBatchSize()`, `-mat_mumps_icntl_27`, `-matmatmult_Bbn`
1322: @*/
1323: PetscErrorCode KSPSetMatSolveBatchSize(KSP ksp, PetscInt bs)
1324: {
1325:   PetscFunctionBegin;
1328:   ksp->nmax = bs;
1329:   PetscFunctionReturn(PETSC_SUCCESS);
1330: }

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

1335:   Input Parameter:
1336: . ksp - iterative context

1338:   Output Parameter:
1339: . bs - batch size

1341:   Level: advanced

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

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

1357:   Collective

1359:   Input Parameter:
1360: . ksp - iterative context obtained from `KSPCreate()`

1362:   Level: beginner

1364: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPSolve()`, `KSPSetFromOptions()`, `KSP`
1365: @*/
1366: PetscErrorCode KSPResetViewers(KSP ksp)
1367: {
1368:   PetscFunctionBegin;
1370:   if (!ksp) PetscFunctionReturn(PETSC_SUCCESS);
1371:   PetscCall(PetscViewerDestroy(&ksp->viewer));
1372:   PetscCall(PetscViewerDestroy(&ksp->viewerPre));
1373:   PetscCall(PetscViewerDestroy(&ksp->viewerRate));
1374:   PetscCall(PetscViewerDestroy(&ksp->viewerMat));
1375:   PetscCall(PetscViewerDestroy(&ksp->viewerPMat));
1376:   PetscCall(PetscViewerDestroy(&ksp->viewerRhs));
1377:   PetscCall(PetscViewerDestroy(&ksp->viewerSol));
1378:   PetscCall(PetscViewerDestroy(&ksp->viewerMatExp));
1379:   PetscCall(PetscViewerDestroy(&ksp->viewerEV));
1380:   PetscCall(PetscViewerDestroy(&ksp->viewerSV));
1381:   PetscCall(PetscViewerDestroy(&ksp->viewerEVExp));
1382:   PetscCall(PetscViewerDestroy(&ksp->viewerFinalRes));
1383:   PetscCall(PetscViewerDestroy(&ksp->viewerPOpExp));
1384:   PetscCall(PetscViewerDestroy(&ksp->viewerDScale));
1385:   ksp->view         = PETSC_FALSE;
1386:   ksp->viewPre      = PETSC_FALSE;
1387:   ksp->viewMat      = PETSC_FALSE;
1388:   ksp->viewPMat     = PETSC_FALSE;
1389:   ksp->viewRhs      = PETSC_FALSE;
1390:   ksp->viewSol      = PETSC_FALSE;
1391:   ksp->viewMatExp   = PETSC_FALSE;
1392:   ksp->viewEV       = PETSC_FALSE;
1393:   ksp->viewSV       = PETSC_FALSE;
1394:   ksp->viewEVExp    = PETSC_FALSE;
1395:   ksp->viewFinalRes = PETSC_FALSE;
1396:   ksp->viewPOpExp   = PETSC_FALSE;
1397:   ksp->viewDScale   = PETSC_FALSE;
1398:   PetscFunctionReturn(PETSC_SUCCESS);
1399: }

1401: /*@
1402:   KSPReset - Resets a `KSP` context to the kspsetupcalled = 0 state and removes any allocated Vecs and Mats

1404:   Collective

1406:   Input Parameter:
1407: . ksp - iterative context obtained from `KSPCreate()`

1409:   Level: beginner

1411: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPSolve()`, `KSP`
1412: @*/
1413: PetscErrorCode KSPReset(KSP ksp)
1414: {
1415:   PetscFunctionBegin;
1417:   if (!ksp) PetscFunctionReturn(PETSC_SUCCESS);
1418:   PetscTryTypeMethod(ksp, reset);
1419:   if (ksp->pc) PetscCall(PCReset(ksp->pc));
1420:   if (ksp->guess) {
1421:     KSPGuess guess = ksp->guess;
1422:     PetscTryTypeMethod(guess, reset);
1423:   }
1424:   PetscCall(VecDestroyVecs(ksp->nwork, &ksp->work));
1425:   PetscCall(VecDestroy(&ksp->vec_rhs));
1426:   PetscCall(VecDestroy(&ksp->vec_sol));
1427:   PetscCall(VecDestroy(&ksp->diagonal));
1428:   PetscCall(VecDestroy(&ksp->truediagonal));

1430:   PetscCall(KSPResetViewers(ksp));

1432:   ksp->setupstage = KSP_SETUP_NEW;
1433:   ksp->nmax       = PETSC_DECIDE;
1434:   PetscFunctionReturn(PETSC_SUCCESS);
1435: }

1437: /*@C
1438:   KSPDestroy - Destroys a `KSP` context.

1440:   Collective

1442:   Input Parameter:
1443: . ksp - iterative context obtained from `KSPCreate()`

1445:   Level: beginner

1447: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPSolve()`, `KSP`
1448: @*/
1449: PetscErrorCode KSPDestroy(KSP *ksp)
1450: {
1451:   PC pc;

1453:   PetscFunctionBegin;
1454:   if (!*ksp) PetscFunctionReturn(PETSC_SUCCESS);
1456:   if (--((PetscObject)(*ksp))->refct > 0) {
1457:     *ksp = NULL;
1458:     PetscFunctionReturn(PETSC_SUCCESS);
1459:   }

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

1463:   /*
1464:    Avoid a cascading call to PCReset(ksp->pc) from the following call:
1465:    PCReset() shouldn't be called from KSPDestroy() as it is unprotected by pc's
1466:    refcount (and may be shared, e.g., by other ksps).
1467:    */
1468:   pc         = (*ksp)->pc;
1469:   (*ksp)->pc = NULL;
1470:   PetscCall(KSPReset((*ksp)));
1471:   (*ksp)->pc = pc;
1472:   PetscTryTypeMethod((*ksp), destroy);

1474:   if ((*ksp)->transpose.use_explicittranspose) {
1475:     PetscCall(MatDestroy(&(*ksp)->transpose.AT));
1476:     PetscCall(MatDestroy(&(*ksp)->transpose.BT));
1477:     (*ksp)->transpose.reuse_transpose = PETSC_FALSE;
1478:   }

1480:   PetscCall(KSPGuessDestroy(&(*ksp)->guess));
1481:   PetscCall(DMDestroy(&(*ksp)->dm));
1482:   PetscCall(PCDestroy(&(*ksp)->pc));
1483:   PetscCall(PetscFree((*ksp)->res_hist_alloc));
1484:   PetscCall(PetscFree((*ksp)->err_hist_alloc));
1485:   if ((*ksp)->convergeddestroy) PetscCall((*(*ksp)->convergeddestroy)((*ksp)->cnvP));
1486:   PetscCall(KSPMonitorCancel((*ksp)));
1487:   PetscCall(KSPConvergedReasonViewCancel((*ksp)));
1488:   PetscCall(PetscHeaderDestroy(ksp));
1489:   PetscFunctionReturn(PETSC_SUCCESS);
1490: }

1492: /*@
1493:   KSPSetPCSide - Sets the preconditioning side.

1495:   Logically Collective

1497:   Input Parameter:
1498: . ksp - iterative context obtained from `KSPCreate()`

1500:   Output Parameter:
1501: . side - the preconditioning side, where side is one of
1502: .vb
1503:       PC_LEFT - left preconditioning (default)
1504:       PC_RIGHT - right preconditioning
1505:       PC_SYMMETRIC - symmetric preconditioning
1506: .ve

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

1511:   Level: intermediate

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

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

1518:   Symmetric preconditioning is currently available only for the `KSPQCG` method. However, note that
1519:   symmetric preconditioning can be emulated by using either right or left
1520:   preconditioning and a pre or post processing step.

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

1524: .seealso: [](ch_ksp), `KSPGetPCSide()`, `KSPSetNormType()`, `KSPGetNormType()`, `KSP`
1525: @*/
1526: PetscErrorCode KSPSetPCSide(KSP ksp, PCSide side)
1527: {
1528:   PetscFunctionBegin;
1531:   ksp->pc_side = ksp->pc_side_set = side;
1532:   PetscFunctionReturn(PETSC_SUCCESS);
1533: }

1535: /*@
1536:   KSPGetPCSide - Gets the preconditioning side.

1538:   Not Collective

1540:   Input Parameter:
1541: . ksp - iterative context obtained from `KSPCreate()`

1543:   Output Parameter:
1544: . side - the preconditioning side, where side is one of
1545: .vb
1546:       PC_LEFT - left preconditioning (default)
1547:       PC_RIGHT - right preconditioning
1548:       PC_SYMMETRIC - symmetric preconditioning
1549: .ve

1551:   Level: intermediate

1553: .seealso: [](ch_ksp), `KSPSetPCSide()`, `KSP`
1554: @*/
1555: PetscErrorCode KSPGetPCSide(KSP ksp, PCSide *side)
1556: {
1557:   PetscFunctionBegin;
1559:   PetscAssertPointer(side, 2);
1560:   PetscCall(KSPSetUpNorms_Private(ksp, PETSC_TRUE, &ksp->normtype, &ksp->pc_side));
1561:   *side = ksp->pc_side;
1562:   PetscFunctionReturn(PETSC_SUCCESS);
1563: }

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

1569:   Not Collective

1571:   Input Parameter:
1572: . ksp - the Krylov subspace context

1574:   Output Parameters:
1575: + rtol   - the relative convergence tolerance
1576: . abstol - the absolute convergence tolerance
1577: . dtol   - the divergence tolerance
1578: - maxits - maximum number of iterations

1580:   Level: intermediate

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

1585: .seealso: [](ch_ksp), `KSPSetTolerances()`, `KSP`, `KSPSetMinimumIterations()`, `KSPGetMinimumIterations()`
1586: @*/
1587: PetscErrorCode KSPGetTolerances(KSP ksp, PetscReal *rtol, PetscReal *abstol, PetscReal *dtol, PetscInt *maxits)
1588: {
1589:   PetscFunctionBegin;
1591:   if (abstol) *abstol = ksp->abstol;
1592:   if (rtol) *rtol = ksp->rtol;
1593:   if (dtol) *dtol = ksp->divtol;
1594:   if (maxits) *maxits = ksp->max_it;
1595:   PetscFunctionReturn(PETSC_SUCCESS);
1596: }

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

1602:   Logically Collective

1604:   Input Parameters:
1605: + ksp    - the Krylov subspace context
1606: . rtol   - the relative convergence tolerance, relative decrease in the (possibly preconditioned) residual norm
1607: . abstol - the absolute convergence tolerance   absolute size of the (possibly preconditioned) residual norm
1608: . dtol   - the divergence tolerance,   amount (possibly preconditioned) residual norm can increase before `KSPConvergedDefault()` concludes that the method is diverging
1609: - maxits - maximum number of iterations to use

1611:   Options Database Keys:
1612: + -ksp_atol <abstol>   - Sets `abstol`
1613: . -ksp_rtol <rtol>     - Sets `rtol`
1614: . -ksp_divtol <dtol>   - Sets `dtol`
1615: - -ksp_max_it <maxits> - Sets `maxits`

1617:   Level: intermediate

1619:   Notes:
1620:   Use `PETSC_DEFAULT` to retain the default value of any of the tolerances.

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

1625: .seealso: [](ch_ksp), `KSPGetTolerances()`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSP`, `KSPSetMinimumIterations()`
1626: @*/
1627: PetscErrorCode KSPSetTolerances(KSP ksp, PetscReal rtol, PetscReal abstol, PetscReal dtol, PetscInt maxits)
1628: {
1629:   PetscFunctionBegin;

1636:   if (rtol != (PetscReal)PETSC_DEFAULT) {
1637:     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);
1638:     ksp->rtol = rtol;
1639:   }
1640:   if (abstol != (PetscReal)PETSC_DEFAULT) {
1641:     PetscCheck(abstol >= 0.0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %g must be non-negative", (double)abstol);
1642:     ksp->abstol = abstol;
1643:   }
1644:   if (dtol != (PetscReal)PETSC_DEFAULT) {
1645:     PetscCheck(dtol >= 0.0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Divergence tolerance %g must be larger than 1.0", (double)dtol);
1646:     ksp->divtol = dtol;
1647:   }
1648:   if (maxits != PETSC_DEFAULT) {
1649:     PetscCheck(maxits >= 0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", maxits);
1650:     ksp->max_it = maxits;
1651:   }
1652:   PetscFunctionReturn(PETSC_SUCCESS);
1653: }

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

1658:   Logically Collective

1660:   Input Parameters:
1661: + ksp   - the Krylov subspace context
1662: - minit - minimum number of iterations to use

1664:   Options Database Key:
1665: . -ksp_min_it <minits> - Sets `minit`

1667:   Level: intermediate

1669:   Notes:
1670:   Use `KSPSetTolerances()` to set a variety of other tolerances

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

1675: .seealso: [](ch_ksp), `KSPGetTolerances()`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSP`, `KSPSetTolerances()`, `KSPGetMinimumIterations()`
1676: @*/
1677: PetscErrorCode KSPSetMinimumIterations(KSP ksp, PetscInt minit)
1678: {
1679:   PetscFunctionBegin;

1683:   PetscCheck(minit >= 0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Minimum number of iterations %" PetscInt_FMT " must be non-negative", minit);
1684:   ksp->min_it = minit;
1685:   PetscFunctionReturn(PETSC_SUCCESS);
1686: }

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

1691:   Not Collective

1693:   Input Parameter:
1694: . ksp - the Krylov subspace context

1696:   Output Parameter:
1697: . minit - minimum number of iterations to use

1699:   Level: intermediate

1701: .seealso: [](ch_ksp), `KSPGetTolerances()`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSP`, `KSPSetTolerances()`, `KSPSetMinimumIterations()`
1702: @*/
1703: PetscErrorCode KSPGetMinimumIterations(KSP ksp, PetscInt *minit)
1704: {
1705:   PetscFunctionBegin;
1707:   PetscAssertPointer(minit, 2);

1709:   *minit = ksp->min_it;
1710:   PetscFunctionReturn(PETSC_SUCCESS);
1711: }

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

1718:   Logically Collective

1720:   Input Parameters:
1721: + ksp - iterative context obtained from `KSPCreate()`
1722: - flg - ``PETSC_TRUE`` indicates the guess is non-zero, `PETSC_FALSE` indicates the guess is zero

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

1727:   Level: beginner

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

1732: .seealso: [](ch_ksp), `KSPGetInitialGuessNonzero()`, `KSPSetGuessType()`, `KSPGuessType`, `KSP`
1733: @*/
1734: PetscErrorCode KSPSetInitialGuessNonzero(KSP ksp, PetscBool flg)
1735: {
1736:   PetscFunctionBegin;
1739:   ksp->guess_zero = (PetscBool) !(int)flg;
1740:   PetscFunctionReturn(PETSC_SUCCESS);
1741: }

1743: /*@
1744:   KSPGetInitialGuessNonzero - Determines whether the `KSP` solver is using
1745:   a zero initial guess.

1747:   Not Collective

1749:   Input Parameter:
1750: . ksp - iterative context obtained from `KSPCreate()`

1752:   Output Parameter:
1753: . flag - `PETSC_TRUE` if guess is nonzero, else `PETSC_FALSE`

1755:   Level: intermediate

1757: .seealso: [](ch_ksp), `KSPSetInitialGuessNonzero()`, `KSP`
1758: @*/
1759: PetscErrorCode KSPGetInitialGuessNonzero(KSP ksp, PetscBool *flag)
1760: {
1761:   PetscFunctionBegin;
1763:   PetscAssertPointer(flag, 2);
1764:   if (ksp->guess_zero) *flag = PETSC_FALSE;
1765:   else *flag = PETSC_TRUE;
1766:   PetscFunctionReturn(PETSC_SUCCESS);
1767: }

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

1772:   Logically Collective

1774:   Input Parameters:
1775: + ksp - iterative context obtained from `KSPCreate()`
1776: - flg - `PETSC_TRUE` indicates you want the error generated

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

1781:   Level: intermediate

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

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

1789: .seealso: [](ch_ksp), `KSPGetErrorIfNotConverged()`, `KSP`
1790: @*/
1791: PetscErrorCode KSPSetErrorIfNotConverged(KSP ksp, PetscBool flg)
1792: {
1793:   PetscFunctionBegin;
1796:   ksp->errorifnotconverged = flg;
1797:   PetscFunctionReturn(PETSC_SUCCESS);
1798: }

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

1803:   Not Collective

1805:   Input Parameter:
1806: . ksp - iterative context obtained from KSPCreate()

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

1811:   Level: intermediate

1813: .seealso: [](ch_ksp), `KSPSetErrorIfNotConverged()`, `KSP`
1814: @*/
1815: PetscErrorCode KSPGetErrorIfNotConverged(KSP ksp, PetscBool *flag)
1816: {
1817:   PetscFunctionBegin;
1819:   PetscAssertPointer(flag, 2);
1820:   *flag = ksp->errorifnotconverged;
1821:   PetscFunctionReturn(PETSC_SUCCESS);
1822: }

1824: /*@
1825:   KSPSetInitialGuessKnoll - Tells the iterative solver to use `PCApply()` to compute the initial guess (The Knoll trick)

1827:   Logically Collective

1829:   Input Parameters:
1830: + ksp - iterative context obtained from `KSPCreate()`
1831: - flg - `PETSC_TRUE` or `PETSC_FALSE`

1833:   Level: advanced

1835:   Developer Note:
1836:   The Knoll trick is not currently implemented using the `KSPGuess` class

1838: .seealso: [](ch_ksp), `KSPGetInitialGuessKnoll()`, `KSPSetInitialGuessNonzero()`, `KSPGetInitialGuessNonzero()`, `KSP`
1839: @*/
1840: PetscErrorCode KSPSetInitialGuessKnoll(KSP ksp, PetscBool flg)
1841: {
1842:   PetscFunctionBegin;
1845:   ksp->guess_knoll = flg;
1846:   PetscFunctionReturn(PETSC_SUCCESS);
1847: }

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

1853:   Not Collective

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

1858:   Output Parameter:
1859: . flag - `PETSC_TRUE` if using Knoll trick, else `PETSC_FALSE`

1861:   Level: advanced

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

1874: /*@
1875:   KSPGetComputeSingularValues - Gets the flag indicating whether the extreme singular
1876:   values will be calculated via a Lanczos or Arnoldi process as the linear
1877:   system is solved.

1879:   Not Collective

1881:   Input Parameter:
1882: . ksp - iterative context obtained from `KSPCreate()`

1884:   Output Parameter:
1885: . flg - `PETSC_TRUE` or `PETSC_FALSE`

1887:   Options Database Key:
1888: . -ksp_monitor_singular_value - Activates `KSPSetComputeSingularValues()`

1890:   Level: advanced

1892:   Notes:
1893:   Currently this option is not valid for all iterative methods.

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

1899: .seealso: [](ch_ksp), `KSPComputeExtremeSingularValues()`, `KSPMonitorSingularValue()`, `KSP`
1900: @*/
1901: PetscErrorCode KSPGetComputeSingularValues(KSP ksp, PetscBool *flg)
1902: {
1903:   PetscFunctionBegin;
1905:   PetscAssertPointer(flg, 2);
1906:   *flg = ksp->calc_sings;
1907:   PetscFunctionReturn(PETSC_SUCCESS);
1908: }

1910: /*@
1911:   KSPSetComputeSingularValues - Sets a flag so that the extreme singular
1912:   values will be calculated via a Lanczos or Arnoldi process as the linear
1913:   system is solved.

1915:   Logically Collective

1917:   Input Parameters:
1918: + ksp - iterative context obtained from `KSPCreate()`
1919: - flg - `PETSC_TRUE` or `PETSC_FALSE`

1921:   Options Database Key:
1922: . -ksp_monitor_singular_value - Activates `KSPSetComputeSingularValues()`

1924:   Level: advanced

1926:   Notes:
1927:   Currently this option is not valid for all iterative methods.

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

1933: .seealso: [](ch_ksp), `KSPComputeExtremeSingularValues()`, `KSPMonitorSingularValue()`, `KSP`, `KSPSetComputeRitz()`
1934: @*/
1935: PetscErrorCode KSPSetComputeSingularValues(KSP ksp, PetscBool flg)
1936: {
1937:   PetscFunctionBegin;
1940:   ksp->calc_sings = flg;
1941:   PetscFunctionReturn(PETSC_SUCCESS);
1942: }

1944: /*@
1945:   KSPGetComputeEigenvalues - Gets the flag indicating that the extreme eigenvalues
1946:   values will be calculated via a Lanczos or Arnoldi process as the linear
1947:   system is solved.

1949:   Not Collective

1951:   Input Parameter:
1952: . ksp - iterative context obtained from `KSPCreate()`

1954:   Output Parameter:
1955: . flg - `PETSC_TRUE` or `PETSC_FALSE`

1957:   Level: advanced

1959:   Note:
1960:   Currently this option is not valid for all iterative methods.

1962: .seealso: [](ch_ksp), `KSPComputeEigenvalues()`, `KSPComputeEigenvaluesExplicitly()`, `KSP`, `KSPSetComputeRitz()`
1963: @*/
1964: PetscErrorCode KSPGetComputeEigenvalues(KSP ksp, PetscBool *flg)
1965: {
1966:   PetscFunctionBegin;
1968:   PetscAssertPointer(flg, 2);
1969:   *flg = ksp->calc_sings;
1970:   PetscFunctionReturn(PETSC_SUCCESS);
1971: }

1973: /*@
1974:   KSPSetComputeEigenvalues - Sets a flag so that the extreme eigenvalues
1975:   values will be calculated via a Lanczos or Arnoldi process as the linear
1976:   system is solved.

1978:   Logically Collective

1980:   Input Parameters:
1981: + ksp - iterative context obtained from `KSPCreate()`
1982: - flg - `PETSC_TRUE` or `PETSC_FALSE`

1984:   Level: advanced

1986:   Note:
1987:   Currently this option is not valid for all iterative methods.

1989: .seealso: [](ch_ksp), `KSPComputeEigenvalues()`, `KSPComputeEigenvaluesExplicitly()`, `KSP`, `KSPSetComputeRitz()`
1990: @*/
1991: PetscErrorCode KSPSetComputeEigenvalues(KSP ksp, PetscBool flg)
1992: {
1993:   PetscFunctionBegin;
1996:   ksp->calc_sings = flg;
1997:   PetscFunctionReturn(PETSC_SUCCESS);
1998: }

2000: /*@
2001:   KSPSetComputeRitz - Sets a flag so that the Ritz or harmonic Ritz pairs
2002:   will be calculated via a Lanczos or Arnoldi process as the linear
2003:   system is solved.

2005:   Logically Collective

2007:   Input Parameters:
2008: + ksp - iterative context obtained from `KSPCreate()`
2009: - flg - `PETSC_TRUE` or `PETSC_FALSE`

2011:   Level: advanced

2013:   Note:
2014:   Currently this option is only valid for the `KSPGMRES` method.

2016: .seealso: [](ch_ksp), `KSPComputeRitz()`, `KSP`
2017: @*/
2018: PetscErrorCode KSPSetComputeRitz(KSP ksp, PetscBool flg)
2019: {
2020:   PetscFunctionBegin;
2023:   ksp->calc_ritz = flg;
2024:   PetscFunctionReturn(PETSC_SUCCESS);
2025: }

2027: /*@
2028:   KSPGetRhs - Gets the right-hand-side vector for the linear system to
2029:   be solved.

2031:   Not Collective

2033:   Input Parameter:
2034: . ksp - iterative context obtained from `KSPCreate()`

2036:   Output Parameter:
2037: . r - right-hand-side vector

2039:   Level: developer

2041: .seealso: [](ch_ksp), `KSPGetSolution()`, `KSPSolve()`, `KSP`
2042: @*/
2043: PetscErrorCode KSPGetRhs(KSP ksp, Vec *r)
2044: {
2045:   PetscFunctionBegin;
2047:   PetscAssertPointer(r, 2);
2048:   *r = ksp->vec_rhs;
2049:   PetscFunctionReturn(PETSC_SUCCESS);
2050: }

2052: /*@
2053:   KSPGetSolution - Gets the location of the solution for the
2054:   linear system to be solved.  Note that this may not be where the solution
2055:   is stored during the iterative process; see `KSPBuildSolution()`.

2057:   Not Collective

2059:   Input Parameter:
2060: . ksp - iterative context obtained from `KSPCreate()`

2062:   Output Parameter:
2063: . v - solution vector

2065:   Level: developer

2067: .seealso: [](ch_ksp), `KSPGetRhs()`, `KSPBuildSolution()`, `KSPSolve()`, `KSP`
2068: @*/
2069: PetscErrorCode KSPGetSolution(KSP ksp, Vec *v)
2070: {
2071:   PetscFunctionBegin;
2073:   PetscAssertPointer(v, 2);
2074:   *v = ksp->vec_sol;
2075:   PetscFunctionReturn(PETSC_SUCCESS);
2076: }

2078: /*@
2079:   KSPSetPC - Sets the preconditioner to be used to calculate the
2080:   application of the preconditioner on a vector.

2082:   Collective

2084:   Input Parameters:
2085: + ksp - iterative context obtained from `KSPCreate()`
2086: - pc  - the preconditioner object (can be `NULL`)

2088:   Level: developer

2090:   Note:
2091:   Use `KSPGetPC()` to retrieve the preconditioner context.

2093: .seealso: [](ch_ksp), `KSPGetPC()`, `KSP`
2094: @*/
2095: PetscErrorCode KSPSetPC(KSP ksp, PC pc)
2096: {
2097:   PetscFunctionBegin;
2099:   if (pc) {
2101:     PetscCheckSameComm(ksp, 1, pc, 2);
2102:   }
2103:   if (ksp->pc != pc && ksp->setupstage) ksp->setupstage = KSP_SETUP_NEWMATRIX;
2104:   PetscCall(PetscObjectReference((PetscObject)pc));
2105:   PetscCall(PCDestroy(&ksp->pc));
2106:   ksp->pc = pc;
2107:   PetscFunctionReturn(PETSC_SUCCESS);
2108: }

2110: PETSC_INTERN PetscErrorCode PCCreate_MPI(PC);

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

2116:    Collective

2118:    Input Parameter:
2119: .  ksp - iterative context obtained from `KSPCreate()`

2121:    Level: developer

2123: .seealso: [](ch_ksp), `KSPSetPC()`, `KSP`, `PCMPIServerBegin()`, `PCMPIServerEnd()`
2124: @*/
2125: PETSC_INTERN PetscErrorCode KSPCheckPCMPI(KSP ksp)
2126: {
2127:   PetscBool isPCMPI;

2129:   PetscFunctionBegin;
2131:   PetscCall(PetscObjectTypeCompare((PetscObject)ksp->pc, PCMPI, &isPCMPI));
2132:   if (PCMPIServerActive && ksp->nestlevel == 0 && !isPCMPI) {
2133:     const char *prefix;
2134:     char       *found = NULL;

2136:     PetscCall(KSPGetOptionsPrefix(ksp, &prefix));
2137:     if (prefix) PetscCall(PetscStrstr(prefix, "mpi_linear_solver_server_", &found));
2138:     if (!found) PetscCall(KSPAppendOptionsPrefix(ksp, "mpi_linear_solver_server_"));
2139:     PetscCall(PCSetType(ksp->pc, PCMPI));
2140:   }
2141:   PetscFunctionReturn(PETSC_SUCCESS);
2142: }

2144: /*@
2145:   KSPGetPC - Returns a pointer to the preconditioner context
2146:   set with `KSPSetPC()`.

2148:   Not Collective

2150:   Input Parameter:
2151: . ksp - iterative context obtained from `KSPCreate()`

2153:   Output Parameter:
2154: . pc - preconditioner context

2156:   Level: developer

2158: .seealso: [](ch_ksp), `KSPSetPC()`, `KSP`, `PC`
2159: @*/
2160: PetscErrorCode KSPGetPC(KSP ksp, PC *pc)
2161: {
2162:   PetscFunctionBegin;
2164:   PetscAssertPointer(pc, 2);
2165:   if (!ksp->pc) {
2166:     PetscCall(PCCreate(PetscObjectComm((PetscObject)ksp), &ksp->pc));
2167:     PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp->pc, (PetscObject)ksp, 0));
2168:     PetscCall(PetscObjectSetOptions((PetscObject)ksp->pc, ((PetscObject)ksp)->options));
2169:     PetscCall(PCSetKSPNestLevel(ksp->pc, ksp->nestlevel));
2170:   }
2171:   PetscCall(KSPCheckPCMPI(ksp));
2172:   *pc = ksp->pc;
2173:   PetscFunctionReturn(PETSC_SUCCESS);
2174: }

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

2179:   Collective

2181:   Input Parameters:
2182: + ksp   - iterative context obtained from `KSPCreate()`
2183: . it    - iteration number
2184: - rnorm - relative norm of the residual

2186:   Level: developer

2188:   Note:
2189:   This routine is called by the `KSP` implementations.
2190:   It does not typically need to be called by the user.

2192: .seealso: [](ch_ksp), `KSPMonitorSet()`
2193: @*/
2194: PetscErrorCode KSPMonitor(KSP ksp, PetscInt it, PetscReal rnorm)
2195: {
2196:   PetscInt i, n = ksp->numbermonitors;

2198:   PetscFunctionBegin;
2199:   for (i = 0; i < n; i++) PetscCall((*ksp->monitor[i])(ksp, it, rnorm, ksp->monitorcontext[i]));
2200:   PetscFunctionReturn(PETSC_SUCCESS);
2201: }

2203: /*@C
2204:   KSPMonitorSet - Sets an ADDITIONAL function to be called at every iteration to monitor
2205:   the residual/error etc.

2207:   Logically Collective

2209:   Input Parameters:
2210: + ksp            - iterative context obtained from `KSPCreate()`
2211: . monitor        - pointer to function (if this is `NULL`, it turns off monitoring
2212: . ctx            - [optional] context for private data for the monitor routine (use `NULL` if no context is needed)
2213: - monitordestroy - [optional] routine that frees monitor context (may be `NULL`)

2215:   Calling sequence of `monitor`:
2216: + ksp   - iterative context obtained from `KSPCreate()`
2217: . it    - iteration number
2218: . rnorm - (estimated) 2-norm of (preconditioned) residual
2219: - ctx   - optional monitoring context, as set by `KSPMonitorSet()`

2221:   Calling sequence of `monitordestroy`:
2222: . ctx - optional monitoring context, as set by `KSPMonitorSet()`

2224:   Options Database Keys:
2225: + -ksp_monitor                             - sets `KSPMonitorResidual()`
2226: . -ksp_monitor draw                        - sets `KSPMonitorResidualDraw()` and plots residual
2227: . -ksp_monitor draw::draw_lg               - sets `KSPMonitorResidualDrawLG()` and plots residual
2228: . -ksp_monitor_pause_final                 - Pauses any graphics when the solve finishes (only works for internal monitors)
2229: . -ksp_monitor_true_residual               - sets `KSPMonitorTrueResidual()`
2230: . -ksp_monitor_true_residual draw::draw_lg - sets `KSPMonitorTrueResidualDrawLG()` and plots residual
2231: . -ksp_monitor_max                         - sets `KSPMonitorTrueResidualMax()`
2232: . -ksp_monitor_singular_value              - sets `KSPMonitorSingularValue()`
2233: - -ksp_monitor_cancel                      - cancels all monitors that have
2234:                           been hardwired into a code by
2235:                           calls to `KSPMonitorSet()`, but
2236:                           does not cancel those set via
2237:                           the options database.

2239:   Level: beginner

2241:   Notes:
2242:   The default is to do nothing.  To print the residual, or preconditioned
2243:   residual if `KSPSetNormType`(ksp,`KSP_NORM_PRECONDITIONED`) was called, use
2244:   `KSPMonitorResidual()` as the monitoring routine, with a `PETSCVIEWERASCII` as the
2245:   context.

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

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

2254: .seealso: [](ch_ksp), `KSPMonitorResidual()`, `KSPMonitorCancel()`, `KSP`
2255: @*/
2256: PetscErrorCode KSPMonitorSet(KSP ksp, PetscErrorCode (*monitor)(KSP ksp, PetscInt it, PetscReal rnorm, void *ctx), void *ctx, PetscErrorCode (*monitordestroy)(void **ctx))
2257: {
2258:   PetscInt  i;
2259:   PetscBool identical;

2261:   PetscFunctionBegin;
2263:   for (i = 0; i < ksp->numbermonitors; i++) {
2264:     PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))monitor, ctx, monitordestroy, (PetscErrorCode(*)(void))ksp->monitor[i], ksp->monitorcontext[i], ksp->monitordestroy[i], &identical));
2265:     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
2266:   }
2267:   PetscCheck(ksp->numbermonitors < MAXKSPMONITORS, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Too many KSP monitors set");
2268:   ksp->monitor[ksp->numbermonitors]          = monitor;
2269:   ksp->monitordestroy[ksp->numbermonitors]   = monitordestroy;
2270:   ksp->monitorcontext[ksp->numbermonitors++] = (void *)ctx;
2271:   PetscFunctionReturn(PETSC_SUCCESS);
2272: }

2274: /*@
2275:   KSPMonitorCancel - Clears all monitors for a `KSP` object.

2277:   Logically Collective

2279:   Input Parameter:
2280: . ksp - iterative context obtained from `KSPCreate()`

2282:   Options Database Key:
2283: . -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.

2285:   Level: intermediate

2287: .seealso: [](ch_ksp), `KSPMonitorResidual()`, `KSPMonitorSet()`, `KSP`
2288: @*/
2289: PetscErrorCode KSPMonitorCancel(KSP ksp)
2290: {
2291:   PetscInt i;

2293:   PetscFunctionBegin;
2295:   for (i = 0; i < ksp->numbermonitors; i++) {
2296:     if (ksp->monitordestroy[i]) PetscCall((*ksp->monitordestroy[i])(&ksp->monitorcontext[i]));
2297:   }
2298:   ksp->numbermonitors = 0;
2299:   PetscFunctionReturn(PETSC_SUCCESS);
2300: }

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

2305:   Not Collective

2307:   Input Parameter:
2308: . ksp - iterative context obtained from `KSPCreate()`

2310:   Output Parameter:
2311: . ctx - monitoring context

2313:   Level: intermediate

2315: .seealso: [](ch_ksp), `KSPMonitorResidual()`, `KSP`
2316: @*/
2317: PetscErrorCode KSPGetMonitorContext(KSP ksp, void *ctx)
2318: {
2319:   PetscFunctionBegin;
2321:   *(void **)ctx = ksp->monitorcontext[0];
2322:   PetscFunctionReturn(PETSC_SUCCESS);
2323: }

2325: /*@
2326:   KSPSetResidualHistory - Sets the array used to hold the residual history.
2327:   If set, this array will contain the residual norms computed at each
2328:   iteration of the solver.

2330:   Not Collective

2332:   Input Parameters:
2333: + ksp   - iterative context obtained from `KSPCreate()`
2334: . a     - array to hold history
2335: . na    - size of `a`
2336: - reset - `PETSC_TRUE` indicates the history counter is reset to zero
2337:            for each new linear solve

2339:   Level: advanced

2341:   Notes:
2342:   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.
2343:   If 'a' is `NULL` then space is allocated for the history. If 'na' `PETSC_DECIDE` or `PETSC_DEFAULT` then a
2344:   default array of length 10000 is allocated.

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

2348: .seealso: [](ch_ksp), `KSPGetResidualHistory()`, `KSP`
2349: @*/
2350: PetscErrorCode KSPSetResidualHistory(KSP ksp, PetscReal a[], PetscInt na, PetscBool reset)
2351: {
2352:   PetscFunctionBegin;

2355:   PetscCall(PetscFree(ksp->res_hist_alloc));
2356:   if (na != PETSC_DECIDE && na != PETSC_DEFAULT && a) {
2357:     ksp->res_hist     = a;
2358:     ksp->res_hist_max = (size_t)na;
2359:   } else {
2360:     if (na != PETSC_DECIDE && na != PETSC_DEFAULT) ksp->res_hist_max = (size_t)na;
2361:     else ksp->res_hist_max = 10000; /* like default ksp->max_it */
2362:     PetscCall(PetscCalloc1(ksp->res_hist_max, &ksp->res_hist_alloc));

2364:     ksp->res_hist = ksp->res_hist_alloc;
2365:   }
2366:   ksp->res_hist_len   = 0;
2367:   ksp->res_hist_reset = reset;
2368:   PetscFunctionReturn(PETSC_SUCCESS);
2369: }

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

2374:   Not Collective

2376:   Input Parameter:
2377: . ksp - iterative context obtained from `KSPCreate()`

2379:   Output Parameters:
2380: + a  - pointer to array to hold history (or `NULL`)
2381: - na - number of used entries in a (or `NULL`)

2383:   Level: advanced

2385:   Note:
2386:   This array is borrowed and should not be freed by the caller.

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

2390:   Fortran Note:
2391:   The Fortran version of this routine has a calling sequence
2392: $   call KSPGetResidualHistory(KSP ksp, integer na, integer ierr)
2393:   note that you have passed a Fortran array into `KSPSetResidualHistory()` and you need
2394:   to access the residual values from this Fortran array you provided. Only the `na` (number of
2395:   residual norms currently held) is set.

2397: .seealso: [](ch_ksp), `KSPSetResidualHistory()`, `KSP`
2398: @*/
2399: PetscErrorCode KSPGetResidualHistory(KSP ksp, const PetscReal *a[], PetscInt *na)
2400: {
2401:   PetscFunctionBegin;
2403:   if (a) *a = ksp->res_hist;
2404:   if (na) *na = (PetscInt)ksp->res_hist_len;
2405:   PetscFunctionReturn(PETSC_SUCCESS);
2406: }

2408: /*@
2409:   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.

2411:   Not Collective

2413:   Input Parameters:
2414: + ksp   - iterative context obtained from `KSPCreate()`
2415: . a     - array to hold history
2416: . na    - size of `a`
2417: - reset - `PETSC_TRUE` indicates the history counter is reset to zero for each new linear solve

2419:   Level: advanced

2421:   Notes:
2422:   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.
2423:   If 'a' is `NULL` then space is allocated for the history. If 'na' is `PETSC_DECIDE` or `PETSC_DEFAULT` then a default array of length 10000 is allocated.

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

2427: .seealso: [](ch_ksp), `KSPGetErrorHistory()`, `KSPSetResidualHistory()`, `KSP`
2428: @*/
2429: PetscErrorCode KSPSetErrorHistory(KSP ksp, PetscReal a[], PetscInt na, PetscBool reset)
2430: {
2431:   PetscFunctionBegin;

2434:   PetscCall(PetscFree(ksp->err_hist_alloc));
2435:   if (na != PETSC_DECIDE && na != PETSC_DEFAULT && a) {
2436:     ksp->err_hist     = a;
2437:     ksp->err_hist_max = (size_t)na;
2438:   } else {
2439:     if (na != PETSC_DECIDE && na != PETSC_DEFAULT) ksp->err_hist_max = (size_t)na;
2440:     else ksp->err_hist_max = 10000; /* like default ksp->max_it */
2441:     PetscCall(PetscCalloc1(ksp->err_hist_max, &ksp->err_hist_alloc));

2443:     ksp->err_hist = ksp->err_hist_alloc;
2444:   }
2445:   ksp->err_hist_len   = 0;
2446:   ksp->err_hist_reset = reset;
2447:   PetscFunctionReturn(PETSC_SUCCESS);
2448: }

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

2453:   Not Collective

2455:   Input Parameter:
2456: . ksp - iterative context obtained from `KSPCreate()`

2458:   Output Parameters:
2459: + a  - pointer to array to hold history (or `NULL`)
2460: - na - number of used entries in a (or `NULL`)

2462:   Level: advanced

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

2468:   Fortran Note:
2469:   The Fortran version of this routine has a calling sequence
2470: $   call KSPGetErrorHistory(KSP ksp, integer na, integer ierr)
2471:   note that you have passed a Fortran array into `KSPSetErrorHistory()` and you need
2472:   to access the residual values from this Fortran array you provided. Only the `na` (number of
2473:   residual norms currently held) is set.

2475: .seealso: [](ch_ksp), `KSPSetErrorHistory()`, `KSPGetResidualHistory()`, `KSP`
2476: @*/
2477: PetscErrorCode KSPGetErrorHistory(KSP ksp, const PetscReal *a[], PetscInt *na)
2478: {
2479:   PetscFunctionBegin;
2481:   if (a) *a = ksp->err_hist;
2482:   if (na) *na = (PetscInt)ksp->err_hist_len;
2483:   PetscFunctionReturn(PETSC_SUCCESS);
2484: }

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

2489:   Not Collective

2491:   Input Parameter:
2492: . ksp - The `KSP`

2494:   Output Parameters:
2495: + cr   - The residual contraction rate
2496: . rRsq - The coefficient of determination, $R^2$, indicating the linearity of the data
2497: . ce   - The error contraction rate
2498: - eRsq - The coefficient of determination, $R^2$, indicating the linearity of the data

2500:   Level: advanced

2502:   Note:
2503:   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,
2504:   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)}$,

2506: .seealso: [](ch_ksp), `KSP`, `KSPConvergedRateView()`
2507: @*/
2508: PetscErrorCode KSPComputeConvergenceRate(KSP ksp, PetscReal *cr, PetscReal *rRsq, PetscReal *ce, PetscReal *eRsq)
2509: {
2510:   PetscReal const *hist;
2511:   PetscReal       *x, *y, slope, intercept, mean = 0.0, var = 0.0, res = 0.0;
2512:   PetscInt         n, k;

2514:   PetscFunctionBegin;
2515:   if (cr || rRsq) {
2516:     PetscCall(KSPGetResidualHistory(ksp, &hist, &n));
2517:     if (!n) {
2518:       if (cr) *cr = 0.0;
2519:       if (rRsq) *rRsq = -1.0;
2520:     } else {
2521:       PetscCall(PetscMalloc2(n, &x, n, &y));
2522:       for (k = 0; k < n; ++k) {
2523:         x[k] = k;
2524:         y[k] = PetscLogReal(hist[k]);
2525:         mean += y[k];
2526:       }
2527:       mean /= n;
2528:       PetscCall(PetscLinearRegression(n, x, y, &slope, &intercept));
2529:       for (k = 0; k < n; ++k) {
2530:         res += PetscSqr(y[k] - (slope * x[k] + intercept));
2531:         var += PetscSqr(y[k] - mean);
2532:       }
2533:       PetscCall(PetscFree2(x, y));
2534:       if (cr) *cr = PetscExpReal(slope);
2535:       if (rRsq) *rRsq = var < PETSC_MACHINE_EPSILON ? 0.0 : 1.0 - (res / var);
2536:     }
2537:   }
2538:   if (ce || eRsq) {
2539:     PetscCall(KSPGetErrorHistory(ksp, &hist, &n));
2540:     if (!n) {
2541:       if (ce) *ce = 0.0;
2542:       if (eRsq) *eRsq = -1.0;
2543:     } else {
2544:       PetscCall(PetscMalloc2(n, &x, n, &y));
2545:       for (k = 0; k < n; ++k) {
2546:         x[k] = k;
2547:         y[k] = PetscLogReal(hist[k]);
2548:         mean += y[k];
2549:       }
2550:       mean /= n;
2551:       PetscCall(PetscLinearRegression(n, x, y, &slope, &intercept));
2552:       for (k = 0; k < n; ++k) {
2553:         res += PetscSqr(y[k] - (slope * x[k] + intercept));
2554:         var += PetscSqr(y[k] - mean);
2555:       }
2556:       PetscCall(PetscFree2(x, y));
2557:       if (ce) *ce = PetscExpReal(slope);
2558:       if (eRsq) *eRsq = var < PETSC_MACHINE_EPSILON ? 0.0 : 1.0 - (res / var);
2559:     }
2560:   }
2561:   PetscFunctionReturn(PETSC_SUCCESS);
2562: }

2564: /*@C
2565:   KSPSetConvergenceTest - Sets the function to be used to determine convergence.

2567:   Logically Collective

2569:   Input Parameters:
2570: + ksp      - iterative context obtained from `KSPCreate()`
2571: . converge - pointer to the function
2572: . ctx      - context for private data for the convergence routine (may be `NULL`)
2573: - destroy  - a routine for destroying the context (may be `NULL`)

2575:   Calling sequence of `converge`:
2576: + ksp    - iterative context obtained from `KSPCreate()`
2577: . it     - iteration number
2578: . rnorm  - (estimated) 2-norm of (preconditioned) residual
2579: . reason - the reason why it has converged or diverged
2580: - ctx    - optional convergence context, as set by `KSPSetConvergenceTest()`

2582:   Calling sequence of `destroy`:
2583: . ctx - the context

2585:   Level: advanced

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

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

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

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

2601: .seealso: [](ch_ksp), `KSP`, `KSPConvergedDefault()`, `KSPGetConvergenceContext()`, `KSPSetTolerances()`, `KSPGetConvergenceTest()`, `KSPGetAndClearConvergenceTest()`
2602: @*/
2603: PetscErrorCode KSPSetConvergenceTest(KSP ksp, PetscErrorCode (*converge)(KSP ksp, PetscInt it, PetscReal rnorm, KSPConvergedReason *reason, void *ctx), void *ctx, PetscErrorCode (*destroy)(void *ctx))
2604: {
2605:   PetscFunctionBegin;
2607:   if (ksp->convergeddestroy) PetscCall((*ksp->convergeddestroy)(ksp->cnvP));
2608:   ksp->converged        = converge;
2609:   ksp->convergeddestroy = destroy;
2610:   ksp->cnvP             = (void *)ctx;
2611:   PetscFunctionReturn(PETSC_SUCCESS);
2612: }

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

2617:   Logically Collective

2619:   Input Parameter:
2620: . ksp - iterative context obtained from `KSPCreate()`

2622:   Output Parameters:
2623: + converge - pointer to convergence test function
2624: . ctx      - context for private data for the convergence routine (may be `NULL`)
2625: - destroy  - a routine for destroying the context (may be `NULL`)

2627:   Calling sequence of `converge`:
2628: + ksp    - iterative context obtained from `KSPCreate()`
2629: . it     - iteration number
2630: . rnorm  - (estimated) 2-norm of (preconditioned) residual
2631: . reason - the reason why it has converged or diverged
2632: - ctx    - optional convergence context, as set by `KSPSetConvergenceTest()`

2634:   Calling sequence of `destroy`:
2635: . ctx - the convergence test context

2637:   Level: advanced

2639: .seealso: [](ch_ksp), `KSP`, `KSPConvergedDefault()`, `KSPGetConvergenceContext()`, `KSPSetTolerances()`, `KSPSetConvergenceTest()`, `KSPGetAndClearConvergenceTest()`
2640: @*/
2641: PetscErrorCode KSPGetConvergenceTest(KSP ksp, PetscErrorCode (**converge)(KSP ksp, PetscInt it, PetscReal rnorm, KSPConvergedReason *reason, void *ctx), void **ctx, PetscErrorCode (**destroy)(void *ctx))
2642: {
2643:   PetscFunctionBegin;
2645:   if (converge) *converge = ksp->converged;
2646:   if (destroy) *destroy = ksp->convergeddestroy;
2647:   if (ctx) *ctx = ksp->cnvP;
2648:   PetscFunctionReturn(PETSC_SUCCESS);
2649: }

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

2654:   Logically Collective

2656:   Input Parameter:
2657: . ksp - iterative context obtained from `KSPCreate()`

2659:   Output Parameters:
2660: + converge - pointer to convergence test function
2661: . ctx      - context for private data for the convergence routine
2662: - destroy  - a routine for destroying the context

2664:   Calling sequence of `converge`:
2665: + ksp    - iterative context obtained from `KSPCreate()`
2666: . it     - iteration number
2667: . rnorm  - (estimated) 2-norm of (preconditioned) residual
2668: . reason - the reason why it has converged or diverged
2669: - ctx    - optional convergence context, as set by `KSPSetConvergenceTest()`

2671:   Calling sequence of `destroy`:
2672: . ctx - the convergence test context

2674:   Level: advanced

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

2682: .seealso: [](ch_ksp), `KSP`, `KSPConvergedDefault()`, `KSPGetConvergenceContext()`, `KSPSetTolerances()`, `KSPSetConvergenceTest()`, `KSPGetConvergenceTest()`
2683: @*/
2684: PetscErrorCode KSPGetAndClearConvergenceTest(KSP ksp, PetscErrorCode (**converge)(KSP ksp, PetscInt it, PetscReal rnorm, KSPConvergedReason *reason, void *ctx), void **ctx, PetscErrorCode (**destroy)(void *ctx))
2685: {
2686:   PetscFunctionBegin;
2688:   *converge             = ksp->converged;
2689:   *destroy              = ksp->convergeddestroy;
2690:   *ctx                  = ksp->cnvP;
2691:   ksp->converged        = NULL;
2692:   ksp->cnvP             = NULL;
2693:   ksp->convergeddestroy = NULL;
2694:   PetscFunctionReturn(PETSC_SUCCESS);
2695: }

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

2700:   Not Collective

2702:   Input Parameter:
2703: . ksp - iterative context obtained from `KSPCreate()`

2705:   Output Parameter:
2706: . ctx - monitoring context

2708:   Level: advanced

2710: .seealso: [](ch_ksp), `KSP`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSPGetConvergenceTest()`
2711: @*/
2712: PetscErrorCode KSPGetConvergenceContext(KSP ksp, void *ctx)
2713: {
2714:   PetscFunctionBegin;
2716:   *(void **)ctx = ksp->cnvP;
2717:   PetscFunctionReturn(PETSC_SUCCESS);
2718: }

2720: /*@C
2721:   KSPBuildSolution - Builds the approximate solution in a vector provided.

2723:   Collective

2725:   Input Parameter:
2726: . ksp - iterative context obtained from `KSPCreate()`

2728:   Output Parameter:
2729:    Provide exactly one of
2730: + v - location to stash solution.
2731: - V - the solution is returned in this location. This vector is created
2732:        internally. This vector should NOT be destroyed by the user with
2733:        `VecDestroy()`.

2735:   Level: developer

2737:   Notes:
2738:   This routine can be used in one of two ways
2739: .vb
2740:       KSPBuildSolution(ksp,NULL,&V);
2741:    or
2742:       KSPBuildSolution(ksp,v,NULL); or KSPBuildSolution(ksp,v,&v);
2743: .ve
2744:   In the first case an internal vector is allocated to store the solution
2745:   (the user cannot destroy this vector). In the second case the solution
2746:   is generated in the vector that the user provides. Note that for certain
2747:   methods, such as `KSPCG`, the second case requires a copy of the solution,
2748:   while in the first case the call is essentially free since it simply
2749:   returns the vector where the solution already is stored. For some methods
2750:   like `KSPGMRES` this is a reasonably expensive operation and should only be
2751:   used in truly needed.

2753: .seealso: [](ch_ksp), `KSPGetSolution()`, `KSPBuildResidual()`, `KSP`
2754: @*/
2755: PetscErrorCode KSPBuildSolution(KSP ksp, Vec v, Vec *V)
2756: {
2757:   PetscFunctionBegin;
2759:   PetscCheck(V || v, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_WRONG, "Must provide either v or V");
2760:   if (!V) V = &v;
2761:   PetscUseTypeMethod(ksp, buildsolution, v, V);
2762:   PetscFunctionReturn(PETSC_SUCCESS);
2763: }

2765: /*@C
2766:   KSPBuildResidual - Builds the residual in a vector provided.

2768:   Collective

2770:   Input Parameter:
2771: . ksp - iterative context obtained from `KSPCreate()`

2773:   Output Parameters:
2774: + v - optional location to stash residual.  If `v` is not provided,
2775:        then a location is generated.
2776: . t - work vector.  If not provided then one is generated.
2777: - V - the residual

2779:   Level: advanced

2781:   Note:
2782:   Regardless of whether or not `v` is provided, the residual is
2783:   returned in `V`.

2785: .seealso: [](ch_ksp), `KSP`, `KSPBuildSolution()`
2786: @*/
2787: PetscErrorCode KSPBuildResidual(KSP ksp, Vec t, Vec v, Vec *V)
2788: {
2789:   PetscBool flag = PETSC_FALSE;
2790:   Vec       w = v, tt = t;

2792:   PetscFunctionBegin;
2794:   if (!w) PetscCall(VecDuplicate(ksp->vec_rhs, &w));
2795:   if (!tt) {
2796:     PetscCall(VecDuplicate(ksp->vec_sol, &tt));
2797:     flag = PETSC_TRUE;
2798:   }
2799:   PetscUseTypeMethod(ksp, buildresidual, tt, w, V);
2800:   if (flag) PetscCall(VecDestroy(&tt));
2801:   PetscFunctionReturn(PETSC_SUCCESS);
2802: }

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

2808:   Logically Collective

2810:   Input Parameters:
2811: + ksp   - the `KSP` context
2812: - scale - `PETSC_TRUE` or `PETSC_FALSE`

2814:   Options Database Keys:
2815: + -ksp_diagonal_scale     - perform a diagonal scaling before the solve
2816: - -ksp_diagonal_scale_fix - scale the matrix back AFTER the solve

2818:   Level: advanced

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

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

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

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

2835: .seealso: [](ch_ksp), `KSPGetDiagonalScale()`, `KSPSetDiagonalScaleFix()`, `KSP`
2836: @*/
2837: PetscErrorCode KSPSetDiagonalScale(KSP ksp, PetscBool scale)
2838: {
2839:   PetscFunctionBegin;
2842:   ksp->dscale = scale;
2843:   PetscFunctionReturn(PETSC_SUCCESS);
2844: }

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

2849:   Not Collective

2851:   Input Parameter:
2852: . ksp - the `KSP` context

2854:   Output Parameter:
2855: . scale - `PETSC_TRUE` or `PETSC_FALSE`

2857:   Level: intermediate

2859: .seealso: [](ch_ksp), `KSP`, `KSPSetDiagonalScale()`, `KSPSetDiagonalScaleFix()`
2860: @*/
2861: PetscErrorCode KSPGetDiagonalScale(KSP ksp, PetscBool *scale)
2862: {
2863:   PetscFunctionBegin;
2865:   PetscAssertPointer(scale, 2);
2866:   *scale = ksp->dscale;
2867:   PetscFunctionReturn(PETSC_SUCCESS);
2868: }

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

2873:   Logically Collective

2875:   Input Parameters:
2876: + ksp - the `KSP` context
2877: - fix - `PETSC_TRUE` to scale back after the system solve, `PETSC_FALSE` to not
2878:          rescale (default)

2880:   Level: intermediate

2882:   Notes:
2883:   Must be called after `KSPSetDiagonalScale()`

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

2890: .seealso: [](ch_ksp), `KSPGetDiagonalScale()`, `KSPSetDiagonalScale()`, `KSPGetDiagonalScaleFix()`, `KSP`
2891: @*/
2892: PetscErrorCode KSPSetDiagonalScaleFix(KSP ksp, PetscBool fix)
2893: {
2894:   PetscFunctionBegin;
2897:   ksp->dscalefix = fix;
2898:   PetscFunctionReturn(PETSC_SUCCESS);
2899: }

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

2904:   Not Collective

2906:   Input Parameter:
2907: . ksp - the `KSP` context

2909:   Output Parameter:
2910: . fix - `PETSC_TRUE` to scale back after the system solve, `PETSC_FALSE` to not
2911:          rescale (default)

2913:   Level: intermediate

2915: .seealso: [](ch_ksp), `KSPGetDiagonalScale()`, `KSPSetDiagonalScale()`, `KSPSetDiagonalScaleFix()`, `KSP`
2916: @*/
2917: PetscErrorCode KSPGetDiagonalScaleFix(KSP ksp, PetscBool *fix)
2918: {
2919:   PetscFunctionBegin;
2921:   PetscAssertPointer(fix, 2);
2922:   *fix = ksp->dscalefix;
2923:   PetscFunctionReturn(PETSC_SUCCESS);
2924: }

2926: /*@C
2927:   KSPSetComputeOperators - set routine to compute the linear operators

2929:   Logically Collective

2931:   Input Parameters:
2932: + ksp  - the `KSP` context
2933: . func - function to compute the operators
2934: - ctx  - optional context

2936:   Calling sequence of `func`:
2937: + ksp - the `KSP` context
2938: . A   - the linear operator
2939: . B   - the matrix from which the preconditioner is built, often `A`
2940: - ctx - optional user-provided context

2942:   Level: beginner

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

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

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

2955: .seealso: [](ch_ksp), `KSP`, `KSPSetOperators()`, `KSPSetComputeRHS()`, `DMKSPSetComputeOperators()`, `KSPSetComputeInitialGuess()`
2956: @*/
2957: PetscErrorCode KSPSetComputeOperators(KSP ksp, PetscErrorCode (*func)(KSP ksp, Mat A, Mat B, void *ctx), void *ctx)
2958: {
2959:   DM dm;

2961:   PetscFunctionBegin;
2963:   PetscCall(KSPGetDM(ksp, &dm));
2964:   PetscCall(DMKSPSetComputeOperators(dm, func, ctx));
2965:   if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX;
2966:   PetscFunctionReturn(PETSC_SUCCESS);
2967: }

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

2972:   Logically Collective

2974:   Input Parameters:
2975: + ksp  - the `KSP` context
2976: . func - function to compute the right hand side
2977: - ctx  - optional context

2979:   Calling sequence of `func`:
2980: + ksp - the `KSP` context
2981: . b   - right hand side of linear system
2982: - ctx - optional user-provided context

2984:   Level: beginner

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

2989: .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `DMKSPSetComputeRHS()`, `KSPSetComputeOperators()`, `KSPSetOperators()`
2990: @*/
2991: PetscErrorCode KSPSetComputeRHS(KSP ksp, PetscErrorCode (*func)(KSP ksp, Vec b, void *ctx), void *ctx)
2992: {
2993:   DM dm;

2995:   PetscFunctionBegin;
2997:   PetscCall(KSPGetDM(ksp, &dm));
2998:   PetscCall(DMKSPSetComputeRHS(dm, func, ctx));
2999:   PetscFunctionReturn(PETSC_SUCCESS);
3000: }

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

3005:   Logically Collective

3007:   Input Parameters:
3008: + ksp  - the `KSP` context
3009: . func - function to compute the initial guess
3010: - ctx  - optional context

3012:   Calling sequence of `func`:
3013: + ksp - the `KSP` context
3014: . x   - solution vector
3015: - ctx - optional user-provided context

3017:   Level: beginner

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

3023: .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `KSPSetComputeRHS()`, `KSPSetComputeOperators()`, `DMKSPSetComputeInitialGuess()`, `KSPSetInitialGuessNonzero()`
3024: @*/
3025: PetscErrorCode KSPSetComputeInitialGuess(KSP ksp, PetscErrorCode (*func)(KSP ksp, Vec x, void *ctx), void *ctx)
3026: {
3027:   DM dm;

3029:   PetscFunctionBegin;
3031:   PetscCall(KSPGetDM(ksp, &dm));
3032:   PetscCall(DMKSPSetComputeInitialGuess(dm, func, ctx));
3033:   PetscFunctionReturn(PETSC_SUCCESS);
3034: }

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

3040:   Logically Collective

3042:   Input Parameter:
3043: . ksp - the `KSP` context

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

3048:   Level: advanced

3050: .seealso: [](ch_ksp), `KSPSolveTranspose()`, `KSP`
3051: @*/
3052: PetscErrorCode KSPSetUseExplicitTranspose(KSP ksp, PetscBool flg)
3053: {
3054:   PetscFunctionBegin;
3057:   ksp->transpose.use_explicittranspose = flg;
3058:   PetscFunctionReturn(PETSC_SUCCESS);
3059: }