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