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