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