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