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