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