xref: /petsc/src/ksp/ksp/interface/itcreate.c (revision f5d0f301a17fda14f5f53985d938d1a5d5001cde)
1 /*
2      The basic KSP routines, Create, View etc. are here.
3 */
4 #include <petsc/private/kspimpl.h> /*I "petscksp.h" I*/
5 
6 /* Logging support */
7 PetscClassId  KSP_CLASSID;
8 PetscClassId  DMKSP_CLASSID;
9 PetscClassId  KSPGUESS_CLASSID;
10 PetscLogEvent KSP_GMRESOrthogonalization, KSP_SetUp, KSP_Solve, KSP_SolveTranspose, KSP_MatSolve, KSP_MatSolveTranspose;
11 
12 /*
13    Contains the list of registered KSP routines
14 */
15 PetscFunctionList KSPList              = NULL;
16 PetscBool         KSPRegisterAllCalled = PETSC_FALSE;
17 
18 /*
19    Contains the list of registered KSP monitors
20 */
21 PetscFunctionList KSPMonitorList              = NULL;
22 PetscFunctionList KSPMonitorCreateList        = NULL;
23 PetscFunctionList KSPMonitorDestroyList       = NULL;
24 PetscBool         KSPMonitorRegisterAllCalled = PETSC_FALSE;
25 
26 /*@
27   KSPLoad - Loads a `KSP` that has been stored in a `PETSCVIEWERBINARY`  with `KSPView()`.
28 
29   Collective
30 
31   Input Parameters:
32 + newdm  - the newly loaded `KSP`, this needs to have been created with `KSPCreate()` or
33            some related function before a call to `KSPLoad()`.
34 - viewer - binary file viewer, obtained from `PetscViewerBinaryOpen()`
35 
36   Level: intermediate
37 
38   Note:
39   The type is determined by the data in the file, any type set into the `KSP` before this call is ignored.
40 
41 .seealso: [](ch_ksp), `KSP`, `PetscViewerBinaryOpen()`, `KSPView()`, `MatLoad()`, `VecLoad()`
42 @*/
43 PetscErrorCode KSPLoad(KSP newdm, PetscViewer viewer)
44 {
45   PetscBool isbinary;
46   PetscInt  classid;
47   char      type[256];
48   PC        pc;
49 
50   PetscFunctionBegin;
51   PetscValidHeaderSpecific(newdm, KSP_CLASSID, 1);
52   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
53   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
54   PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()");
55 
56   PetscCall(PetscViewerBinaryRead(viewer, &classid, 1, NULL, PETSC_INT));
57   PetscCheck(classid == KSP_FILE_CLASSID, PetscObjectComm((PetscObject)newdm), PETSC_ERR_ARG_WRONG, "Not KSP next in file");
58   PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR));
59   PetscCall(KSPSetType(newdm, type));
60   PetscTryTypeMethod(newdm, load, viewer);
61   PetscCall(KSPGetPC(newdm, &pc));
62   PetscCall(PCLoad(pc, viewer));
63   PetscFunctionReturn(PETSC_SUCCESS);
64 }
65 
66 #include <petscdraw.h>
67 #if defined(PETSC_HAVE_SAWS)
68   #include <petscviewersaws.h>
69 #endif
70 /*@
71   KSPView - Prints the various parameters currently set in the `KSP` object. For example, the convergence tolerances and `KSPType`.
72   Also views the `PC` and `Mat` contained by the `KSP` with `PCView()` and `MatView()`.
73 
74   Collective
75 
76   Input Parameters:
77 + ksp    - the Krylov space context
78 - viewer - visualization context
79 
80   Options Database Key:
81 . -ksp_view - print the `KSP` data structure at the end of each `KSPSolve()` call
82 
83   Level: beginner
84 
85   Notes:
86   The available visualization contexts include
87 +     `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
88 -     `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
89   output where only the first processor opens
90   the file.  All other processors send their
91   data to the first processor to print.
92 
93   The available formats include
94 +     `PETSC_VIEWER_DEFAULT` - standard output (default)
95 -     `PETSC_VIEWER_ASCII_INFO_DETAIL` - more verbose output for `PCBJACOBI` and `PCASM`
96 
97   The user can open an alternative visualization context with
98   `PetscViewerASCIIOpen()` - output to a specified file.
99 
100   Use `KSPViewFromOptions()` to allow the user to select many different `PetscViewerType` and formats from the options database.
101 
102   In the debugger you can do call `KSPView(ksp,0)` to display the `KSP`. (The same holds for any PETSc object viewer).
103 
104 .seealso: [](ch_ksp), `KSP`, `PetscViewer`, `PCView()`, `PetscViewerASCIIOpen()`, `KSPViewFromOptions()`
105 @*/
106 PetscErrorCode KSPView(KSP ksp, PetscViewer viewer)
107 {
108   PetscBool isascii, isbinary, isdraw, isstring;
109 #if defined(PETSC_HAVE_SAWS)
110   PetscBool issaws;
111 #endif
112 
113   PetscFunctionBegin;
114   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
115   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ksp), &viewer));
116   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
117   PetscCheckSameComm(ksp, 1, viewer, 2);
118 
119   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
120   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
121   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
122   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
123 #if defined(PETSC_HAVE_SAWS)
124   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSAWS, &issaws));
125 #endif
126   if (isascii) {
127     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)ksp, viewer));
128     PetscCall(PetscViewerASCIIPushTab(viewer));
129     PetscTryTypeMethod(ksp, view, viewer);
130     PetscCall(PetscViewerASCIIPopTab(viewer));
131     if (ksp->guess_zero) {
132       PetscCall(PetscViewerASCIIPrintf(viewer, "  maximum iterations=%" PetscInt_FMT ", initial guess is zero\n", ksp->max_it));
133     } else {
134       PetscCall(PetscViewerASCIIPrintf(viewer, "  maximum iterations=%" PetscInt_FMT ", nonzero initial guess\n", ksp->max_it));
135     }
136     if (ksp->min_it) PetscCall(PetscViewerASCIIPrintf(viewer, "  minimum iterations=%" PetscInt_FMT "\n", ksp->min_it));
137     if (ksp->guess_knoll) PetscCall(PetscViewerASCIIPrintf(viewer, "  using preconditioner applied to right-hand side for initial guess\n"));
138     PetscCall(PetscViewerASCIIPrintf(viewer, "  tolerances: relative=%g, absolute=%g, divergence=%g\n", (double)ksp->rtol, (double)ksp->abstol, (double)ksp->divtol));
139     if (ksp->pc_side == PC_RIGHT) {
140       PetscCall(PetscViewerASCIIPrintf(viewer, "  right preconditioning\n"));
141     } else if (ksp->pc_side == PC_SYMMETRIC) {
142       PetscCall(PetscViewerASCIIPrintf(viewer, "  symmetric preconditioning\n"));
143     } else {
144       PetscCall(PetscViewerASCIIPrintf(viewer, "  left preconditioning\n"));
145     }
146     if (ksp->guess) {
147       PetscCall(PetscViewerASCIIPushTab(viewer));
148       PetscCall(KSPGuessView(ksp->guess, viewer));
149       PetscCall(PetscViewerASCIIPopTab(viewer));
150     }
151     if (ksp->dscale) PetscCall(PetscViewerASCIIPrintf(viewer, "  diagonally scaled system\n"));
152     if (ksp->converged == KSPConvergedSkip || ksp->normtype == KSP_NORM_NONE) PetscCall(PetscViewerASCIIPrintf(viewer, "  not checking for convergence\n"));
153     else PetscCall(PetscViewerASCIIPrintf(viewer, "  using %s norm type for convergence test\n", KSPNormTypes[ksp->normtype]));
154   } else if (isbinary) {
155     PetscInt    classid = KSP_FILE_CLASSID;
156     MPI_Comm    comm;
157     PetscMPIInt rank;
158     char        type[256];
159 
160     PetscCall(PetscObjectGetComm((PetscObject)ksp, &comm));
161     PetscCallMPI(MPI_Comm_rank(comm, &rank));
162     if (rank == 0) {
163       PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT));
164       PetscCall(PetscStrncpy(type, ((PetscObject)ksp)->type_name, 256));
165       PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR));
166     }
167     PetscTryTypeMethod(ksp, view, viewer);
168   } else if (isstring) {
169     const char *type;
170     PetscCall(KSPGetType(ksp, &type));
171     PetscCall(PetscViewerStringSPrintf(viewer, " KSPType: %-7.7s", type));
172     PetscTryTypeMethod(ksp, view, viewer);
173   } else if (isdraw) {
174     PetscDraw draw;
175     char      str[36];
176     PetscReal x, y, bottom, h;
177     PetscBool flg;
178 
179     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
180     PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
181     PetscCall(PetscObjectTypeCompare((PetscObject)ksp, KSPPREONLY, &flg));
182     if (!flg) {
183       PetscCall(PetscStrncpy(str, "KSP: ", sizeof(str)));
184       PetscCall(PetscStrlcat(str, ((PetscObject)ksp)->type_name, sizeof(str)));
185       PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h));
186       bottom = y - h;
187     } else {
188       bottom = y;
189     }
190     PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
191 #if defined(PETSC_HAVE_SAWS)
192   } else if (issaws) {
193     PetscMPIInt rank;
194     const char *name;
195 
196     PetscCall(PetscObjectGetName((PetscObject)ksp, &name));
197     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
198     if (!((PetscObject)ksp)->amsmem && rank == 0) {
199       char dir[1024];
200 
201       PetscCall(PetscObjectViewSAWs((PetscObject)ksp, viewer));
202       PetscCall(PetscSNPrintf(dir, 1024, "/PETSc/Objects/%s/its", name));
203       PetscCallSAWs(SAWs_Register, (dir, &ksp->its, 1, SAWs_READ, SAWs_INT));
204       if (!ksp->res_hist) PetscCall(KSPSetResidualHistory(ksp, NULL, PETSC_DECIDE, PETSC_TRUE));
205       PetscCall(PetscSNPrintf(dir, 1024, "/PETSc/Objects/%s/res_hist", name));
206       PetscCallSAWs(SAWs_Register, (dir, ksp->res_hist, 10, SAWs_READ, SAWs_DOUBLE));
207     }
208 #endif
209   } else PetscTryTypeMethod(ksp, view, viewer);
210   if (ksp->pc) PetscCall(PCView(ksp->pc, viewer));
211   if (isdraw) {
212     PetscDraw draw;
213     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
214     PetscCall(PetscDrawPopCurrentPoint(draw));
215   }
216   PetscFunctionReturn(PETSC_SUCCESS);
217 }
218 
219 /*@
220   KSPViewFromOptions - View (print) a `KSP` object based on values in the options database. Also views the `PC` and `Mat` contained by the `KSP`
221   with `PCView()` and `MatView()`.
222 
223   Collective
224 
225   Input Parameters:
226 + A    - Krylov solver context
227 . obj  - Optional object that provides the options prefix used to query the options database
228 - name - command line option
229 
230   Level: intermediate
231 
232 .seealso: [](ch_ksp), `KSP`, `KSPView()`, `PetscObjectViewFromOptions()`, `KSPCreate()`
233 @*/
234 PetscErrorCode KSPViewFromOptions(KSP A, PetscObject obj, const char name[])
235 {
236   PetscFunctionBegin;
237   PetscValidHeaderSpecific(A, KSP_CLASSID, 1);
238   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
239   PetscFunctionReturn(PETSC_SUCCESS);
240 }
241 
242 /*@
243   KSPSetNormType - Sets the type of residual norm that is used for convergence testing in `KSPSolve()` for the given `KSP` context
244 
245   Logically Collective
246 
247   Input Parameters:
248 + ksp      - Krylov solver context
249 - normtype - one of
250 .vb
251    KSP_NORM_NONE             - skips computing the norm, this should generally only be used if you are using
252                                the Krylov method as a smoother with a fixed small number of iterations.
253                                Implicitly sets `KSPConvergedSkip()` as the `KSP` convergence test.
254                                Note that certain algorithms such as `KSPGMRES` ALWAYS require the norm calculation,
255                                for these methods the norms are still computed, they are just not used in
256                                the convergence test.
257    KSP_NORM_PRECONDITIONED   - the default for left-preconditioned solves, uses the 2-norm
258                                of the preconditioned residual  $B^{-1}(b - A x)$.
259    KSP_NORM_UNPRECONDITIONED - uses the 2-norm of the true $b - Ax$ residual.
260    KSP_NORM_NATURAL          - uses the $A$ norm of the true $b - Ax$ residual; supported by `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS`
261 .ve
262 
263   Options Database Key:
264 . -ksp_norm_type <none,preconditioned,unpreconditioned,natural> - set `KSP` norm type
265 
266   Level: advanced
267 
268   Notes:
269   The norm is always of the equations residual $\| b - A x^n \|$  (or an approximation to that norm), they are never a norm of the error in the equation.
270 
271   Not all combinations of preconditioner side (see `KSPSetPCSide()`) and norm types are supported by all Krylov methods.
272   If only one is set, PETSc tries to automatically change the other to find a compatible pair.  If no such combination
273   is supported, PETSc will generate an error.
274 
275   Developer Note:
276   Supported combinations of norm and preconditioner side are set using `KSPSetSupportedNorm()` for each `KSPType`.
277 
278 .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSPConvergedSkip()`, `KSPSetCheckNormIteration()`, `KSPSetPCSide()`, `KSPGetPCSide()`, `KSPNormType`
279 @*/
280 PetscErrorCode KSPSetNormType(KSP ksp, KSPNormType normtype)
281 {
282   PetscFunctionBegin;
283   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
284   PetscValidLogicalCollectiveEnum(ksp, normtype, 2);
285   ksp->normtype = ksp->normtype_set = normtype;
286   PetscFunctionReturn(PETSC_SUCCESS);
287 }
288 
289 /*@
290   KSPSetCheckNormIteration - Sets the first iteration at which the norm of the residual will be
291   computed and used in the convergence test of `KSPSolve()` for the given `KSP` context
292 
293   Logically Collective
294 
295   Input Parameters:
296 + ksp - Krylov solver context
297 - it  - use -1 to check at all iterations
298 
299   Level: advanced
300 
301   Notes:
302   Currently only works with `KSPCG`, `KSPBCGS` and `KSPIBCGS`
303 
304   Use `KSPSetNormType`(ksp,`KSP_NORM_NONE`) to never check the norm
305 
306   On steps where the norm is not computed, the previous norm is still in the variable, so if you run with, for example,
307   `-ksp_monitor` the residual norm will appear to be unchanged for several iterations (though it is not really unchanged).
308 
309   Certain methods such as `KSPGMRES` always compute the residual norm, this routine will not change that computation, but it will
310   prevent the computed norm from being checked.
311 
312 .seealso: [](ch_ksp), `KSP`, `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSPConvergedSkip()`, `KSPSetNormType()`, `KSPSetLagNorm()`
313 @*/
314 PetscErrorCode KSPSetCheckNormIteration(KSP ksp, PetscInt it)
315 {
316   PetscFunctionBegin;
317   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
318   PetscValidLogicalCollectiveInt(ksp, it, 2);
319   ksp->chknorm = it;
320   PetscFunctionReturn(PETSC_SUCCESS);
321 }
322 
323 /*@
324   KSPSetLagNorm - Lags the residual norm calculation so that it is computed as part of the `MPI_Allreduce()` used for
325   computing the inner products needed for the next iteration.
326 
327   Logically Collective
328 
329   Input Parameters:
330 + ksp - Krylov solver context
331 - flg - `PETSC_TRUE` or `PETSC_FALSE`
332 
333   Options Database Key:
334 . -ksp_lag_norm - lag the calculated residual norm
335 
336   Level: advanced
337 
338   Notes:
339   Currently only works with `KSPIBCGS`.
340 
341   This can reduce communication costs at the expense of doing
342   one additional iteration because the norm used in the convergence test of `KSPSolve()` is one iteration behind the actual
343   current residual norm (which has not yet been computed due to the lag).
344 
345   Use `KSPSetNormType`(ksp,`KSP_NORM_NONE`) to never check the norm
346 
347   If you lag the norm and run with, for example, `-ksp_monitor`, the residual norm reported will be the lagged one.
348 
349   `KSPSetCheckNormIteration()` is an alternative way of avoiding the expense of computing the residual norm at each iteration.
350 
351 .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSPConvergedSkip()`, `KSPSetNormType()`, `KSPSetCheckNormIteration()`
352 @*/
353 PetscErrorCode KSPSetLagNorm(KSP ksp, PetscBool flg)
354 {
355   PetscFunctionBegin;
356   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
357   PetscValidLogicalCollectiveBool(ksp, flg, 2);
358   ksp->lagnorm = flg;
359   PetscFunctionReturn(PETSC_SUCCESS);
360 }
361 
362 /*@
363   KSPSetSupportedNorm - Sets a norm and preconditioner side supported by a `KSPType`
364 
365   Logically Collective
366 
367   Input Parameters:
368 + ksp      - Krylov method
369 . normtype - supported norm type of the type `KSPNormType`
370 . pcside   - preconditioner side, of the type `PCSide` that can be used with this `KSPNormType`
371 - priority - positive integer preference for this combination; larger values have higher priority
372 
373   Level: developer
374 
375   Notes:
376   This function should be called from the implementation files `KSPCreate_XXX()` to declare
377   which norms and preconditioner sides are supported. Users should not call this
378   function.
379 
380   This function can be called multiple times for each combination of `KSPNormType` and `PCSide`
381   the `KSPType` supports
382 
383 .seealso: [](ch_ksp), `KSP`, `KSPNormType`, `PCSide`, `KSPSetNormType()`, `KSPSetPCSide()`
384 @*/
385 PetscErrorCode KSPSetSupportedNorm(KSP ksp, KSPNormType normtype, PCSide pcside, PetscInt priority)
386 {
387   PetscFunctionBegin;
388   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
389   ksp->normsupporttable[normtype][pcside] = priority;
390   PetscFunctionReturn(PETSC_SUCCESS);
391 }
392 
393 static PetscErrorCode KSPNormSupportTableReset_Private(KSP ksp)
394 {
395   PetscFunctionBegin;
396   PetscCall(PetscMemzero(ksp->normsupporttable, sizeof(ksp->normsupporttable)));
397   ksp->pc_side  = ksp->pc_side_set;
398   ksp->normtype = ksp->normtype_set;
399   PetscFunctionReturn(PETSC_SUCCESS);
400 }
401 
402 PetscErrorCode KSPSetUpNorms_Private(KSP ksp, PetscBool errorifnotsupported, KSPNormType *normtype, PCSide *pcside)
403 {
404   PetscInt i, j, best, ibest = 0, jbest = 0;
405 
406   PetscFunctionBegin;
407   best = 0;
408   for (i = 0; i < KSP_NORM_MAX; i++) {
409     for (j = 0; j < PC_SIDE_MAX; j++) {
410       if ((ksp->normtype == KSP_NORM_DEFAULT || ksp->normtype == i) && (ksp->pc_side == PC_SIDE_DEFAULT || ksp->pc_side == j) && ksp->normsupporttable[i][j] > best) {
411         best  = ksp->normsupporttable[i][j];
412         ibest = i;
413         jbest = j;
414       }
415     }
416   }
417   if (best < 1 && errorifnotsupported) {
418     PetscCheck(ksp->normtype != KSP_NORM_DEFAULT || ksp->pc_side != PC_SIDE_DEFAULT, PetscObjectComm((PetscObject)ksp), PETSC_ERR_PLIB, "The %s KSP implementation did not call KSPSetSupportedNorm()", ((PetscObject)ksp)->type_name);
419     PetscCheck(ksp->normtype != KSP_NORM_DEFAULT, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "KSP %s does not support preconditioner side %s", ((PetscObject)ksp)->type_name, PCSides[ksp->pc_side]);
420     PetscCheck(ksp->pc_side != PC_SIDE_DEFAULT, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "KSP %s does not support norm type %s", ((PetscObject)ksp)->type_name, KSPNormTypes[ksp->normtype]);
421     SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "KSP %s does not support norm type %s with preconditioner side %s", ((PetscObject)ksp)->type_name, KSPNormTypes[ksp->normtype], PCSides[ksp->pc_side]);
422   }
423   if (normtype) *normtype = (KSPNormType)ibest;
424   if (pcside) *pcside = (PCSide)jbest;
425   PetscFunctionReturn(PETSC_SUCCESS);
426 }
427 
428 /*@
429   KSPGetNormType - Gets the `KSPNormType` that is used for convergence testing during `KSPSolve()` for this `KSP` context
430 
431   Not Collective
432 
433   Input Parameter:
434 . ksp - Krylov solver context
435 
436   Output Parameter:
437 . normtype - the `KSPNormType` that is used for convergence testing
438 
439   Level: advanced
440 
441 .seealso: [](ch_ksp), `KSPNormType`, `KSPSetNormType()`, `KSPConvergedSkip()`
442 @*/
443 PetscErrorCode KSPGetNormType(KSP ksp, KSPNormType *normtype)
444 {
445   PetscFunctionBegin;
446   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
447   PetscAssertPointer(normtype, 2);
448   PetscCall(KSPSetUpNorms_Private(ksp, PETSC_TRUE, &ksp->normtype, &ksp->pc_side));
449   *normtype = ksp->normtype;
450   PetscFunctionReturn(PETSC_SUCCESS);
451 }
452 
453 #if defined(PETSC_HAVE_SAWS)
454   #include <petscviewersaws.h>
455 #endif
456 
457 /*@
458   KSPSetOperators - Sets the matrix associated with the linear system
459   and a (possibly) different one from which the preconditioner will be built into the `KSP` context. The matrix will then be used during `KSPSolve()`
460 
461   Collective
462 
463   Input Parameters:
464 + ksp  - the `KSP` context
465 . Amat - the matrix that defines the linear system
466 - Pmat - the matrix to be used in constructing the preconditioner, usually the same as `Amat`.
467 
468   Level: beginner
469 
470   Notes:
471   If you know the operator `Amat` has a null space you can use `MatSetNullSpace()` and `MatSetTransposeNullSpace()` to supply the null
472   space to `Amat` and the `KSP` solvers will automatically use that null space as needed during the solution process.
473 
474   All future calls to `KSPSetOperators()` must use the same size matrices, unless `KSPReset()` is called!
475 
476   Passing a `NULL` for `Amat` or `Pmat` removes the matrix that is currently being used from the `KSP` context.
477 
478   If you wish to replace either `Amat` or `Pmat` but leave the other one untouched then
479   first call `KSPGetOperators()` to get the one you wish to keep, call `PetscObjectReference()`
480   on it and then pass it back in your call to `KSPSetOperators()`.
481 
482   Developer Notes:
483   If the operators have NOT been set with `KSPSetOperators()` then the operators
484   are created in the `PC` and returned to the user. In this case, if both operators
485   mat and pmat are requested, two DIFFERENT operators will be returned. If
486   only one is requested both operators in the `PC` will be the same (i.e. as
487   if one had called `KSPSetOperators()` with the same argument for both `Mat`s).
488   The user must set the sizes of the returned matrices and their type etc just
489   as if the user created them with `MatCreate()`. For example,
490 
491 .vb
492          KSPGetOperators(ksp/pc,&mat,NULL); is equivalent to
493            set size, type, etc of mat
494 
495          MatCreate(comm,&mat);
496          KSP/PCSetOperators(ksp/pc,mat,mat);
497          PetscObjectDereference((PetscObject)mat);
498            set size, type, etc of mat
499 
500      and
501 
502          KSP/PCGetOperators(ksp/pc,&mat,&pmat); is equivalent to
503            set size, type, etc of mat and pmat
504 
505          MatCreate(comm,&mat);
506          MatCreate(comm,&pmat);
507          KSP/PCSetOperators(ksp/pc,mat,pmat);
508          PetscObjectDereference((PetscObject)mat);
509          PetscObjectDereference((PetscObject)pmat);
510            set size, type, etc of mat and pmat
511 .ve
512 
513   The rationale for this support is so that when creating a `TS`, `SNES`, or `KSP` the hierarchy
514   of underlying objects (i.e. `SNES`, `KSP`, `PC`, `Mat`) and their lifespans can be completely
515   managed by the top most level object (i.e. the `TS`, `SNES`, or `KSP`). Another way to look
516   at this is when you create a `SNES` you do not NEED to create a `KSP` and attach it to
517   the `SNES` object (the `SNES` object manages it for you). Similarly when you create a `KSP`
518   you do not need to attach a `PC` to it (the `KSP` object manages the `PC` object for you).
519   Thus, why should YOU have to create the `Mat` and attach it to the `SNES`/`KSP`/`PC`, when
520   it can be created for you?
521 
522 .seealso: [](ch_ksp), `KSP`, `Mat`, `KSPSolve()`, `KSPGetPC()`, `PCGetOperators()`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetComputeOperators()`, `KSPSetComputeInitialGuess()`, `KSPSetComputeRHS()`
523 @*/
524 PetscErrorCode KSPSetOperators(KSP ksp, Mat Amat, Mat Pmat)
525 {
526   PetscFunctionBegin;
527   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
528   if (Amat) PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2);
529   if (Pmat) PetscValidHeaderSpecific(Pmat, MAT_CLASSID, 3);
530   if (Amat) PetscCheckSameComm(ksp, 1, Amat, 2);
531   if (Pmat) PetscCheckSameComm(ksp, 1, Pmat, 3);
532   if (!ksp->pc) PetscCall(KSPGetPC(ksp, &ksp->pc));
533   PetscCall(PCSetOperators(ksp->pc, Amat, Pmat));
534   if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX; /* so that next solve call will call PCSetUp() on new matrix */
535   PetscFunctionReturn(PETSC_SUCCESS);
536 }
537 
538 /*@
539   KSPGetOperators - Gets the matrix associated with the linear system
540   and a (possibly) different one used to construct the preconditioner from the `KSP` context
541 
542   Collective
543 
544   Input Parameter:
545 . ksp - the `KSP` context
546 
547   Output Parameters:
548 + Amat - the matrix that defines the linear system
549 - Pmat - the matrix to be used in constructing the preconditioner, usually the same as `Amat`.
550 
551   Level: intermediate
552 
553   Notes:
554   If `KSPSetOperators()` has not been called then the `KSP` object will attempt to automatically create the matrix `Amat` and return it
555 
556   Use `KSPGetOperatorsSet()` to determine if matrices have been provided.
557 
558   DOES NOT increase the reference counts of the matrix, so you should NOT destroy them.
559 
560 .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `KSPGetPC()`, `PCSetOperators()`, `KSPSetOperators()`, `KSPGetOperatorsSet()`
561 @*/
562 PetscErrorCode KSPGetOperators(KSP ksp, Mat *Amat, Mat *Pmat)
563 {
564   PetscFunctionBegin;
565   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
566   if (!ksp->pc) PetscCall(KSPGetPC(ksp, &ksp->pc));
567   PetscCall(PCGetOperators(ksp->pc, Amat, Pmat));
568   PetscFunctionReturn(PETSC_SUCCESS);
569 }
570 
571 /*@
572   KSPGetOperatorsSet - Determines if the matrix associated with the linear system and
573   possibly a different one from which the preconditioner will be built have been set in the `KSP` with `KSPSetOperators()`
574 
575   Not Collective, though the results on all processes will be the same
576 
577   Input Parameter:
578 . ksp - the `KSP` context
579 
580   Output Parameters:
581 + mat  - the matrix associated with the linear system was set
582 - pmat - matrix from which the preconditioner will be built, usually the same as `mat` was set
583 
584   Level: intermediate
585 
586   Note:
587   This routine exists because if you call `KSPGetOperators()` on a `KSP` that does not yet have operators they are
588   automatically created in the call.
589 
590 .seealso: [](ch_ksp), `KSP`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperators()`, `PCGetOperatorsSet()`
591 @*/
592 PetscErrorCode KSPGetOperatorsSet(KSP ksp, PetscBool *mat, PetscBool *pmat)
593 {
594   PetscFunctionBegin;
595   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
596   if (!ksp->pc) PetscCall(KSPGetPC(ksp, &ksp->pc));
597   PetscCall(PCGetOperatorsSet(ksp->pc, mat, pmat));
598   PetscFunctionReturn(PETSC_SUCCESS);
599 }
600 
601 /*@C
602   KSPSetPreSolve - Sets a function that is called at the beginning of each `KSPSolve()`. Used in conjunction with `KSPSetPostSolve()`.
603 
604   Logically Collective
605 
606   Input Parameters:
607 + ksp      - the solver object
608 . presolve - the function to call before the solve, see` KSPPSolveFn`
609 - ctx      - an optional context needed by the function
610 
611   Level: developer
612 
613   Notes:
614   The function provided here `presolve` is used to modify the right hand side, and possibly the matrix, of the linear system to be solved.
615   The function provided with `KSPSetPostSolve()` then modifies the resulting solution of that linear system to obtain the correct solution
616   to the initial linear system.
617 
618   The functions `PCPreSolve()` and `PCPostSolve()` provide a similar functionality and are used, for example with `PCEISENSTAT`.
619 
620 .seealso: [](ch_ksp), `KSPPSolveFn`, `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPSetPostSolve()`, `PCEISENSTAT`, `PCPreSolve()`, `PCPostSolve()`
621 @*/
622 PetscErrorCode KSPSetPreSolve(KSP ksp, KSPPSolveFn *presolve, void *ctx)
623 {
624   PetscFunctionBegin;
625   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
626   ksp->presolve = presolve;
627   ksp->prectx   = ctx;
628   PetscFunctionReturn(PETSC_SUCCESS);
629 }
630 
631 /*@C
632   KSPSetPostSolve - Sets a function that is called at the end of each `KSPSolve()` (whether it converges or not). Used in conjunction with `KSPSetPreSolve()`.
633 
634   Logically Collective
635 
636   Input Parameters:
637 + ksp       - the solver object
638 . postsolve - the function to call after the solve, see` KSPPSolveFn`
639 - ctx       - an optional context needed by the function
640 
641   Level: developer
642 
643 .seealso: [](ch_ksp), `KSPPSolveFn`, `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPSetPreSolve()`, `PCEISENSTAT`
644 @*/
645 PetscErrorCode KSPSetPostSolve(KSP ksp, KSPPSolveFn *postsolve, void *ctx)
646 {
647   PetscFunctionBegin;
648   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
649   ksp->postsolve = postsolve;
650   ksp->postctx   = ctx;
651   PetscFunctionReturn(PETSC_SUCCESS);
652 }
653 
654 /*@
655   KSPSetNestLevel - sets the amount of nesting the `KSP` has. That is the number of levels of `KSP` above this `KSP` in a linear solve.
656 
657   Collective
658 
659   Input Parameters:
660 + ksp   - the `KSP`
661 - level - the nest level
662 
663   Level: developer
664 
665   Note:
666   For example, the `KSP` in each block of a `KSPBJACOBI` has a level of 1, while the outer `KSP` has a level of 0.
667 
668 .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPGetNestLevel()`, `PCSetKSPNestLevel()`, `PCGetKSPNestLevel()`
669 @*/
670 PetscErrorCode KSPSetNestLevel(KSP ksp, PetscInt level)
671 {
672   PetscFunctionBegin;
673   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
674   PetscValidLogicalCollectiveInt(ksp, level, 2);
675   ksp->nestlevel = level;
676   PetscFunctionReturn(PETSC_SUCCESS);
677 }
678 
679 /*@
680   KSPGetNestLevel - gets the amount of nesting the `KSP` has
681 
682   Not Collective
683 
684   Input Parameter:
685 . ksp - the `KSP`
686 
687   Output Parameter:
688 . level - the nest level
689 
690   Level: developer
691 
692 .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPSetNestLevel()`, `PCSetKSPNestLevel()`, `PCGetKSPNestLevel()`
693 @*/
694 PetscErrorCode KSPGetNestLevel(KSP ksp, PetscInt *level)
695 {
696   PetscFunctionBegin;
697   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
698   PetscAssertPointer(level, 2);
699   *level = ksp->nestlevel;
700   PetscFunctionReturn(PETSC_SUCCESS);
701 }
702 
703 /*@
704   KSPCreate - Creates the `KSP` context. This `KSP` context is used in PETSc to solve linear systems with `KSPSolve()`
705 
706   Collective
707 
708   Input Parameter:
709 . comm - MPI communicator
710 
711   Output Parameter:
712 . inksp - location to put the `KSP` context
713 
714   Level: beginner
715 
716   Note:
717   The default `KSPType` is `KSPGMRES` with a restart of 30, using modified Gram-Schmidt orthogonalization. The `KSPType` may be
718   changed with `KSPSetType()`
719 
720 .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPSetType()`
721 @*/
722 PetscErrorCode KSPCreate(MPI_Comm comm, KSP *inksp)
723 {
724   KSP   ksp;
725   void *ctx;
726 
727   PetscFunctionBegin;
728   PetscAssertPointer(inksp, 2);
729   PetscCall(KSPInitializePackage());
730 
731   PetscCall(PetscHeaderCreate(ksp, KSP_CLASSID, "KSP", "Krylov Method", "KSP", comm, KSPDestroy, KSPView));
732   ksp->default_max_it = ksp->max_it = 10000;
733   ksp->pc_side = ksp->pc_side_set = PC_SIDE_DEFAULT;
734 
735   ksp->default_rtol = ksp->rtol = 1.e-5;
736   ksp->default_abstol = ksp->abstol = PetscDefined(USE_REAL_SINGLE) ? 1.e-25 : 1.e-50;
737   ksp->default_divtol = ksp->divtol = 1.e4;
738 
739   ksp->chknorm  = -1;
740   ksp->normtype = ksp->normtype_set = KSP_NORM_DEFAULT;
741   ksp->rnorm                        = 0.0;
742   ksp->its                          = 0;
743   ksp->guess_zero                   = PETSC_TRUE;
744   ksp->calc_sings                   = PETSC_FALSE;
745   ksp->res_hist                     = NULL;
746   ksp->res_hist_alloc               = NULL;
747   ksp->res_hist_len                 = 0;
748   ksp->res_hist_max                 = 0;
749   ksp->res_hist_reset               = PETSC_TRUE;
750   ksp->err_hist                     = NULL;
751   ksp->err_hist_alloc               = NULL;
752   ksp->err_hist_len                 = 0;
753   ksp->err_hist_max                 = 0;
754   ksp->err_hist_reset               = PETSC_TRUE;
755   ksp->numbermonitors               = 0;
756   ksp->numberreasonviews            = 0;
757   ksp->setfromoptionscalled         = 0;
758   ksp->nmax                         = PETSC_DECIDE;
759 
760   PetscCall(KSPConvergedDefaultCreate(&ctx));
761   PetscCall(KSPSetConvergenceTest(ksp, KSPConvergedDefault, ctx, KSPConvergedDefaultDestroy));
762   ksp->ops->buildsolution = KSPBuildSolutionDefault;
763   ksp->ops->buildresidual = KSPBuildResidualDefault;
764 
765   ksp->vec_sol    = NULL;
766   ksp->vec_rhs    = NULL;
767   ksp->pc         = NULL;
768   ksp->data       = NULL;
769   ksp->nwork      = 0;
770   ksp->work       = NULL;
771   ksp->reason     = KSP_CONVERGED_ITERATING;
772   ksp->setupstage = KSP_SETUP_NEW;
773 
774   PetscCall(KSPNormSupportTableReset_Private(ksp));
775 
776   *inksp = ksp;
777   PetscFunctionReturn(PETSC_SUCCESS);
778 }
779 
780 /*@
781   KSPSetType - Sets the algorithm/method to be used to solve the linear system with the given `KSP`
782 
783   Logically Collective
784 
785   Input Parameters:
786 + ksp  - the Krylov space context
787 - type - a known method
788 
789   Options Database Key:
790 . -ksp_type  <method> - Sets the method; see `KSPGType` or use `-help` for a list  of available methods (for instance, cg or gmres)
791 
792   Level: intermediate
793 
794   Notes:
795   See `KSPType` for available methods (for instance, `KSPCG` or `KSPGMRES`).
796 
797   Normally, it is best to use the `KSPSetFromOptions()` command and
798   then set the `KSP` type from the options database rather than by using
799   this routine.  Using the options database provides the user with
800   maximum flexibility in evaluating the many different Krylov methods.
801   The `KSPSetType()` routine is provided for those situations where it
802   is necessary to set the iterative solver independently of the command
803   line or options database.  This might be the case, for example, when
804   the choice of iterative solver changes during the execution of the
805   program, and the user's application is taking responsibility for
806   choosing the appropriate method.  In other words, this routine is
807   not for beginners.
808 
809   Developer Note:
810   `KSPRegister()` is used to add Krylov types to `KSPList` from which they are accessed by `KSPSetType()`.
811 
812 .seealso: [](ch_ksp), `PCSetType()`, `KSPType`, `KSPRegister()`, `KSPCreate()`, `KSP`
813 @*/
814 PetscErrorCode KSPSetType(KSP ksp, KSPType type)
815 {
816   PetscBool match;
817   PetscErrorCode (*r)(KSP);
818 
819   PetscFunctionBegin;
820   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
821   PetscAssertPointer(type, 2);
822 
823   PetscCall(PetscObjectTypeCompare((PetscObject)ksp, type, &match));
824   if (match) PetscFunctionReturn(PETSC_SUCCESS);
825 
826   PetscCall(PetscFunctionListFind(KSPList, type, &r));
827   PetscCheck(r, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested KSP type %s", type);
828   /* Destroy the previous private KSP context */
829   PetscTryTypeMethod(ksp, destroy);
830 
831   /* Reinitialize function pointers in KSPOps structure */
832   PetscCall(PetscMemzero(ksp->ops, sizeof(struct _KSPOps)));
833   ksp->ops->buildsolution = KSPBuildSolutionDefault;
834   ksp->ops->buildresidual = KSPBuildResidualDefault;
835   PetscCall(KSPNormSupportTableReset_Private(ksp));
836   ksp->converged_neg_curve = PETSC_FALSE; // restore default
837   ksp->setupnewmatrix      = PETSC_FALSE; // restore default (setup not called in case of new matrix)
838   /* Call the KSPCreate_XXX routine for this particular Krylov solver */
839   ksp->setupstage     = KSP_SETUP_NEW;
840   ksp->guess_not_read = PETSC_FALSE; // restore default
841   PetscCall((*r)(ksp));
842   PetscCall(PetscObjectChangeTypeName((PetscObject)ksp, type));
843   PetscFunctionReturn(PETSC_SUCCESS);
844 }
845 
846 /*@
847   KSPGetType - Gets the `KSP` type as a string from the `KSP` object.
848 
849   Not Collective
850 
851   Input Parameter:
852 . ksp - Krylov context
853 
854   Output Parameter:
855 . type - name of the `KSP` method
856 
857   Level: intermediate
858 
859 .seealso: [](ch_ksp), `KSPType`, `KSP`, `KSPSetType()`
860 @*/
861 PetscErrorCode KSPGetType(KSP ksp, KSPType *type)
862 {
863   PetscFunctionBegin;
864   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
865   PetscAssertPointer(type, 2);
866   *type = ((PetscObject)ksp)->type_name;
867   PetscFunctionReturn(PETSC_SUCCESS);
868 }
869 
870 /*@C
871   KSPRegister -  Adds a method, `KSPType`, to the Krylov subspace solver package.
872 
873   Not Collective, No Fortran Support
874 
875   Input Parameters:
876 + sname    - name of a new user-defined solver
877 - function - routine to create method
878 
879   Level: advanced
880 
881   Note:
882   `KSPRegister()` may be called multiple times to add several user-defined solvers.
883 
884   Example Usage:
885 .vb
886    KSPRegister("my_solver", MySolverCreate);
887 .ve
888 
889   Then, your solver can be chosen with the procedural interface via
890 .vb
891   KSPSetType(ksp, "my_solver")
892 .ve
893   or at runtime via the option `-ksp_type my_solver`
894 
895 .seealso: [](ch_ksp), `KSP`, `KSPType`, `KSPSetType`, `KSPRegisterAll()`
896 @*/
897 PetscErrorCode KSPRegister(const char sname[], PetscErrorCode (*function)(KSP))
898 {
899   PetscFunctionBegin;
900   PetscCall(KSPInitializePackage());
901   PetscCall(PetscFunctionListAdd(&KSPList, sname, function));
902   PetscFunctionReturn(PETSC_SUCCESS);
903 }
904 
905 PetscErrorCode KSPMonitorMakeKey_Internal(const char name[], PetscViewerType vtype, PetscViewerFormat format, char key[])
906 {
907   PetscFunctionBegin;
908   PetscCall(PetscStrncpy(key, name, PETSC_MAX_PATH_LEN));
909   PetscCall(PetscStrlcat(key, ":", PETSC_MAX_PATH_LEN));
910   PetscCall(PetscStrlcat(key, vtype, PETSC_MAX_PATH_LEN));
911   PetscCall(PetscStrlcat(key, ":", PETSC_MAX_PATH_LEN));
912   PetscCall(PetscStrlcat(key, PetscViewerFormats[format], PETSC_MAX_PATH_LEN));
913   PetscFunctionReturn(PETSC_SUCCESS);
914 }
915 
916 /*@C
917   KSPMonitorRegister -  Registers a Krylov subspace solver monitor routine that may be accessed with `KSPMonitorSetFromOptions()`
918 
919   Not Collective
920 
921   Input Parameters:
922 + name    - name of a new monitor type
923 . vtype   - A `PetscViewerType` for the output
924 . format  - A `PetscViewerFormat` for the output
925 . monitor - Monitor routine, see `KSPMonitorRegisterFn`
926 . create  - Creation routine, or `NULL`
927 - destroy - Destruction routine, or `NULL`
928 
929   Level: advanced
930 
931   Notes:
932   `KSPMonitorRegister()` may be called multiple times to add several user-defined monitors.
933 
934   The calling sequence for the given function matches the calling sequence used by `KSPMonitorFn` functions passed to `KSPMonitorSet()` with the additional
935   requirement that its final argument be a `PetscViewerAndFormat`.
936 
937   Example Usage:
938 .vb
939   KSPMonitorRegister("my_monitor", PETSCVIEWERASCII, PETSC_VIEWER_ASCII_INFO_DETAIL, MyMonitor, NULL, NULL);
940 .ve
941 
942   Then, your monitor can be chosen with the procedural interface via
943 .vb
944   KSPMonitorSetFromOptions(ksp, "-ksp_monitor_my_monitor", "my_monitor", NULL)
945 .ve
946   or at runtime via the option `-ksp_monitor_my_monitor`
947 
948 .seealso: [](ch_ksp), `KSP`, `KSPMonitorSet()`, `KSPMonitorRegisterAll()`, `KSPMonitorSetFromOptions()`
949 @*/
950 PetscErrorCode KSPMonitorRegister(const char name[], PetscViewerType vtype, PetscViewerFormat format, KSPMonitorRegisterFn *monitor, KSPMonitorRegisterCreateFn *create, KSPMonitorRegisterDestroyFn *destroy)
951 {
952   char key[PETSC_MAX_PATH_LEN];
953 
954   PetscFunctionBegin;
955   PetscCall(KSPInitializePackage());
956   PetscCall(KSPMonitorMakeKey_Internal(name, vtype, format, key));
957   PetscCall(PetscFunctionListAdd(&KSPMonitorList, key, monitor));
958   if (create) PetscCall(PetscFunctionListAdd(&KSPMonitorCreateList, key, create));
959   if (destroy) PetscCall(PetscFunctionListAdd(&KSPMonitorDestroyList, key, destroy));
960   PetscFunctionReturn(PETSC_SUCCESS);
961 }
962