xref: /petsc/src/ksp/ksp/interface/itcreate.c (revision 6d02ad95930ea9a86b3c6dfb8d8dcdefc1c73438)
1 /*$Id: itcreate.c,v 1.206 2001/08/06 21:16:38 bsmith Exp $*/
2 /*
3      The basic KSP routines, Create, View etc. are here.
4 */
5 #include "src/ksp/ksp/kspimpl.h"      /*I "petscksp.h" I*/
6 #include "petscsys.h"
7 
8 /* Logging support */
9 int KSP_COOKIE = 0;
10 int KSP_GMRESOrthogonalization = 0, KSP_SetUp = 0, KSP_Solve = 0;
11 
12 
13 PetscTruth KSPRegisterAllCalled = PETSC_FALSE;
14 
15 #undef __FUNCT__
16 #define __FUNCT__ "KSPView"
17 /*@C
18    KSPView - Prints the KSP data structure.
19 
20    Collective on KSP
21 
22    Input Parameters:
23 +  ksp - the Krylov space context
24 -  viewer - visualization context
25 
26    Note:
27    The available visualization contexts include
28 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
29 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
30          output where only the first processor opens
31          the file.  All other processors send their
32          data to the first processor to print.
33 
34    The user can open an alternative visualization context with
35    PetscViewerASCIIOpen() - output to a specified file.
36 
37    Level: developer
38 
39 .keywords: KSP, view
40 
41 .seealso: PCView(), PetscViewerASCIIOpen()
42 @*/
43 int KSPView(KSP ksp,PetscViewer viewer)
44 {
45   char        *type;
46   int         ierr;
47   PetscTruth  isascii;
48 
49   PetscFunctionBegin;
50   PetscValidHeaderSpecific(ksp,KSP_COOKIE);
51   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(ksp->comm);
52   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE);
53   PetscCheckSameComm(ksp,viewer);
54 
55   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
56   if (isascii) {
57     ierr = KSPGetType(ksp,&type);CHKERRQ(ierr);
58     if (ksp->prefix) {
59       ierr = PetscViewerASCIIPrintf(viewer,"KSP Object:(%s)\n",ksp->prefix);CHKERRQ(ierr);
60     } else {
61       ierr = PetscViewerASCIIPrintf(viewer,"KSP Object:\n");CHKERRQ(ierr);
62     }
63     if (type) {
64       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
65     } else {
66       ierr = PetscViewerASCIIPrintf(viewer,"  type: not yet set\n");CHKERRQ(ierr);
67     }
68     if (ksp->ops->view) {
69       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
70       ierr = (*ksp->ops->view)(ksp,viewer);CHKERRQ(ierr);
71       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
72     }
73     if (ksp->guess_zero) {ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%d, initial guess is zero\n",ksp->max_it);CHKERRQ(ierr);}
74     else                 {ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%d\n", ksp->max_it);CHKERRQ(ierr);}
75     if (ksp->guess_knoll) {ierr = PetscViewerASCIIPrintf(viewer,"  using preconditioner applied to right hand side for initial guess\n");CHKERRQ(ierr);}
76     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances:  relative=%g, absolute=%g, divergence=%g\n",ksp->rtol,ksp->atol,ksp->divtol);CHKERRQ(ierr);
77     if (ksp->pc_side == PC_RIGHT)          {ierr = PetscViewerASCIIPrintf(viewer,"  right preconditioning\n");CHKERRQ(ierr);}
78     else if (ksp->pc_side == PC_SYMMETRIC) {ierr = PetscViewerASCIIPrintf(viewer,"  symmetric preconditioning\n");CHKERRQ(ierr);}
79     else                                   {ierr = PetscViewerASCIIPrintf(viewer,"  left preconditioning\n");CHKERRQ(ierr);}
80   } else {
81     if (ksp->ops->view) {
82       ierr = (*ksp->ops->view)(ksp,viewer);CHKERRQ(ierr);
83     }
84   }
85   ierr = PCView(ksp->B,viewer);CHKERRQ(ierr);
86   PetscFunctionReturn(0);
87 }
88 
89 /*
90    Contains the list of registered KSP routines
91 */
92 PetscFList KSPList = 0;
93 
94 #undef __FUNCT__
95 #define __FUNCT__ "KSPSetNormType"
96 /*@C
97    KSPSetNormType - Sets the norm that is used for convergence testing.
98 
99    Collective on KSP
100 
101    Input Parameter:
102 +  ksp - Krylov solver context
103 -  normtype - one of
104 $   KSP_NO_NORM - skips computing the norm, this should only be used if you are using
105 $                 the Krylov method as a smoother with a fixed small number of iterations.
106 $                 You must also call KSPSetConvergenceTest(ksp,KSPSkipConverged,PETSC_NULL);
107 $                 supported only by CG, Richardson, Bi-CG-stab, CR, and CGS methods.
108 $   KSP_PRECONDITIONED_NORM - the default for left preconditioned solves, uses the l2 norm
109 $                 of the preconditioned residual
110 $   KSP_UNPRECONDITIONED_NORM - uses the l2 norm of the true b - Ax residual, supported only by
111 $                 CG, CHEBYCHEV, and RICHARDSON
112 $   KSP_NATURAL_NORM - supported  by cg, cr, and cgs
113 
114 
115    Options Database Key:
116 .   -ksp_norm_type <none,preconditioned,unpreconditioned,natural>
117 
118    Notes:
119    Currently only works with the CG, Richardson, Bi-CG-stab, CR, and CGS methods.
120 
121    Level: advanced
122 
123 .keywords: KSP, create, context, norms
124 
125 .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPSkipConverged()
126 @*/
127 int KSPSetNormType(KSP ksp,KSPNormType normtype)
128 {
129 
130   PetscFunctionBegin;
131   PetscValidHeaderSpecific(ksp,KSP_COOKIE);
132   ksp->normtype = normtype;
133   if (normtype == KSP_NO_NORM) {
134     PetscLogInfo(ksp,"KSPSetNormType:Warning seting KSPNormType to skip computing the norm\n\
135   make sure you set the KSP convergence test to KSPSkipConvergence\n");
136   }
137   PetscFunctionReturn(0);
138 }
139 
140 #undef __FUNCT__
141 #define __FUNCT__ "KSPPublish_Petsc"
142 static int KSPPublish_Petsc(PetscObject obj)
143 {
144 #if defined(PETSC_HAVE_AMS)
145   KSP          v = (KSP) obj;
146   int          ierr;
147 #endif
148 
149   PetscFunctionBegin;
150 
151 #if defined(PETSC_HAVE_AMS)
152   /* if it is already published then return */
153   if (v->amem >=0) PetscFunctionReturn(0);
154 
155   ierr = PetscObjectPublishBaseBegin(obj);CHKERRQ(ierr);
156   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Iteration",&v->its,1,AMS_INT,AMS_READ,
157                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
158   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Residual",&v->rnorm,1,AMS_DOUBLE,AMS_READ,
159                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
160 
161   if (v->res_hist_max > 0) {
162     ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"ResidualNormsCount",&v->res_hist_len,1,AMS_INT,AMS_READ,
163                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
164     ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"ResidualNormsCountMax",&v->res_hist_max,1,AMS_INT,AMS_READ,
165                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
166     ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"ResidualNorms",v->res_hist,v->res_hist_max,AMS_DOUBLE,AMS_READ,
167                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
168   }
169 
170   ierr = PetscObjectPublishBaseEnd(obj);CHKERRQ(ierr);
171 #endif
172   PetscFunctionReturn(0);
173 }
174 
175 #undef __FUNCT__
176 #define __FUNCT__ "KSPSetOperators"
177 /*@
178    KSPSetOperators - Sets the matrix associated with the linear system
179    and a (possibly) different one associated with the preconditioner.
180 
181    Collective on KSP and Mat
182 
183    Input Parameters:
184 +  ksp - the KSP context
185 .  Amat - the matrix associated with the linear system
186 .  Pmat - the matrix to be used in constructing the preconditioner, usually the
187           same as Amat.
188 -  flag - flag indicating information about the preconditioner matrix structure
189    during successive linear solves.  This flag is ignored the first time a
190    linear system is solved, and thus is irrelevant when solving just one linear
191    system.
192 
193    Notes:
194    The flag can be used to eliminate unnecessary work in the preconditioner
195    during the repeated solution of linear systems of the same size.  The
196    available options are
197 $    SAME_PRECONDITIONER -
198 $      Pmat is identical during successive linear solves.
199 $      This option is intended for folks who are using
200 $      different Amat and Pmat matrices and want to reuse the
201 $      same preconditioner matrix.  For example, this option
202 $      saves work by not recomputing incomplete factorization
203 $      for ILU/ICC preconditioners.
204 $    SAME_NONZERO_PATTERN -
205 $      Pmat has the same nonzero structure during
206 $      successive linear solves.
207 $    DIFFERENT_NONZERO_PATTERN -
208 $      Pmat does not have the same nonzero structure.
209 
210     Caution:
211     If you specify SAME_NONZERO_PATTERN, PETSc believes your assertion
212     and does not check the structure of the matrix.  If you erroneously
213     claim that the structure is the same when it actually is not, the new
214     preconditioner will not function correctly.  Thus, use this optimization
215     feature carefully!
216 
217     If in doubt about whether your preconditioner matrix has changed
218     structure or not, use the flag DIFFERENT_NONZERO_PATTERN.
219 
220     Level: beginner
221 
222 .keywords: KSP, set, operators, matrix, preconditioner, linear system
223 
224 .seealso: KSPSolve(), KSPGetPC(), PCGetOperators(), PCSetOperators()
225 @*/
226 int KSPSetOperators(KSP ksp,Mat Amat,Mat Pmat,MatStructure flag)
227 {
228   PetscFunctionBegin;
229   PetscValidHeaderSpecific(ksp,KSP_COOKIE);
230   PetscValidHeaderSpecific(Amat,MAT_COOKIE);
231   PetscValidHeaderSpecific(Pmat,MAT_COOKIE);
232   PCSetOperators(ksp->B,Amat,Pmat,flag);
233   ksp->setupcalled = 0;  /* so that next solve call will call setup */
234   PetscFunctionReturn(0);
235 }
236 
237 #undef __FUNCT__
238 #define __FUNCT__ "KSPCreate"
239 /*@C
240    KSPCreate - Creates the default KSP context.
241 
242    Collective on MPI_Comm
243 
244    Input Parameter:
245 .  comm - MPI communicator
246 
247    Output Parameter:
248 .  ksp - location to put the KSP context
249 
250    Notes:
251    The default KSP type is GMRES with a restart of 30, using modified Gram-Schmidt
252    orthogonalization.
253 
254    Level: developer
255 
256 .keywords: KSP, create, context
257 
258 .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP
259 @*/
260 int KSPCreate(MPI_Comm comm,KSP *inksp)
261 {
262   KSP ksp;
263   int ierr;
264 
265   PetscFunctionBegin;
266   PetscValidPointer(inksp);
267   *inksp = 0;
268 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
269   ierr = KSPInitializePackage(PETSC_NULL);                                                               CHKERRQ(ierr);
270 #endif
271 
272   PetscHeaderCreate(ksp,_p_KSP,struct _KSPOps,KSP_COOKIE,-1,"KSP",comm,KSPDestroy,KSPView);
273   PetscLogObjectCreate(ksp);
274   *inksp             = ksp;
275   ksp->bops->publish = KSPPublish_Petsc;
276 
277   ksp->type          = -1;
278   ksp->max_it        = 10000;
279   ksp->pc_side       = PC_LEFT;
280   ksp->rtol          = 1.e-5;
281   ksp->atol          = 1.e-50;
282   ksp->divtol        = 1.e4;
283 
284   ksp->normtype            = KSP_PRECONDITIONED_NORM;
285   ksp->rnorm               = 0.0;
286   ksp->its                 = 0;
287   ksp->guess_zero          = PETSC_TRUE;
288   ksp->calc_sings          = PETSC_FALSE;
289   ksp->res_hist            = PETSC_NULL;
290   ksp->res_hist_len        = 0;
291   ksp->res_hist_max        = 0;
292   ksp->res_hist_reset      = PETSC_TRUE;
293   ksp->numbermonitors      = 0;
294   ksp->converged           = KSPDefaultConverged;
295   ksp->ops->buildsolution  = KSPDefaultBuildSolution;
296   ksp->ops->buildresidual  = KSPDefaultBuildResidual;
297 
298   ksp->ops->setfromoptions = 0;
299 
300   ksp->vec_sol         = 0;
301   ksp->vec_rhs         = 0;
302   ksp->B               = 0;
303 
304   ksp->ops->solve      = 0;
305   ksp->ops->setup      = 0;
306   ksp->ops->destroy    = 0;
307 
308   ksp->data            = 0;
309   ksp->nwork           = 0;
310   ksp->work            = 0;
311 
312   ksp->cnvP            = 0;
313 
314   ksp->reason          = KSP_CONVERGED_ITERATING;
315 
316   ksp->setupcalled     = 0;
317   ierr = PetscPublishAll(ksp);CHKERRQ(ierr);
318   ierr = PCCreate(comm,&ksp->B);CHKERRQ(ierr);
319   PetscFunctionReturn(0);
320 }
321 
322 #undef __FUNCT__
323 #define __FUNCT__ "KSPSetType"
324 /*@C
325    KSPSetType - Builds KSP for a particular solver.
326 
327    Collective on KSP
328 
329    Input Parameters:
330 +  ksp      - the Krylov space context
331 -  type - a known method
332 
333    Options Database Key:
334 .  -ksp_type  <method> - Sets the method; use -help for a list
335     of available methods (for instance, cg or gmres)
336 
337    Notes:
338    See "petsc/include/petscksp.h" for available methods (for instance,
339    KSPCG or KSPGMRES).
340 
341   Normally, it is best to use the KSPSetFromOptions() command and
342   then set the KSP type from the options database rather than by using
343   this routine.  Using the options database provides the user with
344   maximum flexibility in evaluating the many different Krylov methods.
345   The KSPSetType() routine is provided for those situations where it
346   is necessary to set the iterative solver independently of the command
347   line or options database.  This might be the case, for example, when
348   the choice of iterative solver changes during the execution of the
349   program, and the user's application is taking responsibility for
350   choosing the appropriate method.  In other words, this routine is
351   not for beginners.
352 
353   Level: intermediate
354 
355 .keywords: KSP, set, method
356 
357 .seealso: PCSetType(), KSPType
358 
359 @*/
360 int KSPSetType(KSP ksp,const KSPType type)
361 {
362   int        ierr,(*r)(KSP);
363   PetscTruth match;
364 
365   PetscFunctionBegin;
366   PetscValidHeaderSpecific(ksp,KSP_COOKIE);
367   PetscValidCharPointer(type);
368 
369   ierr = PetscTypeCompare((PetscObject)ksp,type,&match);CHKERRQ(ierr);
370   if (match) PetscFunctionReturn(0);
371 
372   if (ksp->data) {
373     /* destroy the old private KSP context */
374     ierr = (*ksp->ops->destroy)(ksp);CHKERRQ(ierr);
375     ksp->data = 0;
376   }
377   /* Get the function pointers for the iterative method requested */
378   if (!KSPRegisterAllCalled) {ierr = KSPRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
379 
380   ierr =  PetscFListFind(ksp->comm,KSPList,type,(void (**)(void)) &r);CHKERRQ(ierr);
381 
382   if (!r) SETERRQ1(1,"Unknown KSP type given: %s",type);
383 
384   ksp->setupcalled = 0;
385   ierr = (*r)(ksp);CHKERRQ(ierr);
386 
387   ierr = PetscObjectChangeTypeName((PetscObject)ksp,type);CHKERRQ(ierr);
388   PetscFunctionReturn(0);
389 }
390 
391 #undef __FUNCT__
392 #define __FUNCT__ "KSPRegisterDestroy"
393 /*@C
394    KSPRegisterDestroy - Frees the list of KSP methods that were
395    registered by KSPRegisterDynamic().
396 
397    Not Collective
398 
399    Level: advanced
400 
401 .keywords: KSP, register, destroy
402 
403 .seealso: KSPRegisterDynamic(), KSPRegisterAll()
404 @*/
405 int KSPRegisterDestroy(void)
406 {
407   int ierr;
408 
409   PetscFunctionBegin;
410   if (KSPList) {
411     ierr = PetscFListDestroy(&KSPList);CHKERRQ(ierr);
412     KSPList = 0;
413   }
414   KSPRegisterAllCalled = PETSC_FALSE;
415   PetscFunctionReturn(0);
416 }
417 
418 #undef __FUNCT__
419 #define __FUNCT__ "KSPGetType"
420 /*@C
421    KSPGetType - Gets the KSP type as a string from the KSP object.
422 
423    Not Collective
424 
425    Input Parameter:
426 .  ksp - Krylov context
427 
428    Output Parameter:
429 .  name - name of KSP method
430 
431    Level: intermediate
432 
433 .keywords: KSP, get, method, name
434 
435 .seealso: KSPSetType()
436 @*/
437 int KSPGetType(KSP ksp,KSPType *type)
438 {
439   PetscFunctionBegin;
440   PetscValidHeaderSpecific(ksp,KSP_COOKIE);
441   *type = ksp->type_name;
442   PetscFunctionReturn(0);
443 }
444 
445 #undef __FUNCT__
446 #define __FUNCT__ "KSPSetFromOptions"
447 /*@
448    KSPSetFromOptions - Sets KSP options from the options database.
449    This routine must be called before KSPSetUp() if the user is to be
450    allowed to set the Krylov type.
451 
452    Collective on KSP
453 
454    Input Parameters:
455 .  ksp - the Krylov space context
456 
457    Options Database Keys:
458 +   -ksp_max_it - maximum number of linear iterations
459 .   -ksp_rtol rtol - relative tolerance used in default determination of convergence, i.e.
460                 if residual norm decreases by this factor than convergence is declared
461 .   -ksp_atol atol - absolute tolerance used in default convergence test, i.e. if residual
462                 norm is less than this then convergence is declared
463 .   -ksp_divtol tol - if residual norm increases by this factor than divergence is declared
464 .   -ksp_norm_type - none - skip norms used in convergence tests (useful only when not using
465 $                       convergence test (say you always want to run with 5 iterations) to
466 $                       save on communication overhead
467 $                    preconditioned - default for left preconditioning
468 $                    unpreconditioned - see KSPSetNormType()
469 $                    natural - see KSPSetNormType()
470 .   -ksp_knoll - compute initial guess by applying the preconditioner to the right hand side
471 .   -ksp_cancelmonitors - cancel all previous convergene monitor routines set
472 .   -ksp_monitor - print residual norm at each iteration
473 .   -ksp_xmonitor - plot residual norm at each iteration
474 .   -ksp_vecmonitor - plot solution at each iteration
475 -   -ksp_singmonitor - monitor extremem singular values at each iteration
476 
477    Notes:
478    To see all options, run your program with the -help option
479    or consult the users manual.
480 
481    Level: developer
482 
483 .keywords: KSP, set, from, options, database
484 
485 .seealso:
486 @*/
487 int KSPSetFromOptions(KSP ksp)
488 {
489   int        ierr,indx;
490   char       type[256],*stype[] = {"none","preconditioned","unpreconditioned","natural"};
491   PetscTruth flg;
492 
493   PetscFunctionBegin;
494   PetscValidHeaderSpecific(ksp,KSP_COOKIE);
495   ierr = PCSetFromOptions(ksp->B);CHKERRQ(ierr);
496 
497   if (!KSPRegisterAllCalled) {ierr = KSPRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
498   ierr = PetscOptionsBegin(ksp->comm,ksp->prefix,"Krylov Method (KSP) Options","KSP");CHKERRQ(ierr);
499     ierr = PetscOptionsList("-ksp_type","Krylov method","KSPSetType",KSPList,(char*)(ksp->type_name?ksp->type_name:KSPGMRES),type,256,&flg);CHKERRQ(ierr);
500     if (flg) {
501       ierr = KSPSetType(ksp,type);CHKERRQ(ierr);
502     }
503     /*
504       Set the type if it was never set.
505     */
506     if (!ksp->type_name) {
507       ierr = KSPSetType(ksp,KSPGMRES);CHKERRQ(ierr);
508     }
509 
510     ierr = PetscOptionsInt("-ksp_max_it","Maximum number of iterations","KSPSetTolerances",ksp->max_it,&ksp->max_it,PETSC_NULL);CHKERRQ(ierr);
511     ierr = PetscOptionsReal("-ksp_rtol","Relative decrease in residual norm","KSPSetTolerances",ksp->rtol,&ksp->rtol,PETSC_NULL);CHKERRQ(ierr);
512     ierr = PetscOptionsReal("-ksp_atol","Absolute value of residual norm","KSPSetTolerances",ksp->atol,&ksp->atol,PETSC_NULL);CHKERRQ(ierr);
513     ierr = PetscOptionsReal("-ksp_divtol","Residual norm increase cause divergence","KSPSetTolerances",ksp->divtol,&ksp->divtol,PETSC_NULL);CHKERRQ(ierr);
514     ierr = PetscOptionsLogical("-ksp_knoll","Use preconditioner applied to b for initial guess","KSPSetInitialGuessKnoll",ksp->guess_knoll,
515                                   &ksp->guess_knoll,PETSC_NULL);CHKERRQ(ierr);
516 
517     ierr = PetscOptionsEList("-ksp_norm_type","KSP Norm type","KSPSetNormType",stype,4,"preconditioned",&indx,&flg);CHKERRQ(ierr);
518     if (flg) {
519       switch (indx) {
520       case 0:
521         ierr = KSPSetNormType(ksp,KSP_NO_NORM);CHKERRQ(ierr);
522         ierr = KSPSetConvergenceTest(ksp,KSPSkipConverged,0);CHKERRQ(ierr);
523         break;
524       case 1:
525         ierr = KSPSetNormType(ksp,KSP_PRECONDITIONED_NORM);CHKERRQ(ierr);
526         break;
527       case 2:
528         ierr = KSPSetNormType(ksp,KSP_UNPRECONDITIONED_NORM);CHKERRQ(ierr);
529         break;
530       case 3:
531         ierr = KSPSetNormType(ksp,KSP_NATURAL_NORM);CHKERRQ(ierr);
532         break;
533       }
534     }
535 
536     ierr = PetscOptionsName("-ksp_diagonal_scale","Diagonal scale matrix before building preconditioner","KSPSetDiagonalScale",&flg);CHKERRQ(ierr);
537     if (flg) {
538       ierr = KSPSetDiagonalScale(ksp,PETSC_TRUE);CHKERRQ(ierr);
539     }
540     ierr = PetscOptionsName("-ksp_diagonal_scale_fix","Fix diagonaled scaled matrix after solve","KSPSetDiagonalScaleFix",&flg);CHKERRQ(ierr);
541     if (flg) {
542       ierr = KSPSetDiagonalScaleFix(ksp,PETSC_TRUE);CHKERRQ(ierr);
543     }
544 
545     /*
546       Prints reason for convergence or divergence of each linear solve
547     */
548     ierr = PetscOptionsName("-ksp_converged_reason","Print reason for converged or diverged","KSPSolve",&flg);CHKERRQ(ierr);
549     if (flg) {
550       ksp->printreason = PETSC_TRUE;
551     }
552 
553     ierr = PetscOptionsName("-ksp_cancelmonitors","Remove any hardwired monitor routines","KSPClearMonitor",&flg);CHKERRQ(ierr);
554     /* -----------------------------------------------------------------------*/
555     /*
556       Cancels all monitors hardwired into code before call to KSPSetFromOptions()
557     */
558     if (flg) {
559       ierr = KSPClearMonitor(ksp);CHKERRQ(ierr);
560     }
561     /*
562       Prints preconditioned residual norm at each iteration
563     */
564     ierr = PetscOptionsName("-ksp_monitor","Monitor preconditioned residual norm","KSPSetMonitor",&flg);CHKERRQ(ierr);
565     if (flg) {
566       ierr = KSPSetMonitor(ksp,KSPDefaultMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
567     }
568     /*
569       Plots the vector solution
570     */
571     ierr = PetscOptionsName("-ksp_vecmonitor","Monitor solution graphically","KSPSetMonitor",&flg);CHKERRQ(ierr);
572     if (flg) {
573       ierr = KSPSetMonitor(ksp,KSPVecViewMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
574     }
575     /*
576       Prints preconditioned and true residual norm at each iteration
577     */
578     ierr = PetscOptionsName("-ksp_truemonitor","Monitor true (unpreconditioned) residual norm","KSPSetMonitor",&flg);CHKERRQ(ierr);
579     if (flg) {
580       ierr = KSPSetMonitor(ksp,KSPTrueMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
581     }
582     /*
583       Prints extreme eigenvalue estimates at each iteration
584     */
585     ierr = PetscOptionsName("-ksp_singmonitor","Monitor singular values","KSPSetMonitor",&flg);CHKERRQ(ierr);
586     if (flg) {
587       ierr = KSPSetComputeSingularValues(ksp,PETSC_TRUE);CHKERRQ(ierr);
588       ierr = KSPSetMonitor(ksp,KSPSingularValueMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
589     }
590     /*
591       Prints preconditioned residual norm with fewer digits
592     */
593     ierr = PetscOptionsName("-ksp_smonitor","Monitor preconditioned residual norm with fewer digitis","KSPSetMonitor",&flg);CHKERRQ(ierr);
594     if (flg) {
595       ierr = KSPSetMonitor(ksp,KSPDefaultSMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
596     }
597     /*
598       Graphically plots preconditioned residual norm
599     */
600     ierr = PetscOptionsName("-ksp_xmonitor","Monitor graphically preconditioned residual norm","KSPSetMonitor",&flg);CHKERRQ(ierr);
601     if (flg) {
602       ierr = KSPSetMonitor(ksp,KSPLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
603     }
604     /*
605       Graphically plots preconditioned and true residual norm
606     */
607     ierr = PetscOptionsName("-ksp_xtruemonitor","Monitor graphically true residual norm","KSPSetMonitor",&flg);CHKERRQ(ierr);
608     if (flg){
609       ierr = KSPSetMonitor(ksp,KSPLGTrueMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
610     }
611 
612     /* -----------------------------------------------------------------------*/
613 
614     ierr = PetscOptionsLogicalGroupBegin("-ksp_left_pc","Use left preconditioning","KSPSetPreconditionerSide",&flg);CHKERRQ(ierr);
615     if (flg) { ierr = KSPSetPreconditionerSide(ksp,PC_LEFT);CHKERRQ(ierr); }
616     ierr = PetscOptionsLogicalGroup("-ksp_right_pc","Use right preconditioning","KSPSetPreconditionerSide",&flg);CHKERRQ(ierr);
617     if (flg) { ierr = KSPSetPreconditionerSide(ksp,PC_RIGHT);CHKERRQ(ierr);}
618     ierr = PetscOptionsLogicalGroupEnd("-ksp_symmetric_pc","Use symmetric (factorized) preconditioning","KSPSetPreconditionerSide",&flg);CHKERRQ(ierr);
619     if (flg) { ierr = KSPSetPreconditionerSide(ksp,PC_SYMMETRIC);CHKERRQ(ierr);}
620 
621     ierr = PetscOptionsName("-ksp_compute_singularvalues","Compute singular values of preconditioned operator","KSPSetComputeSingularValues",&flg);CHKERRQ(ierr);
622     if (flg) { ierr = KSPSetComputeSingularValues(ksp,PETSC_TRUE);CHKERRQ(ierr); }
623     ierr = PetscOptionsName("-ksp_compute_eigenvalues","Compute eigenvalues of preconditioned operator","KSPSetComputeSingularValues",&flg);CHKERRQ(ierr);
624     if (flg) { ierr = KSPSetComputeSingularValues(ksp,PETSC_TRUE);CHKERRQ(ierr); }
625     ierr = PetscOptionsName("-ksp_plot_eigenvalues","Scatter plot extreme eigenvalues","KSPSetComputeSingularValues",&flg);CHKERRQ(ierr);
626     if (flg) { ierr = KSPSetComputeSingularValues(ksp,PETSC_TRUE);CHKERRQ(ierr); }
627 
628     if (ksp->ops->setfromoptions) {
629       ierr = (*ksp->ops->setfromoptions)(ksp);CHKERRQ(ierr);
630     }
631     ierr = PetscOptionsName("-ksp_view","View linear solver parameters","KSPView",&flg);CHKERRQ(ierr);
632     if (flg) {
633       ierr = KSPView(ksp,PETSC_VIEWER_STDOUT_(ksp->comm));CHKERRQ(ierr);
634     }
635   ierr = PetscOptionsEnd();CHKERRQ(ierr);
636   PetscFunctionReturn(0);
637 }
638 
639 
640 #undef __FUNCT__
641 #define __FUNCT__ "KSPRegister"
642 /*@C
643   KSPRegister - See KSPRegisterDynamic()
644 
645   Level: advanced
646 @*/
647 int KSPRegister(const char sname[],const char path[],const char name[],int (*function)(KSP))
648 {
649   int  ierr;
650   char fullname[256];
651 
652   PetscFunctionBegin;
653   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
654   ierr = PetscFListAdd(&KSPList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
655   PetscFunctionReturn(0);
656 }
657