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: `KSPConvergedReason` view routine (use `NULL` if no context is desired)
498: - reasonviewdestroy - [optional] routine that frees `vctx` (may be `NULL`)
500: Options Database Keys:
501: + -ksp_converged_reason - sets a default `KSPConvergedReasonView()`
502: - -ksp_converged_reason_view_cancel - cancels all converged reason viewers that have been hardwired into a code by
503: calls to `KSPConvergedReasonViewSet()`, but does not cancel those set via the options database.
505: Level: intermediate
507: Note:
508: Several different converged reason view routines may be set by calling
509: `KSPConvergedReasonViewSet()` multiple times; all will be called in the
510: order in which they were set.
512: Developer Note:
513: Should be named KSPConvergedReasonViewAdd().
515: .seealso: [](ch_ksp), `KSPConvergedReasonView()`, `KSPConvergedReasonViewCancel()`
516: @*/
517: PetscErrorCode KSPConvergedReasonViewSet(KSP ksp, PetscErrorCode (*f)(KSP, void *), void *vctx, PetscErrorCode (*reasonviewdestroy)(void **))
518: {
519: PetscInt i;
520: PetscBool identical;
522: PetscFunctionBegin;
524: for (i = 0; i < ksp->numberreasonviews; i++) {
525: PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))f, vctx, reasonviewdestroy, (PetscErrorCode(*)(void))ksp->reasonview[i], ksp->reasonviewcontext[i], ksp->reasonviewdestroy[i], &identical));
526: if (identical) PetscFunctionReturn(PETSC_SUCCESS);
527: }
528: PetscCheck(ksp->numberreasonviews < MAXKSPREASONVIEWS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many KSP reasonview set");
529: ksp->reasonview[ksp->numberreasonviews] = f;
530: ksp->reasonviewdestroy[ksp->numberreasonviews] = reasonviewdestroy;
531: ksp->reasonviewcontext[ksp->numberreasonviews++] = (void *)vctx;
532: PetscFunctionReturn(PETSC_SUCCESS);
533: }
535: /*@
536: KSPConvergedReasonViewCancel - Clears all the reasonview functions for a `KSP` object set with `KSPConvergedReasonViewSet()`.
538: Collective
540: Input Parameter:
541: . ksp - iterative context obtained from `KSPCreate()`
543: Level: intermediate
545: .seealso: [](ch_ksp), `KSPCreate()`, `KSPDestroy()`, `KSPReset()`, `KSPConvergedReasonViewSet()`
546: @*/
547: PetscErrorCode KSPConvergedReasonViewCancel(KSP ksp)
548: {
549: PetscInt i;
551: PetscFunctionBegin;
553: for (i = 0; i < ksp->numberreasonviews; i++) {
554: if (ksp->reasonviewdestroy[i]) PetscCall((*ksp->reasonviewdestroy[i])(&ksp->reasonviewcontext[i]));
555: }
556: ksp->numberreasonviews = 0;
557: PetscFunctionReturn(PETSC_SUCCESS);
558: }
560: /*@
561: KSPConvergedReasonViewFromOptions - Processes command line options to determine if/how a `KSPReason` is to be viewed.
563: Collective
565: Input Parameter:
566: . ksp - the `KSP` object
568: Level: intermediate
570: .seealso: [](ch_ksp), `KSPConvergedReasonView()`, `KSPConvergedReasonViewSet()`
571: @*/
572: PetscErrorCode KSPConvergedReasonViewFromOptions(KSP ksp)
573: {
574: PetscViewer viewer;
575: PetscBool flg;
576: PetscViewerFormat format;
577: PetscInt i;
579: PetscFunctionBegin;
580: /* Call all user-provided reason review routines */
581: for (i = 0; i < ksp->numberreasonviews; i++) PetscCall((*ksp->reasonview[i])(ksp, ksp->reasonviewcontext[i]));
583: /* Call the default PETSc routine */
584: PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp), ((PetscObject)ksp)->options, ((PetscObject)ksp)->prefix, "-ksp_converged_reason", &viewer, &format, &flg));
585: if (flg) {
586: PetscCall(PetscViewerPushFormat(viewer, format));
587: PetscCall(KSPConvergedReasonView(ksp, viewer));
588: PetscCall(PetscViewerPopFormat(viewer));
589: PetscCall(PetscOptionsRestoreViewer(&viewer));
590: }
591: PetscFunctionReturn(PETSC_SUCCESS);
592: }
594: /*@C
595: KSPConvergedRateView - Displays the convergence rate <https://en.wikipedia.org/wiki/Coefficient_of_determination> of `KSPSolve()` to a viewer
597: Collective
599: Input Parameters:
600: + ksp - iterative context obtained from `KSPCreate()`
601: - viewer - the viewer to display the reason
603: Options Database Key:
604: . -ksp_converged_rate - print reason for convergence or divergence and the convergence rate (or 0.0 for divergence)
606: Level: intermediate
608: Notes:
609: To change the format of the output, call `PetscViewerPushFormat`(`viewer`,`format`) before this call.
611: 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,
612: 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)}$,
614: .seealso: [](ch_ksp), `KSPConvergedReasonView()`, `KSPGetConvergedRate()`, `KSPSetTolerances()`, `KSPConvergedDefault()`
615: @*/
616: PetscErrorCode KSPConvergedRateView(KSP ksp, PetscViewer viewer)
617: {
618: PetscViewerFormat format;
619: PetscBool isAscii;
620: PetscReal rrate, rRsq, erate = 0.0, eRsq = 0.0;
621: PetscInt its;
622: const char *prefix, *reason = KSPConvergedReasons[ksp->reason];
624: PetscFunctionBegin;
625: PetscCall(KSPGetOptionsPrefix(ksp, &prefix));
626: PetscCall(KSPGetIterationNumber(ksp, &its));
627: PetscCall(KSPComputeConvergenceRate(ksp, &rrate, &rRsq, &erate, &eRsq));
628: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)ksp));
629: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isAscii));
630: if (isAscii) {
631: PetscCall(PetscViewerGetFormat(viewer, &format));
632: PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ksp)->tablevel));
633: if (ksp->reason > 0) {
634: if (prefix) PetscCall(PetscViewerASCIIPrintf(viewer, "Linear %s solve converged due to %s iterations %" PetscInt_FMT, prefix, reason, its));
635: else PetscCall(PetscViewerASCIIPrintf(viewer, "Linear solve converged due to %s iterations %" PetscInt_FMT, reason, its));
636: PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
637: if (rRsq >= 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, " res rate %g R^2 %g", (double)rrate, (double)rRsq));
638: if (eRsq >= 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, " error rate %g R^2 %g", (double)erate, (double)eRsq));
639: PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
640: PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
641: } else if (ksp->reason <= 0) {
642: if (prefix) PetscCall(PetscViewerASCIIPrintf(viewer, "Linear %s solve did not converge due to %s iterations %" PetscInt_FMT, prefix, reason, its));
643: else PetscCall(PetscViewerASCIIPrintf(viewer, "Linear solve did not converge due to %s iterations %" PetscInt_FMT, reason, its));
644: PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
645: if (rRsq >= 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, " res rate %g R^2 %g", (double)rrate, (double)rRsq));
646: if (eRsq >= 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, " error rate %g R^2 %g", (double)erate, (double)eRsq));
647: PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
648: PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
649: if (ksp->reason == KSP_DIVERGED_PC_FAILED) {
650: PCFailedReason reason;
651: PetscCall(PCGetFailedReason(ksp->pc, &reason));
652: PetscCall(PetscViewerASCIIPrintf(viewer, " PC failed due to %s \n", PCFailedReasons[reason]));
653: }
654: }
655: PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ksp)->tablevel));
656: }
657: PetscFunctionReturn(PETSC_SUCCESS);
658: }
660: #include <petscdraw.h>
662: static PetscErrorCode KSPViewEigenvalues_Internal(KSP ksp, PetscBool isExplicit, PetscViewer viewer, PetscViewerFormat format)
663: {
664: PetscReal *r, *c;
665: PetscInt n, i, neig;
666: PetscBool isascii, isdraw;
667: PetscMPIInt rank;
669: PetscFunctionBegin;
670: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)ksp), &rank));
671: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
672: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
673: if (isExplicit) {
674: PetscCall(VecGetSize(ksp->vec_sol, &n));
675: PetscCall(PetscMalloc2(n, &r, n, &c));
676: PetscCall(KSPComputeEigenvaluesExplicitly(ksp, n, r, c));
677: neig = n;
678: } else {
679: PetscInt nits;
681: PetscCall(KSPGetIterationNumber(ksp, &nits));
682: n = nits + 2;
683: if (!nits) {
684: PetscCall(PetscViewerASCIIPrintf(viewer, "Zero iterations in solver, cannot approximate any eigenvalues\n"));
685: PetscFunctionReturn(PETSC_SUCCESS);
686: }
687: PetscCall(PetscMalloc2(n, &r, n, &c));
688: PetscCall(KSPComputeEigenvalues(ksp, n, r, c, &neig));
689: }
690: if (isascii) {
691: PetscCall(PetscViewerASCIIPrintf(viewer, "%s computed eigenvalues\n", isExplicit ? "Explicitly" : "Iteratively"));
692: for (i = 0; i < neig; ++i) {
693: if (c[i] >= 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, "%g + %gi\n", (double)r[i], (double)c[i]));
694: else PetscCall(PetscViewerASCIIPrintf(viewer, "%g - %gi\n", (double)r[i], -(double)c[i]));
695: }
696: } else if (isdraw && rank == 0) {
697: PetscDraw draw;
698: PetscDrawSP drawsp;
700: if (format == PETSC_VIEWER_DRAW_CONTOUR) {
701: PetscCall(KSPPlotEigenContours_Private(ksp, neig, r, c));
702: } else {
703: PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
704: PetscCall(PetscDrawSPCreate(draw, 1, &drawsp));
705: PetscCall(PetscDrawSPReset(drawsp));
706: for (i = 0; i < neig; ++i) PetscCall(PetscDrawSPAddPoint(drawsp, r + i, c + i));
707: PetscCall(PetscDrawSPDraw(drawsp, PETSC_TRUE));
708: PetscCall(PetscDrawSPSave(drawsp));
709: PetscCall(PetscDrawSPDestroy(&drawsp));
710: }
711: }
712: PetscCall(PetscFree2(r, c));
713: PetscFunctionReturn(PETSC_SUCCESS);
714: }
716: static PetscErrorCode KSPViewSingularvalues_Internal(KSP ksp, PetscViewer viewer, PetscViewerFormat format)
717: {
718: PetscReal smax, smin;
719: PetscInt nits;
720: PetscBool isascii;
722: PetscFunctionBegin;
723: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
724: PetscCall(KSPGetIterationNumber(ksp, &nits));
725: if (!nits) {
726: PetscCall(PetscViewerASCIIPrintf(viewer, "Zero iterations in solver, cannot approximate any singular values\n"));
727: PetscFunctionReturn(PETSC_SUCCESS);
728: }
729: PetscCall(KSPComputeExtremeSingularValues(ksp, &smax, &smin));
730: if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, "Iteratively computed extreme %svalues: max %g min %g max/min %g\n", smin < 0 ? "eigen" : "singular ", (double)smax, (double)smin, (double)(smax / smin)));
731: PetscFunctionReturn(PETSC_SUCCESS);
732: }
734: static PetscErrorCode KSPViewFinalResidual_Internal(KSP ksp, PetscViewer viewer, PetscViewerFormat format)
735: {
736: PetscBool isascii;
738: PetscFunctionBegin;
739: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
740: 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");
741: if (isascii) {
742: Mat A;
743: Vec t;
744: PetscReal norm;
746: PetscCall(PCGetOperators(ksp->pc, &A, NULL));
747: PetscCall(VecDuplicate(ksp->vec_rhs, &t));
748: PetscCall(KSP_MatMult(ksp, A, ksp->vec_sol, t));
749: PetscCall(VecAYPX(t, -1.0, ksp->vec_rhs));
750: PetscCall(VecViewFromOptions(t, (PetscObject)ksp, "-ksp_view_final_residual_vec"));
751: PetscCall(VecNorm(t, NORM_2, &norm));
752: PetscCall(VecDestroy(&t));
753: PetscCall(PetscViewerASCIIPrintf(viewer, "KSP final norm of residual %g\n", (double)norm));
754: }
755: PetscFunctionReturn(PETSC_SUCCESS);
756: }
758: static PetscErrorCode KSPMonitorPauseFinal_Internal(KSP ksp)
759: {
760: PetscInt i;
762: PetscFunctionBegin;
763: if (!ksp->pauseFinal) PetscFunctionReturn(PETSC_SUCCESS);
764: for (i = 0; i < ksp->numbermonitors; ++i) {
765: PetscViewerAndFormat *vf = (PetscViewerAndFormat *)ksp->monitorcontext[i];
766: PetscDraw draw;
767: PetscReal lpause;
769: if (!vf) continue;
770: if (vf->lg) {
771: if (!PetscCheckPointer(vf->lg, PETSC_OBJECT)) continue;
772: if (((PetscObject)vf->lg)->classid != PETSC_DRAWLG_CLASSID) continue;
773: PetscCall(PetscDrawLGGetDraw(vf->lg, &draw));
774: PetscCall(PetscDrawGetPause(draw, &lpause));
775: PetscCall(PetscDrawSetPause(draw, -1.0));
776: PetscCall(PetscDrawPause(draw));
777: PetscCall(PetscDrawSetPause(draw, lpause));
778: } else {
779: PetscBool isdraw;
781: if (!PetscCheckPointer(vf->viewer, PETSC_OBJECT)) continue;
782: if (((PetscObject)vf->viewer)->classid != PETSC_VIEWER_CLASSID) continue;
783: PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERDRAW, &isdraw));
784: if (!isdraw) continue;
785: PetscCall(PetscViewerDrawGetDraw(vf->viewer, 0, &draw));
786: PetscCall(PetscDrawGetPause(draw, &lpause));
787: PetscCall(PetscDrawSetPause(draw, -1.0));
788: PetscCall(PetscDrawPause(draw));
789: PetscCall(PetscDrawSetPause(draw, lpause));
790: }
791: }
792: PetscFunctionReturn(PETSC_SUCCESS);
793: }
795: static PetscErrorCode KSPSolve_Private(KSP ksp, Vec b, Vec x)
796: {
797: PetscBool flg = PETSC_FALSE, inXisinB = PETSC_FALSE, guess_zero;
798: Mat mat, pmat;
799: MPI_Comm comm;
800: MatNullSpace nullsp;
801: Vec btmp, vec_rhs = NULL;
803: PetscFunctionBegin;
804: level++;
805: comm = PetscObjectComm((PetscObject)ksp);
806: if (x && x == b) {
807: PetscCheck(ksp->guess_zero, comm, PETSC_ERR_ARG_INCOMP, "Cannot use x == b with nonzero initial guess");
808: PetscCall(VecDuplicate(b, &x));
809: inXisinB = PETSC_TRUE;
810: }
811: if (b) {
812: PetscCall(PetscObjectReference((PetscObject)b));
813: PetscCall(VecDestroy(&ksp->vec_rhs));
814: ksp->vec_rhs = b;
815: }
816: if (x) {
817: PetscCall(PetscObjectReference((PetscObject)x));
818: PetscCall(VecDestroy(&ksp->vec_sol));
819: ksp->vec_sol = x;
820: }
822: if (ksp->viewPre) PetscCall(ObjectView((PetscObject)ksp, ksp->viewerPre, ksp->formatPre));
824: if (ksp->presolve) PetscCall((*ksp->presolve)(ksp, ksp->vec_rhs, ksp->vec_sol, ksp->prectx));
826: /* reset the residual history list if requested */
827: if (ksp->res_hist_reset) ksp->res_hist_len = 0;
828: if (ksp->err_hist_reset) ksp->err_hist_len = 0;
830: /* KSPSetUp() scales the matrix if needed */
831: PetscCall(KSPSetUp(ksp));
832: PetscCall(KSPSetUpOnBlocks(ksp));
834: if (ksp->guess) {
835: PetscObjectState ostate, state;
837: PetscCall(KSPGuessSetUp(ksp->guess));
838: PetscCall(PetscObjectStateGet((PetscObject)ksp->vec_sol, &ostate));
839: PetscCall(KSPGuessFormGuess(ksp->guess, ksp->vec_rhs, ksp->vec_sol));
840: PetscCall(PetscObjectStateGet((PetscObject)ksp->vec_sol, &state));
841: if (state != ostate) {
842: ksp->guess_zero = PETSC_FALSE;
843: } else {
844: PetscCall(PetscInfo(ksp, "Using zero initial guess since the KSPGuess object did not change the vector\n"));
845: ksp->guess_zero = PETSC_TRUE;
846: }
847: }
849: PetscCall(VecSetErrorIfLocked(ksp->vec_sol, 3));
851: PetscCall(PetscLogEventBegin(!ksp->transpose_solve ? KSP_Solve : KSP_SolveTranspose, ksp, ksp->vec_rhs, ksp->vec_sol, 0));
852: PetscCall(PCGetOperators(ksp->pc, &mat, &pmat));
853: /* diagonal scale RHS if called for */
854: if (ksp->dscale) {
855: PetscCall(VecPointwiseMult(ksp->vec_rhs, ksp->vec_rhs, ksp->diagonal));
856: /* second time in, but matrix was scaled back to original */
857: if (ksp->dscalefix && ksp->dscalefix2) {
858: Mat mat, pmat;
860: PetscCall(PCGetOperators(ksp->pc, &mat, &pmat));
861: PetscCall(MatDiagonalScale(pmat, ksp->diagonal, ksp->diagonal));
862: if (mat != pmat) PetscCall(MatDiagonalScale(mat, ksp->diagonal, ksp->diagonal));
863: }
865: /* scale initial guess */
866: if (!ksp->guess_zero) {
867: if (!ksp->truediagonal) {
868: PetscCall(VecDuplicate(ksp->diagonal, &ksp->truediagonal));
869: PetscCall(VecCopy(ksp->diagonal, ksp->truediagonal));
870: PetscCall(VecReciprocal(ksp->truediagonal));
871: }
872: PetscCall(VecPointwiseMult(ksp->vec_sol, ksp->vec_sol, ksp->truediagonal));
873: }
874: }
875: PetscCall(PCPreSolve(ksp->pc, ksp));
877: if (ksp->guess_zero && !ksp->guess_not_read) PetscCall(VecSet(ksp->vec_sol, 0.0));
878: if (ksp->guess_knoll) { /* The Knoll trick is independent on the KSPGuess specified */
879: PetscCall(PCApply(ksp->pc, ksp->vec_rhs, ksp->vec_sol));
880: PetscCall(KSP_RemoveNullSpace(ksp, ksp->vec_sol));
881: ksp->guess_zero = PETSC_FALSE;
882: }
884: /* can we mark the initial guess as zero for this solve? */
885: guess_zero = ksp->guess_zero;
886: if (!ksp->guess_zero) {
887: PetscReal norm;
889: PetscCall(VecNormAvailable(ksp->vec_sol, NORM_2, &flg, &norm));
890: if (flg && !norm) ksp->guess_zero = PETSC_TRUE;
891: }
892: if (ksp->transpose_solve) {
893: PetscCall(MatGetNullSpace(pmat, &nullsp));
894: } else {
895: PetscCall(MatGetTransposeNullSpace(pmat, &nullsp));
896: }
897: if (nullsp) {
898: PetscCall(VecDuplicate(ksp->vec_rhs, &btmp));
899: PetscCall(VecCopy(ksp->vec_rhs, btmp));
900: PetscCall(MatNullSpaceRemove(nullsp, btmp));
901: vec_rhs = ksp->vec_rhs;
902: ksp->vec_rhs = btmp;
903: }
904: PetscCall(VecLockReadPush(ksp->vec_rhs));
905: PetscUseTypeMethod(ksp, solve);
906: PetscCall(KSPMonitorPauseFinal_Internal(ksp));
908: PetscCall(VecLockReadPop(ksp->vec_rhs));
909: if (nullsp) {
910: ksp->vec_rhs = vec_rhs;
911: PetscCall(VecDestroy(&btmp));
912: }
914: ksp->guess_zero = guess_zero;
916: PetscCheck(ksp->reason, comm, PETSC_ERR_PLIB, "Internal error, solver returned without setting converged reason");
917: ksp->totalits += ksp->its;
919: PetscCall(KSPConvergedReasonViewFromOptions(ksp));
921: if (ksp->viewRate) {
922: PetscCall(PetscViewerPushFormat(ksp->viewerRate, ksp->formatRate));
923: PetscCall(KSPConvergedRateView(ksp, ksp->viewerRate));
924: PetscCall(PetscViewerPopFormat(ksp->viewerRate));
925: }
926: PetscCall(PCPostSolve(ksp->pc, ksp));
928: /* diagonal scale solution if called for */
929: if (ksp->dscale) {
930: PetscCall(VecPointwiseMult(ksp->vec_sol, ksp->vec_sol, ksp->diagonal));
931: /* unscale right hand side and matrix */
932: if (ksp->dscalefix) {
933: Mat mat, pmat;
935: PetscCall(VecReciprocal(ksp->diagonal));
936: PetscCall(VecPointwiseMult(ksp->vec_rhs, ksp->vec_rhs, ksp->diagonal));
937: PetscCall(PCGetOperators(ksp->pc, &mat, &pmat));
938: PetscCall(MatDiagonalScale(pmat, ksp->diagonal, ksp->diagonal));
939: if (mat != pmat) PetscCall(MatDiagonalScale(mat, ksp->diagonal, ksp->diagonal));
940: PetscCall(VecReciprocal(ksp->diagonal));
941: ksp->dscalefix2 = PETSC_TRUE;
942: }
943: }
944: PetscCall(PetscLogEventEnd(!ksp->transpose_solve ? KSP_Solve : KSP_SolveTranspose, ksp, ksp->vec_rhs, ksp->vec_sol, 0));
945: if (ksp->guess) PetscCall(KSPGuessUpdate(ksp->guess, ksp->vec_rhs, ksp->vec_sol));
946: if (ksp->postsolve) PetscCall((*ksp->postsolve)(ksp, ksp->vec_rhs, ksp->vec_sol, ksp->postctx));
948: PetscCall(PCGetOperators(ksp->pc, &mat, &pmat));
949: if (ksp->viewEV) PetscCall(KSPViewEigenvalues_Internal(ksp, PETSC_FALSE, ksp->viewerEV, ksp->formatEV));
950: if (ksp->viewEVExp) PetscCall(KSPViewEigenvalues_Internal(ksp, PETSC_TRUE, ksp->viewerEVExp, ksp->formatEVExp));
951: if (ksp->viewSV) PetscCall(KSPViewSingularvalues_Internal(ksp, ksp->viewerSV, ksp->formatSV));
952: if (ksp->viewFinalRes) PetscCall(KSPViewFinalResidual_Internal(ksp, ksp->viewerFinalRes, ksp->formatFinalRes));
953: if (ksp->viewMat) PetscCall(ObjectView((PetscObject)mat, ksp->viewerMat, ksp->formatMat));
954: if (ksp->viewPMat) PetscCall(ObjectView((PetscObject)pmat, ksp->viewerPMat, ksp->formatPMat));
955: if (ksp->viewRhs) PetscCall(ObjectView((PetscObject)ksp->vec_rhs, ksp->viewerRhs, ksp->formatRhs));
956: if (ksp->viewSol) PetscCall(ObjectView((PetscObject)ksp->vec_sol, ksp->viewerSol, ksp->formatSol));
957: if (ksp->view) PetscCall(ObjectView((PetscObject)ksp, ksp->viewer, ksp->format));
958: if (ksp->viewDScale) PetscCall(ObjectView((PetscObject)ksp->diagonal, ksp->viewerDScale, ksp->formatDScale));
959: if (ksp->viewMatExp) {
960: Mat A, B;
962: PetscCall(PCGetOperators(ksp->pc, &A, NULL));
963: if (ksp->transpose_solve) {
964: Mat AT;
966: PetscCall(MatCreateTranspose(A, &AT));
967: PetscCall(MatComputeOperator(AT, MATAIJ, &B));
968: PetscCall(MatDestroy(&AT));
969: } else {
970: PetscCall(MatComputeOperator(A, MATAIJ, &B));
971: }
972: PetscCall(ObjectView((PetscObject)B, ksp->viewerMatExp, ksp->formatMatExp));
973: PetscCall(MatDestroy(&B));
974: }
975: if (ksp->viewPOpExp) {
976: Mat B;
978: PetscCall(KSPComputeOperator(ksp, MATAIJ, &B));
979: PetscCall(ObjectView((PetscObject)B, ksp->viewerPOpExp, ksp->formatPOpExp));
980: PetscCall(MatDestroy(&B));
981: }
983: if (inXisinB) {
984: PetscCall(VecCopy(x, b));
985: PetscCall(VecDestroy(&x));
986: }
987: PetscCall(PetscObjectSAWsBlock((PetscObject)ksp));
988: if (ksp->errorifnotconverged && ksp->reason < 0 && ((level == 1) || (ksp->reason != KSP_DIVERGED_ITS))) {
989: PCFailedReason reason;
991: 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]);
992: PetscCall(PCGetFailedReason(ksp->pc, &reason));
993: 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]);
994: }
995: level--;
996: PetscFunctionReturn(PETSC_SUCCESS);
997: }
999: /*@
1000: KSPSolve - Solves linear system.
1002: Collective
1004: Input Parameters:
1005: + ksp - iterative context obtained from `KSPCreate()`
1006: . b - the right hand side vector
1007: - x - the solution (this may be the same vector as `b`, then `b` will be overwritten with answer)
1009: Options Database Keys:
1010: + -ksp_view_eigenvalues - compute preconditioned operators eigenvalues
1011: . -ksp_view_eigenvalues_explicit - compute the eigenvalues by forming the dense operator and using LAPACK
1012: . -ksp_view_mat binary - save matrix to the default binary viewer
1013: . -ksp_view_pmat binary - save matrix used to build preconditioner to the default binary viewer
1014: . -ksp_view_rhs binary - save right hand side vector to the default binary viewer
1015: . -ksp_view_solution binary - save computed solution vector to the default binary viewer
1016: (can be read later with src/ksp/tutorials/ex10.c for testing solvers)
1017: . -ksp_view_mat_explicit - for matrix-free operators, computes the matrix entries and views them
1018: . -ksp_view_preconditioned_operator_explicit - computes the product of the preconditioner and matrix as an explicit matrix and views it
1019: . -ksp_converged_reason - print reason for converged or diverged, also prints number of iterations
1020: . -ksp_view_final_residual - print 2-norm of true linear system residual at the end of the solution process
1021: . -ksp_error_if_not_converged - stop the program as soon as an error is detected in a `KSPSolve()`
1022: . -ksp_view_pre - print the ksp data structure before the system solution
1023: - -ksp_view - print the ksp data structure at the end of the system solution
1025: Level: beginner
1027: Notes:
1028: If one uses `KSPSetDM()` then `x` or `b` need not be passed. Use `KSPGetSolution()` to access the solution in this case.
1030: The operator is specified with `KSPSetOperators()`.
1032: `KSPSolve()` will normally return without generating an error regardless of whether the linear system was solved or if constructing the preconditioner failed.
1033: Call `KSPGetConvergedReason()` to determine if the solver converged or failed and why. The option -ksp_error_if_not_converged or function `KSPSetErrorIfNotConverged()`
1034: 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
1035: it may be fine that inner solvers in the preconditioner do not converge during the solution process.
1037: The number of iterations can be obtained from `KSPGetIterationNumber()`.
1039: If you provide a matrix that has a `MatSetNullSpace()` and `MatSetTransposeNullSpace()` this will use that information to solve singular systems
1040: in the least squares sense with a norm minimizing solution.
1042: $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()`).
1044: `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
1045: 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
1046: direction thus the solution which is a linear combination of the search directions has no component in the nullspace(A).
1048: We recommend always using `KSPGMRES` for such singular systems.
1049: If nullspace(A) = nullspace(A') (note symmetric matrices always satisfy this property) then both left and right preconditioning will work
1050: If nullspace(A) != nullspace(A') then left preconditioning will work but right preconditioning may not work (or it may).
1052: Developer Notes:
1053: The reason we cannot always solve nullspace(A) != nullspace(A') systems with right preconditioning is because we need to remove at each iteration
1054: 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
1055: such as diagonal scaling we cannot apply the inverse of the preconditioner to a vector and thus cannot compute the nullspace(AB).
1057: If using a direct method (e.g., via the `KSP` solver
1058: `KSPPREONLY` and a preconditioner such as `PCLU` or `PCILU`,
1059: then its=1. See `KSPSetTolerances()` and `KSPConvergedDefault()`
1060: for more details.
1062: Understanding Convergence\:
1063: The routines `KSPMonitorSet()`, `KSPComputeEigenvalues()`, and
1064: `KSPComputeEigenvaluesExplicitly()` provide information on additional
1065: options to monitor convergence and print eigenvalue information.
1067: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPDestroy()`, `KSPSetTolerances()`, `KSPConvergedDefault()`,
1068: `KSPSolveTranspose()`, `KSPGetIterationNumber()`, `MatNullSpaceCreate()`, `MatSetNullSpace()`, `MatSetTransposeNullSpace()`, `KSP`,
1069: `KSPConvergedReasonView()`, `KSPCheckSolve()`, `KSPSetErrorIfNotConverged()`
1070: @*/
1071: PetscErrorCode KSPSolve(KSP ksp, Vec b, Vec x)
1072: {
1073: PetscFunctionBegin;
1077: ksp->transpose_solve = PETSC_FALSE;
1078: PetscCall(KSPSolve_Private(ksp, b, x));
1079: PetscFunctionReturn(PETSC_SUCCESS);
1080: }
1082: /*@
1083: KSPSolveTranspose - Solves a linear system with the transpose of the matrix, $ A^T x = b$.
1085: Collective
1087: Input Parameters:
1088: + ksp - iterative context obtained from `KSPCreate()`
1089: . b - right hand side vector
1090: - x - solution vector
1092: Level: developer
1094: Note:
1095: For complex numbers this solve the non-Hermitian transpose system.
1097: Developer Note:
1098: We need to implement a `KSPSolveHermitianTranspose()`
1100: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPDestroy()`, `KSPSetTolerances()`, `KSPConvergedDefault()`,
1101: `KSPSolve()`, `KSP`
1102: @*/
1103: PetscErrorCode KSPSolveTranspose(KSP ksp, Vec b, Vec x)
1104: {
1105: PetscFunctionBegin;
1109: if (ksp->transpose.use_explicittranspose) {
1110: Mat J, Jpre;
1111: PetscCall(KSPGetOperators(ksp, &J, &Jpre));
1112: if (!ksp->transpose.reuse_transpose) {
1113: PetscCall(MatTranspose(J, MAT_INITIAL_MATRIX, &ksp->transpose.AT));
1114: if (J != Jpre) PetscCall(MatTranspose(Jpre, MAT_INITIAL_MATRIX, &ksp->transpose.BT));
1115: ksp->transpose.reuse_transpose = PETSC_TRUE;
1116: } else {
1117: PetscCall(MatTranspose(J, MAT_REUSE_MATRIX, &ksp->transpose.AT));
1118: if (J != Jpre) PetscCall(MatTranspose(Jpre, MAT_REUSE_MATRIX, &ksp->transpose.BT));
1119: }
1120: if (J == Jpre && ksp->transpose.BT != ksp->transpose.AT) {
1121: PetscCall(PetscObjectReference((PetscObject)ksp->transpose.AT));
1122: ksp->transpose.BT = ksp->transpose.AT;
1123: }
1124: PetscCall(KSPSetOperators(ksp, ksp->transpose.AT, ksp->transpose.BT));
1125: } else {
1126: ksp->transpose_solve = PETSC_TRUE;
1127: }
1128: PetscCall(KSPSolve_Private(ksp, b, x));
1129: PetscFunctionReturn(PETSC_SUCCESS);
1130: }
1132: static PetscErrorCode KSPViewFinalMatResidual_Internal(KSP ksp, Mat B, Mat X, PetscViewer viewer, PetscViewerFormat format, PetscInt shift)
1133: {
1134: Mat A, R;
1135: PetscReal *norms;
1136: PetscInt i, N;
1137: PetscBool flg;
1139: PetscFunctionBegin;
1140: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &flg));
1141: if (flg) {
1142: PetscCall(PCGetOperators(ksp->pc, &A, NULL));
1143: if (!ksp->transpose_solve) PetscCall(MatMatMult(A, X, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &R));
1144: else PetscCall(MatTransposeMatMult(A, X, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &R));
1145: PetscCall(MatAYPX(R, -1.0, B, SAME_NONZERO_PATTERN));
1146: PetscCall(MatGetSize(R, NULL, &N));
1147: PetscCall(PetscMalloc1(N, &norms));
1148: PetscCall(MatGetColumnNorms(R, NORM_2, norms));
1149: PetscCall(MatDestroy(&R));
1150: 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]));
1151: PetscCall(PetscFree(norms));
1152: }
1153: PetscFunctionReturn(PETSC_SUCCESS);
1154: }
1156: static PetscErrorCode KSPMatSolve_Private(KSP ksp, Mat B, Mat X)
1157: {
1158: Mat A, P, vB, vX;
1159: Vec cb, cx;
1160: PetscInt n1, N1, n2, N2, Bbn = PETSC_DECIDE;
1161: PetscBool match;
1163: PetscFunctionBegin;
1167: PetscCheckSameComm(ksp, 1, B, 2);
1168: PetscCheckSameComm(ksp, 1, X, 3);
1169: PetscCheckSameType(B, 2, X, 3);
1170: PetscCheck(B->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
1171: MatCheckPreallocated(X, 3);
1172: if (!X->assembled) {
1173: PetscCall(MatSetOption(X, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
1174: PetscCall(MatAssemblyBegin(X, MAT_FINAL_ASSEMBLY));
1175: PetscCall(MatAssemblyEnd(X, MAT_FINAL_ASSEMBLY));
1176: }
1177: PetscCheck(B != X, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_IDN, "B and X must be different matrices");
1178: PetscCheck(!ksp->transpose_solve || !ksp->transpose.use_explicittranspose, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "KSPMatSolveTranspose() does not support -ksp_use_explicittranspose");
1179: PetscCall(KSPGetOperators(ksp, &A, &P));
1180: PetscCall(MatGetLocalSize(B, NULL, &n2));
1181: PetscCall(MatGetLocalSize(X, NULL, &n1));
1182: PetscCall(MatGetSize(B, NULL, &N2));
1183: PetscCall(MatGetSize(X, NULL, &N1));
1184: 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);
1185: PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)B, &match, MATSEQDENSE, MATMPIDENSE, ""));
1186: PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of right-hand sides not stored in a dense Mat");
1187: PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)X, &match, MATSEQDENSE, MATMPIDENSE, ""));
1188: PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of solutions not stored in a dense Mat");
1189: PetscCall(KSPSetUp(ksp));
1190: PetscCall(KSPSetUpOnBlocks(ksp));
1191: if (ksp->ops->matsolve) {
1192: level++;
1193: if (ksp->guess_zero) PetscCall(MatZeroEntries(X));
1194: PetscCall(PetscLogEventBegin(!ksp->transpose_solve ? KSP_MatSolve : KSP_MatSolveTranspose, ksp, B, X, 0));
1195: PetscCall(KSPGetMatSolveBatchSize(ksp, &Bbn));
1196: /* by default, do a single solve with all columns */
1197: if (Bbn == PETSC_DECIDE) Bbn = N2;
1198: else PetscCheck(Bbn >= 1, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "KSPMatSolve() batch size %" PetscInt_FMT " must be positive", Bbn);
1199: PetscCall(PetscInfo(ksp, "KSP type %s solving using batches of width at most %" PetscInt_FMT "\n", ((PetscObject)ksp)->type_name, Bbn));
1200: /* if -ksp_matsolve_batch_size is greater than the actual number of columns, do a single solve with all columns */
1201: if (Bbn >= N2) {
1202: PetscUseTypeMethod(ksp, matsolve, B, X);
1203: if (ksp->viewFinalRes) PetscCall(KSPViewFinalMatResidual_Internal(ksp, B, X, ksp->viewerFinalRes, ksp->formatFinalRes, 0));
1205: PetscCall(KSPConvergedReasonViewFromOptions(ksp));
1207: if (ksp->viewRate) {
1208: PetscCall(PetscViewerPushFormat(ksp->viewerRate, PETSC_VIEWER_DEFAULT));
1209: PetscCall(KSPConvergedRateView(ksp, ksp->viewerRate));
1210: PetscCall(PetscViewerPopFormat(ksp->viewerRate));
1211: }
1212: } else {
1213: for (n2 = 0; n2 < N2; n2 += Bbn) {
1214: PetscCall(MatDenseGetSubMatrix(B, PETSC_DECIDE, PETSC_DECIDE, n2, PetscMin(n2 + Bbn, N2), &vB));
1215: PetscCall(MatDenseGetSubMatrix(X, PETSC_DECIDE, PETSC_DECIDE, n2, PetscMin(n2 + Bbn, N2), &vX));
1216: PetscUseTypeMethod(ksp, matsolve, vB, vX);
1217: if (ksp->viewFinalRes) PetscCall(KSPViewFinalMatResidual_Internal(ksp, vB, vX, ksp->viewerFinalRes, ksp->formatFinalRes, n2));
1219: PetscCall(KSPConvergedReasonViewFromOptions(ksp));
1221: if (ksp->viewRate) {
1222: PetscCall(PetscViewerPushFormat(ksp->viewerRate, PETSC_VIEWER_DEFAULT));
1223: PetscCall(KSPConvergedRateView(ksp, ksp->viewerRate));
1224: PetscCall(PetscViewerPopFormat(ksp->viewerRate));
1225: }
1226: PetscCall(MatDenseRestoreSubMatrix(B, &vB));
1227: PetscCall(MatDenseRestoreSubMatrix(X, &vX));
1228: }
1229: }
1230: if (ksp->viewMat) PetscCall(ObjectView((PetscObject)A, ksp->viewerMat, ksp->formatMat));
1231: if (ksp->viewPMat) PetscCall(ObjectView((PetscObject)P, ksp->viewerPMat, ksp->formatPMat));
1232: if (ksp->viewRhs) PetscCall(ObjectView((PetscObject)B, ksp->viewerRhs, ksp->formatRhs));
1233: if (ksp->viewSol) PetscCall(ObjectView((PetscObject)X, ksp->viewerSol, ksp->formatSol));
1234: if (ksp->view) PetscCall(KSPView(ksp, ksp->viewer));
1235: PetscCall(PetscLogEventEnd(!ksp->transpose_solve ? KSP_MatSolve : KSP_MatSolveTranspose, ksp, B, X, 0));
1236: if (ksp->errorifnotconverged && ksp->reason < 0 && (level == 1 || ksp->reason != KSP_DIVERGED_ITS)) {
1237: PCFailedReason reason;
1239: 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]);
1240: PetscCall(PCGetFailedReason(ksp->pc, &reason));
1241: 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]);
1242: }
1243: level--;
1244: } else {
1245: PetscCall(PetscInfo(ksp, "KSP type %s solving column by column\n", ((PetscObject)ksp)->type_name));
1246: for (n2 = 0; n2 < N2; ++n2) {
1247: PetscCall(MatDenseGetColumnVecRead(B, n2, &cb));
1248: PetscCall(MatDenseGetColumnVecWrite(X, n2, &cx));
1249: PetscCall(KSPSolve_Private(ksp, cb, cx));
1250: PetscCall(MatDenseRestoreColumnVecWrite(X, n2, &cx));
1251: PetscCall(MatDenseRestoreColumnVecRead(B, n2, &cb));
1252: }
1253: }
1254: PetscFunctionReturn(PETSC_SUCCESS);
1255: }
1257: /*@
1258: KSPMatSolve - Solves a linear system with multiple right-hand sides stored as a `MATDENSE`. Unlike `KSPSolve()`, `B` and `X` must be different matrices.
1260: Input Parameters:
1261: + ksp - iterative context
1262: - B - block of right-hand sides
1264: Output Parameter:
1265: . X - block of solutions
1267: Level: intermediate
1269: Note:
1270: This is a stripped-down version of `KSPSolve()`, which only handles `-ksp_view`, `-ksp_converged_reason`, `-ksp_converged_rate`, and `-ksp_view_final_residual`.
1272: .seealso: [](ch_ksp), `KSPSolve()`, `MatMatSolve()`, `KSPMatSolveTranspose()`, `MATDENSE`, `KSPHPDDM`, `PCBJACOBI`, `PCASM`
1273: @*/
1274: PetscErrorCode KSPMatSolve(KSP ksp, Mat B, Mat X)
1275: {
1276: PetscFunctionBegin;
1277: ksp->transpose_solve = PETSC_FALSE;
1278: PetscCall(KSPMatSolve_Private(ksp, B, X));
1279: PetscFunctionReturn(PETSC_SUCCESS);
1280: }
1282: /*@
1283: KSPMatSolveTranspose - Solves a linear system with the transposed matrix with multiple right-hand sides stored as a `MATDENSE`. Unlike `KSPSolveTranspose()`,
1284: `B` and `X` must be different matrices and the transposed matrix cannot be assembled explicitly for the user.
1286: Input Parameters:
1287: + ksp - iterative context
1288: - B - block of right-hand sides
1290: Output Parameter:
1291: . X - block of solutions
1293: Level: intermediate
1295: Note:
1296: This is a stripped-down version of `KSPSolveTranspose()`, which only handles `-ksp_view`, `-ksp_converged_reason`, `-ksp_converged_rate`, and `-ksp_view_final_residual`.
1298: .seealso: [](ch_ksp), `KSPSolveTranspose()`, `MatMatTransposeSolve()`, `KSPMatSolve()`, `MATDENSE`, `KSPHPDDM`, `PCBJACOBI`, `PCASM`
1299: @*/
1300: PetscErrorCode KSPMatSolveTranspose(KSP ksp, Mat B, Mat X)
1301: {
1302: PetscFunctionBegin;
1303: ksp->transpose_solve = PETSC_TRUE;
1304: PetscCall(KSPMatSolve_Private(ksp, B, X));
1305: PetscFunctionReturn(PETSC_SUCCESS);
1306: }
1308: /*@
1309: KSPSetMatSolveBatchSize - Sets the maximum number of columns treated simultaneously in `KSPMatSolve()`.
1311: Logically Collective
1313: Input Parameters:
1314: + ksp - iterative context
1315: - bs - batch size
1317: Level: advanced
1319: .seealso: [](ch_ksp), `KSPMatSolve()`, `KSPGetMatSolveBatchSize()`, `-mat_mumps_icntl_27`, `-matmatmult_Bbn`
1320: @*/
1321: PetscErrorCode KSPSetMatSolveBatchSize(KSP ksp, PetscInt bs)
1322: {
1323: PetscFunctionBegin;
1326: ksp->nmax = bs;
1327: PetscFunctionReturn(PETSC_SUCCESS);
1328: }
1330: /*@
1331: KSPGetMatSolveBatchSize - Gets the maximum number of columns treated simultaneously in `KSPMatSolve()`.
1333: Input Parameter:
1334: . ksp - iterative context
1336: Output Parameter:
1337: . bs - batch size
1339: Level: advanced
1341: .seealso: [](ch_ksp), `KSPMatSolve()`, `KSPSetMatSolveBatchSize()`, `-mat_mumps_icntl_27`, `-matmatmult_Bbn`
1342: @*/
1343: PetscErrorCode KSPGetMatSolveBatchSize(KSP ksp, PetscInt *bs)
1344: {
1345: PetscFunctionBegin;
1347: PetscAssertPointer(bs, 2);
1348: *bs = ksp->nmax;
1349: PetscFunctionReturn(PETSC_SUCCESS);
1350: }
1352: /*@
1353: KSPResetViewers - Resets all the viewers set from the options database during `KSPSetFromOptions()`
1355: Collective
1357: Input Parameter:
1358: . ksp - iterative context obtained from `KSPCreate()`
1360: Level: beginner
1362: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPSolve()`, `KSPSetFromOptions()`, `KSP`
1363: @*/
1364: PetscErrorCode KSPResetViewers(KSP ksp)
1365: {
1366: PetscFunctionBegin;
1368: if (!ksp) PetscFunctionReturn(PETSC_SUCCESS);
1369: PetscCall(PetscOptionsRestoreViewer(&ksp->viewer));
1370: PetscCall(PetscOptionsRestoreViewer(&ksp->viewerPre));
1371: PetscCall(PetscOptionsRestoreViewer(&ksp->viewerRate));
1372: PetscCall(PetscOptionsRestoreViewer(&ksp->viewerMat));
1373: PetscCall(PetscOptionsRestoreViewer(&ksp->viewerPMat));
1374: PetscCall(PetscOptionsRestoreViewer(&ksp->viewerRhs));
1375: PetscCall(PetscOptionsRestoreViewer(&ksp->viewerSol));
1376: PetscCall(PetscOptionsRestoreViewer(&ksp->viewerMatExp));
1377: PetscCall(PetscOptionsRestoreViewer(&ksp->viewerEV));
1378: PetscCall(PetscOptionsRestoreViewer(&ksp->viewerSV));
1379: PetscCall(PetscOptionsRestoreViewer(&ksp->viewerEVExp));
1380: PetscCall(PetscOptionsRestoreViewer(&ksp->viewerFinalRes));
1381: PetscCall(PetscOptionsRestoreViewer(&ksp->viewerPOpExp));
1382: PetscCall(PetscOptionsRestoreViewer(&ksp->viewerDScale));
1383: ksp->view = PETSC_FALSE;
1384: ksp->viewPre = PETSC_FALSE;
1385: ksp->viewMat = PETSC_FALSE;
1386: ksp->viewPMat = PETSC_FALSE;
1387: ksp->viewRhs = PETSC_FALSE;
1388: ksp->viewSol = PETSC_FALSE;
1389: ksp->viewMatExp = PETSC_FALSE;
1390: ksp->viewEV = PETSC_FALSE;
1391: ksp->viewSV = PETSC_FALSE;
1392: ksp->viewEVExp = PETSC_FALSE;
1393: ksp->viewFinalRes = PETSC_FALSE;
1394: ksp->viewPOpExp = PETSC_FALSE;
1395: ksp->viewDScale = PETSC_FALSE;
1396: PetscFunctionReturn(PETSC_SUCCESS);
1397: }
1399: /*@
1400: KSPReset - Resets a `KSP` context to the kspsetupcalled = 0 state and removes any allocated Vecs and Mats
1402: Collective
1404: Input Parameter:
1405: . ksp - iterative context obtained from `KSPCreate()`
1407: Level: beginner
1409: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPSolve()`, `KSP`
1410: @*/
1411: PetscErrorCode KSPReset(KSP ksp)
1412: {
1413: PetscFunctionBegin;
1415: if (!ksp) PetscFunctionReturn(PETSC_SUCCESS);
1416: PetscTryTypeMethod(ksp, reset);
1417: if (ksp->pc) PetscCall(PCReset(ksp->pc));
1418: if (ksp->guess) {
1419: KSPGuess guess = ksp->guess;
1420: PetscTryTypeMethod(guess, reset);
1421: }
1422: PetscCall(VecDestroyVecs(ksp->nwork, &ksp->work));
1423: PetscCall(VecDestroy(&ksp->vec_rhs));
1424: PetscCall(VecDestroy(&ksp->vec_sol));
1425: PetscCall(VecDestroy(&ksp->diagonal));
1426: PetscCall(VecDestroy(&ksp->truediagonal));
1428: ksp->setupstage = KSP_SETUP_NEW;
1429: ksp->nmax = PETSC_DECIDE;
1430: PetscFunctionReturn(PETSC_SUCCESS);
1431: }
1433: /*@C
1434: KSPDestroy - Destroys a `KSP` context.
1436: Collective
1438: Input Parameter:
1439: . ksp - iterative context obtained from `KSPCreate()`
1441: Level: beginner
1443: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPSolve()`, `KSP`
1444: @*/
1445: PetscErrorCode KSPDestroy(KSP *ksp)
1446: {
1447: PC pc;
1449: PetscFunctionBegin;
1450: if (!*ksp) PetscFunctionReturn(PETSC_SUCCESS);
1452: if (--((PetscObject)*ksp)->refct > 0) {
1453: *ksp = NULL;
1454: PetscFunctionReturn(PETSC_SUCCESS);
1455: }
1457: PetscCall(PetscObjectSAWsViewOff((PetscObject)*ksp));
1459: /*
1460: Avoid a cascading call to PCReset(ksp->pc) from the following call:
1461: PCReset() shouldn't be called from KSPDestroy() as it is unprotected by pc's
1462: refcount (and may be shared, e.g., by other ksps).
1463: */
1464: pc = (*ksp)->pc;
1465: (*ksp)->pc = NULL;
1466: PetscCall(KSPReset(*ksp));
1467: PetscCall(KSPResetViewers(*ksp));
1468: (*ksp)->pc = pc;
1469: PetscTryTypeMethod(*ksp, destroy);
1471: if ((*ksp)->transpose.use_explicittranspose) {
1472: PetscCall(MatDestroy(&(*ksp)->transpose.AT));
1473: PetscCall(MatDestroy(&(*ksp)->transpose.BT));
1474: (*ksp)->transpose.reuse_transpose = PETSC_FALSE;
1475: }
1477: PetscCall(KSPGuessDestroy(&(*ksp)->guess));
1478: PetscCall(DMDestroy(&(*ksp)->dm));
1479: PetscCall(PCDestroy(&(*ksp)->pc));
1480: PetscCall(PetscFree((*ksp)->res_hist_alloc));
1481: PetscCall(PetscFree((*ksp)->err_hist_alloc));
1482: if ((*ksp)->convergeddestroy) PetscCall((*(*ksp)->convergeddestroy)((*ksp)->cnvP));
1483: PetscCall(KSPMonitorCancel(*ksp));
1484: PetscCall(KSPConvergedReasonViewCancel(*ksp));
1485: PetscCall(PetscHeaderDestroy(ksp));
1486: PetscFunctionReturn(PETSC_SUCCESS);
1487: }
1489: /*@
1490: KSPSetPCSide - Sets the preconditioning side.
1492: Logically Collective
1494: Input Parameter:
1495: . ksp - iterative context obtained from `KSPCreate()`
1497: Output Parameter:
1498: . side - the preconditioning side, where side is one of
1499: .vb
1500: PC_LEFT - left preconditioning (default)
1501: PC_RIGHT - right preconditioning
1502: PC_SYMMETRIC - symmetric preconditioning
1503: .ve
1505: Options Database Key:
1506: . -ksp_pc_side <right,left,symmetric> - `KSP` preconditioner side
1508: Level: intermediate
1510: Notes:
1511: Left preconditioning is used by default for most Krylov methods except `KSPFGMRES` which only supports right preconditioning.
1513: For methods changing the side of the preconditioner changes the norm type that is used, see `KSPSetNormType()`.
1515: Symmetric preconditioning is currently available only for the `KSPQCG` method. However, note that
1516: symmetric preconditioning can be emulated by using either right or left
1517: preconditioning and a pre or post processing step.
1519: Setting the `PCSide` often affects the default norm type. See `KSPSetNormType()` for details.
1521: .seealso: [](ch_ksp), `KSPGetPCSide()`, `KSPSetNormType()`, `KSPGetNormType()`, `KSP`
1522: @*/
1523: PetscErrorCode KSPSetPCSide(KSP ksp, PCSide side)
1524: {
1525: PetscFunctionBegin;
1528: ksp->pc_side = ksp->pc_side_set = side;
1529: PetscFunctionReturn(PETSC_SUCCESS);
1530: }
1532: /*@
1533: KSPGetPCSide - Gets the preconditioning side.
1535: Not Collective
1537: Input Parameter:
1538: . ksp - iterative context obtained from `KSPCreate()`
1540: Output Parameter:
1541: . side - the preconditioning side, where side is one of
1542: .vb
1543: PC_LEFT - left preconditioning (default)
1544: PC_RIGHT - right preconditioning
1545: PC_SYMMETRIC - symmetric preconditioning
1546: .ve
1548: Level: intermediate
1550: .seealso: [](ch_ksp), `KSPSetPCSide()`, `KSP`
1551: @*/
1552: PetscErrorCode KSPGetPCSide(KSP ksp, PCSide *side)
1553: {
1554: PetscFunctionBegin;
1556: PetscAssertPointer(side, 2);
1557: PetscCall(KSPSetUpNorms_Private(ksp, PETSC_TRUE, &ksp->normtype, &ksp->pc_side));
1558: *side = ksp->pc_side;
1559: PetscFunctionReturn(PETSC_SUCCESS);
1560: }
1562: /*@
1563: KSPGetTolerances - Gets the relative, absolute, divergence, and maximum
1564: iteration tolerances used by the default `KSP` convergence tests.
1566: Not Collective
1568: Input Parameter:
1569: . ksp - the Krylov subspace context
1571: Output Parameters:
1572: + rtol - the relative convergence tolerance
1573: . abstol - the absolute convergence tolerance
1574: . dtol - the divergence tolerance
1575: - maxits - maximum number of iterations
1577: Level: intermediate
1579: Note:
1580: The user can specify `NULL` for any parameter that is not needed.
1582: .seealso: [](ch_ksp), `KSPSetTolerances()`, `KSP`, `KSPSetMinimumIterations()`, `KSPGetMinimumIterations()`
1583: @*/
1584: PetscErrorCode KSPGetTolerances(KSP ksp, PetscReal *rtol, PetscReal *abstol, PetscReal *dtol, PetscInt *maxits)
1585: {
1586: PetscFunctionBegin;
1588: if (abstol) *abstol = ksp->abstol;
1589: if (rtol) *rtol = ksp->rtol;
1590: if (dtol) *dtol = ksp->divtol;
1591: if (maxits) *maxits = ksp->max_it;
1592: PetscFunctionReturn(PETSC_SUCCESS);
1593: }
1595: /*@
1596: KSPSetTolerances - Sets the relative, absolute, divergence, and maximum
1597: iteration tolerances used by the default `KSP` convergence testers.
1599: Logically Collective
1601: Input Parameters:
1602: + ksp - the Krylov subspace context
1603: . rtol - the relative convergence tolerance, relative decrease in the (possibly preconditioned) residual norm
1604: . abstol - the absolute convergence tolerance absolute size of the (possibly preconditioned) residual norm
1605: . dtol - the divergence tolerance, amount (possibly preconditioned) residual norm can increase before `KSPConvergedDefault()` concludes that the method is diverging
1606: - maxits - maximum number of iterations to use
1608: Options Database Keys:
1609: + -ksp_atol <abstol> - Sets `abstol`
1610: . -ksp_rtol <rtol> - Sets `rtol`
1611: . -ksp_divtol <dtol> - Sets `dtol`
1612: - -ksp_max_it <maxits> - Sets `maxits`
1614: Level: intermediate
1616: Notes:
1617: Use `PETSC_DEFAULT` to retain the default value of any of the tolerances.
1619: See `KSPConvergedDefault()` for details how these parameters are used in the default convergence test. See also `KSPSetConvergenceTest()`
1620: for setting user-defined stopping criteria.
1622: .seealso: [](ch_ksp), `KSPGetTolerances()`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSP`, `KSPSetMinimumIterations()`
1623: @*/
1624: PetscErrorCode KSPSetTolerances(KSP ksp, PetscReal rtol, PetscReal abstol, PetscReal dtol, PetscInt maxits)
1625: {
1626: PetscFunctionBegin;
1633: if (rtol != (PetscReal)PETSC_DEFAULT) {
1634: 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);
1635: ksp->rtol = rtol;
1636: }
1637: if (abstol != (PetscReal)PETSC_DEFAULT) {
1638: PetscCheck(abstol >= 0.0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %g must be non-negative", (double)abstol);
1639: ksp->abstol = abstol;
1640: }
1641: if (dtol != (PetscReal)PETSC_DEFAULT) {
1642: PetscCheck(dtol >= 0.0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Divergence tolerance %g must be larger than 1.0", (double)dtol);
1643: ksp->divtol = dtol;
1644: }
1645: if (maxits != PETSC_DEFAULT) {
1646: PetscCheck(maxits >= 0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", maxits);
1647: ksp->max_it = maxits;
1648: }
1649: PetscFunctionReturn(PETSC_SUCCESS);
1650: }
1652: /*@
1653: KSPSetMinimumIterations - Sets the minimum number of iterations to use, regardless of the tolerances
1655: Logically Collective
1657: Input Parameters:
1658: + ksp - the Krylov subspace context
1659: - minit - minimum number of iterations to use
1661: Options Database Key:
1662: . -ksp_min_it <minits> - Sets `minit`
1664: Level: intermediate
1666: Notes:
1667: Use `KSPSetTolerances()` to set a variety of other tolerances
1669: See `KSPConvergedDefault()` for details on how these parameters are used in the default convergence test. See also `KSPSetConvergenceTest()`
1670: for setting user-defined stopping criteria.
1672: .seealso: [](ch_ksp), `KSPGetTolerances()`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSP`, `KSPSetTolerances()`, `KSPGetMinimumIterations()`
1673: @*/
1674: PetscErrorCode KSPSetMinimumIterations(KSP ksp, PetscInt minit)
1675: {
1676: PetscFunctionBegin;
1680: PetscCheck(minit >= 0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Minimum number of iterations %" PetscInt_FMT " must be non-negative", minit);
1681: ksp->min_it = minit;
1682: PetscFunctionReturn(PETSC_SUCCESS);
1683: }
1685: /*@
1686: KSPGetMinimumIterations - Gets the minimum number of iterations to use, regardless of the tolerances, that was set with `KSPSetMinimumIterations()` or `-ksp_min_it`
1688: Not Collective
1690: Input Parameter:
1691: . ksp - the Krylov subspace context
1693: Output Parameter:
1694: . minit - minimum number of iterations to use
1696: Level: intermediate
1698: .seealso: [](ch_ksp), `KSPGetTolerances()`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSP`, `KSPSetTolerances()`, `KSPSetMinimumIterations()`
1699: @*/
1700: PetscErrorCode KSPGetMinimumIterations(KSP ksp, PetscInt *minit)
1701: {
1702: PetscFunctionBegin;
1704: PetscAssertPointer(minit, 2);
1706: *minit = ksp->min_it;
1707: PetscFunctionReturn(PETSC_SUCCESS);
1708: }
1710: /*@
1711: KSPSetInitialGuessNonzero - Tells the iterative solver that the
1712: initial guess is nonzero; otherwise `KSP` assumes the initial guess
1713: is to be zero (and thus zeros it out before solving).
1715: Logically Collective
1717: Input Parameters:
1718: + ksp - iterative context obtained from `KSPCreate()`
1719: - flg - ``PETSC_TRUE`` indicates the guess is non-zero, `PETSC_FALSE` indicates the guess is zero
1721: Options Database Key:
1722: . -ksp_initial_guess_nonzero <true,false> - use nonzero initial guess
1724: Level: beginner
1726: Note:
1727: If this is not called the X vector is zeroed in the call to `KSPSolve()`.
1729: .seealso: [](ch_ksp), `KSPGetInitialGuessNonzero()`, `KSPGuessSetType()`, `KSPGuessType`, `KSP`
1730: @*/
1731: PetscErrorCode KSPSetInitialGuessNonzero(KSP ksp, PetscBool flg)
1732: {
1733: PetscFunctionBegin;
1736: ksp->guess_zero = (PetscBool) !(int)flg;
1737: PetscFunctionReturn(PETSC_SUCCESS);
1738: }
1740: /*@
1741: KSPGetInitialGuessNonzero - Determines whether the `KSP` solver is using
1742: a zero initial guess.
1744: Not Collective
1746: Input Parameter:
1747: . ksp - iterative context obtained from `KSPCreate()`
1749: Output Parameter:
1750: . flag - `PETSC_TRUE` if guess is nonzero, else `PETSC_FALSE`
1752: Level: intermediate
1754: .seealso: [](ch_ksp), `KSPSetInitialGuessNonzero()`, `KSP`
1755: @*/
1756: PetscErrorCode KSPGetInitialGuessNonzero(KSP ksp, PetscBool *flag)
1757: {
1758: PetscFunctionBegin;
1760: PetscAssertPointer(flag, 2);
1761: if (ksp->guess_zero) *flag = PETSC_FALSE;
1762: else *flag = PETSC_TRUE;
1763: PetscFunctionReturn(PETSC_SUCCESS);
1764: }
1766: /*@
1767: KSPSetErrorIfNotConverged - Causes `KSPSolve()` to generate an error if the solver has not converged as soon as the error is detected.
1769: Logically Collective
1771: Input Parameters:
1772: + ksp - iterative context obtained from `KSPCreate()`
1773: - flg - `PETSC_TRUE` indicates you want the error generated
1775: Options Database Key:
1776: . -ksp_error_if_not_converged <true,false> - generate an error and stop the program
1778: Level: intermediate
1780: Notes:
1781: Normally PETSc continues if a linear solver fails to converge, you can call `KSPGetConvergedReason()` after a `KSPSolve()`
1782: to determine if it has converged.
1784: A `KSP_DIVERGED_ITS` will not generate an error in a `KSPSolve()` inside a nested linear solver
1786: .seealso: [](ch_ksp), `KSPGetErrorIfNotConverged()`, `KSP`
1787: @*/
1788: PetscErrorCode KSPSetErrorIfNotConverged(KSP ksp, PetscBool flg)
1789: {
1790: PetscFunctionBegin;
1793: ksp->errorifnotconverged = flg;
1794: PetscFunctionReturn(PETSC_SUCCESS);
1795: }
1797: /*@
1798: KSPGetErrorIfNotConverged - Will `KSPSolve()` generate an error if the solver does not converge?
1800: Not Collective
1802: Input Parameter:
1803: . ksp - iterative context obtained from KSPCreate()
1805: Output Parameter:
1806: . flag - `PETSC_TRUE` if it will generate an error, else `PETSC_FALSE`
1808: Level: intermediate
1810: .seealso: [](ch_ksp), `KSPSetErrorIfNotConverged()`, `KSP`
1811: @*/
1812: PetscErrorCode KSPGetErrorIfNotConverged(KSP ksp, PetscBool *flag)
1813: {
1814: PetscFunctionBegin;
1816: PetscAssertPointer(flag, 2);
1817: *flag = ksp->errorifnotconverged;
1818: PetscFunctionReturn(PETSC_SUCCESS);
1819: }
1821: /*@
1822: KSPSetInitialGuessKnoll - Tells the iterative solver to use `PCApply()` to compute the initial guess (The Knoll trick)
1824: Logically Collective
1826: Input Parameters:
1827: + ksp - iterative context obtained from `KSPCreate()`
1828: - flg - `PETSC_TRUE` or `PETSC_FALSE`
1830: Level: advanced
1832: Developer Note:
1833: The Knoll trick is not currently implemented using the `KSPGuess` class
1835: .seealso: [](ch_ksp), `KSPGetInitialGuessKnoll()`, `KSPSetInitialGuessNonzero()`, `KSPGetInitialGuessNonzero()`, `KSP`
1836: @*/
1837: PetscErrorCode KSPSetInitialGuessKnoll(KSP ksp, PetscBool flg)
1838: {
1839: PetscFunctionBegin;
1842: ksp->guess_knoll = flg;
1843: PetscFunctionReturn(PETSC_SUCCESS);
1844: }
1846: /*@
1847: KSPGetInitialGuessKnoll - Determines whether the `KSP` solver is using the Knoll trick (using PCApply(pc,b,...) to compute
1848: the initial guess
1850: Not Collective
1852: Input Parameter:
1853: . ksp - iterative context obtained from `KSPCreate()`
1855: Output Parameter:
1856: . flag - `PETSC_TRUE` if using Knoll trick, else `PETSC_FALSE`
1858: Level: advanced
1860: .seealso: [](ch_ksp), `KSPSetInitialGuessKnoll()`, `KSPSetInitialGuessNonzero()`, `KSPGetInitialGuessNonzero()`, `KSP`
1861: @*/
1862: PetscErrorCode KSPGetInitialGuessKnoll(KSP ksp, PetscBool *flag)
1863: {
1864: PetscFunctionBegin;
1866: PetscAssertPointer(flag, 2);
1867: *flag = ksp->guess_knoll;
1868: PetscFunctionReturn(PETSC_SUCCESS);
1869: }
1871: /*@
1872: KSPGetComputeSingularValues - Gets the flag indicating whether the extreme singular
1873: values will be calculated via a Lanczos or Arnoldi process as the linear
1874: system is solved.
1876: Not Collective
1878: Input Parameter:
1879: . ksp - iterative context obtained from `KSPCreate()`
1881: Output Parameter:
1882: . flg - `PETSC_TRUE` or `PETSC_FALSE`
1884: Options Database Key:
1885: . -ksp_monitor_singular_value - Activates `KSPSetComputeSingularValues()`
1887: Level: advanced
1889: Notes:
1890: Currently this option is not valid for all iterative methods.
1892: Many users may just want to use the monitoring routine
1893: `KSPMonitorSingularValue()` (which can be set with option -ksp_monitor_singular_value)
1894: to print the singular values at each iteration of the linear solve.
1896: .seealso: [](ch_ksp), `KSPComputeExtremeSingularValues()`, `KSPMonitorSingularValue()`, `KSP`
1897: @*/
1898: PetscErrorCode KSPGetComputeSingularValues(KSP ksp, PetscBool *flg)
1899: {
1900: PetscFunctionBegin;
1902: PetscAssertPointer(flg, 2);
1903: *flg = ksp->calc_sings;
1904: PetscFunctionReturn(PETSC_SUCCESS);
1905: }
1907: /*@
1908: KSPSetComputeSingularValues - Sets a flag so that the extreme singular
1909: values will be calculated via a Lanczos or Arnoldi process as the linear
1910: system is solved.
1912: Logically Collective
1914: Input Parameters:
1915: + ksp - iterative context obtained from `KSPCreate()`
1916: - flg - `PETSC_TRUE` or `PETSC_FALSE`
1918: Options Database Key:
1919: . -ksp_monitor_singular_value - Activates `KSPSetComputeSingularValues()`
1921: Level: advanced
1923: Notes:
1924: Currently this option is not valid for all iterative methods.
1926: Many users may just want to use the monitoring routine
1927: `KSPMonitorSingularValue()` (which can be set with option -ksp_monitor_singular_value)
1928: to print the singular values at each iteration of the linear solve.
1930: .seealso: [](ch_ksp), `KSPComputeExtremeSingularValues()`, `KSPMonitorSingularValue()`, `KSP`, `KSPSetComputeRitz()`
1931: @*/
1932: PetscErrorCode KSPSetComputeSingularValues(KSP ksp, PetscBool flg)
1933: {
1934: PetscFunctionBegin;
1937: ksp->calc_sings = flg;
1938: PetscFunctionReturn(PETSC_SUCCESS);
1939: }
1941: /*@
1942: KSPGetComputeEigenvalues - Gets the flag indicating that the extreme eigenvalues
1943: values will be calculated via a Lanczos or Arnoldi process as the linear
1944: system is solved.
1946: Not Collective
1948: Input Parameter:
1949: . ksp - iterative context obtained from `KSPCreate()`
1951: Output Parameter:
1952: . flg - `PETSC_TRUE` or `PETSC_FALSE`
1954: Level: advanced
1956: Note:
1957: Currently this option is not valid for all iterative methods.
1959: .seealso: [](ch_ksp), `KSPComputeEigenvalues()`, `KSPComputeEigenvaluesExplicitly()`, `KSP`, `KSPSetComputeRitz()`
1960: @*/
1961: PetscErrorCode KSPGetComputeEigenvalues(KSP ksp, PetscBool *flg)
1962: {
1963: PetscFunctionBegin;
1965: PetscAssertPointer(flg, 2);
1966: *flg = ksp->calc_sings;
1967: PetscFunctionReturn(PETSC_SUCCESS);
1968: }
1970: /*@
1971: KSPSetComputeEigenvalues - Sets a flag so that the extreme eigenvalues
1972: values will be calculated via a Lanczos or Arnoldi process as the linear
1973: system is solved.
1975: Logically Collective
1977: Input Parameters:
1978: + ksp - iterative context obtained from `KSPCreate()`
1979: - flg - `PETSC_TRUE` or `PETSC_FALSE`
1981: Level: advanced
1983: Note:
1984: Currently this option is not valid for all iterative methods.
1986: .seealso: [](ch_ksp), `KSPComputeEigenvalues()`, `KSPComputeEigenvaluesExplicitly()`, `KSP`, `KSPSetComputeRitz()`
1987: @*/
1988: PetscErrorCode KSPSetComputeEigenvalues(KSP ksp, PetscBool flg)
1989: {
1990: PetscFunctionBegin;
1993: ksp->calc_sings = flg;
1994: PetscFunctionReturn(PETSC_SUCCESS);
1995: }
1997: /*@
1998: KSPSetComputeRitz - Sets a flag so that the Ritz or harmonic Ritz pairs
1999: will be calculated via a Lanczos or Arnoldi process as the linear
2000: system is solved.
2002: Logically Collective
2004: Input Parameters:
2005: + ksp - iterative context obtained from `KSPCreate()`
2006: - flg - `PETSC_TRUE` or `PETSC_FALSE`
2008: Level: advanced
2010: Note:
2011: Currently this option is only valid for the `KSPGMRES` method.
2013: .seealso: [](ch_ksp), `KSPComputeRitz()`, `KSP`
2014: @*/
2015: PetscErrorCode KSPSetComputeRitz(KSP ksp, PetscBool flg)
2016: {
2017: PetscFunctionBegin;
2020: ksp->calc_ritz = flg;
2021: PetscFunctionReturn(PETSC_SUCCESS);
2022: }
2024: /*@
2025: KSPGetRhs - Gets the right-hand-side vector for the linear system to
2026: be solved.
2028: Not Collective
2030: Input Parameter:
2031: . ksp - iterative context obtained from `KSPCreate()`
2033: Output Parameter:
2034: . r - right-hand-side vector
2036: Level: developer
2038: .seealso: [](ch_ksp), `KSPGetSolution()`, `KSPSolve()`, `KSP`
2039: @*/
2040: PetscErrorCode KSPGetRhs(KSP ksp, Vec *r)
2041: {
2042: PetscFunctionBegin;
2044: PetscAssertPointer(r, 2);
2045: *r = ksp->vec_rhs;
2046: PetscFunctionReturn(PETSC_SUCCESS);
2047: }
2049: /*@
2050: KSPGetSolution - Gets the location of the solution for the
2051: linear system to be solved. Note that this may not be where the solution
2052: is stored during the iterative process; see `KSPBuildSolution()`.
2054: Not Collective
2056: Input Parameter:
2057: . ksp - iterative context obtained from `KSPCreate()`
2059: Output Parameter:
2060: . v - solution vector
2062: Level: developer
2064: .seealso: [](ch_ksp), `KSPGetRhs()`, `KSPBuildSolution()`, `KSPSolve()`, `KSP`
2065: @*/
2066: PetscErrorCode KSPGetSolution(KSP ksp, Vec *v)
2067: {
2068: PetscFunctionBegin;
2070: PetscAssertPointer(v, 2);
2071: *v = ksp->vec_sol;
2072: PetscFunctionReturn(PETSC_SUCCESS);
2073: }
2075: /*@
2076: KSPSetPC - Sets the preconditioner to be used to calculate the
2077: application of the preconditioner on a vector.
2079: Collective
2081: Input Parameters:
2082: + ksp - iterative context obtained from `KSPCreate()`
2083: - pc - the preconditioner object (can be `NULL`)
2085: Level: developer
2087: Note:
2088: Use `KSPGetPC()` to retrieve the preconditioner context.
2090: .seealso: [](ch_ksp), `KSPGetPC()`, `KSP`
2091: @*/
2092: PetscErrorCode KSPSetPC(KSP ksp, PC pc)
2093: {
2094: PetscFunctionBegin;
2096: if (pc) {
2098: PetscCheckSameComm(ksp, 1, pc, 2);
2099: }
2100: if (ksp->pc != pc && ksp->setupstage) ksp->setupstage = KSP_SETUP_NEWMATRIX;
2101: PetscCall(PetscObjectReference((PetscObject)pc));
2102: PetscCall(PCDestroy(&ksp->pc));
2103: ksp->pc = pc;
2104: PetscFunctionReturn(PETSC_SUCCESS);
2105: }
2107: PETSC_INTERN PetscErrorCode PCCreate_MPI(PC);
2109: // PetscClangLinter pragma disable: -fdoc-internal-linkage
2110: /*@C
2111: KSPCheckPCMPI - Checks if `-mpi_linear_solver_server` is active and the `PC` should be changed to `PCMPI`
2113: Collective
2115: Input Parameter:
2116: . ksp - iterative context obtained from `KSPCreate()`
2118: Level: developer
2120: .seealso: [](ch_ksp), `KSPSetPC()`, `KSP`, `PCMPIServerBegin()`, `PCMPIServerEnd()`
2121: @*/
2122: PETSC_INTERN PetscErrorCode KSPCheckPCMPI(KSP ksp)
2123: {
2124: PetscBool isPCMPI;
2126: PetscFunctionBegin;
2128: PetscCall(PetscObjectTypeCompare((PetscObject)ksp->pc, PCMPI, &isPCMPI));
2129: if (PCMPIServerActive && ksp->nestlevel == 0 && !isPCMPI) {
2130: const char *prefix;
2131: char *found = NULL;
2133: PetscCall(KSPGetOptionsPrefix(ksp, &prefix));
2134: if (prefix) PetscCall(PetscStrstr(prefix, "mpi_linear_solver_server_", &found));
2135: if (!found) PetscCall(KSPAppendOptionsPrefix(ksp, "mpi_linear_solver_server_"));
2136: PetscCall(PCSetType(ksp->pc, PCMPI));
2137: }
2138: PetscFunctionReturn(PETSC_SUCCESS);
2139: }
2141: /*@
2142: KSPGetPC - Returns a pointer to the preconditioner context
2143: set with `KSPSetPC()`.
2145: Not Collective
2147: Input Parameter:
2148: . ksp - iterative context obtained from `KSPCreate()`
2150: Output Parameter:
2151: . pc - preconditioner context
2153: Level: developer
2155: .seealso: [](ch_ksp), `KSPSetPC()`, `KSP`, `PC`
2156: @*/
2157: PetscErrorCode KSPGetPC(KSP ksp, PC *pc)
2158: {
2159: PetscFunctionBegin;
2161: PetscAssertPointer(pc, 2);
2162: if (!ksp->pc) {
2163: PetscCall(PCCreate(PetscObjectComm((PetscObject)ksp), &ksp->pc));
2164: PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp->pc, (PetscObject)ksp, 0));
2165: PetscCall(PetscObjectSetOptions((PetscObject)ksp->pc, ((PetscObject)ksp)->options));
2166: PetscCall(PCSetKSPNestLevel(ksp->pc, ksp->nestlevel));
2167: }
2168: PetscCall(KSPCheckPCMPI(ksp));
2169: *pc = ksp->pc;
2170: PetscFunctionReturn(PETSC_SUCCESS);
2171: }
2173: /*@
2174: KSPMonitor - runs the user provided monitor routines, if they exist
2176: Collective
2178: Input Parameters:
2179: + ksp - iterative context obtained from `KSPCreate()`
2180: . it - iteration number
2181: - rnorm - relative norm of the residual
2183: Level: developer
2185: Note:
2186: This routine is called by the `KSP` implementations.
2187: It does not typically need to be called by the user.
2189: .seealso: [](ch_ksp), `KSPMonitorSet()`
2190: @*/
2191: PetscErrorCode KSPMonitor(KSP ksp, PetscInt it, PetscReal rnorm)
2192: {
2193: PetscInt i, n = ksp->numbermonitors;
2195: PetscFunctionBegin;
2196: for (i = 0; i < n; i++) PetscCall((*ksp->monitor[i])(ksp, it, rnorm, ksp->monitorcontext[i]));
2197: PetscFunctionReturn(PETSC_SUCCESS);
2198: }
2200: /*@C
2201: KSPMonitorSet - Sets an ADDITIONAL function to be called at every iteration to monitor
2202: the residual/error etc.
2204: Logically Collective
2206: Input Parameters:
2207: + ksp - iterative context obtained from `KSPCreate()`
2208: . monitor - pointer to function (if this is `NULL`, it turns off monitoring
2209: . ctx - [optional] context for private data for the monitor routine (use `NULL` if no context is needed)
2210: - monitordestroy - [optional] routine that frees monitor context (may be `NULL`)
2212: Calling sequence of `monitor`:
2213: + ksp - iterative context obtained from `KSPCreate()`
2214: . it - iteration number
2215: . rnorm - (estimated) 2-norm of (preconditioned) residual
2216: - ctx - optional monitoring context, as set by `KSPMonitorSet()`
2218: Calling sequence of `monitordestroy`:
2219: . ctx - optional monitoring context, as set by `KSPMonitorSet()`
2221: Options Database Keys:
2222: + -ksp_monitor - sets `KSPMonitorResidual()`
2223: . -ksp_monitor draw - sets `KSPMonitorResidualDraw()` and plots residual
2224: . -ksp_monitor draw::draw_lg - sets `KSPMonitorResidualDrawLG()` and plots residual
2225: . -ksp_monitor_pause_final - Pauses any graphics when the solve finishes (only works for internal monitors)
2226: . -ksp_monitor_true_residual - sets `KSPMonitorTrueResidual()`
2227: . -ksp_monitor_true_residual draw::draw_lg - sets `KSPMonitorTrueResidualDrawLG()` and plots residual
2228: . -ksp_monitor_max - sets `KSPMonitorTrueResidualMax()`
2229: . -ksp_monitor_singular_value - sets `KSPMonitorSingularValue()`
2230: - -ksp_monitor_cancel - cancels all monitors that have been hardwired into a code by calls to `KSPMonitorSet()`, but
2231: does not cancel those set via the options database.
2233: Level: beginner
2235: Notes:
2236: The default is to do nothing. To print the residual, or preconditioned
2237: residual if `KSPSetNormType`(ksp,`KSP_NORM_PRECONDITIONED`) was called, use
2238: `KSPMonitorResidual()` as the monitoring routine, with a `PETSCVIEWERASCII` as the
2239: context.
2241: Several different monitoring routines may be set by calling
2242: `KSPMonitorSet()` multiple times; all will be called in the
2243: order in which they were set.
2245: Fortran Note:
2246: Only a single monitor function can be set for each `KSP` object
2248: .seealso: [](ch_ksp), `KSPMonitorResidual()`, `KSPMonitorCancel()`, `KSP`
2249: @*/
2250: PetscErrorCode KSPMonitorSet(KSP ksp, PetscErrorCode (*monitor)(KSP ksp, PetscInt it, PetscReal rnorm, void *ctx), void *ctx, PetscErrorCode (*monitordestroy)(void **ctx))
2251: {
2252: PetscInt i;
2253: PetscBool identical;
2255: PetscFunctionBegin;
2257: for (i = 0; i < ksp->numbermonitors; i++) {
2258: PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))monitor, ctx, monitordestroy, (PetscErrorCode(*)(void))ksp->monitor[i], ksp->monitorcontext[i], ksp->monitordestroy[i], &identical));
2259: if (identical) PetscFunctionReturn(PETSC_SUCCESS);
2260: }
2261: PetscCheck(ksp->numbermonitors < MAXKSPMONITORS, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Too many KSP monitors set");
2262: ksp->monitor[ksp->numbermonitors] = monitor;
2263: ksp->monitordestroy[ksp->numbermonitors] = monitordestroy;
2264: ksp->monitorcontext[ksp->numbermonitors++] = (void *)ctx;
2265: PetscFunctionReturn(PETSC_SUCCESS);
2266: }
2268: /*@
2269: KSPMonitorCancel - Clears all monitors for a `KSP` object.
2271: Logically Collective
2273: Input Parameter:
2274: . ksp - iterative context obtained from `KSPCreate()`
2276: Options Database Key:
2277: . -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.
2279: Level: intermediate
2281: .seealso: [](ch_ksp), `KSPMonitorResidual()`, `KSPMonitorSet()`, `KSP`
2282: @*/
2283: PetscErrorCode KSPMonitorCancel(KSP ksp)
2284: {
2285: PetscInt i;
2287: PetscFunctionBegin;
2289: for (i = 0; i < ksp->numbermonitors; i++) {
2290: if (ksp->monitordestroy[i]) PetscCall((*ksp->monitordestroy[i])(&ksp->monitorcontext[i]));
2291: }
2292: ksp->numbermonitors = 0;
2293: PetscFunctionReturn(PETSC_SUCCESS);
2294: }
2296: /*@C
2297: KSPGetMonitorContext - Gets the monitoring context, as set by `KSPMonitorSet()` for the FIRST monitor only.
2299: Not Collective
2301: Input Parameter:
2302: . ksp - iterative context obtained from `KSPCreate()`
2304: Output Parameter:
2305: . ctx - monitoring context
2307: Level: intermediate
2309: .seealso: [](ch_ksp), `KSPMonitorResidual()`, `KSP`
2310: @*/
2311: PetscErrorCode KSPGetMonitorContext(KSP ksp, void *ctx)
2312: {
2313: PetscFunctionBegin;
2315: *(void **)ctx = ksp->monitorcontext[0];
2316: PetscFunctionReturn(PETSC_SUCCESS);
2317: }
2319: /*@
2320: KSPSetResidualHistory - Sets the array used to hold the residual history.
2321: If set, this array will contain the residual norms computed at each
2322: iteration of the solver.
2324: Not Collective
2326: Input Parameters:
2327: + ksp - iterative context obtained from `KSPCreate()`
2328: . a - array to hold history
2329: . na - size of `a`
2330: - reset - `PETSC_TRUE` indicates the history counter is reset to zero
2331: for each new linear solve
2333: Level: advanced
2335: Notes:
2336: 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.
2337: If 'a' is `NULL` then space is allocated for the history. If 'na' `PETSC_DECIDE` or `PETSC_DEFAULT` then a
2338: default array of length 10000 is allocated.
2340: If the array is not long enough then once the iterations is longer than the array length `KSPSolve()` stops recording the history
2342: .seealso: [](ch_ksp), `KSPGetResidualHistory()`, `KSP`
2343: @*/
2344: PetscErrorCode KSPSetResidualHistory(KSP ksp, PetscReal a[], PetscInt na, PetscBool reset)
2345: {
2346: PetscFunctionBegin;
2349: PetscCall(PetscFree(ksp->res_hist_alloc));
2350: if (na != PETSC_DECIDE && na != PETSC_DEFAULT && a) {
2351: ksp->res_hist = a;
2352: ksp->res_hist_max = (size_t)na;
2353: } else {
2354: if (na != PETSC_DECIDE && na != PETSC_DEFAULT) ksp->res_hist_max = (size_t)na;
2355: else ksp->res_hist_max = 10000; /* like default ksp->max_it */
2356: PetscCall(PetscCalloc1(ksp->res_hist_max, &ksp->res_hist_alloc));
2358: ksp->res_hist = ksp->res_hist_alloc;
2359: }
2360: ksp->res_hist_len = 0;
2361: ksp->res_hist_reset = reset;
2362: PetscFunctionReturn(PETSC_SUCCESS);
2363: }
2365: /*@C
2366: KSPGetResidualHistory - Gets the array used to hold the residual history and the number of residuals it contains.
2368: Not Collective
2370: Input Parameter:
2371: . ksp - iterative context obtained from `KSPCreate()`
2373: Output Parameters:
2374: + a - pointer to array to hold history (or `NULL`)
2375: - na - number of used entries in a (or `NULL`). Note this has different meanings depending on the `reset` argument to `KSPSetResidualHistory()`
2377: Level: advanced
2379: Note:
2380: This array is borrowed and should not be freed by the caller.
2382: Can only be called after a `KSPSetResidualHistory()` otherwise `a` and `na` are set to `NULL` and zero
2384: When `reset` was `PETSC_TRUE` since a residual is computed before the first iteration, the value of `na` is generally one more than the value
2385: returned with `KSPGetIterationNumber()`.
2387: Some Krylov methods may not compute the final residual norm when convergence is declared because the maximum number of iterations allowed has been reached.
2388: In this situation, when `reset` was `PETSC_TRUE`, `na` will then equal the number of iterations reported with `KSPGetIterationNumber()`
2390: Some Krylov methods (such as `KSPSTCG`), under certain circumstances, do not compute the final residual norm. In this situation, when `reset` was `PETSC_TRUE`,
2391: `na` will then equal the number of iterations reported with `KSPGetIterationNumber()`
2393: `KSPBCGSL` does not record the residual norms for the "subiterations" hence the results from `KSPGetResidualHistory()` and `KSPGetIterationNumber()` will be different
2395: Fortran Note:
2396: The Fortran version of this routine has a calling sequence
2397: .vb
2398: call KSPGetResidualHistory(KSP ksp, integer na, integer ierr)
2399: .ve
2400: note that you have passed a Fortran array into `KSPSetResidualHistory()` and you need
2401: to access the residual values from this Fortran array you provided. Only the `na` (number of
2402: residual norms currently held) is set.
2404: .seealso: [](ch_ksp), `KSPSetResidualHistory()`, `KSP`, `KSPGetIterationNumber()`, `KSPSTCG`, `KSPBCGSL`
2405: @*/
2406: PetscErrorCode KSPGetResidualHistory(KSP ksp, const PetscReal *a[], PetscInt *na)
2407: {
2408: PetscFunctionBegin;
2410: if (a) *a = ksp->res_hist;
2411: if (na) *na = (PetscInt)ksp->res_hist_len;
2412: PetscFunctionReturn(PETSC_SUCCESS);
2413: }
2415: /*@
2416: 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.
2418: Not Collective
2420: Input Parameters:
2421: + ksp - iterative context obtained from `KSPCreate()`
2422: . a - array to hold history
2423: . na - size of `a`
2424: - reset - `PETSC_TRUE` indicates the history counter is reset to zero for each new linear solve
2426: Level: advanced
2428: Notes:
2429: 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.
2430: 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.
2432: If the array is not long enough then once the iterations is longer than the array length `KSPSolve()` stops recording the history
2434: .seealso: [](ch_ksp), `KSPGetErrorHistory()`, `KSPSetResidualHistory()`, `KSP`
2435: @*/
2436: PetscErrorCode KSPSetErrorHistory(KSP ksp, PetscReal a[], PetscInt na, PetscBool reset)
2437: {
2438: PetscFunctionBegin;
2441: PetscCall(PetscFree(ksp->err_hist_alloc));
2442: if (na != PETSC_DECIDE && na != PETSC_DEFAULT && a) {
2443: ksp->err_hist = a;
2444: ksp->err_hist_max = (size_t)na;
2445: } else {
2446: if (na != PETSC_DECIDE && na != PETSC_DEFAULT) ksp->err_hist_max = (size_t)na;
2447: else ksp->err_hist_max = 10000; /* like default ksp->max_it */
2448: PetscCall(PetscCalloc1(ksp->err_hist_max, &ksp->err_hist_alloc));
2450: ksp->err_hist = ksp->err_hist_alloc;
2451: }
2452: ksp->err_hist_len = 0;
2453: ksp->err_hist_reset = reset;
2454: PetscFunctionReturn(PETSC_SUCCESS);
2455: }
2457: /*@C
2458: KSPGetErrorHistory - Gets the array used to hold the error history and the number of residuals it contains.
2460: Not Collective
2462: Input Parameter:
2463: . ksp - iterative context obtained from `KSPCreate()`
2465: Output Parameters:
2466: + a - pointer to array to hold history (or `NULL`)
2467: - na - number of used entries in a (or `NULL`)
2469: Level: advanced
2471: Note:
2472: This array is borrowed and should not be freed by the caller.
2473: Can only be called after a `KSPSetErrorHistory()` otherwise `a` and `na` are set to `NULL` and zero
2475: Fortran Note:
2476: The Fortran version of this routine has a calling sequence
2477: .vb
2478: call KSPGetErrorHistory(KSP ksp, integer na, integer ierr)
2479: .ve
2480: note that you have passed a Fortran array into `KSPSetErrorHistory()` and you need
2481: to access the residual values from this Fortran array you provided. Only the `na` (number of
2482: residual norms currently held) is set.
2484: .seealso: [](ch_ksp), `KSPSetErrorHistory()`, `KSPGetResidualHistory()`, `KSP`
2485: @*/
2486: PetscErrorCode KSPGetErrorHistory(KSP ksp, const PetscReal *a[], PetscInt *na)
2487: {
2488: PetscFunctionBegin;
2490: if (a) *a = ksp->err_hist;
2491: if (na) *na = (PetscInt)ksp->err_hist_len;
2492: PetscFunctionReturn(PETSC_SUCCESS);
2493: }
2495: /*@
2496: KSPComputeConvergenceRate - Compute the convergence rate for the iteration <https:/en.wikipedia.org/wiki/Coefficient_of_determination>
2498: Not Collective
2500: Input Parameter:
2501: . ksp - The `KSP`
2503: Output Parameters:
2504: + cr - The residual contraction rate
2505: . rRsq - The coefficient of determination, $R^2$, indicating the linearity of the data
2506: . ce - The error contraction rate
2507: - eRsq - The coefficient of determination, $R^2$, indicating the linearity of the data
2509: Level: advanced
2511: Note:
2512: 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,
2513: 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)}$,
2515: .seealso: [](ch_ksp), `KSP`, `KSPConvergedRateView()`
2516: @*/
2517: PetscErrorCode KSPComputeConvergenceRate(KSP ksp, PetscReal *cr, PetscReal *rRsq, PetscReal *ce, PetscReal *eRsq)
2518: {
2519: PetscReal const *hist;
2520: PetscReal *x, *y, slope, intercept, mean = 0.0, var = 0.0, res = 0.0;
2521: PetscInt n, k;
2523: PetscFunctionBegin;
2524: if (cr || rRsq) {
2525: PetscCall(KSPGetResidualHistory(ksp, &hist, &n));
2526: if (!n) {
2527: if (cr) *cr = 0.0;
2528: if (rRsq) *rRsq = -1.0;
2529: } else {
2530: PetscCall(PetscMalloc2(n, &x, n, &y));
2531: for (k = 0; k < n; ++k) {
2532: x[k] = k;
2533: y[k] = PetscLogReal(hist[k]);
2534: mean += y[k];
2535: }
2536: mean /= n;
2537: PetscCall(PetscLinearRegression(n, x, y, &slope, &intercept));
2538: for (k = 0; k < n; ++k) {
2539: res += PetscSqr(y[k] - (slope * x[k] + intercept));
2540: var += PetscSqr(y[k] - mean);
2541: }
2542: PetscCall(PetscFree2(x, y));
2543: if (cr) *cr = PetscExpReal(slope);
2544: if (rRsq) *rRsq = var < PETSC_MACHINE_EPSILON ? 0.0 : 1.0 - (res / var);
2545: }
2546: }
2547: if (ce || eRsq) {
2548: PetscCall(KSPGetErrorHistory(ksp, &hist, &n));
2549: if (!n) {
2550: if (ce) *ce = 0.0;
2551: if (eRsq) *eRsq = -1.0;
2552: } else {
2553: PetscCall(PetscMalloc2(n, &x, n, &y));
2554: for (k = 0; k < n; ++k) {
2555: x[k] = k;
2556: y[k] = PetscLogReal(hist[k]);
2557: mean += y[k];
2558: }
2559: mean /= n;
2560: PetscCall(PetscLinearRegression(n, x, y, &slope, &intercept));
2561: for (k = 0; k < n; ++k) {
2562: res += PetscSqr(y[k] - (slope * x[k] + intercept));
2563: var += PetscSqr(y[k] - mean);
2564: }
2565: PetscCall(PetscFree2(x, y));
2566: if (ce) *ce = PetscExpReal(slope);
2567: if (eRsq) *eRsq = var < PETSC_MACHINE_EPSILON ? 0.0 : 1.0 - (res / var);
2568: }
2569: }
2570: PetscFunctionReturn(PETSC_SUCCESS);
2571: }
2573: /*@C
2574: KSPSetConvergenceTest - Sets the function to be used to determine convergence.
2576: Logically Collective
2578: Input Parameters:
2579: + ksp - iterative context obtained from `KSPCreate()`
2580: . converge - pointer to the function
2581: . ctx - context for private data for the convergence routine (may be `NULL`)
2582: - destroy - a routine for destroying the context (may be `NULL`)
2584: Calling sequence of `converge`:
2585: + ksp - iterative context obtained from `KSPCreate()`
2586: . it - iteration number
2587: . rnorm - (estimated) 2-norm of (preconditioned) residual
2588: . reason - the reason why it has converged or diverged
2589: - ctx - optional convergence context, as set by `KSPSetConvergenceTest()`
2591: Calling sequence of `destroy`:
2592: . ctx - the context
2594: Level: advanced
2596: Notes:
2597: Must be called after the `KSP` type has been set so put this after
2598: a call to `KSPSetType()`, or `KSPSetFromOptions()`.
2600: The default convergence test, `KSPConvergedDefault()`, aborts if the
2601: residual grows to more than 10000 times the initial residual.
2603: The default is a combination of relative and absolute tolerances.
2604: The residual value that is tested may be an approximation; routines
2605: that need exact values should compute them.
2607: In the default PETSc convergence test, the precise values of reason
2608: are macros such as `KSP_CONVERGED_RTOL`, which are defined in petscksp.h.
2610: .seealso: [](ch_ksp), `KSP`, `KSPConvergedDefault()`, `KSPGetConvergenceContext()`, `KSPSetTolerances()`, `KSPGetConvergenceTest()`, `KSPGetAndClearConvergenceTest()`
2611: @*/
2612: PetscErrorCode KSPSetConvergenceTest(KSP ksp, PetscErrorCode (*converge)(KSP ksp, PetscInt it, PetscReal rnorm, KSPConvergedReason *reason, void *ctx), void *ctx, PetscErrorCode (*destroy)(void *ctx))
2613: {
2614: PetscFunctionBegin;
2616: if (ksp->convergeddestroy) PetscCall((*ksp->convergeddestroy)(ksp->cnvP));
2617: ksp->converged = converge;
2618: ksp->convergeddestroy = destroy;
2619: ksp->cnvP = (void *)ctx;
2620: PetscFunctionReturn(PETSC_SUCCESS);
2621: }
2623: /*@C
2624: KSPGetConvergenceTest - Gets the function to be used to determine convergence.
2626: Logically Collective
2628: Input Parameter:
2629: . ksp - iterative context obtained from `KSPCreate()`
2631: Output Parameters:
2632: + converge - pointer to convergence test function
2633: . ctx - context for private data for the convergence routine (may be `NULL`)
2634: - destroy - a routine for destroying the context (may be `NULL`)
2636: Calling sequence of `converge`:
2637: + ksp - iterative context obtained from `KSPCreate()`
2638: . it - iteration number
2639: . rnorm - (estimated) 2-norm of (preconditioned) residual
2640: . reason - the reason why it has converged or diverged
2641: - ctx - optional convergence context, as set by `KSPSetConvergenceTest()`
2643: Calling sequence of `destroy`:
2644: . ctx - the convergence test context
2646: Level: advanced
2648: .seealso: [](ch_ksp), `KSP`, `KSPConvergedDefault()`, `KSPGetConvergenceContext()`, `KSPSetTolerances()`, `KSPSetConvergenceTest()`, `KSPGetAndClearConvergenceTest()`
2649: @*/
2650: PetscErrorCode KSPGetConvergenceTest(KSP ksp, PetscErrorCode (**converge)(KSP ksp, PetscInt it, PetscReal rnorm, KSPConvergedReason *reason, void *ctx), void **ctx, PetscErrorCode (**destroy)(void *ctx))
2651: {
2652: PetscFunctionBegin;
2654: if (converge) *converge = ksp->converged;
2655: if (destroy) *destroy = ksp->convergeddestroy;
2656: if (ctx) *ctx = ksp->cnvP;
2657: PetscFunctionReturn(PETSC_SUCCESS);
2658: }
2660: /*@C
2661: KSPGetAndClearConvergenceTest - Gets the function to be used to determine convergence. Removes the current test without calling destroy on the test context
2663: Logically Collective
2665: Input Parameter:
2666: . ksp - iterative context obtained from `KSPCreate()`
2668: Output Parameters:
2669: + converge - pointer to convergence test function
2670: . ctx - context for private data for the convergence routine
2671: - destroy - a routine for destroying the context
2673: Calling sequence of `converge`:
2674: + ksp - iterative context obtained from `KSPCreate()`
2675: . it - iteration number
2676: . rnorm - (estimated) 2-norm of (preconditioned) residual
2677: . reason - the reason why it has converged or diverged
2678: - ctx - optional convergence context, as set by `KSPSetConvergenceTest()`
2680: Calling sequence of `destroy`:
2681: . ctx - the convergence test context
2683: Level: advanced
2685: Note:
2686: This is intended to be used to allow transferring the convergence test (and its context) to another testing object (for example another `KSP`)
2687: and then calling `KSPSetConvergenceTest()` on this original `KSP`. If you just called `KSPGetConvergenceTest()` followed
2688: by `KSPSetConvergenceTest()` the original context information
2689: would be destroyed and hence the transferred context would be invalid and trigger a crash on use
2691: .seealso: [](ch_ksp), `KSP`, `KSPConvergedDefault()`, `KSPGetConvergenceContext()`, `KSPSetTolerances()`, `KSPSetConvergenceTest()`, `KSPGetConvergenceTest()`
2692: @*/
2693: PetscErrorCode KSPGetAndClearConvergenceTest(KSP ksp, PetscErrorCode (**converge)(KSP ksp, PetscInt it, PetscReal rnorm, KSPConvergedReason *reason, void *ctx), void **ctx, PetscErrorCode (**destroy)(void *ctx))
2694: {
2695: PetscFunctionBegin;
2697: *converge = ksp->converged;
2698: *destroy = ksp->convergeddestroy;
2699: *ctx = ksp->cnvP;
2700: ksp->converged = NULL;
2701: ksp->cnvP = NULL;
2702: ksp->convergeddestroy = NULL;
2703: PetscFunctionReturn(PETSC_SUCCESS);
2704: }
2706: /*@C
2707: KSPGetConvergenceContext - Gets the convergence context set with `KSPSetConvergenceTest()`.
2709: Not Collective
2711: Input Parameter:
2712: . ksp - iterative context obtained from `KSPCreate()`
2714: Output Parameter:
2715: . ctx - monitoring context
2717: Level: advanced
2719: .seealso: [](ch_ksp), `KSP`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSPGetConvergenceTest()`
2720: @*/
2721: PetscErrorCode KSPGetConvergenceContext(KSP ksp, void *ctx)
2722: {
2723: PetscFunctionBegin;
2725: *(void **)ctx = ksp->cnvP;
2726: PetscFunctionReturn(PETSC_SUCCESS);
2727: }
2729: /*@C
2730: KSPBuildSolution - Builds the approximate solution in a vector provided.
2732: Collective
2734: Input Parameter:
2735: . ksp - iterative context obtained from `KSPCreate()`
2737: Output Parameter:
2738: Provide exactly one of
2739: + v - location to stash solution.
2740: - V - the solution is returned in this location. This vector is created internally. This vector should NOT be destroyed by the user with `VecDestroy()`.
2742: Level: developer
2744: Notes:
2745: This routine can be used in one of two ways
2746: .vb
2747: KSPBuildSolution(ksp,NULL,&V);
2748: or
2749: KSPBuildSolution(ksp,v,NULL); or KSPBuildSolution(ksp,v,&v);
2750: .ve
2751: In the first case an internal vector is allocated to store the solution
2752: (the user cannot destroy this vector). In the second case the solution
2753: is generated in the vector that the user provides. Note that for certain
2754: methods, such as `KSPCG`, the second case requires a copy of the solution,
2755: while in the first case the call is essentially free since it simply
2756: returns the vector where the solution already is stored. For some methods
2757: like `KSPGMRES` this is a reasonably expensive operation and should only be
2758: used in truly needed.
2760: .seealso: [](ch_ksp), `KSPGetSolution()`, `KSPBuildResidual()`, `KSP`
2761: @*/
2762: PetscErrorCode KSPBuildSolution(KSP ksp, Vec v, Vec *V)
2763: {
2764: PetscFunctionBegin;
2766: PetscCheck(V || v, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_WRONG, "Must provide either v or V");
2767: if (!V) V = &v;
2768: PetscUseTypeMethod(ksp, buildsolution, v, V);
2769: PetscFunctionReturn(PETSC_SUCCESS);
2770: }
2772: /*@C
2773: KSPBuildResidual - Builds the residual in a vector provided.
2775: Collective
2777: Input Parameter:
2778: . ksp - iterative context obtained from `KSPCreate()`
2780: Output Parameters:
2781: + v - optional location to stash residual. If `v` is not provided,
2782: then a location is generated.
2783: . t - work vector. If not provided then one is generated.
2784: - V - the residual
2786: Level: advanced
2788: Note:
2789: Regardless of whether or not `v` is provided, the residual is
2790: returned in `V`.
2792: .seealso: [](ch_ksp), `KSP`, `KSPBuildSolution()`
2793: @*/
2794: PetscErrorCode KSPBuildResidual(KSP ksp, Vec t, Vec v, Vec *V)
2795: {
2796: PetscBool flag = PETSC_FALSE;
2797: Vec w = v, tt = t;
2799: PetscFunctionBegin;
2801: if (!w) PetscCall(VecDuplicate(ksp->vec_rhs, &w));
2802: if (!tt) {
2803: PetscCall(VecDuplicate(ksp->vec_sol, &tt));
2804: flag = PETSC_TRUE;
2805: }
2806: PetscUseTypeMethod(ksp, buildresidual, tt, w, V);
2807: if (flag) PetscCall(VecDestroy(&tt));
2808: PetscFunctionReturn(PETSC_SUCCESS);
2809: }
2811: /*@
2812: KSPSetDiagonalScale - Tells `KSP` to symmetrically diagonally scale the system
2813: before solving. This actually CHANGES the matrix (and right hand side).
2815: Logically Collective
2817: Input Parameters:
2818: + ksp - the `KSP` context
2819: - scale - `PETSC_TRUE` or `PETSC_FALSE`
2821: Options Database Keys:
2822: + -ksp_diagonal_scale - perform a diagonal scaling before the solve
2823: - -ksp_diagonal_scale_fix - scale the matrix back AFTER the solve
2825: Level: advanced
2827: Notes:
2828: Scales the matrix by $D^{-1/2} A D^{-1/2} [D^{1/2} x ] = D^{-1/2} b $
2829: where $D_{ii}$ is $1/abs(A_{ii}) $ unless $A_{ii}$ is zero and then it is 1.
2831: BE CAREFUL with this routine: it actually scales the matrix and right
2832: hand side that define the system. After the system is solved the matrix
2833: and right hand side remain scaled unless you use `KSPSetDiagonalScaleFix()`
2835: This should NOT be used within the `SNES` solves if you are using a line
2836: search.
2838: If you use this with the `PCType` `PCEISENSTAT` preconditioner than you can
2839: use the `PCEisenstatSetNoDiagonalScaling()` option, or `-pc_eisenstat_no_diagonal_scaling`
2840: to save some unneeded, redundant flops.
2842: .seealso: [](ch_ksp), `KSPGetDiagonalScale()`, `KSPSetDiagonalScaleFix()`, `KSP`
2843: @*/
2844: PetscErrorCode KSPSetDiagonalScale(KSP ksp, PetscBool scale)
2845: {
2846: PetscFunctionBegin;
2849: ksp->dscale = scale;
2850: PetscFunctionReturn(PETSC_SUCCESS);
2851: }
2853: /*@
2854: KSPGetDiagonalScale - Checks if `KSP` solver scales the matrix and right hand side, that is if `KSPSetDiagonalScale()` has been called
2856: Not Collective
2858: Input Parameter:
2859: . ksp - the `KSP` context
2861: Output Parameter:
2862: . scale - `PETSC_TRUE` or `PETSC_FALSE`
2864: Level: intermediate
2866: .seealso: [](ch_ksp), `KSP`, `KSPSetDiagonalScale()`, `KSPSetDiagonalScaleFix()`
2867: @*/
2868: PetscErrorCode KSPGetDiagonalScale(KSP ksp, PetscBool *scale)
2869: {
2870: PetscFunctionBegin;
2872: PetscAssertPointer(scale, 2);
2873: *scale = ksp->dscale;
2874: PetscFunctionReturn(PETSC_SUCCESS);
2875: }
2877: /*@
2878: KSPSetDiagonalScaleFix - Tells `KSP` to diagonally scale the system back after solving.
2880: Logically Collective
2882: Input Parameters:
2883: + ksp - the `KSP` context
2884: - fix - `PETSC_TRUE` to scale back after the system solve, `PETSC_FALSE` to not
2885: rescale (default)
2887: Level: intermediate
2889: Notes:
2890: Must be called after `KSPSetDiagonalScale()`
2892: Using this will slow things down, because it rescales the matrix before and
2893: after each linear solve. This is intended mainly for testing to allow one
2894: to easily get back the original system to make sure the solution computed is
2895: accurate enough.
2897: .seealso: [](ch_ksp), `KSPGetDiagonalScale()`, `KSPSetDiagonalScale()`, `KSPGetDiagonalScaleFix()`, `KSP`
2898: @*/
2899: PetscErrorCode KSPSetDiagonalScaleFix(KSP ksp, PetscBool fix)
2900: {
2901: PetscFunctionBegin;
2904: ksp->dscalefix = fix;
2905: PetscFunctionReturn(PETSC_SUCCESS);
2906: }
2908: /*@
2909: KSPGetDiagonalScaleFix - Determines if `KSP` diagonally scales the system back after solving. That is `KSPSetDiagonalScaleFix()` has been called
2911: Not Collective
2913: Input Parameter:
2914: . ksp - the `KSP` context
2916: Output Parameter:
2917: . fix - `PETSC_TRUE` to scale back after the system solve, `PETSC_FALSE` to not
2918: rescale (default)
2920: Level: intermediate
2922: .seealso: [](ch_ksp), `KSPGetDiagonalScale()`, `KSPSetDiagonalScale()`, `KSPSetDiagonalScaleFix()`, `KSP`
2923: @*/
2924: PetscErrorCode KSPGetDiagonalScaleFix(KSP ksp, PetscBool *fix)
2925: {
2926: PetscFunctionBegin;
2928: PetscAssertPointer(fix, 2);
2929: *fix = ksp->dscalefix;
2930: PetscFunctionReturn(PETSC_SUCCESS);
2931: }
2933: /*@C
2934: KSPSetComputeOperators - set routine to compute the linear operators
2936: Logically Collective
2938: Input Parameters:
2939: + ksp - the `KSP` context
2940: . func - function to compute the operators, see `KSPComputeOperatorsFn` for the calling sequence
2941: - ctx - optional context
2943: Level: beginner
2945: Notes:
2946: `func()` will be called automatically at the very next call to `KSPSolve()`. It will NOT be called at future `KSPSolve()` calls
2947: unless either `KSPSetComputeOperators()` or `KSPSetOperators()` is called before that `KSPSolve()` is called. This allows the same system to be solved several times
2948: with different right hand side functions but is a confusing API since one might expect it to be called for each `KSPSolve()`
2950: To reuse the same preconditioner for the next `KSPSolve()` and not compute a new one based on the most recently computed matrix call `KSPSetReusePreconditioner()`
2952: Developer Note:
2953: Perhaps this routine and `KSPSetComputeRHS()` could be combined into a new API that makes clear when new matrices are computing without requiring call this
2954: routine to indicate when the new matrix should be computed.
2956: .seealso: [](ch_ksp), `KSP`, `KSPSetOperators()`, `KSPSetComputeRHS()`, `DMKSPSetComputeOperators()`, `KSPSetComputeInitialGuess()`, `KSPComputeOperatorsFn`
2957: @*/
2958: PetscErrorCode KSPSetComputeOperators(KSP ksp, KSPComputeOperatorsFn *func, void *ctx)
2959: {
2960: DM dm;
2962: PetscFunctionBegin;
2964: PetscCall(KSPGetDM(ksp, &dm));
2965: PetscCall(DMKSPSetComputeOperators(dm, func, ctx));
2966: if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX;
2967: PetscFunctionReturn(PETSC_SUCCESS);
2968: }
2970: /*@C
2971: KSPSetComputeRHS - set routine to compute the right hand side of the linear system
2973: Logically Collective
2975: Input Parameters:
2976: + ksp - the `KSP` context
2977: . func - function to compute the right hand side, see `KSPComputeRHSFn` for the calling squence
2978: - ctx - optional context
2980: Level: beginner
2982: Note:
2983: The routine you provide will be called EACH you call `KSPSolve()` to prepare the new right hand side for that solve
2985: .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `DMKSPSetComputeRHS()`, `KSPSetComputeOperators()`, `KSPSetOperators()`, `KSPComputeRHSFn`
2986: @*/
2987: PetscErrorCode KSPSetComputeRHS(KSP ksp, KSPComputeRHSFn *func, void *ctx)
2988: {
2989: DM dm;
2991: PetscFunctionBegin;
2993: PetscCall(KSPGetDM(ksp, &dm));
2994: PetscCall(DMKSPSetComputeRHS(dm, func, ctx));
2995: PetscFunctionReturn(PETSC_SUCCESS);
2996: }
2998: /*@C
2999: KSPSetComputeInitialGuess - set routine to compute the initial guess of the linear system
3001: Logically Collective
3003: Input Parameters:
3004: + ksp - the `KSP` context
3005: . func - function to compute the initial guess, see `KSPComputeInitialGuessFn` for calling sequence
3006: - ctx - optional context
3008: Level: beginner
3010: Note:
3011: This should only be used in conjunction with `KSPSetComputeRHS()` and `KSPSetComputeOperators()`, otherwise
3012: call `KSPSetInitialGuessNonzero()` and set the initial guess values in the solution vector passed to `KSPSolve()` before calling the solver
3014: .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `KSPSetComputeRHS()`, `KSPSetComputeOperators()`, `DMKSPSetComputeInitialGuess()`, `KSPSetInitialGuessNonzero()`,
3015: `KSPComputeInitialGuessFn`
3016: @*/
3017: PetscErrorCode KSPSetComputeInitialGuess(KSP ksp, KSPComputeInitialGuessFn *func, void *ctx)
3018: {
3019: DM dm;
3021: PetscFunctionBegin;
3023: PetscCall(KSPGetDM(ksp, &dm));
3024: PetscCall(DMKSPSetComputeInitialGuess(dm, func, ctx));
3025: PetscFunctionReturn(PETSC_SUCCESS);
3026: }
3028: /*@
3029: KSPSetUseExplicitTranspose - Determines the explicit transpose of the operator is formed in `KSPSolveTranspose()`. In some configurations (like GPUs) it may
3030: be explicitly formed since the solve is much more efficient.
3032: Logically Collective
3034: Input Parameter:
3035: . ksp - the `KSP` context
3037: Output Parameter:
3038: . flg - `PETSC_TRUE` to transpose the system in `KSPSolveTranspose()`, `PETSC_FALSE` to not transpose (default)
3040: Level: advanced
3042: .seealso: [](ch_ksp), `KSPSolveTranspose()`, `KSP`
3043: @*/
3044: PetscErrorCode KSPSetUseExplicitTranspose(KSP ksp, PetscBool flg)
3045: {
3046: PetscFunctionBegin;
3049: ksp->transpose.use_explicittranspose = flg;
3050: PetscFunctionReturn(PETSC_SUCCESS);
3051: }