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;
581: /* Call all user-provided reason review routines */
582: for (i = 0; i < ksp->numberreasonviews; i++) PetscCall((*ksp->reasonview[i])(ksp, ksp->reasonviewcontext[i]));
584: /* Call the default PETSc routine */
585: PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp), ((PetscObject)ksp)->options, ((PetscObject)ksp)->prefix, "-ksp_converged_reason", &viewer, &format, &flg));
586: if (flg) {
587: PetscCall(PetscViewerPushFormat(viewer, format));
588: PetscCall(KSPConvergedReasonView(ksp, viewer));
589: PetscCall(PetscViewerPopFormat(viewer));
590: PetscCall(PetscViewerDestroy(&viewer));
591: }
592: PetscFunctionReturn(PETSC_SUCCESS);
593: }
595: /*@C
596: KSPConvergedRateView - Displays the convergence rate <https://en.wikipedia.org/wiki/Coefficient_of_determination> of `KSPSolve()` to a viewer
598: Collective
600: Input Parameters:
601: + ksp - iterative context obtained from `KSPCreate()`
602: - viewer - the viewer to display the reason
604: Options Database Key:
605: . -ksp_converged_rate - print reason for convergence or divergence and the convergence rate (or 0.0 for divergence)
607: Level: intermediate
609: Notes:
610: To change the format of the output, call `PetscViewerPushFormat`(`viewer`,`format`) before this call.
612: 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,
613: 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)}$,
615: .seealso: [](ch_ksp), `KSPConvergedReasonView()`, `KSPGetConvergedRate()`, `KSPSetTolerances()`, `KSPConvergedDefault()`
616: @*/
617: PetscErrorCode KSPConvergedRateView(KSP ksp, PetscViewer viewer)
618: {
619: PetscViewerFormat format;
620: PetscBool isAscii;
621: PetscReal rrate, rRsq, erate = 0.0, eRsq = 0.0;
622: PetscInt its;
623: const char *prefix, *reason = KSPConvergedReasons[ksp->reason];
625: PetscFunctionBegin;
626: PetscCall(KSPGetOptionsPrefix(ksp, &prefix));
627: PetscCall(KSPGetIterationNumber(ksp, &its));
628: PetscCall(KSPComputeConvergenceRate(ksp, &rrate, &rRsq, &erate, &eRsq));
629: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)ksp));
630: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isAscii));
631: if (isAscii) {
632: PetscCall(PetscViewerGetFormat(viewer, &format));
633: PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ksp)->tablevel));
634: if (ksp->reason > 0) {
635: if (prefix) PetscCall(PetscViewerASCIIPrintf(viewer, "Linear %s solve converged due to %s iterations %" PetscInt_FMT, prefix, reason, its));
636: else PetscCall(PetscViewerASCIIPrintf(viewer, "Linear solve converged due to %s iterations %" PetscInt_FMT, reason, its));
637: PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
638: if (rRsq >= 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, " res rate %g R^2 %g", (double)rrate, (double)rRsq));
639: if (eRsq >= 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, " error rate %g R^2 %g", (double)erate, (double)eRsq));
640: PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
641: PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
642: } else if (ksp->reason <= 0) {
643: if (prefix) PetscCall(PetscViewerASCIIPrintf(viewer, "Linear %s solve did not converge due to %s iterations %" PetscInt_FMT, prefix, reason, its));
644: else PetscCall(PetscViewerASCIIPrintf(viewer, "Linear solve did not converge due to %s iterations %" PetscInt_FMT, reason, its));
645: PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
646: if (rRsq >= 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, " res rate %g R^2 %g", (double)rrate, (double)rRsq));
647: if (eRsq >= 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, " error rate %g R^2 %g", (double)erate, (double)eRsq));
648: PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
649: PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
650: if (ksp->reason == KSP_DIVERGED_PC_FAILED) {
651: PCFailedReason reason;
652: PetscCall(PCGetFailedReason(ksp->pc, &reason));
653: PetscCall(PetscViewerASCIIPrintf(viewer, " PC failed due to %s \n", PCFailedReasons[reason]));
654: }
655: }
656: PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ksp)->tablevel));
657: }
658: PetscFunctionReturn(PETSC_SUCCESS);
659: }
661: #include <petscdraw.h>
663: static PetscErrorCode KSPViewEigenvalues_Internal(KSP ksp, PetscBool isExplicit, PetscViewer viewer, PetscViewerFormat format)
664: {
665: PetscReal *r, *c;
666: PetscInt n, i, neig;
667: PetscBool isascii, isdraw;
668: PetscMPIInt rank;
670: PetscFunctionBegin;
671: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)ksp), &rank));
672: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
673: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
674: if (isExplicit) {
675: PetscCall(VecGetSize(ksp->vec_sol, &n));
676: PetscCall(PetscMalloc2(n, &r, n, &c));
677: PetscCall(KSPComputeEigenvaluesExplicitly(ksp, n, r, c));
678: neig = n;
679: } else {
680: PetscInt nits;
682: PetscCall(KSPGetIterationNumber(ksp, &nits));
683: n = nits + 2;
684: if (!nits) {
685: PetscCall(PetscViewerASCIIPrintf(viewer, "Zero iterations in solver, cannot approximate any eigenvalues\n"));
686: PetscFunctionReturn(PETSC_SUCCESS);
687: }
688: PetscCall(PetscMalloc2(n, &r, n, &c));
689: PetscCall(KSPComputeEigenvalues(ksp, n, r, c, &neig));
690: }
691: if (isascii) {
692: PetscCall(PetscViewerASCIIPrintf(viewer, "%s computed eigenvalues\n", isExplicit ? "Explicitly" : "Iteratively"));
693: for (i = 0; i < neig; ++i) {
694: if (c[i] >= 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, "%g + %gi\n", (double)r[i], (double)c[i]));
695: else PetscCall(PetscViewerASCIIPrintf(viewer, "%g - %gi\n", (double)r[i], -(double)c[i]));
696: }
697: } else if (isdraw && rank == 0) {
698: PetscDraw draw;
699: PetscDrawSP drawsp;
701: if (format == PETSC_VIEWER_DRAW_CONTOUR) {
702: PetscCall(KSPPlotEigenContours_Private(ksp, neig, r, c));
703: } else {
704: PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
705: PetscCall(PetscDrawSPCreate(draw, 1, &drawsp));
706: PetscCall(PetscDrawSPReset(drawsp));
707: for (i = 0; i < neig; ++i) PetscCall(PetscDrawSPAddPoint(drawsp, r + i, c + i));
708: PetscCall(PetscDrawSPDraw(drawsp, PETSC_TRUE));
709: PetscCall(PetscDrawSPSave(drawsp));
710: PetscCall(PetscDrawSPDestroy(&drawsp));
711: }
712: }
713: PetscCall(PetscFree2(r, c));
714: PetscFunctionReturn(PETSC_SUCCESS);
715: }
717: static PetscErrorCode KSPViewSingularvalues_Internal(KSP ksp, PetscViewer viewer, PetscViewerFormat format)
718: {
719: PetscReal smax, smin;
720: PetscInt nits;
721: PetscBool isascii;
723: PetscFunctionBegin;
724: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
725: PetscCall(KSPGetIterationNumber(ksp, &nits));
726: if (!nits) {
727: PetscCall(PetscViewerASCIIPrintf(viewer, "Zero iterations in solver, cannot approximate any singular values\n"));
728: PetscFunctionReturn(PETSC_SUCCESS);
729: }
730: PetscCall(KSPComputeExtremeSingularValues(ksp, &smax, &smin));
731: if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, "Iteratively computed extreme singular values: max %g min %g max/min %g\n", (double)smax, (double)smin, (double)(smax / smin)));
732: PetscFunctionReturn(PETSC_SUCCESS);
733: }
735: static PetscErrorCode KSPViewFinalResidual_Internal(KSP ksp, PetscViewer viewer, PetscViewerFormat format)
736: {
737: PetscBool isascii;
739: PetscFunctionBegin;
740: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
741: 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");
742: if (isascii) {
743: Mat A;
744: Vec t;
745: PetscReal norm;
747: PetscCall(PCGetOperators(ksp->pc, &A, NULL));
748: PetscCall(VecDuplicate(ksp->vec_rhs, &t));
749: PetscCall(KSP_MatMult(ksp, A, ksp->vec_sol, t));
750: PetscCall(VecAYPX(t, -1.0, ksp->vec_rhs));
751: PetscCall(VecViewFromOptions(t, (PetscObject)ksp, "-ksp_view_final_residual_vec"));
752: PetscCall(VecNorm(t, NORM_2, &norm));
753: PetscCall(VecDestroy(&t));
754: PetscCall(PetscViewerASCIIPrintf(viewer, "KSP final norm of residual %g\n", (double)norm));
755: }
756: PetscFunctionReturn(PETSC_SUCCESS);
757: }
759: static PetscErrorCode KSPMonitorPauseFinal_Internal(KSP ksp)
760: {
761: PetscInt i;
763: PetscFunctionBegin;
764: if (!ksp->pauseFinal) PetscFunctionReturn(PETSC_SUCCESS);
765: for (i = 0; i < ksp->numbermonitors; ++i) {
766: PetscViewerAndFormat *vf = (PetscViewerAndFormat *)ksp->monitorcontext[i];
767: PetscDraw draw;
768: PetscReal lpause;
770: if (!vf) continue;
771: if (vf->lg) {
772: if (!PetscCheckPointer(vf->lg, PETSC_OBJECT)) continue;
773: if (((PetscObject)vf->lg)->classid != PETSC_DRAWLG_CLASSID) continue;
774: PetscCall(PetscDrawLGGetDraw(vf->lg, &draw));
775: PetscCall(PetscDrawGetPause(draw, &lpause));
776: PetscCall(PetscDrawSetPause(draw, -1.0));
777: PetscCall(PetscDrawPause(draw));
778: PetscCall(PetscDrawSetPause(draw, lpause));
779: } else {
780: PetscBool isdraw;
782: if (!PetscCheckPointer(vf->viewer, PETSC_OBJECT)) continue;
783: if (((PetscObject)vf->viewer)->classid != PETSC_VIEWER_CLASSID) continue;
784: PetscCall(PetscObjectTypeCompare((PetscObject)vf->viewer, PETSCVIEWERDRAW, &isdraw));
785: if (!isdraw) continue;
786: PetscCall(PetscViewerDrawGetDraw(vf->viewer, 0, &draw));
787: PetscCall(PetscDrawGetPause(draw, &lpause));
788: PetscCall(PetscDrawSetPause(draw, -1.0));
789: PetscCall(PetscDrawPause(draw));
790: PetscCall(PetscDrawSetPause(draw, lpause));
791: }
792: }
793: PetscFunctionReturn(PETSC_SUCCESS);
794: }
796: static PetscErrorCode KSPSolve_Private(KSP ksp, Vec b, Vec x)
797: {
798: PetscBool flg = PETSC_FALSE, inXisinB = PETSC_FALSE, guess_zero;
799: Mat mat, pmat;
800: MPI_Comm comm;
801: MatNullSpace nullsp;
802: Vec btmp, vec_rhs = NULL;
804: PetscFunctionBegin;
805: level++;
806: comm = PetscObjectComm((PetscObject)ksp);
807: if (x && x == b) {
808: PetscCheck(ksp->guess_zero, comm, PETSC_ERR_ARG_INCOMP, "Cannot use x == b with nonzero initial guess");
809: PetscCall(VecDuplicate(b, &x));
810: inXisinB = PETSC_TRUE;
811: }
812: if (b) {
813: PetscCall(PetscObjectReference((PetscObject)b));
814: PetscCall(VecDestroy(&ksp->vec_rhs));
815: ksp->vec_rhs = b;
816: }
817: if (x) {
818: PetscCall(PetscObjectReference((PetscObject)x));
819: PetscCall(VecDestroy(&ksp->vec_sol));
820: ksp->vec_sol = x;
821: }
823: if (ksp->viewPre) PetscCall(ObjectView((PetscObject)ksp, ksp->viewerPre, ksp->formatPre));
825: if (ksp->presolve) PetscCall((*ksp->presolve)(ksp, ksp->vec_rhs, ksp->vec_sol, ksp->prectx));
827: /* reset the residual history list if requested */
828: if (ksp->res_hist_reset) ksp->res_hist_len = 0;
829: if (ksp->err_hist_reset) ksp->err_hist_len = 0;
831: /* KSPSetUp() scales the matrix if needed */
832: PetscCall(KSPSetUp(ksp));
833: PetscCall(KSPSetUpOnBlocks(ksp));
835: if (ksp->guess) {
836: PetscObjectState ostate, state;
838: PetscCall(KSPGuessSetUp(ksp->guess));
839: PetscCall(PetscObjectStateGet((PetscObject)ksp->vec_sol, &ostate));
840: PetscCall(KSPGuessFormGuess(ksp->guess, ksp->vec_rhs, ksp->vec_sol));
841: PetscCall(PetscObjectStateGet((PetscObject)ksp->vec_sol, &state));
842: if (state != ostate) {
843: ksp->guess_zero = PETSC_FALSE;
844: } else {
845: PetscCall(PetscInfo(ksp, "Using zero initial guess since the KSPGuess object did not change the vector\n"));
846: ksp->guess_zero = PETSC_TRUE;
847: }
848: }
850: PetscCall(VecSetErrorIfLocked(ksp->vec_sol, 3));
852: PetscCall(PetscLogEventBegin(!ksp->transpose_solve ? KSP_Solve : KSP_SolveTranspose, ksp, ksp->vec_rhs, ksp->vec_sol, 0));
853: PetscCall(PCGetOperators(ksp->pc, &mat, &pmat));
854: /* diagonal scale RHS if called for */
855: if (ksp->dscale) {
856: PetscCall(VecPointwiseMult(ksp->vec_rhs, ksp->vec_rhs, ksp->diagonal));
857: /* second time in, but matrix was scaled back to original */
858: if (ksp->dscalefix && ksp->dscalefix2) {
859: Mat mat, pmat;
861: PetscCall(PCGetOperators(ksp->pc, &mat, &pmat));
862: PetscCall(MatDiagonalScale(pmat, ksp->diagonal, ksp->diagonal));
863: if (mat != pmat) PetscCall(MatDiagonalScale(mat, ksp->diagonal, ksp->diagonal));
864: }
866: /* scale initial guess */
867: if (!ksp->guess_zero) {
868: if (!ksp->truediagonal) {
869: PetscCall(VecDuplicate(ksp->diagonal, &ksp->truediagonal));
870: PetscCall(VecCopy(ksp->diagonal, ksp->truediagonal));
871: PetscCall(VecReciprocal(ksp->truediagonal));
872: }
873: PetscCall(VecPointwiseMult(ksp->vec_sol, ksp->vec_sol, ksp->truediagonal));
874: }
875: }
876: PetscCall(PCPreSolve(ksp->pc, ksp));
878: if (ksp->guess_zero && !ksp->guess_not_read) PetscCall(VecSet(ksp->vec_sol, 0.0));
879: if (ksp->guess_knoll) { /* The Knoll trick is independent on the KSPGuess specified */
880: PetscCall(PCApply(ksp->pc, ksp->vec_rhs, ksp->vec_sol));
881: PetscCall(KSP_RemoveNullSpace(ksp, ksp->vec_sol));
882: ksp->guess_zero = PETSC_FALSE;
883: }
885: /* can we mark the initial guess as zero for this solve? */
886: guess_zero = ksp->guess_zero;
887: if (!ksp->guess_zero) {
888: PetscReal norm;
890: PetscCall(VecNormAvailable(ksp->vec_sol, NORM_2, &flg, &norm));
891: if (flg && !norm) ksp->guess_zero = PETSC_TRUE;
892: }
893: if (ksp->transpose_solve) {
894: PetscCall(MatGetNullSpace(pmat, &nullsp));
895: } else {
896: PetscCall(MatGetTransposeNullSpace(pmat, &nullsp));
897: }
898: if (nullsp) {
899: PetscCall(VecDuplicate(ksp->vec_rhs, &btmp));
900: PetscCall(VecCopy(ksp->vec_rhs, btmp));
901: PetscCall(MatNullSpaceRemove(nullsp, btmp));
902: vec_rhs = ksp->vec_rhs;
903: ksp->vec_rhs = btmp;
904: }
905: PetscCall(VecLockReadPush(ksp->vec_rhs));
906: PetscUseTypeMethod(ksp, solve);
907: PetscCall(KSPMonitorPauseFinal_Internal(ksp));
909: PetscCall(VecLockReadPop(ksp->vec_rhs));
910: if (nullsp) {
911: ksp->vec_rhs = vec_rhs;
912: PetscCall(VecDestroy(&btmp));
913: }
915: ksp->guess_zero = guess_zero;
917: PetscCheck(ksp->reason, comm, PETSC_ERR_PLIB, "Internal error, solver returned without setting converged reason");
918: ksp->totalits += ksp->its;
920: PetscCall(KSPConvergedReasonViewFromOptions(ksp));
922: if (ksp->viewRate) {
923: PetscCall(PetscViewerPushFormat(ksp->viewerRate, ksp->formatRate));
924: PetscCall(KSPConvergedRateView(ksp, ksp->viewerRate));
925: PetscCall(PetscViewerPopFormat(ksp->viewerRate));
926: }
927: PetscCall(PCPostSolve(ksp->pc, ksp));
929: /* diagonal scale solution if called for */
930: if (ksp->dscale) {
931: PetscCall(VecPointwiseMult(ksp->vec_sol, ksp->vec_sol, ksp->diagonal));
932: /* unscale right hand side and matrix */
933: if (ksp->dscalefix) {
934: Mat mat, pmat;
936: PetscCall(VecReciprocal(ksp->diagonal));
937: PetscCall(VecPointwiseMult(ksp->vec_rhs, ksp->vec_rhs, ksp->diagonal));
938: PetscCall(PCGetOperators(ksp->pc, &mat, &pmat));
939: PetscCall(MatDiagonalScale(pmat, ksp->diagonal, ksp->diagonal));
940: if (mat != pmat) PetscCall(MatDiagonalScale(mat, ksp->diagonal, ksp->diagonal));
941: PetscCall(VecReciprocal(ksp->diagonal));
942: ksp->dscalefix2 = PETSC_TRUE;
943: }
944: }
945: PetscCall(PetscLogEventEnd(!ksp->transpose_solve ? KSP_Solve : KSP_SolveTranspose, ksp, ksp->vec_rhs, ksp->vec_sol, 0));
946: if (ksp->guess) PetscCall(KSPGuessUpdate(ksp->guess, ksp->vec_rhs, ksp->vec_sol));
947: if (ksp->postsolve) PetscCall((*ksp->postsolve)(ksp, ksp->vec_rhs, ksp->vec_sol, ksp->postctx));
949: PetscCall(PCGetOperators(ksp->pc, &mat, &pmat));
950: if (ksp->viewEV) PetscCall(KSPViewEigenvalues_Internal(ksp, PETSC_FALSE, ksp->viewerEV, ksp->formatEV));
951: if (ksp->viewEVExp) PetscCall(KSPViewEigenvalues_Internal(ksp, PETSC_TRUE, ksp->viewerEVExp, ksp->formatEVExp));
952: if (ksp->viewSV) PetscCall(KSPViewSingularvalues_Internal(ksp, ksp->viewerSV, ksp->formatSV));
953: if (ksp->viewFinalRes) PetscCall(KSPViewFinalResidual_Internal(ksp, ksp->viewerFinalRes, ksp->formatFinalRes));
954: if (ksp->viewMat) PetscCall(ObjectView((PetscObject)mat, ksp->viewerMat, ksp->formatMat));
955: if (ksp->viewPMat) PetscCall(ObjectView((PetscObject)pmat, ksp->viewerPMat, ksp->formatPMat));
956: if (ksp->viewRhs) PetscCall(ObjectView((PetscObject)ksp->vec_rhs, ksp->viewerRhs, ksp->formatRhs));
957: if (ksp->viewSol) PetscCall(ObjectView((PetscObject)ksp->vec_sol, ksp->viewerSol, ksp->formatSol));
958: if (ksp->view) PetscCall(ObjectView((PetscObject)ksp, ksp->viewer, ksp->format));
959: if (ksp->viewDScale) PetscCall(ObjectView((PetscObject)ksp->diagonal, ksp->viewerDScale, ksp->formatDScale));
960: if (ksp->viewMatExp) {
961: Mat A, B;
963: PetscCall(PCGetOperators(ksp->pc, &A, NULL));
964: if (ksp->transpose_solve) {
965: Mat AT;
967: PetscCall(MatCreateTranspose(A, &AT));
968: PetscCall(MatComputeOperator(AT, MATAIJ, &B));
969: PetscCall(MatDestroy(&AT));
970: } else {
971: PetscCall(MatComputeOperator(A, MATAIJ, &B));
972: }
973: PetscCall(ObjectView((PetscObject)B, ksp->viewerMatExp, ksp->formatMatExp));
974: PetscCall(MatDestroy(&B));
975: }
976: if (ksp->viewPOpExp) {
977: Mat B;
979: PetscCall(KSPComputeOperator(ksp, MATAIJ, &B));
980: PetscCall(ObjectView((PetscObject)B, ksp->viewerPOpExp, ksp->formatPOpExp));
981: PetscCall(MatDestroy(&B));
982: }
984: if (inXisinB) {
985: PetscCall(VecCopy(x, b));
986: PetscCall(VecDestroy(&x));
987: }
988: PetscCall(PetscObjectSAWsBlock((PetscObject)ksp));
989: if (ksp->errorifnotconverged && ksp->reason < 0 && ((level == 1) || (ksp->reason != KSP_DIVERGED_ITS))) {
990: PCFailedReason reason;
992: 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]);
993: PetscCall(PCGetFailedReason(ksp->pc, &reason));
994: 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]);
995: }
996: level--;
997: PetscFunctionReturn(PETSC_SUCCESS);
998: }
1000: /*@
1001: KSPSolve - Solves linear system.
1003: Collective
1005: Input Parameters:
1006: + ksp - iterative context obtained from `KSPCreate()`
1007: . b - the right hand side vector
1008: - x - the solution (this may be the same vector as `b`, then `b` will be overwritten with answer)
1010: Options Database Keys:
1011: + -ksp_view_eigenvalues - compute preconditioned operators eigenvalues
1012: . -ksp_view_eigenvalues_explicit - compute the eigenvalues by forming the dense operator and using LAPACK
1013: . -ksp_view_mat binary - save matrix to the default binary viewer
1014: . -ksp_view_pmat binary - save matrix used to build preconditioner to the default binary viewer
1015: . -ksp_view_rhs binary - save right hand side vector to the default binary viewer
1016: . -ksp_view_solution binary - save computed solution vector to the default binary viewer
1017: (can be read later with src/ksp/tutorials/ex10.c for testing solvers)
1018: . -ksp_view_mat_explicit - for matrix-free operators, computes the matrix entries and views them
1019: . -ksp_view_preconditioned_operator_explicit - computes the product of the preconditioner and matrix as an explicit matrix and views it
1020: . -ksp_converged_reason - print reason for converged or diverged, also prints number of iterations
1021: . -ksp_view_final_residual - print 2-norm of true linear system residual at the end of the solution process
1022: . -ksp_error_if_not_converged - stop the program as soon as an error is detected in a `KSPSolve()`
1023: . -ksp_view_pre - print the ksp data structure before the system solution
1024: - -ksp_view - print the ksp data structure at the end of the system solution
1026: Level: beginner
1028: Notes:
1029: If one uses `KSPSetDM()` then `x` or `b` need not be passed. Use `KSPGetSolution()` to access the solution in this case.
1031: The operator is specified with `KSPSetOperators()`.
1033: `KSPSolve()` will normally return without generating an error regardless of whether the linear system was solved or if constructing the preconditioner failed.
1034: Call `KSPGetConvergedReason()` to determine if the solver converged or failed and why. The option -ksp_error_if_not_converged or function `KSPSetErrorIfNotConverged()`
1035: 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
1036: it may be fine that inner solvers in the preconditioner do not converge during the solution process.
1038: The number of iterations can be obtained from `KSPGetIterationNumber()`.
1040: If you provide a matrix that has a `MatSetNullSpace()` and `MatSetTransposeNullSpace()` this will use that information to solve singular systems
1041: in the least squares sense with a norm minimizing solution.
1043: $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()`).
1045: `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
1046: 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
1047: direction thus the solution which is a linear combination of the search directions has no component in the nullspace(A).
1049: We recommend always using `KSPGMRES` for such singular systems.
1050: If nullspace(A) = nullspace(A') (note symmetric matrices always satisfy this property) then both left and right preconditioning will work
1051: If nullspace(A) != nullspace(A') then left preconditioning will work but right preconditioning may not work (or it may).
1053: Developer Notes:
1054: The reason we cannot always solve nullspace(A) != nullspace(A') systems with right preconditioning is because we need to remove at each iteration
1055: 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
1056: such as diagonal scaling we cannot apply the inverse of the preconditioner to a vector and thus cannot compute the nullspace(AB).
1058: If using a direct method (e.g., via the `KSP` solver
1059: `KSPPREONLY` and a preconditioner such as `PCLU` or `PCILU`,
1060: then its=1. See `KSPSetTolerances()` and `KSPConvergedDefault()`
1061: for more details.
1063: Understanding Convergence\:
1064: The routines `KSPMonitorSet()`, `KSPComputeEigenvalues()`, and
1065: `KSPComputeEigenvaluesExplicitly()` provide information on additional
1066: options to monitor convergence and print eigenvalue information.
1068: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPDestroy()`, `KSPSetTolerances()`, `KSPConvergedDefault()`,
1069: `KSPSolveTranspose()`, `KSPGetIterationNumber()`, `MatNullSpaceCreate()`, `MatSetNullSpace()`, `MatSetTransposeNullSpace()`, `KSP`,
1070: `KSPConvergedReasonView()`, `KSPCheckSolve()`, `KSPSetErrorIfNotConverged()`
1071: @*/
1072: PetscErrorCode KSPSolve(KSP ksp, Vec b, Vec x)
1073: {
1074: PetscFunctionBegin;
1078: ksp->transpose_solve = PETSC_FALSE;
1079: PetscCall(KSPSolve_Private(ksp, b, x));
1080: PetscFunctionReturn(PETSC_SUCCESS);
1081: }
1083: /*@
1084: KSPSolveTranspose - Solves a linear system with the transpose of the matrix, $ A^T x = b$.
1086: Collective
1088: Input Parameters:
1089: + ksp - iterative context obtained from `KSPCreate()`
1090: . b - right hand side vector
1091: - x - solution vector
1093: Level: developer
1095: Note:
1096: For complex numbers this solve the non-Hermitian transpose system.
1098: Developer Note:
1099: We need to implement a `KSPSolveHermitianTranspose()`
1101: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPDestroy()`, `KSPSetTolerances()`, `KSPConvergedDefault()`,
1102: `KSPSolve()`, `KSP`
1103: @*/
1104: PetscErrorCode KSPSolveTranspose(KSP ksp, Vec b, Vec x)
1105: {
1106: PetscFunctionBegin;
1110: if (ksp->transpose.use_explicittranspose) {
1111: Mat J, Jpre;
1112: PetscCall(KSPGetOperators(ksp, &J, &Jpre));
1113: if (!ksp->transpose.reuse_transpose) {
1114: PetscCall(MatTranspose(J, MAT_INITIAL_MATRIX, &ksp->transpose.AT));
1115: if (J != Jpre) PetscCall(MatTranspose(Jpre, MAT_INITIAL_MATRIX, &ksp->transpose.BT));
1116: ksp->transpose.reuse_transpose = PETSC_TRUE;
1117: } else {
1118: PetscCall(MatTranspose(J, MAT_REUSE_MATRIX, &ksp->transpose.AT));
1119: if (J != Jpre) PetscCall(MatTranspose(Jpre, MAT_REUSE_MATRIX, &ksp->transpose.BT));
1120: }
1121: if (J == Jpre && ksp->transpose.BT != ksp->transpose.AT) {
1122: PetscCall(PetscObjectReference((PetscObject)ksp->transpose.AT));
1123: ksp->transpose.BT = ksp->transpose.AT;
1124: }
1125: PetscCall(KSPSetOperators(ksp, ksp->transpose.AT, ksp->transpose.BT));
1126: } else {
1127: ksp->transpose_solve = PETSC_TRUE;
1128: }
1129: PetscCall(KSPSolve_Private(ksp, b, x));
1130: PetscFunctionReturn(PETSC_SUCCESS);
1131: }
1133: static PetscErrorCode KSPViewFinalMatResidual_Internal(KSP ksp, Mat B, Mat X, PetscViewer viewer, PetscViewerFormat format, PetscInt shift)
1134: {
1135: Mat A, R;
1136: PetscReal *norms;
1137: PetscInt i, N;
1138: PetscBool flg;
1140: PetscFunctionBegin;
1141: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &flg));
1142: if (flg) {
1143: PetscCall(PCGetOperators(ksp->pc, &A, NULL));
1144: if (!ksp->transpose_solve) PetscCall(MatMatMult(A, X, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &R));
1145: else PetscCall(MatTransposeMatMult(A, X, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &R));
1146: PetscCall(MatAYPX(R, -1.0, B, SAME_NONZERO_PATTERN));
1147: PetscCall(MatGetSize(R, NULL, &N));
1148: PetscCall(PetscMalloc1(N, &norms));
1149: PetscCall(MatGetColumnNorms(R, NORM_2, norms));
1150: PetscCall(MatDestroy(&R));
1151: 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]));
1152: PetscCall(PetscFree(norms));
1153: }
1154: PetscFunctionReturn(PETSC_SUCCESS);
1155: }
1157: static PetscErrorCode KSPMatSolve_Private(KSP ksp, Mat B, Mat X)
1158: {
1159: Mat A, P, vB, vX;
1160: Vec cb, cx;
1161: PetscInt n1, N1, n2, N2, Bbn = PETSC_DECIDE;
1162: PetscBool match;
1164: PetscFunctionBegin;
1168: PetscCheckSameComm(ksp, 1, B, 2);
1169: PetscCheckSameComm(ksp, 1, X, 3);
1170: PetscCheckSameType(B, 2, X, 3);
1171: PetscCheck(B->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
1172: MatCheckPreallocated(X, 3);
1173: if (!X->assembled) {
1174: PetscCall(MatSetOption(X, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
1175: PetscCall(MatAssemblyBegin(X, MAT_FINAL_ASSEMBLY));
1176: PetscCall(MatAssemblyEnd(X, MAT_FINAL_ASSEMBLY));
1177: }
1178: PetscCheck(B != X, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_IDN, "B and X must be different matrices");
1179: PetscCheck(!ksp->transpose_solve || !ksp->transpose.use_explicittranspose, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "KSPMatSolveTranspose() does not support -ksp_use_explicittranspose");
1180: PetscCall(KSPGetOperators(ksp, &A, &P));
1181: PetscCall(MatGetLocalSize(B, NULL, &n2));
1182: PetscCall(MatGetLocalSize(X, NULL, &n1));
1183: PetscCall(MatGetSize(B, NULL, &N2));
1184: PetscCall(MatGetSize(X, NULL, &N1));
1185: 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);
1186: PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)B, &match, MATSEQDENSE, MATMPIDENSE, ""));
1187: PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of right-hand sides not stored in a dense Mat");
1188: PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)X, &match, MATSEQDENSE, MATMPIDENSE, ""));
1189: PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of solutions not stored in a dense Mat");
1190: PetscCall(KSPSetUp(ksp));
1191: PetscCall(KSPSetUpOnBlocks(ksp));
1192: if (ksp->ops->matsolve) {
1193: level++;
1194: if (ksp->guess_zero) PetscCall(MatZeroEntries(X));
1195: PetscCall(PetscLogEventBegin(!ksp->transpose_solve ? KSP_MatSolve : KSP_MatSolveTranspose, ksp, B, X, 0));
1196: PetscCall(KSPGetMatSolveBatchSize(ksp, &Bbn));
1197: /* by default, do a single solve with all columns */
1198: if (Bbn == PETSC_DECIDE) Bbn = N2;
1199: else PetscCheck(Bbn >= 1, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "KSPMatSolve() batch size %" PetscInt_FMT " must be positive", Bbn);
1200: PetscCall(PetscInfo(ksp, "KSP type %s solving using batches of width at most %" PetscInt_FMT "\n", ((PetscObject)ksp)->type_name, Bbn));
1201: /* if -ksp_matsolve_batch_size is greater than the actual number of columns, do a single solve with all columns */
1202: if (Bbn >= N2) {
1203: PetscUseTypeMethod(ksp, matsolve, B, X);
1204: if (ksp->viewFinalRes) PetscCall(KSPViewFinalMatResidual_Internal(ksp, B, X, ksp->viewerFinalRes, ksp->formatFinalRes, 0));
1206: PetscCall(KSPConvergedReasonViewFromOptions(ksp));
1208: if (ksp->viewRate) {
1209: PetscCall(PetscViewerPushFormat(ksp->viewerRate, PETSC_VIEWER_DEFAULT));
1210: PetscCall(KSPConvergedRateView(ksp, ksp->viewerRate));
1211: PetscCall(PetscViewerPopFormat(ksp->viewerRate));
1212: }
1213: } else {
1214: for (n2 = 0; n2 < N2; n2 += Bbn) {
1215: PetscCall(MatDenseGetSubMatrix(B, PETSC_DECIDE, PETSC_DECIDE, n2, PetscMin(n2 + Bbn, N2), &vB));
1216: PetscCall(MatDenseGetSubMatrix(X, PETSC_DECIDE, PETSC_DECIDE, n2, PetscMin(n2 + Bbn, N2), &vX));
1217: PetscUseTypeMethod(ksp, matsolve, vB, vX);
1218: if (ksp->viewFinalRes) PetscCall(KSPViewFinalMatResidual_Internal(ksp, vB, vX, ksp->viewerFinalRes, ksp->formatFinalRes, n2));
1220: PetscCall(KSPConvergedReasonViewFromOptions(ksp));
1222: if (ksp->viewRate) {
1223: PetscCall(PetscViewerPushFormat(ksp->viewerRate, PETSC_VIEWER_DEFAULT));
1224: PetscCall(KSPConvergedRateView(ksp, ksp->viewerRate));
1225: PetscCall(PetscViewerPopFormat(ksp->viewerRate));
1226: }
1227: PetscCall(MatDenseRestoreSubMatrix(B, &vB));
1228: PetscCall(MatDenseRestoreSubMatrix(X, &vX));
1229: }
1230: }
1231: if (ksp->viewMat) PetscCall(ObjectView((PetscObject)A, ksp->viewerMat, ksp->formatMat));
1232: if (ksp->viewPMat) PetscCall(ObjectView((PetscObject)P, ksp->viewerPMat, ksp->formatPMat));
1233: if (ksp->viewRhs) PetscCall(ObjectView((PetscObject)B, ksp->viewerRhs, ksp->formatRhs));
1234: if (ksp->viewSol) PetscCall(ObjectView((PetscObject)X, ksp->viewerSol, ksp->formatSol));
1235: if (ksp->view) PetscCall(KSPView(ksp, ksp->viewer));
1236: PetscCall(PetscLogEventEnd(!ksp->transpose_solve ? KSP_MatSolve : KSP_MatSolveTranspose, ksp, B, X, 0));
1237: if (ksp->errorifnotconverged && ksp->reason < 0 && (level == 1 || ksp->reason != KSP_DIVERGED_ITS)) {
1238: PCFailedReason reason;
1240: 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]);
1241: PetscCall(PCGetFailedReason(ksp->pc, &reason));
1242: 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]);
1243: }
1244: level--;
1245: } else {
1246: PetscCall(PetscInfo(ksp, "KSP type %s solving column by column\n", ((PetscObject)ksp)->type_name));
1247: for (n2 = 0; n2 < N2; ++n2) {
1248: PetscCall(MatDenseGetColumnVecRead(B, n2, &cb));
1249: PetscCall(MatDenseGetColumnVecWrite(X, n2, &cx));
1250: PetscCall(KSPSolve_Private(ksp, cb, cx));
1251: PetscCall(MatDenseRestoreColumnVecWrite(X, n2, &cx));
1252: PetscCall(MatDenseRestoreColumnVecRead(B, n2, &cb));
1253: }
1254: }
1255: PetscFunctionReturn(PETSC_SUCCESS);
1256: }
1258: /*@
1259: KSPMatSolve - Solves a linear system with multiple right-hand sides stored as a `MATDENSE`. Unlike `KSPSolve()`, `B` and `X` must be different matrices.
1261: Input Parameters:
1262: + ksp - iterative context
1263: - B - block of right-hand sides
1265: Output Parameter:
1266: . X - block of solutions
1268: Level: intermediate
1270: Note:
1271: This is a stripped-down version of `KSPSolve()`, which only handles `-ksp_view`, `-ksp_converged_reason`, `-ksp_converged_rate`, and `-ksp_view_final_residual`.
1273: .seealso: [](ch_ksp), `KSPSolve()`, `MatMatSolve()`, `KSPMatSolveTranspose()`, `MATDENSE`, `KSPHPDDM`, `PCBJACOBI`, `PCASM`
1274: @*/
1275: PetscErrorCode KSPMatSolve(KSP ksp, Mat B, Mat X)
1276: {
1277: PetscFunctionBegin;
1278: ksp->transpose_solve = PETSC_FALSE;
1279: PetscCall(KSPMatSolve_Private(ksp, B, X));
1280: PetscFunctionReturn(PETSC_SUCCESS);
1281: }
1283: /*@
1284: KSPMatSolveTranspose - Solves a linear system with the transposed matrix with multiple right-hand sides stored as a `MATDENSE`. Unlike `KSPSolveTranspose()`,
1285: `B` and `X` must be different matrices and the transposed matrix cannot be assembled explicitly for the user.
1287: Input Parameters:
1288: + ksp - iterative context
1289: - B - block of right-hand sides
1291: Output Parameter:
1292: . X - block of solutions
1294: Level: intermediate
1296: Note:
1297: This is a stripped-down version of `KSPSolveTranspose()`, which only handles `-ksp_view`, `-ksp_converged_reason`, `-ksp_converged_rate`, and `-ksp_view_final_residual`.
1299: .seealso: [](ch_ksp), `KSPSolveTranspose()`, `MatMatTransposeSolve()`, `KSPMatSolve()`, `MATDENSE`, `KSPHPDDM`, `PCBJACOBI`, `PCASM`
1300: @*/
1301: PetscErrorCode KSPMatSolveTranspose(KSP ksp, Mat B, Mat X)
1302: {
1303: PetscFunctionBegin;
1304: ksp->transpose_solve = PETSC_TRUE;
1305: PetscCall(KSPMatSolve_Private(ksp, B, X));
1306: PetscFunctionReturn(PETSC_SUCCESS);
1307: }
1309: /*@
1310: KSPSetMatSolveBatchSize - Sets the maximum number of columns treated simultaneously in `KSPMatSolve()`.
1312: Logically Collective
1314: Input Parameters:
1315: + ksp - iterative context
1316: - bs - batch size
1318: Level: advanced
1320: .seealso: [](ch_ksp), `KSPMatSolve()`, `KSPGetMatSolveBatchSize()`, `-mat_mumps_icntl_27`, `-matmatmult_Bbn`
1321: @*/
1322: PetscErrorCode KSPSetMatSolveBatchSize(KSP ksp, PetscInt bs)
1323: {
1324: PetscFunctionBegin;
1327: ksp->nmax = bs;
1328: PetscFunctionReturn(PETSC_SUCCESS);
1329: }
1331: /*@
1332: KSPGetMatSolveBatchSize - Gets the maximum number of columns treated simultaneously in `KSPMatSolve()`.
1334: Input Parameter:
1335: . ksp - iterative context
1337: Output Parameter:
1338: . bs - batch size
1340: Level: advanced
1342: .seealso: [](ch_ksp), `KSPMatSolve()`, `KSPSetMatSolveBatchSize()`, `-mat_mumps_icntl_27`, `-matmatmult_Bbn`
1343: @*/
1344: PetscErrorCode KSPGetMatSolveBatchSize(KSP ksp, PetscInt *bs)
1345: {
1346: PetscFunctionBegin;
1348: PetscAssertPointer(bs, 2);
1349: *bs = ksp->nmax;
1350: PetscFunctionReturn(PETSC_SUCCESS);
1351: }
1353: /*@
1354: KSPResetViewers - Resets all the viewers set from the options database during `KSPSetFromOptions()`
1356: Collective
1358: Input Parameter:
1359: . ksp - iterative context obtained from `KSPCreate()`
1361: Level: beginner
1363: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPSolve()`, `KSPSetFromOptions()`, `KSP`
1364: @*/
1365: PetscErrorCode KSPResetViewers(KSP ksp)
1366: {
1367: PetscFunctionBegin;
1369: if (!ksp) PetscFunctionReturn(PETSC_SUCCESS);
1370: PetscCall(PetscViewerDestroy(&ksp->viewer));
1371: PetscCall(PetscViewerDestroy(&ksp->viewerPre));
1372: PetscCall(PetscViewerDestroy(&ksp->viewerRate));
1373: PetscCall(PetscViewerDestroy(&ksp->viewerMat));
1374: PetscCall(PetscViewerDestroy(&ksp->viewerPMat));
1375: PetscCall(PetscViewerDestroy(&ksp->viewerRhs));
1376: PetscCall(PetscViewerDestroy(&ksp->viewerSol));
1377: PetscCall(PetscViewerDestroy(&ksp->viewerMatExp));
1378: PetscCall(PetscViewerDestroy(&ksp->viewerEV));
1379: PetscCall(PetscViewerDestroy(&ksp->viewerSV));
1380: PetscCall(PetscViewerDestroy(&ksp->viewerEVExp));
1381: PetscCall(PetscViewerDestroy(&ksp->viewerFinalRes));
1382: PetscCall(PetscViewerDestroy(&ksp->viewerPOpExp));
1383: PetscCall(PetscViewerDestroy(&ksp->viewerDScale));
1384: ksp->view = PETSC_FALSE;
1385: ksp->viewPre = PETSC_FALSE;
1386: ksp->viewMat = PETSC_FALSE;
1387: ksp->viewPMat = PETSC_FALSE;
1388: ksp->viewRhs = PETSC_FALSE;
1389: ksp->viewSol = PETSC_FALSE;
1390: ksp->viewMatExp = PETSC_FALSE;
1391: ksp->viewEV = PETSC_FALSE;
1392: ksp->viewSV = PETSC_FALSE;
1393: ksp->viewEVExp = PETSC_FALSE;
1394: ksp->viewFinalRes = PETSC_FALSE;
1395: ksp->viewPOpExp = PETSC_FALSE;
1396: ksp->viewDScale = PETSC_FALSE;
1397: PetscFunctionReturn(PETSC_SUCCESS);
1398: }
1400: /*@
1401: KSPReset - Resets a `KSP` context to the kspsetupcalled = 0 state and removes any allocated Vecs and Mats
1403: Collective
1405: Input Parameter:
1406: . ksp - iterative context obtained from `KSPCreate()`
1408: Level: beginner
1410: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPSolve()`, `KSP`
1411: @*/
1412: PetscErrorCode KSPReset(KSP ksp)
1413: {
1414: PetscFunctionBegin;
1416: if (!ksp) PetscFunctionReturn(PETSC_SUCCESS);
1417: PetscTryTypeMethod(ksp, reset);
1418: if (ksp->pc) PetscCall(PCReset(ksp->pc));
1419: if (ksp->guess) {
1420: KSPGuess guess = ksp->guess;
1421: PetscTryTypeMethod(guess, reset);
1422: }
1423: PetscCall(VecDestroyVecs(ksp->nwork, &ksp->work));
1424: PetscCall(VecDestroy(&ksp->vec_rhs));
1425: PetscCall(VecDestroy(&ksp->vec_sol));
1426: PetscCall(VecDestroy(&ksp->diagonal));
1427: PetscCall(VecDestroy(&ksp->truediagonal));
1429: PetscCall(KSPResetViewers(ksp));
1431: ksp->setupstage = KSP_SETUP_NEW;
1432: ksp->nmax = PETSC_DECIDE;
1433: PetscFunctionReturn(PETSC_SUCCESS);
1434: }
1436: /*@C
1437: KSPDestroy - Destroys a `KSP` context.
1439: Collective
1441: Input Parameter:
1442: . ksp - iterative context obtained from `KSPCreate()`
1444: Level: beginner
1446: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetUp()`, `KSPSolve()`, `KSP`
1447: @*/
1448: PetscErrorCode KSPDestroy(KSP *ksp)
1449: {
1450: PC pc;
1452: PetscFunctionBegin;
1453: if (!*ksp) PetscFunctionReturn(PETSC_SUCCESS);
1455: if (--((PetscObject)(*ksp))->refct > 0) {
1456: *ksp = NULL;
1457: PetscFunctionReturn(PETSC_SUCCESS);
1458: }
1460: PetscCall(PetscObjectSAWsViewOff((PetscObject)*ksp));
1462: /*
1463: Avoid a cascading call to PCReset(ksp->pc) from the following call:
1464: PCReset() shouldn't be called from KSPDestroy() as it is unprotected by pc's
1465: refcount (and may be shared, e.g., by other ksps).
1466: */
1467: pc = (*ksp)->pc;
1468: (*ksp)->pc = NULL;
1469: PetscCall(KSPReset((*ksp)));
1470: (*ksp)->pc = pc;
1471: PetscTryTypeMethod((*ksp), destroy);
1473: if ((*ksp)->transpose.use_explicittranspose) {
1474: PetscCall(MatDestroy(&(*ksp)->transpose.AT));
1475: PetscCall(MatDestroy(&(*ksp)->transpose.BT));
1476: (*ksp)->transpose.reuse_transpose = PETSC_FALSE;
1477: }
1479: PetscCall(KSPGuessDestroy(&(*ksp)->guess));
1480: PetscCall(DMDestroy(&(*ksp)->dm));
1481: PetscCall(PCDestroy(&(*ksp)->pc));
1482: PetscCall(PetscFree((*ksp)->res_hist_alloc));
1483: PetscCall(PetscFree((*ksp)->err_hist_alloc));
1484: if ((*ksp)->convergeddestroy) PetscCall((*(*ksp)->convergeddestroy)((*ksp)->cnvP));
1485: PetscCall(KSPMonitorCancel((*ksp)));
1486: PetscCall(KSPConvergedReasonViewCancel((*ksp)));
1487: PetscCall(PetscHeaderDestroy(ksp));
1488: PetscFunctionReturn(PETSC_SUCCESS);
1489: }
1491: /*@
1492: KSPSetPCSide - Sets the preconditioning side.
1494: Logically Collective
1496: Input Parameter:
1497: . ksp - iterative context obtained from `KSPCreate()`
1499: Output Parameter:
1500: . side - the preconditioning side, where side is one of
1501: .vb
1502: PC_LEFT - left preconditioning (default)
1503: PC_RIGHT - right preconditioning
1504: PC_SYMMETRIC - symmetric preconditioning
1505: .ve
1507: Options Database Key:
1508: . -ksp_pc_side <right,left,symmetric> - `KSP` preconditioner side
1510: Level: intermediate
1512: Notes:
1513: Left preconditioning is used by default for most Krylov methods except `KSPFGMRES` which only supports right preconditioning.
1515: For methods changing the side of the preconditioner changes the norm type that is used, see `KSPSetNormType()`.
1517: Symmetric preconditioning is currently available only for the `KSPQCG` method. However, note that
1518: symmetric preconditioning can be emulated by using either right or left
1519: preconditioning and a pre or post processing step.
1521: Setting the `PCSide` often affects the default norm type. See `KSPSetNormType()` for details.
1523: .seealso: [](ch_ksp), `KSPGetPCSide()`, `KSPSetNormType()`, `KSPGetNormType()`, `KSP`
1524: @*/
1525: PetscErrorCode KSPSetPCSide(KSP ksp, PCSide side)
1526: {
1527: PetscFunctionBegin;
1530: ksp->pc_side = ksp->pc_side_set = side;
1531: PetscFunctionReturn(PETSC_SUCCESS);
1532: }
1534: /*@
1535: KSPGetPCSide - Gets the preconditioning side.
1537: Not Collective
1539: Input Parameter:
1540: . ksp - iterative context obtained from `KSPCreate()`
1542: Output Parameter:
1543: . side - the preconditioning side, where side is one of
1544: .vb
1545: PC_LEFT - left preconditioning (default)
1546: PC_RIGHT - right preconditioning
1547: PC_SYMMETRIC - symmetric preconditioning
1548: .ve
1550: Level: intermediate
1552: .seealso: [](ch_ksp), `KSPSetPCSide()`, `KSP`
1553: @*/
1554: PetscErrorCode KSPGetPCSide(KSP ksp, PCSide *side)
1555: {
1556: PetscFunctionBegin;
1558: PetscAssertPointer(side, 2);
1559: PetscCall(KSPSetUpNorms_Private(ksp, PETSC_TRUE, &ksp->normtype, &ksp->pc_side));
1560: *side = ksp->pc_side;
1561: PetscFunctionReturn(PETSC_SUCCESS);
1562: }
1564: /*@
1565: KSPGetTolerances - Gets the relative, absolute, divergence, and maximum
1566: iteration tolerances used by the default `KSP` convergence tests.
1568: Not Collective
1570: Input Parameter:
1571: . ksp - the Krylov subspace context
1573: Output Parameters:
1574: + rtol - the relative convergence tolerance
1575: . abstol - the absolute convergence tolerance
1576: . dtol - the divergence tolerance
1577: - maxits - maximum number of iterations
1579: Level: intermediate
1581: Note:
1582: The user can specify `NULL` for any parameter that is not needed.
1584: .seealso: [](ch_ksp), `KSPSetTolerances()`, `KSP`, `KSPSetMinimumIterations()`, `KSPGetMinimumIterations()`
1585: @*/
1586: PetscErrorCode KSPGetTolerances(KSP ksp, PetscReal *rtol, PetscReal *abstol, PetscReal *dtol, PetscInt *maxits)
1587: {
1588: PetscFunctionBegin;
1590: if (abstol) *abstol = ksp->abstol;
1591: if (rtol) *rtol = ksp->rtol;
1592: if (dtol) *dtol = ksp->divtol;
1593: if (maxits) *maxits = ksp->max_it;
1594: PetscFunctionReturn(PETSC_SUCCESS);
1595: }
1597: /*@
1598: KSPSetTolerances - Sets the relative, absolute, divergence, and maximum
1599: iteration tolerances used by the default `KSP` convergence testers.
1601: Logically Collective
1603: Input Parameters:
1604: + ksp - the Krylov subspace context
1605: . rtol - the relative convergence tolerance, relative decrease in the (possibly preconditioned) residual norm
1606: . abstol - the absolute convergence tolerance absolute size of the (possibly preconditioned) residual norm
1607: . dtol - the divergence tolerance, amount (possibly preconditioned) residual norm can increase before `KSPConvergedDefault()` concludes that the method is diverging
1608: - maxits - maximum number of iterations to use
1610: Options Database Keys:
1611: + -ksp_atol <abstol> - Sets `abstol`
1612: . -ksp_rtol <rtol> - Sets `rtol`
1613: . -ksp_divtol <dtol> - Sets `dtol`
1614: - -ksp_max_it <maxits> - Sets `maxits`
1616: Level: intermediate
1618: Notes:
1619: Use `PETSC_DEFAULT` to retain the default value of any of the tolerances.
1621: See `KSPConvergedDefault()` for details how these parameters are used in the default convergence test. See also `KSPSetConvergenceTest()`
1622: for setting user-defined stopping criteria.
1624: .seealso: [](ch_ksp), `KSPGetTolerances()`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSP`, `KSPSetMinimumIterations()`
1625: @*/
1626: PetscErrorCode KSPSetTolerances(KSP ksp, PetscReal rtol, PetscReal abstol, PetscReal dtol, PetscInt maxits)
1627: {
1628: PetscFunctionBegin;
1635: if (rtol != (PetscReal)PETSC_DEFAULT) {
1636: 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);
1637: ksp->rtol = rtol;
1638: }
1639: if (abstol != (PetscReal)PETSC_DEFAULT) {
1640: PetscCheck(abstol >= 0.0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %g must be non-negative", (double)abstol);
1641: ksp->abstol = abstol;
1642: }
1643: if (dtol != (PetscReal)PETSC_DEFAULT) {
1644: PetscCheck(dtol >= 0.0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Divergence tolerance %g must be larger than 1.0", (double)dtol);
1645: ksp->divtol = dtol;
1646: }
1647: if (maxits != PETSC_DEFAULT) {
1648: PetscCheck(maxits >= 0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", maxits);
1649: ksp->max_it = maxits;
1650: }
1651: PetscFunctionReturn(PETSC_SUCCESS);
1652: }
1654: /*@
1655: KSPSetMinimumIterations - Sets the minimum number of iterations to use, regardless of the tolerances
1657: Logically Collective
1659: Input Parameters:
1660: + ksp - the Krylov subspace context
1661: - minit - minimum number of iterations to use
1663: Options Database Key:
1664: . -ksp_min_it <minits> - Sets `minit`
1666: Level: intermediate
1668: Notes:
1669: Use `KSPSetTolerances()` to set a variety of other tolerances
1671: See `KSPConvergedDefault()` for details on how these parameters are used in the default convergence test. See also `KSPSetConvergenceTest()`
1672: for setting user-defined stopping criteria.
1674: .seealso: [](ch_ksp), `KSPGetTolerances()`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSP`, `KSPSetTolerances()`, `KSPGetMinimumIterations()`
1675: @*/
1676: PetscErrorCode KSPSetMinimumIterations(KSP ksp, PetscInt minit)
1677: {
1678: PetscFunctionBegin;
1682: PetscCheck(minit >= 0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Minimum number of iterations %" PetscInt_FMT " must be non-negative", minit);
1683: ksp->min_it = minit;
1684: PetscFunctionReturn(PETSC_SUCCESS);
1685: }
1687: /*@
1688: KSPGetMinimumIterations - Gets the minimum number of iterations to use, regardless of the tolerances, that was set with `KSPSetMinimumIterations()` or `-ksp_min_it`
1690: Not Collective
1692: Input Parameter:
1693: . ksp - the Krylov subspace context
1695: Output Parameter:
1696: . minit - minimum number of iterations to use
1698: Level: intermediate
1700: .seealso: [](ch_ksp), `KSPGetTolerances()`, `KSPConvergedDefault()`, `KSPSetConvergenceTest()`, `KSP`, `KSPSetTolerances()`, `KSPSetMinimumIterations()`
1701: @*/
1702: PetscErrorCode KSPGetMinimumIterations(KSP ksp, PetscInt *minit)
1703: {
1704: PetscFunctionBegin;
1706: PetscAssertPointer(minit, 2);
1708: *minit = ksp->min_it;
1709: PetscFunctionReturn(PETSC_SUCCESS);
1710: }
1712: /*@
1713: KSPSetInitialGuessNonzero - Tells the iterative solver that the
1714: initial guess is nonzero; otherwise `KSP` assumes the initial guess
1715: is to be zero (and thus zeros it out before solving).
1717: Logically Collective
1719: Input Parameters:
1720: + ksp - iterative context obtained from `KSPCreate()`
1721: - flg - ``PETSC_TRUE`` indicates the guess is non-zero, `PETSC_FALSE` indicates the guess is zero
1723: Options Database Key:
1724: . -ksp_initial_guess_nonzero <true,false> - use nonzero initial guess
1726: Level: beginner
1728: Note:
1729: If this is not called the X vector is zeroed in the call to `KSPSolve()`.
1731: .seealso: [](ch_ksp), `KSPGetInitialGuessNonzero()`, `KSPGuessSetType()`, `KSPGuessType`, `KSP`
1732: @*/
1733: PetscErrorCode KSPSetInitialGuessNonzero(KSP ksp, PetscBool flg)
1734: {
1735: PetscFunctionBegin;
1738: ksp->guess_zero = (PetscBool) !(int)flg;
1739: PetscFunctionReturn(PETSC_SUCCESS);
1740: }
1742: /*@
1743: KSPGetInitialGuessNonzero - Determines whether the `KSP` solver is using
1744: a zero initial guess.
1746: Not Collective
1748: Input Parameter:
1749: . ksp - iterative context obtained from `KSPCreate()`
1751: Output Parameter:
1752: . flag - `PETSC_TRUE` if guess is nonzero, else `PETSC_FALSE`
1754: Level: intermediate
1756: .seealso: [](ch_ksp), `KSPSetInitialGuessNonzero()`, `KSP`
1757: @*/
1758: PetscErrorCode KSPGetInitialGuessNonzero(KSP ksp, PetscBool *flag)
1759: {
1760: PetscFunctionBegin;
1762: PetscAssertPointer(flag, 2);
1763: if (ksp->guess_zero) *flag = PETSC_FALSE;
1764: else *flag = PETSC_TRUE;
1765: PetscFunctionReturn(PETSC_SUCCESS);
1766: }
1768: /*@
1769: KSPSetErrorIfNotConverged - Causes `KSPSolve()` to generate an error if the solver has not converged as soon as the error is detected.
1771: Logically Collective
1773: Input Parameters:
1774: + ksp - iterative context obtained from `KSPCreate()`
1775: - flg - `PETSC_TRUE` indicates you want the error generated
1777: Options Database Key:
1778: . -ksp_error_if_not_converged <true,false> - generate an error and stop the program
1780: Level: intermediate
1782: Notes:
1783: Normally PETSc continues if a linear solver fails to converge, you can call `KSPGetConvergedReason()` after a `KSPSolve()`
1784: to determine if it has converged.
1786: A `KSP_DIVERGED_ITS` will not generate an error in a `KSPSolve()` inside a nested linear solver
1788: .seealso: [](ch_ksp), `KSPGetErrorIfNotConverged()`, `KSP`
1789: @*/
1790: PetscErrorCode KSPSetErrorIfNotConverged(KSP ksp, PetscBool flg)
1791: {
1792: PetscFunctionBegin;
1795: ksp->errorifnotconverged = flg;
1796: PetscFunctionReturn(PETSC_SUCCESS);
1797: }
1799: /*@
1800: KSPGetErrorIfNotConverged - Will `KSPSolve()` generate an error if the solver does not converge?
1802: Not Collective
1804: Input Parameter:
1805: . ksp - iterative context obtained from KSPCreate()
1807: Output Parameter:
1808: . flag - `PETSC_TRUE` if it will generate an error, else `PETSC_FALSE`
1810: Level: intermediate
1812: .seealso: [](ch_ksp), `KSPSetErrorIfNotConverged()`, `KSP`
1813: @*/
1814: PetscErrorCode KSPGetErrorIfNotConverged(KSP ksp, PetscBool *flag)
1815: {
1816: PetscFunctionBegin;
1818: PetscAssertPointer(flag, 2);
1819: *flag = ksp->errorifnotconverged;
1820: PetscFunctionReturn(PETSC_SUCCESS);
1821: }
1823: /*@
1824: KSPSetInitialGuessKnoll - Tells the iterative solver to use `PCApply()` to compute the initial guess (The Knoll trick)
1826: Logically Collective
1828: Input Parameters:
1829: + ksp - iterative context obtained from `KSPCreate()`
1830: - flg - `PETSC_TRUE` or `PETSC_FALSE`
1832: Level: advanced
1834: Developer Note:
1835: The Knoll trick is not currently implemented using the `KSPGuess` class
1837: .seealso: [](ch_ksp), `KSPGetInitialGuessKnoll()`, `KSPSetInitialGuessNonzero()`, `KSPGetInitialGuessNonzero()`, `KSP`
1838: @*/
1839: PetscErrorCode KSPSetInitialGuessKnoll(KSP ksp, PetscBool flg)
1840: {
1841: PetscFunctionBegin;
1844: ksp->guess_knoll = flg;
1845: PetscFunctionReturn(PETSC_SUCCESS);
1846: }
1848: /*@
1849: KSPGetInitialGuessKnoll - Determines whether the `KSP` solver is using the Knoll trick (using PCApply(pc,b,...) to compute
1850: the initial guess
1852: Not Collective
1854: Input Parameter:
1855: . ksp - iterative context obtained from `KSPCreate()`
1857: Output Parameter:
1858: . flag - `PETSC_TRUE` if using Knoll trick, else `PETSC_FALSE`
1860: Level: advanced
1862: .seealso: [](ch_ksp), `KSPSetInitialGuessKnoll()`, `KSPSetInitialGuessNonzero()`, `KSPGetInitialGuessNonzero()`, `KSP`
1863: @*/
1864: PetscErrorCode KSPGetInitialGuessKnoll(KSP ksp, PetscBool *flag)
1865: {
1866: PetscFunctionBegin;
1868: PetscAssertPointer(flag, 2);
1869: *flag = ksp->guess_knoll;
1870: PetscFunctionReturn(PETSC_SUCCESS);
1871: }
1873: /*@
1874: KSPGetComputeSingularValues - Gets the flag indicating whether the extreme singular
1875: values will be calculated via a Lanczos or Arnoldi process as the linear
1876: system is solved.
1878: Not Collective
1880: Input Parameter:
1881: . ksp - iterative context obtained from `KSPCreate()`
1883: Output Parameter:
1884: . flg - `PETSC_TRUE` or `PETSC_FALSE`
1886: Options Database Key:
1887: . -ksp_monitor_singular_value - Activates `KSPSetComputeSingularValues()`
1889: Level: advanced
1891: Notes:
1892: Currently this option is not valid for all iterative methods.
1894: Many users may just want to use the monitoring routine
1895: `KSPMonitorSingularValue()` (which can be set with option -ksp_monitor_singular_value)
1896: to print the singular values at each iteration of the linear solve.
1898: .seealso: [](ch_ksp), `KSPComputeExtremeSingularValues()`, `KSPMonitorSingularValue()`, `KSP`
1899: @*/
1900: PetscErrorCode KSPGetComputeSingularValues(KSP ksp, PetscBool *flg)
1901: {
1902: PetscFunctionBegin;
1904: PetscAssertPointer(flg, 2);
1905: *flg = ksp->calc_sings;
1906: PetscFunctionReturn(PETSC_SUCCESS);
1907: }
1909: /*@
1910: KSPSetComputeSingularValues - Sets a flag so that the extreme singular
1911: values will be calculated via a Lanczos or Arnoldi process as the linear
1912: system is solved.
1914: Logically Collective
1916: Input Parameters:
1917: + ksp - iterative context obtained from `KSPCreate()`
1918: - flg - `PETSC_TRUE` or `PETSC_FALSE`
1920: Options Database Key:
1921: . -ksp_monitor_singular_value - Activates `KSPSetComputeSingularValues()`
1923: Level: advanced
1925: Notes:
1926: Currently this option is not valid for all iterative methods.
1928: Many users may just want to use the monitoring routine
1929: `KSPMonitorSingularValue()` (which can be set with option -ksp_monitor_singular_value)
1930: to print the singular values at each iteration of the linear solve.
1932: .seealso: [](ch_ksp), `KSPComputeExtremeSingularValues()`, `KSPMonitorSingularValue()`, `KSP`, `KSPSetComputeRitz()`
1933: @*/
1934: PetscErrorCode KSPSetComputeSingularValues(KSP ksp, PetscBool flg)
1935: {
1936: PetscFunctionBegin;
1939: ksp->calc_sings = flg;
1940: PetscFunctionReturn(PETSC_SUCCESS);
1941: }
1943: /*@
1944: KSPGetComputeEigenvalues - Gets the flag indicating that the extreme eigenvalues
1945: values will be calculated via a Lanczos or Arnoldi process as the linear
1946: system is solved.
1948: Not Collective
1950: Input Parameter:
1951: . ksp - iterative context obtained from `KSPCreate()`
1953: Output Parameter:
1954: . flg - `PETSC_TRUE` or `PETSC_FALSE`
1956: Level: advanced
1958: Note:
1959: Currently this option is not valid for all iterative methods.
1961: .seealso: [](ch_ksp), `KSPComputeEigenvalues()`, `KSPComputeEigenvaluesExplicitly()`, `KSP`, `KSPSetComputeRitz()`
1962: @*/
1963: PetscErrorCode KSPGetComputeEigenvalues(KSP ksp, PetscBool *flg)
1964: {
1965: PetscFunctionBegin;
1967: PetscAssertPointer(flg, 2);
1968: *flg = ksp->calc_sings;
1969: PetscFunctionReturn(PETSC_SUCCESS);
1970: }
1972: /*@
1973: KSPSetComputeEigenvalues - Sets a flag so that the extreme eigenvalues
1974: values will be calculated via a Lanczos or Arnoldi process as the linear
1975: system is solved.
1977: Logically Collective
1979: Input Parameters:
1980: + ksp - iterative context obtained from `KSPCreate()`
1981: - flg - `PETSC_TRUE` or `PETSC_FALSE`
1983: Level: advanced
1985: Note:
1986: Currently this option is not valid for all iterative methods.
1988: .seealso: [](ch_ksp), `KSPComputeEigenvalues()`, `KSPComputeEigenvaluesExplicitly()`, `KSP`, `KSPSetComputeRitz()`
1989: @*/
1990: PetscErrorCode KSPSetComputeEigenvalues(KSP ksp, PetscBool flg)
1991: {
1992: PetscFunctionBegin;
1995: ksp->calc_sings = flg;
1996: PetscFunctionReturn(PETSC_SUCCESS);
1997: }
1999: /*@
2000: KSPSetComputeRitz - Sets a flag so that the Ritz or harmonic Ritz pairs
2001: will be calculated via a Lanczos or Arnoldi process as the linear
2002: system is solved.
2004: Logically Collective
2006: Input Parameters:
2007: + ksp - iterative context obtained from `KSPCreate()`
2008: - flg - `PETSC_TRUE` or `PETSC_FALSE`
2010: Level: advanced
2012: Note:
2013: Currently this option is only valid for the `KSPGMRES` method.
2015: .seealso: [](ch_ksp), `KSPComputeRitz()`, `KSP`
2016: @*/
2017: PetscErrorCode KSPSetComputeRitz(KSP ksp, PetscBool flg)
2018: {
2019: PetscFunctionBegin;
2022: ksp->calc_ritz = flg;
2023: PetscFunctionReturn(PETSC_SUCCESS);
2024: }
2026: /*@
2027: KSPGetRhs - Gets the right-hand-side vector for the linear system to
2028: be solved.
2030: Not Collective
2032: Input Parameter:
2033: . ksp - iterative context obtained from `KSPCreate()`
2035: Output Parameter:
2036: . r - right-hand-side vector
2038: Level: developer
2040: .seealso: [](ch_ksp), `KSPGetSolution()`, `KSPSolve()`, `KSP`
2041: @*/
2042: PetscErrorCode KSPGetRhs(KSP ksp, Vec *r)
2043: {
2044: PetscFunctionBegin;
2046: PetscAssertPointer(r, 2);
2047: *r = ksp->vec_rhs;
2048: PetscFunctionReturn(PETSC_SUCCESS);
2049: }
2051: /*@
2052: KSPGetSolution - Gets the location of the solution for the
2053: linear system to be solved. Note that this may not be where the solution
2054: is stored during the iterative process; see `KSPBuildSolution()`.
2056: Not Collective
2058: Input Parameter:
2059: . ksp - iterative context obtained from `KSPCreate()`
2061: Output Parameter:
2062: . v - solution vector
2064: Level: developer
2066: .seealso: [](ch_ksp), `KSPGetRhs()`, `KSPBuildSolution()`, `KSPSolve()`, `KSP`
2067: @*/
2068: PetscErrorCode KSPGetSolution(KSP ksp, Vec *v)
2069: {
2070: PetscFunctionBegin;
2072: PetscAssertPointer(v, 2);
2073: *v = ksp->vec_sol;
2074: PetscFunctionReturn(PETSC_SUCCESS);
2075: }
2077: /*@
2078: KSPSetPC - Sets the preconditioner to be used to calculate the
2079: application of the preconditioner on a vector.
2081: Collective
2083: Input Parameters:
2084: + ksp - iterative context obtained from `KSPCreate()`
2085: - pc - the preconditioner object (can be `NULL`)
2087: Level: developer
2089: Note:
2090: Use `KSPGetPC()` to retrieve the preconditioner context.
2092: .seealso: [](ch_ksp), `KSPGetPC()`, `KSP`
2093: @*/
2094: PetscErrorCode KSPSetPC(KSP ksp, PC pc)
2095: {
2096: PetscFunctionBegin;
2098: if (pc) {
2100: PetscCheckSameComm(ksp, 1, pc, 2);
2101: }
2102: if (ksp->pc != pc && ksp->setupstage) ksp->setupstage = KSP_SETUP_NEWMATRIX;
2103: PetscCall(PetscObjectReference((PetscObject)pc));
2104: PetscCall(PCDestroy(&ksp->pc));
2105: ksp->pc = pc;
2106: PetscFunctionReturn(PETSC_SUCCESS);
2107: }
2109: PETSC_INTERN PetscErrorCode PCCreate_MPI(PC);
2111: // PetscClangLinter pragma disable: -fdoc-internal-linkage
2112: /*@C
2113: KSPCheckPCMPI - Checks if `-mpi_linear_solver_server` is active and the `PC` should be changed to `PCMPI`
2115: Collective
2117: Input Parameter:
2118: . ksp - iterative context obtained from `KSPCreate()`
2120: Level: developer
2122: .seealso: [](ch_ksp), `KSPSetPC()`, `KSP`, `PCMPIServerBegin()`, `PCMPIServerEnd()`
2123: @*/
2124: PETSC_INTERN PetscErrorCode KSPCheckPCMPI(KSP ksp)
2125: {
2126: PetscBool isPCMPI;
2128: PetscFunctionBegin;
2130: PetscCall(PetscObjectTypeCompare((PetscObject)ksp->pc, PCMPI, &isPCMPI));
2131: if (PCMPIServerActive && ksp->nestlevel == 0 && !isPCMPI) {
2132: const char *prefix;
2133: char *found = NULL;
2135: PetscCall(KSPGetOptionsPrefix(ksp, &prefix));
2136: if (prefix) PetscCall(PetscStrstr(prefix, "mpi_linear_solver_server_", &found));
2137: if (!found) PetscCall(KSPAppendOptionsPrefix(ksp, "mpi_linear_solver_server_"));
2138: PetscCall(PCSetType(ksp->pc, PCMPI));
2139: }
2140: PetscFunctionReturn(PETSC_SUCCESS);
2141: }
2143: /*@
2144: KSPGetPC - Returns a pointer to the preconditioner context
2145: set with `KSPSetPC()`.
2147: Not Collective
2149: Input Parameter:
2150: . ksp - iterative context obtained from `KSPCreate()`
2152: Output Parameter:
2153: . pc - preconditioner context
2155: Level: developer
2157: .seealso: [](ch_ksp), `KSPSetPC()`, `KSP`, `PC`
2158: @*/
2159: PetscErrorCode KSPGetPC(KSP ksp, PC *pc)
2160: {
2161: PetscFunctionBegin;
2163: PetscAssertPointer(pc, 2);
2164: if (!ksp->pc) {
2165: PetscCall(PCCreate(PetscObjectComm((PetscObject)ksp), &ksp->pc));
2166: PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp->pc, (PetscObject)ksp, 0));
2167: PetscCall(PetscObjectSetOptions((PetscObject)ksp->pc, ((PetscObject)ksp)->options));
2168: PetscCall(PCSetKSPNestLevel(ksp->pc, ksp->nestlevel));
2169: }
2170: PetscCall(KSPCheckPCMPI(ksp));
2171: *pc = ksp->pc;
2172: PetscFunctionReturn(PETSC_SUCCESS);
2173: }
2175: /*@
2176: KSPMonitor - runs the user provided monitor routines, if they exist
2178: Collective
2180: Input Parameters:
2181: + ksp - iterative context obtained from `KSPCreate()`
2182: . it - iteration number
2183: - rnorm - relative norm of the residual
2185: Level: developer
2187: Note:
2188: This routine is called by the `KSP` implementations.
2189: It does not typically need to be called by the user.
2191: .seealso: [](ch_ksp), `KSPMonitorSet()`
2192: @*/
2193: PetscErrorCode KSPMonitor(KSP ksp, PetscInt it, PetscReal rnorm)
2194: {
2195: PetscInt i, n = ksp->numbermonitors;
2197: PetscFunctionBegin;
2198: for (i = 0; i < n; i++) PetscCall((*ksp->monitor[i])(ksp, it, rnorm, ksp->monitorcontext[i]));
2199: PetscFunctionReturn(PETSC_SUCCESS);
2200: }
2202: /*@C
2203: KSPMonitorSet - Sets an ADDITIONAL function to be called at every iteration to monitor
2204: the residual/error etc.
2206: Logically Collective
2208: Input Parameters:
2209: + ksp - iterative context obtained from `KSPCreate()`
2210: . monitor - pointer to function (if this is `NULL`, it turns off monitoring
2211: . ctx - [optional] context for private data for the monitor routine (use `NULL` if no context is needed)
2212: - monitordestroy - [optional] routine that frees monitor context (may be `NULL`)
2214: Calling sequence of `monitor`:
2215: + ksp - iterative context obtained from `KSPCreate()`
2216: . it - iteration number
2217: . rnorm - (estimated) 2-norm of (preconditioned) residual
2218: - ctx - optional monitoring context, as set by `KSPMonitorSet()`
2220: Calling sequence of `monitordestroy`:
2221: . ctx - optional monitoring context, as set by `KSPMonitorSet()`
2223: Options Database Keys:
2224: + -ksp_monitor - sets `KSPMonitorResidual()`
2225: . -ksp_monitor draw - sets `KSPMonitorResidualDraw()` and plots residual
2226: . -ksp_monitor draw::draw_lg - sets `KSPMonitorResidualDrawLG()` and plots residual
2227: . -ksp_monitor_pause_final - Pauses any graphics when the solve finishes (only works for internal monitors)
2228: . -ksp_monitor_true_residual - sets `KSPMonitorTrueResidual()`
2229: . -ksp_monitor_true_residual draw::draw_lg - sets `KSPMonitorTrueResidualDrawLG()` and plots residual
2230: . -ksp_monitor_max - sets `KSPMonitorTrueResidualMax()`
2231: . -ksp_monitor_singular_value - sets `KSPMonitorSingularValue()`
2232: - -ksp_monitor_cancel - cancels all monitors that have been hardwired into a code by calls to `KSPMonitorSet()`, but
2233: does not cancel those set via the options database.
2235: Level: beginner
2237: Notes:
2238: The default is to do nothing. To print the residual, or preconditioned
2239: residual if `KSPSetNormType`(ksp,`KSP_NORM_PRECONDITIONED`) was called, use
2240: `KSPMonitorResidual()` as the monitoring routine, with a `PETSCVIEWERASCII` as the
2241: context.
2243: Several different monitoring routines may be set by calling
2244: `KSPMonitorSet()` multiple times; all will be called in the
2245: order in which they were set.
2247: Fortran Note:
2248: Only a single monitor function can be set for each `KSP` object
2250: .seealso: [](ch_ksp), `KSPMonitorResidual()`, `KSPMonitorCancel()`, `KSP`
2251: @*/
2252: PetscErrorCode KSPMonitorSet(KSP ksp, PetscErrorCode (*monitor)(KSP ksp, PetscInt it, PetscReal rnorm, void *ctx), void *ctx, PetscErrorCode (*monitordestroy)(void **ctx))
2253: {
2254: PetscInt i;
2255: PetscBool identical;
2257: PetscFunctionBegin;
2259: for (i = 0; i < ksp->numbermonitors; i++) {
2260: PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))monitor, ctx, monitordestroy, (PetscErrorCode(*)(void))ksp->monitor[i], ksp->monitorcontext[i], ksp->monitordestroy[i], &identical));
2261: if (identical) PetscFunctionReturn(PETSC_SUCCESS);
2262: }
2263: PetscCheck(ksp->numbermonitors < MAXKSPMONITORS, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "Too many KSP monitors set");
2264: ksp->monitor[ksp->numbermonitors] = monitor;
2265: ksp->monitordestroy[ksp->numbermonitors] = monitordestroy;
2266: ksp->monitorcontext[ksp->numbermonitors++] = (void *)ctx;
2267: PetscFunctionReturn(PETSC_SUCCESS);
2268: }
2270: /*@
2271: KSPMonitorCancel - Clears all monitors for a `KSP` object.
2273: Logically Collective
2275: Input Parameter:
2276: . ksp - iterative context obtained from `KSPCreate()`
2278: Options Database Key:
2279: . -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.
2281: Level: intermediate
2283: .seealso: [](ch_ksp), `KSPMonitorResidual()`, `KSPMonitorSet()`, `KSP`
2284: @*/
2285: PetscErrorCode KSPMonitorCancel(KSP ksp)
2286: {
2287: PetscInt i;
2289: PetscFunctionBegin;
2291: for (i = 0; i < ksp->numbermonitors; i++) {
2292: if (ksp->monitordestroy[i]) PetscCall((*ksp->monitordestroy[i])(&ksp->monitorcontext[i]));
2293: }
2294: ksp->numbermonitors = 0;
2295: PetscFunctionReturn(PETSC_SUCCESS);
2296: }
2298: /*@C
2299: KSPGetMonitorContext - Gets the monitoring context, as set by `KSPMonitorSet()` for the FIRST monitor only.
2301: Not Collective
2303: Input Parameter:
2304: . ksp - iterative context obtained from `KSPCreate()`
2306: Output Parameter:
2307: . ctx - monitoring context
2309: Level: intermediate
2311: .seealso: [](ch_ksp), `KSPMonitorResidual()`, `KSP`
2312: @*/
2313: PetscErrorCode KSPGetMonitorContext(KSP ksp, void *ctx)
2314: {
2315: PetscFunctionBegin;
2317: *(void **)ctx = ksp->monitorcontext[0];
2318: PetscFunctionReturn(PETSC_SUCCESS);
2319: }
2321: /*@
2322: KSPSetResidualHistory - Sets the array used to hold the residual history.
2323: If set, this array will contain the residual norms computed at each
2324: iteration of the solver.
2326: Not Collective
2328: Input Parameters:
2329: + ksp - iterative context obtained from `KSPCreate()`
2330: . a - array to hold history
2331: . na - size of `a`
2332: - reset - `PETSC_TRUE` indicates the history counter is reset to zero
2333: for each new linear solve
2335: Level: advanced
2337: Notes:
2338: 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.
2339: If 'a' is `NULL` then space is allocated for the history. If 'na' `PETSC_DECIDE` or `PETSC_DEFAULT` then a
2340: default array of length 10000 is allocated.
2342: If the array is not long enough then once the iterations is longer than the array length `KSPSolve()` stops recording the history
2344: .seealso: [](ch_ksp), `KSPGetResidualHistory()`, `KSP`
2345: @*/
2346: PetscErrorCode KSPSetResidualHistory(KSP ksp, PetscReal a[], PetscInt na, PetscBool reset)
2347: {
2348: PetscFunctionBegin;
2351: PetscCall(PetscFree(ksp->res_hist_alloc));
2352: if (na != PETSC_DECIDE && na != PETSC_DEFAULT && a) {
2353: ksp->res_hist = a;
2354: ksp->res_hist_max = (size_t)na;
2355: } else {
2356: if (na != PETSC_DECIDE && na != PETSC_DEFAULT) ksp->res_hist_max = (size_t)na;
2357: else ksp->res_hist_max = 10000; /* like default ksp->max_it */
2358: PetscCall(PetscCalloc1(ksp->res_hist_max, &ksp->res_hist_alloc));
2360: ksp->res_hist = ksp->res_hist_alloc;
2361: }
2362: ksp->res_hist_len = 0;
2363: ksp->res_hist_reset = reset;
2364: PetscFunctionReturn(PETSC_SUCCESS);
2365: }
2367: /*@C
2368: KSPGetResidualHistory - Gets the array used to hold the residual history and the number of residuals it contains.
2370: Not Collective
2372: Input Parameter:
2373: . ksp - iterative context obtained from `KSPCreate()`
2375: Output Parameters:
2376: + a - pointer to array to hold history (or `NULL`)
2377: - na - number of used entries in a (or `NULL`). Note this has different meanings depending on the `reset` argument to `KSPSetResidualHistory()`
2379: Level: advanced
2381: Note:
2382: This array is borrowed and should not be freed by the caller.
2384: Can only be called after a `KSPSetResidualHistory()` otherwise `a` and `na` are set to `NULL` and zero
2386: 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
2387: returned with `KSPGetIterationNumber()`.
2389: Some Krylov methods may not compute the final residual norm when convergence is declared because the maximum number of iterations allowed has been reached.
2390: In this situation, when `reset` was `PETSC_TRUE`, `na` will then equal the number of iterations reported with `KSPGetIterationNumber()`
2392: Some Krylov methods (such as `KSPSTCG`), under certain circumstances, do not compute the final residual norm. In this situation, when `reset` was `PETSC_TRUE`,
2393: `na` will then equal the number of iterations reported with `KSPGetIterationNumber()`
2395: `KSPBCGSL` does not record the residual norms for the "subiterations" hence the results from `KSPGetResidualHistory()` and `KSPGetIterationNumber()` will be different
2397: Fortran Note:
2398: The Fortran version of this routine has a calling sequence
2399: .vb
2400: call KSPGetResidualHistory(KSP ksp, integer na, integer ierr)
2401: .ve
2402: note that you have passed a Fortran array into `KSPSetResidualHistory()` and you need
2403: to access the residual values from this Fortran array you provided. Only the `na` (number of
2404: residual norms currently held) is set.
2406: .seealso: [](ch_ksp), `KSPSetResidualHistory()`, `KSP`, `KSPGetIterationNumber()`, `KSPSTCG`, `KSPBCGSL`
2407: @*/
2408: PetscErrorCode KSPGetResidualHistory(KSP ksp, const PetscReal *a[], PetscInt *na)
2409: {
2410: PetscFunctionBegin;
2412: if (a) *a = ksp->res_hist;
2413: if (na) *na = (PetscInt)ksp->res_hist_len;
2414: PetscFunctionReturn(PETSC_SUCCESS);
2415: }
2417: /*@
2418: 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.
2420: Not Collective
2422: Input Parameters:
2423: + ksp - iterative context obtained from `KSPCreate()`
2424: . a - array to hold history
2425: . na - size of `a`
2426: - reset - `PETSC_TRUE` indicates the history counter is reset to zero for each new linear solve
2428: Level: advanced
2430: Notes:
2431: 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.
2432: 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.
2434: If the array is not long enough then once the iterations is longer than the array length `KSPSolve()` stops recording the history
2436: .seealso: [](ch_ksp), `KSPGetErrorHistory()`, `KSPSetResidualHistory()`, `KSP`
2437: @*/
2438: PetscErrorCode KSPSetErrorHistory(KSP ksp, PetscReal a[], PetscInt na, PetscBool reset)
2439: {
2440: PetscFunctionBegin;
2443: PetscCall(PetscFree(ksp->err_hist_alloc));
2444: if (na != PETSC_DECIDE && na != PETSC_DEFAULT && a) {
2445: ksp->err_hist = a;
2446: ksp->err_hist_max = (size_t)na;
2447: } else {
2448: if (na != PETSC_DECIDE && na != PETSC_DEFAULT) ksp->err_hist_max = (size_t)na;
2449: else ksp->err_hist_max = 10000; /* like default ksp->max_it */
2450: PetscCall(PetscCalloc1(ksp->err_hist_max, &ksp->err_hist_alloc));
2452: ksp->err_hist = ksp->err_hist_alloc;
2453: }
2454: ksp->err_hist_len = 0;
2455: ksp->err_hist_reset = reset;
2456: PetscFunctionReturn(PETSC_SUCCESS);
2457: }
2459: /*@C
2460: KSPGetErrorHistory - Gets the array used to hold the error history and the number of residuals it contains.
2462: Not Collective
2464: Input Parameter:
2465: . ksp - iterative context obtained from `KSPCreate()`
2467: Output Parameters:
2468: + a - pointer to array to hold history (or `NULL`)
2469: - na - number of used entries in a (or `NULL`)
2471: Level: advanced
2473: Note:
2474: This array is borrowed and should not be freed by the caller.
2475: Can only be called after a `KSPSetErrorHistory()` otherwise `a` and `na` are set to `NULL` and zero
2477: Fortran Note:
2478: The Fortran version of this routine has a calling sequence
2479: $ call KSPGetErrorHistory(KSP ksp, integer na, integer ierr)
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
2941: - ctx - optional context
2943: Calling sequence of `func`:
2944: + ksp - the `KSP` context
2945: . A - the linear operator
2946: . B - the matrix from which the preconditioner is built, often `A`
2947: - ctx - optional user-provided context
2949: Level: beginner
2951: Notes:
2952: `func()` will be called automatically at the very next call to `KSPSolve()`. It will NOT be called at future `KSPSolve()` calls
2953: unless either `KSPSetComputeOperators()` or `KSPSetOperators()` is called before that `KSPSolve()` is called. This allows the same system to be solved several times
2954: with different right hand side functions but is a confusing API since one might expect it to be called for each `KSPSolve()`
2956: To reuse the same preconditioner for the next `KSPSolve()` and not compute a new one based on the most recently computed matrix call `KSPSetReusePreconditioner()`
2958: Developer Note:
2959: Perhaps this routine and `KSPSetComputeRHS()` could be combined into a new API that makes clear when new matrices are computing without requiring call this
2960: routine to indicate when the new matrix should be computed.
2962: .seealso: [](ch_ksp), `KSP`, `KSPSetOperators()`, `KSPSetComputeRHS()`, `DMKSPSetComputeOperators()`, `KSPSetComputeInitialGuess()`
2963: @*/
2964: PetscErrorCode KSPSetComputeOperators(KSP ksp, PetscErrorCode (*func)(KSP ksp, Mat A, Mat B, void *ctx), void *ctx)
2965: {
2966: DM dm;
2968: PetscFunctionBegin;
2970: PetscCall(KSPGetDM(ksp, &dm));
2971: PetscCall(DMKSPSetComputeOperators(dm, func, ctx));
2972: if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX;
2973: PetscFunctionReturn(PETSC_SUCCESS);
2974: }
2976: /*@C
2977: KSPSetComputeRHS - set routine to compute the right hand side of the linear system
2979: Logically Collective
2981: Input Parameters:
2982: + ksp - the `KSP` context
2983: . func - function to compute the right hand side
2984: - ctx - optional context
2986: Calling sequence of `func`:
2987: + ksp - the `KSP` context
2988: . b - right hand side of linear system
2989: - ctx - optional user-provided context
2991: Level: beginner
2993: Note:
2994: The routine you provide will be called EACH you call `KSPSolve()` to prepare the new right hand side for that solve
2996: .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `DMKSPSetComputeRHS()`, `KSPSetComputeOperators()`, `KSPSetOperators()`
2997: @*/
2998: PetscErrorCode KSPSetComputeRHS(KSP ksp, PetscErrorCode (*func)(KSP ksp, Vec b, void *ctx), void *ctx)
2999: {
3000: DM dm;
3002: PetscFunctionBegin;
3004: PetscCall(KSPGetDM(ksp, &dm));
3005: PetscCall(DMKSPSetComputeRHS(dm, func, ctx));
3006: PetscFunctionReturn(PETSC_SUCCESS);
3007: }
3009: /*@C
3010: KSPSetComputeInitialGuess - set routine to compute the initial guess of the linear system
3012: Logically Collective
3014: Input Parameters:
3015: + ksp - the `KSP` context
3016: . func - function to compute the initial guess
3017: - ctx - optional context
3019: Calling sequence of `func`:
3020: + ksp - the `KSP` context
3021: . x - solution vector
3022: - ctx - optional user-provided context
3024: Level: beginner
3026: Note:
3027: This should only be used in conjunction with `KSPSetComputeRHS()` and `KSPSetComputeOperators()`, otherwise
3028: call `KSPSetInitialGuessNonzero()` and set the initial guess values in the solution vector passed to `KSPSolve()` before calling the solver
3030: .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `KSPSetComputeRHS()`, `KSPSetComputeOperators()`, `DMKSPSetComputeInitialGuess()`, `KSPSetInitialGuessNonzero()`
3031: @*/
3032: PetscErrorCode KSPSetComputeInitialGuess(KSP ksp, PetscErrorCode (*func)(KSP ksp, Vec x, void *ctx), void *ctx)
3033: {
3034: DM dm;
3036: PetscFunctionBegin;
3038: PetscCall(KSPGetDM(ksp, &dm));
3039: PetscCall(DMKSPSetComputeInitialGuess(dm, func, ctx));
3040: PetscFunctionReturn(PETSC_SUCCESS);
3041: }
3043: /*@
3044: KSPSetUseExplicitTranspose - Determines the explicit transpose of the operator is formed in `KSPSolveTranspose()`. In some configurations (like GPUs) it may
3045: be explicitly formed since the solve is much more efficient.
3047: Logically Collective
3049: Input Parameter:
3050: . ksp - the `KSP` context
3052: Output Parameter:
3053: . flg - `PETSC_TRUE` to transpose the system in `KSPSolveTranspose()`, `PETSC_FALSE` to not transpose (default)
3055: Level: advanced
3057: .seealso: [](ch_ksp), `KSPSolveTranspose()`, `KSP`
3058: @*/
3059: PetscErrorCode KSPSetUseExplicitTranspose(KSP ksp, PetscBool flg)
3060: {
3061: PetscFunctionBegin;
3064: ksp->transpose.use_explicittranspose = flg;
3065: PetscFunctionReturn(PETSC_SUCCESS);
3066: }