xref: /petsc/src/ksp/ksp/interface/itcreate.c (revision 292f8084fb157dadf9a2ae26c5bd14368ed7ffcb)
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 #if defined(PETSC_HAVE_AMS)
147   KSP          v = (KSP) obj;
148   PetscErrorCode ierr;
149 #endif
150 
151   PetscFunctionBegin;
152 
153 #if defined(PETSC_HAVE_AMS)
154   /* if it is already published then return */
155   if (v->amem >=0) PetscFunctionReturn(0);
156 
157   ierr = PetscObjectPublishBaseBegin(obj);CHKERRQ(ierr);
158   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Iteration",&v->its,1,AMS_INT,AMS_READ,
159                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
160   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Residual",&v->rnorm,1,AMS_DOUBLE,AMS_READ,
161                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
162 
163   if (v->res_hist_max > 0) {
164     ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"ResidualNormsCount",&v->res_hist_len,1,AMS_INT,AMS_READ,
165                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
166     ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"ResidualNormsCountMax",&v->res_hist_max,1,AMS_INT,AMS_READ,
167                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
168     ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"ResidualNorms",v->res_hist,v->res_hist_max,AMS_DOUBLE,AMS_READ,
169                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
170   }
171 
172   ierr = PetscObjectPublishBaseEnd(obj);CHKERRQ(ierr);
173 #endif
174   PetscFunctionReturn(0);
175 }
176 
177 #undef __FUNCT__
178 #define __FUNCT__ "KSPSetOperators"
179 /*@
180    KSPSetOperators - Sets the matrix associated with the linear system
181    and a (possibly) different one associated with the preconditioner.
182 
183    Collective on KSP and Mat
184 
185    Input Parameters:
186 +  ksp - the KSP context
187 .  Amat - the matrix associated with the linear system
188 .  Pmat - the matrix to be used in constructing the preconditioner, usually the
189           same as Amat.
190 -  flag - flag indicating information about the preconditioner matrix structure
191    during successive linear solves.  This flag is ignored the first time a
192    linear system is solved, and thus is irrelevant when solving just one linear
193    system.
194 
195    Notes:
196    The flag can be used to eliminate unnecessary work in the preconditioner
197    during the repeated solution of linear systems of the same size.  The
198    available options are
199 $    SAME_PRECONDITIONER -
200 $      Pmat is identical during successive linear solves.
201 $      This option is intended for folks who are using
202 $      different Amat and Pmat matrices and want to reuse the
203 $      same preconditioner matrix.  For example, this option
204 $      saves work by not recomputing incomplete factorization
205 $      for ILU/ICC preconditioners.
206 $    SAME_NONZERO_PATTERN -
207 $      Pmat has the same nonzero structure during
208 $      successive linear solves.
209 $    DIFFERENT_NONZERO_PATTERN -
210 $      Pmat does not have the same nonzero structure.
211 
212     Caution:
213     If you specify SAME_NONZERO_PATTERN, PETSc believes your assertion
214     and does not check the structure of the matrix.  If you erroneously
215     claim that the structure is the same when it actually is not, the new
216     preconditioner will not function correctly.  Thus, use this optimization
217     feature carefully!
218 
219     If in doubt about whether your preconditioner matrix has changed
220     structure or not, use the flag DIFFERENT_NONZERO_PATTERN.
221 
222     Level: beginner
223 
224 .keywords: KSP, set, operators, matrix, preconditioner, linear system
225 
226 .seealso: KSPSolve(), KSPGetPC(), PCGetOperators(), PCSetOperators()
227 @*/
228 PetscErrorCode KSPSetOperators(KSP ksp,Mat Amat,Mat Pmat,MatStructure flag)
229 {
230   PetscErrorCode ierr;
231 
232   PetscFunctionBegin;
233   PetscValidHeaderSpecific(ksp,KSP_COOKIE,1);
234   PetscValidHeaderSpecific(Amat,MAT_COOKIE,2);
235   PetscValidHeaderSpecific(Pmat,MAT_COOKIE,3);
236   ierr = PCSetOperators(ksp->pc,Amat,Pmat,flag);CHKERRQ(ierr);
237   if (ksp->setupcalled > 1) ksp->setupcalled = 1;  /* so that next solve call will call setup */
238   PetscFunctionReturn(0);
239 }
240 
241 #undef __FUNCT__
242 #define __FUNCT__ "KSPCreate"
243 /*@C
244    KSPCreate - Creates the default KSP context.
245 
246    Collective on MPI_Comm
247 
248    Input Parameter:
249 .  comm - MPI communicator
250 
251    Output Parameter:
252 .  ksp - location to put the KSP context
253 
254    Notes:
255    The default KSP type is GMRES with a restart of 30, using modified Gram-Schmidt
256    orthogonalization.
257 
258    Level: beginner
259 
260 .keywords: KSP, create, context
261 
262 .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP
263 @*/
264 PetscErrorCode KSPCreate(MPI_Comm comm,KSP *inksp)
265 {
266   KSP ksp;
267   PetscErrorCode ierr;
268 
269   PetscFunctionBegin;
270   PetscValidPointer(inksp,2);
271   *inksp = 0;
272 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
273   ierr = KSPInitializePackage(PETSC_NULL);CHKERRQ(ierr);
274 #endif
275 
276   PetscHeaderCreate(ksp,_p_KSP,struct _KSPOps,KSP_COOKIE,-1,"KSP",comm,KSPDestroy,KSPView);
277   PetscLogObjectCreate(ksp);
278   *inksp             = ksp;
279   ksp->bops->publish = KSPPublish_Petsc;
280 
281   ksp->type          = -1;
282   ksp->max_it        = 10000;
283   ksp->pc_side       = PC_LEFT;
284   ksp->rtol          = 1.e-5;
285   ksp->atol          = 1.e-50;
286   ksp->divtol        = 1.e4;
287 
288   ksp->normtype            = KSP_PRECONDITIONED_NORM;
289   ksp->rnorm               = 0.0;
290   ksp->its                 = 0;
291   ksp->guess_zero          = PETSC_TRUE;
292   ksp->calc_sings          = PETSC_FALSE;
293   ksp->res_hist            = PETSC_NULL;
294   ksp->res_hist_len        = 0;
295   ksp->res_hist_max        = 0;
296   ksp->res_hist_reset      = PETSC_TRUE;
297   ksp->numbermonitors      = 0;
298   ksp->converged           = KSPDefaultConverged;
299   ksp->ops->buildsolution  = KSPDefaultBuildSolution;
300   ksp->ops->buildresidual  = KSPDefaultBuildResidual;
301 
302   ksp->ops->setfromoptions = 0;
303 
304   ksp->vec_sol         = 0;
305   ksp->vec_rhs         = 0;
306   ksp->pc              = 0;
307 
308   ksp->ops->solve      = 0;
309   ksp->ops->setup      = 0;
310   ksp->ops->destroy    = 0;
311 
312   ksp->data            = 0;
313   ksp->nwork           = 0;
314   ksp->work            = 0;
315 
316   ksp->cnvP            = 0;
317 
318   ksp->reason          = KSP_CONVERGED_ITERATING;
319 
320   ksp->setupcalled     = 0;
321   ierr = PetscPublishAll(ksp);CHKERRQ(ierr);
322   ierr = PCCreate(comm,&ksp->pc);CHKERRQ(ierr);
323   PetscFunctionReturn(0);
324 }
325 
326 #undef __FUNCT__
327 #define __FUNCT__ "KSPSetType"
328 /*@C
329    KSPSetType - Builds KSP for a particular solver.
330 
331    Collective on KSP
332 
333    Input Parameters:
334 +  ksp      - the Krylov space context
335 -  type - a known method
336 
337    Options Database Key:
338 .  -ksp_type  <method> - Sets the method; use -help for a list
339     of available methods (for instance, cg or gmres)
340 
341    Notes:
342    See "petsc/include/petscksp.h" for available methods (for instance,
343    KSPCG or KSPGMRES).
344 
345   Normally, it is best to use the KSPSetFromOptions() command and
346   then set the KSP type from the options database rather than by using
347   this routine.  Using the options database provides the user with
348   maximum flexibility in evaluating the many different Krylov methods.
349   The KSPSetType() routine is provided for those situations where it
350   is necessary to set the iterative solver independently of the command
351   line or options database.  This might be the case, for example, when
352   the choice of iterative solver changes during the execution of the
353   program, and the user's application is taking responsibility for
354   choosing the appropriate method.  In other words, this routine is
355   not for beginners.
356 
357   Level: intermediate
358 
359 .keywords: KSP, set, method
360 
361 .seealso: PCSetType(), KSPType
362 
363 @*/
364 PetscErrorCode KSPSetType(KSP ksp,const KSPType type)
365 {
366   PetscErrorCode ierr,(*r)(KSP);
367   PetscTruth match;
368 
369   PetscFunctionBegin;
370   PetscValidHeaderSpecific(ksp,KSP_COOKIE,1);
371   PetscValidCharPointer(type,2);
372 
373   ierr = PetscTypeCompare((PetscObject)ksp,type,&match);CHKERRQ(ierr);
374   if (match) PetscFunctionReturn(0);
375 
376   if (ksp->data) {
377     /* destroy the old private KSP context */
378     ierr = (*ksp->ops->destroy)(ksp);CHKERRQ(ierr);
379     ksp->data = 0;
380   }
381   /* Get the function pointers for the iterative method requested */
382   if (!KSPRegisterAllCalled) {ierr = KSPRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
383   ierr =  PetscFListFind(ksp->comm,KSPList,type,(void (**)(void)) &r);CHKERRQ(ierr);
384   if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown KSP type given: %s",type);
385   ksp->setupcalled = 0;
386   ierr = (*r)(ksp);CHKERRQ(ierr);
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 PetscErrorCode KSPRegisterDestroy(void)
406 {
407   PetscErrorCode 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 PetscErrorCode KSPGetType(KSP ksp,KSPType *type)
438 {
439   PetscFunctionBegin;
440   PetscValidHeaderSpecific(ksp,KSP_COOKIE,1);
441   PetscValidPointer(type,2);
442   *type = ksp->type_name;
443   PetscFunctionReturn(0);
444 }
445 
446 #undef __FUNCT__
447 #define __FUNCT__ "KSPSetFromOptions"
448 /*@
449    KSPSetFromOptions - Sets KSP options from the options database.
450    This routine must be called before KSPSetUp() if the user is to be
451    allowed to set the Krylov type.
452 
453    Collective on KSP
454 
455    Input Parameters:
456 .  ksp - the Krylov space context
457 
458    Options Database Keys:
459 +   -ksp_max_it - maximum number of linear iterations
460 .   -ksp_rtol rtol - relative tolerance used in default determination of convergence, i.e.
461                 if residual norm decreases by this factor than convergence is declared
462 .   -ksp_atol atol - absolute tolerance used in default convergence test, i.e. if residual
463                 norm is less than this then convergence is declared
464 .   -ksp_divtol tol - if residual norm increases by this factor than divergence is declared
465 .   -ksp_norm_type - none - skip norms used in convergence tests (useful only when not using
466 $                       convergence test (say you always want to run with 5 iterations) to
467 $                       save on communication overhead
468 $                    preconditioned - default for left preconditioning
469 $                    unpreconditioned - see KSPSetNormType()
470 $                    natural - see KSPSetNormType()
471 .    -ksp_constant_null_space - assume the operator (matrix) has the constant vector in its null space
472 .    -ksp_test_null_space - tests the null space set with KSPSetNullSpace() to see if it truly is a null space
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: beginner
485 
486 .keywords: KSP, set, from, options, database
487 
488 .seealso:
489 @*/
490 PetscErrorCode KSPSetFromOptions(KSP ksp)
491 {
492   PetscErrorCode ierr;
493   int indx;
494   char       type[256];
495   const char *stype[] = {"none","preconditioned","unpreconditioned","natural"};
496   PetscTruth flg;
497 
498   PetscFunctionBegin;
499   PetscValidHeaderSpecific(ksp,KSP_COOKIE,1);
500   ierr = PCSetFromOptions(ksp->pc);CHKERRQ(ierr);
501 
502   if (!KSPRegisterAllCalled) {ierr = KSPRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
503   ierr = PetscOptionsBegin(ksp->comm,ksp->prefix,"Krylov Method (KSP) Options","KSP");CHKERRQ(ierr);
504     ierr = PetscOptionsList("-ksp_type","Krylov method","KSPSetType",KSPList,(char*)(ksp->type_name?ksp->type_name:KSPGMRES),type,256,&flg);CHKERRQ(ierr);
505     if (flg) {
506       ierr = KSPSetType(ksp,type);CHKERRQ(ierr);
507     }
508     /*
509       Set the type if it was never set.
510     */
511     if (!ksp->type_name) {
512       ierr = KSPSetType(ksp,KSPGMRES);CHKERRQ(ierr);
513     }
514 
515     ierr = PetscOptionsInt("-ksp_max_it","Maximum number of iterations","KSPSetTolerances",ksp->max_it,&ksp->max_it,PETSC_NULL);CHKERRQ(ierr);
516     ierr = PetscOptionsReal("-ksp_rtol","Relative decrease in residual norm","KSPSetTolerances",ksp->rtol,&ksp->rtol,PETSC_NULL);CHKERRQ(ierr);
517     ierr = PetscOptionsReal("-ksp_atol","Absolute value of residual norm","KSPSetTolerances",ksp->atol,&ksp->atol,PETSC_NULL);CHKERRQ(ierr);
518     ierr = PetscOptionsReal("-ksp_divtol","Residual norm increase cause divergence","KSPSetTolerances",ksp->divtol,&ksp->divtol,PETSC_NULL);CHKERRQ(ierr);
519     ierr = PetscOptionsLogical("-ksp_knoll","Use preconditioner applied to b for initial guess","KSPSetInitialGuessKnoll",ksp->guess_knoll,
520                                   &ksp->guess_knoll,PETSC_NULL);CHKERRQ(ierr);
521 
522     ierr = PetscOptionsEList("-ksp_norm_type","KSP Norm type","KSPSetNormType",stype,4,"preconditioned",&indx,&flg);CHKERRQ(ierr);
523     if (flg) {
524       switch (indx) {
525       case 0:
526         ierr = KSPSetNormType(ksp,KSP_NO_NORM);CHKERRQ(ierr);
527         ierr = KSPSetConvergenceTest(ksp,KSPSkipConverged,0);CHKERRQ(ierr);
528         break;
529       case 1:
530         ierr = KSPSetNormType(ksp,KSP_PRECONDITIONED_NORM);CHKERRQ(ierr);
531         break;
532       case 2:
533         ierr = KSPSetNormType(ksp,KSP_UNPRECONDITIONED_NORM);CHKERRQ(ierr);
534         break;
535       case 3:
536         ierr = KSPSetNormType(ksp,KSP_NATURAL_NORM);CHKERRQ(ierr);
537         break;
538       }
539     }
540 
541     ierr = PetscOptionsName("-ksp_diagonal_scale","Diagonal scale matrix before building preconditioner","KSPSetDiagonalScale",&flg);CHKERRQ(ierr);
542     if (flg) {
543       ierr = KSPSetDiagonalScale(ksp,PETSC_TRUE);CHKERRQ(ierr);
544     }
545     ierr = PetscOptionsName("-ksp_diagonal_scale_fix","Fix diagonaled scaled matrix after solve","KSPSetDiagonalScaleFix",&flg);CHKERRQ(ierr);
546     if (flg) {
547       ierr = KSPSetDiagonalScaleFix(ksp,PETSC_TRUE);CHKERRQ(ierr);
548     }
549 
550 
551     ierr = PetscOptionsName("-ksp_constant_null_space","Add constant null space to Krylov solver","KSPSetNullSpace",&flg);CHKERRQ(ierr);
552     if (flg) {
553       MatNullSpace nsp;
554 
555       ierr = MatNullSpaceCreate(ksp->comm,1,0,0,&nsp);CHKERRQ(ierr);
556       ierr = KSPSetNullSpace(ksp,nsp);CHKERRQ(ierr);
557       ierr = MatNullSpaceDestroy(nsp);CHKERRQ(ierr);
558     }
559 
560     /* option is actually checked in KSPSetUp() */
561     if (ksp->nullsp) {
562       ierr = PetscOptionsName("-ksp_test_null_space","Is provided null space correct","None",&flg);CHKERRQ(ierr);
563     }
564 
565     /*
566       Prints reason for convergence or divergence of each linear solve
567     */
568     ierr = PetscOptionsName("-ksp_converged_reason","Print reason for converged or diverged","KSPSolve",&flg);CHKERRQ(ierr);
569     if (flg) {
570       ksp->printreason = PETSC_TRUE;
571     }
572 
573     ierr = PetscOptionsName("-ksp_cancelmonitors","Remove any hardwired monitor routines","KSPClearMonitor",&flg);CHKERRQ(ierr);
574     /* -----------------------------------------------------------------------*/
575     /*
576       Cancels all monitors hardwired into code before call to KSPSetFromOptions()
577     */
578     if (flg) {
579       ierr = KSPClearMonitor(ksp);CHKERRQ(ierr);
580     }
581     /*
582       Prints preconditioned residual norm at each iteration
583     */
584     ierr = PetscOptionsName("-ksp_monitor","Monitor preconditioned residual norm","KSPSetMonitor",&flg);CHKERRQ(ierr);
585     if (flg) {
586       ierr = KSPSetMonitor(ksp,KSPDefaultMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
587     }
588     /*
589       Plots the vector solution
590     */
591     ierr = PetscOptionsName("-ksp_vecmonitor","Monitor solution graphically","KSPSetMonitor",&flg);CHKERRQ(ierr);
592     if (flg) {
593       ierr = KSPSetMonitor(ksp,KSPVecViewMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
594     }
595     /*
596       Prints preconditioned and true residual norm at each iteration
597     */
598     ierr = PetscOptionsName("-ksp_truemonitor","Monitor true (unpreconditioned) residual norm","KSPSetMonitor",&flg);CHKERRQ(ierr);
599     if (flg) {
600       ierr = KSPSetMonitor(ksp,KSPTrueMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
601     }
602     /*
603       Prints extreme eigenvalue estimates at each iteration
604     */
605     ierr = PetscOptionsName("-ksp_singmonitor","Monitor singular values","KSPSetMonitor",&flg);CHKERRQ(ierr);
606     if (flg) {
607       ierr = KSPSetComputeSingularValues(ksp,PETSC_TRUE);CHKERRQ(ierr);
608       ierr = KSPSetMonitor(ksp,KSPSingularValueMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
609     }
610     /*
611       Prints preconditioned residual norm with fewer digits
612     */
613     ierr = PetscOptionsName("-ksp_smonitor","Monitor preconditioned residual norm with fewer digitis","KSPSetMonitor",&flg);CHKERRQ(ierr);
614     if (flg) {
615       ierr = KSPSetMonitor(ksp,KSPDefaultSMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
616     }
617     /*
618       Graphically plots preconditioned residual norm
619     */
620     ierr = PetscOptionsName("-ksp_xmonitor","Monitor graphically preconditioned residual norm","KSPSetMonitor",&flg);CHKERRQ(ierr);
621     if (flg) {
622       ierr = KSPSetMonitor(ksp,KSPLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
623     }
624     /*
625       Graphically plots preconditioned and true residual norm
626     */
627     ierr = PetscOptionsName("-ksp_xtruemonitor","Monitor graphically true residual norm","KSPSetMonitor",&flg);CHKERRQ(ierr);
628     if (flg){
629       ierr = KSPSetMonitor(ksp,KSPLGTrueMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
630     }
631 
632     /* -----------------------------------------------------------------------*/
633 
634     ierr = PetscOptionsLogicalGroupBegin("-ksp_left_pc","Use left preconditioning","KSPSetPreconditionerSide",&flg);CHKERRQ(ierr);
635     if (flg) { ierr = KSPSetPreconditionerSide(ksp,PC_LEFT);CHKERRQ(ierr); }
636     ierr = PetscOptionsLogicalGroup("-ksp_right_pc","Use right preconditioning","KSPSetPreconditionerSide",&flg);CHKERRQ(ierr);
637     if (flg) { ierr = KSPSetPreconditionerSide(ksp,PC_RIGHT);CHKERRQ(ierr);}
638     ierr = PetscOptionsLogicalGroupEnd("-ksp_symmetric_pc","Use symmetric (factorized) preconditioning","KSPSetPreconditionerSide",&flg);CHKERRQ(ierr);
639     if (flg) { ierr = KSPSetPreconditionerSide(ksp,PC_SYMMETRIC);CHKERRQ(ierr);}
640 
641     ierr = PetscOptionsName("-ksp_compute_singularvalues","Compute singular values of preconditioned operator","KSPSetComputeSingularValues",&flg);CHKERRQ(ierr);
642     if (flg) { ierr = KSPSetComputeSingularValues(ksp,PETSC_TRUE);CHKERRQ(ierr); }
643     ierr = PetscOptionsName("-ksp_compute_eigenvalues","Compute eigenvalues of preconditioned operator","KSPSetComputeSingularValues",&flg);CHKERRQ(ierr);
644     if (flg) { ierr = KSPSetComputeSingularValues(ksp,PETSC_TRUE);CHKERRQ(ierr); }
645     ierr = PetscOptionsName("-ksp_plot_eigenvalues","Scatter plot extreme eigenvalues","KSPSetComputeSingularValues",&flg);CHKERRQ(ierr);
646     if (flg) { ierr = KSPSetComputeSingularValues(ksp,PETSC_TRUE);CHKERRQ(ierr); }
647 
648     if (ksp->ops->setfromoptions) {
649       ierr = (*ksp->ops->setfromoptions)(ksp);CHKERRQ(ierr);
650     }
651     ierr = PetscOptionsName("-ksp_view","View linear solver parameters","KSPView",&flg);CHKERRQ(ierr);
652   ierr = PetscOptionsEnd();CHKERRQ(ierr);
653   PetscFunctionReturn(0);
654 }
655 
656 #undef __FUNCT__
657 #define __FUNCT__ "KSPRegister"
658 /*@C
659   KSPRegister - See KSPRegisterDynamic()
660 
661   Level: advanced
662 @*/
663 PetscErrorCode KSPRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(KSP))
664 {
665   PetscErrorCode ierr;
666   char fullname[PETSC_MAX_PATH_LEN];
667 
668   PetscFunctionBegin;
669   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
670   ierr = PetscFListAdd(&KSPList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
671   PetscFunctionReturn(0);
672 }
673 
674 #undef __FUNCT__
675 #define __FUNCT__ "KSPSetNullSpace"
676 /*@C
677   KSPSetNullSpace - Sets the null space of the operator
678 
679   Collective on KSP
680 
681   Input Parameters:
682 +  ksp - the Krylov space object
683 -  nullsp - the null space of the operator
684 
685   Level: advanced
686 
687 .seealso: KSPSetOperators(), MatNullSpaceCreate(), KSPGetNullSpace()
688 @*/
689 PetscErrorCode KSPSetNullSpace(KSP ksp,MatNullSpace nullsp)
690 {
691   PetscErrorCode ierr;
692 
693   PetscFunctionBegin;
694   if (ksp->nullsp) {
695     ierr = MatNullSpaceDestroy(ksp->nullsp);CHKERRQ(ierr);
696   }
697   ksp->nullsp = nullsp;
698   ierr = PetscObjectReference((PetscObject)ksp->nullsp);CHKERRQ(ierr);
699   PetscFunctionReturn(0);
700 }
701 
702 #undef __FUNCT__
703 #define __FUNCT__ "KSPGetNullSpace"
704 /*@C
705   KSPGetNullSpace - Gets the null space of the operator
706 
707   Collective on KSP
708 
709   Input Parameters:
710 +  ksp - the Krylov space object
711 -  nullsp - the null space of the operator
712 
713   Level: advanced
714 
715 .seealso: KSPSetOperators(), MatNullSpaceCreate(), KSPSetNullSpace()
716 @*/
717 PetscErrorCode KSPGetNullSpace(KSP ksp,MatNullSpace *nullsp)
718 {
719   PetscFunctionBegin;
720   *nullsp = ksp->nullsp;
721   PetscFunctionReturn(0);
722 }
723 
724