xref: /petsc/src/ksp/ksp/interface/itcreate.c (revision 4bf112e79d3f87ddad6dbf7b300a6f889c1e7d34)
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()
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__ "KSPCreate"
216 /*@C
217    KSPCreate - Creates the default KSP context.
218 
219    Collective on MPI_Comm
220 
221    Input Parameter:
222 .  comm - MPI communicator
223 
224    Output Parameter:
225 .  ksp - location to put the KSP context
226 
227    Notes:
228    The default KSP type is GMRES with a restart of 30, using modified Gram-Schmidt
229    orthogonalization.
230 
231    Level: beginner
232 
233 .keywords: KSP, create, context
234 
235 .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP
236 @*/
237 PetscErrorCode KSPCreate(MPI_Comm comm,KSP *inksp)
238 {
239   KSP ksp;
240   PetscErrorCode ierr;
241 
242   PetscFunctionBegin;
243   PetscValidPointer(inksp,2);
244   *inksp = 0;
245 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
246   ierr = KSPInitializePackage(PETSC_NULL);CHKERRQ(ierr);
247 #endif
248 
249   PetscHeaderCreate(ksp,_p_KSP,struct _KSPOps,KSP_COOKIE,-1,"KSP",comm,KSPDestroy,KSPView);
250   PetscLogObjectCreate(ksp);
251   *inksp             = ksp;
252   ksp->bops->publish = KSPPublish_Petsc;
253 
254   ksp->type          = -1;
255   ksp->max_it        = 10000;
256   ksp->pc_side       = PC_LEFT;
257   ksp->rtol          = 1.e-5;
258   ksp->atol          = 1.e-50;
259   ksp->divtol        = 1.e4;
260 
261   ksp->normtype            = KSP_PRECONDITIONED_NORM;
262   ksp->rnorm               = 0.0;
263   ksp->its                 = 0;
264   ksp->guess_zero          = PETSC_TRUE;
265   ksp->calc_sings          = PETSC_FALSE;
266   ksp->res_hist            = PETSC_NULL;
267   ksp->res_hist_len        = 0;
268   ksp->res_hist_max        = 0;
269   ksp->res_hist_reset      = PETSC_TRUE;
270   ksp->numbermonitors      = 0;
271   ksp->converged           = KSPDefaultConverged;
272   ksp->ops->buildsolution  = KSPDefaultBuildSolution;
273   ksp->ops->buildresidual  = KSPDefaultBuildResidual;
274 
275   ksp->ops->setfromoptions = 0;
276 
277   ksp->vec_sol         = 0;
278   ksp->vec_rhs         = 0;
279   ksp->pc              = 0;
280 
281   ksp->ops->solve      = 0;
282   ksp->ops->setup      = 0;
283   ksp->ops->destroy    = 0;
284 
285   ksp->data            = 0;
286   ksp->nwork           = 0;
287   ksp->work            = 0;
288 
289   ksp->cnvP            = 0;
290 
291   ksp->reason          = KSP_CONVERGED_ITERATING;
292 
293   ksp->setupcalled     = 0;
294   ierr = PetscPublishAll(ksp);CHKERRQ(ierr);
295   ierr = PCCreate(comm,&ksp->pc);CHKERRQ(ierr);
296   PetscFunctionReturn(0);
297 }
298 
299 #undef __FUNCT__
300 #define __FUNCT__ "KSPSetType"
301 /*@C
302    KSPSetType - Builds KSP for a particular solver.
303 
304    Collective on KSP
305 
306    Input Parameters:
307 +  ksp      - the Krylov space context
308 -  type - a known method
309 
310    Options Database Key:
311 .  -ksp_type  <method> - Sets the method; use -help for a list
312     of available methods (for instance, cg or gmres)
313 
314    Notes:
315    See "petsc/include/petscksp.h" for available methods (for instance,
316    KSPCG or KSPGMRES).
317 
318   Normally, it is best to use the KSPSetFromOptions() command and
319   then set the KSP type from the options database rather than by using
320   this routine.  Using the options database provides the user with
321   maximum flexibility in evaluating the many different Krylov methods.
322   The KSPSetType() routine is provided for those situations where it
323   is necessary to set the iterative solver independently of the command
324   line or options database.  This might be the case, for example, when
325   the choice of iterative solver changes during the execution of the
326   program, and the user's application is taking responsibility for
327   choosing the appropriate method.  In other words, this routine is
328   not for beginners.
329 
330   Level: intermediate
331 
332 .keywords: KSP, set, method
333 
334 .seealso: PCSetType(), KSPType
335 
336 @*/
337 PetscErrorCode KSPSetType(KSP ksp,const KSPType type)
338 {
339   PetscErrorCode ierr,(*r)(KSP);
340   PetscTruth match;
341 
342   PetscFunctionBegin;
343   PetscValidHeaderSpecific(ksp,KSP_COOKIE,1);
344   PetscValidCharPointer(type,2);
345 
346   ierr = PetscTypeCompare((PetscObject)ksp,type,&match);CHKERRQ(ierr);
347   if (match) PetscFunctionReturn(0);
348 
349   if (ksp->data) {
350     /* destroy the old private KSP context */
351     ierr = (*ksp->ops->destroy)(ksp);CHKERRQ(ierr);
352     ksp->data = 0;
353   }
354   /* Get the function pointers for the iterative method requested */
355   if (!KSPRegisterAllCalled) {ierr = KSPRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
356   ierr =  PetscFListFind(ksp->comm,KSPList,type,(void (**)(void)) &r);CHKERRQ(ierr);
357   if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown KSP type given: %s",type);
358   ksp->setupcalled = 0;
359   ierr = (*r)(ksp);CHKERRQ(ierr);
360   ierr = PetscObjectChangeTypeName((PetscObject)ksp,type);CHKERRQ(ierr);
361   PetscFunctionReturn(0);
362 }
363 
364 #undef __FUNCT__
365 #define __FUNCT__ "KSPRegisterDestroy"
366 /*@C
367    KSPRegisterDestroy - Frees the list of KSP methods that were
368    registered by KSPRegisterDynamic().
369 
370    Not Collective
371 
372    Level: advanced
373 
374 .keywords: KSP, register, destroy
375 
376 .seealso: KSPRegisterDynamic(), KSPRegisterAll()
377 @*/
378 PetscErrorCode KSPRegisterDestroy(void)
379 {
380   PetscErrorCode ierr;
381 
382   PetscFunctionBegin;
383   if (KSPList) {
384     ierr = PetscFListDestroy(&KSPList);CHKERRQ(ierr);
385     KSPList = 0;
386   }
387   KSPRegisterAllCalled = PETSC_FALSE;
388   PetscFunctionReturn(0);
389 }
390 
391 #undef __FUNCT__
392 #define __FUNCT__ "KSPGetType"
393 /*@C
394    KSPGetType - Gets the KSP type as a string from the KSP object.
395 
396    Not Collective
397 
398    Input Parameter:
399 .  ksp - Krylov context
400 
401    Output Parameter:
402 .  name - name of KSP method
403 
404    Level: intermediate
405 
406 .keywords: KSP, get, method, name
407 
408 .seealso: KSPSetType()
409 @*/
410 PetscErrorCode KSPGetType(KSP ksp,KSPType *type)
411 {
412   PetscFunctionBegin;
413   PetscValidHeaderSpecific(ksp,KSP_COOKIE,1);
414   PetscValidPointer(type,2);
415   *type = ksp->type_name;
416   PetscFunctionReturn(0);
417 }
418 
419 #undef __FUNCT__
420 #define __FUNCT__ "KSPSetFromOptions"
421 /*@
422    KSPSetFromOptions - Sets KSP options from the options database.
423    This routine must be called before KSPSetUp() if the user is to be
424    allowed to set the Krylov type.
425 
426    Collective on KSP
427 
428    Input Parameters:
429 .  ksp - the Krylov space context
430 
431    Options Database Keys:
432 +   -ksp_max_it - maximum number of linear iterations
433 .   -ksp_rtol rtol - relative tolerance used in default determination of convergence, i.e.
434                 if residual norm decreases by this factor than convergence is declared
435 .   -ksp_atol atol - absolute tolerance used in default convergence test, i.e. if residual
436                 norm is less than this then convergence is declared
437 .   -ksp_divtol tol - if residual norm increases by this factor than divergence is declared
438 .   -ksp_norm_type - none - skip norms used in convergence tests (useful only when not using
439 $                       convergence test (say you always want to run with 5 iterations) to
440 $                       save on communication overhead
441 $                    preconditioned - default for left preconditioning
442 $                    unpreconditioned - see KSPSetNormType()
443 $                    natural - see KSPSetNormType()
444 .    -ksp_constant_null_space - assume the operator (matrix) has the constant vector in its null space
445 .    -ksp_test_null_space - tests the null space set with KSPSetNullSpace() to see if it truly is a null space
446 .   -ksp_knoll - compute initial guess by applying the preconditioner to the right hand side
447 .   -ksp_cancelmonitors - cancel all previous convergene monitor routines set
448 .   -ksp_monitor - print residual norm at each iteration
449 .   -ksp_xmonitor - plot residual norm at each iteration
450 .   -ksp_vecmonitor - plot solution at each iteration
451 -   -ksp_singmonitor - monitor extremem singular values at each iteration
452 
453    Notes:
454    To see all options, run your program with the -help option
455    or consult the users manual.
456 
457    Level: beginner
458 
459 .keywords: KSP, set, from, options, database
460 
461 .seealso:
462 @*/
463 PetscErrorCode KSPSetFromOptions(KSP ksp)
464 {
465   PetscErrorCode ierr;
466   int indx;
467   char       type[256];
468   const char *stype[] = {"none","preconditioned","unpreconditioned","natural"};
469   PetscTruth flg;
470 
471   PetscFunctionBegin;
472   PetscValidHeaderSpecific(ksp,KSP_COOKIE,1);
473   ierr = PCSetFromOptions(ksp->pc);CHKERRQ(ierr);
474 
475   if (!KSPRegisterAllCalled) {ierr = KSPRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
476   ierr = PetscOptionsBegin(ksp->comm,ksp->prefix,"Krylov Method (KSP) Options","KSP");CHKERRQ(ierr);
477     ierr = PetscOptionsList("-ksp_type","Krylov method","KSPSetType",KSPList,(char*)(ksp->type_name?ksp->type_name:KSPGMRES),type,256,&flg);CHKERRQ(ierr);
478     if (flg) {
479       ierr = KSPSetType(ksp,type);CHKERRQ(ierr);
480     }
481     /*
482       Set the type if it was never set.
483     */
484     if (!ksp->type_name) {
485       ierr = KSPSetType(ksp,KSPGMRES);CHKERRQ(ierr);
486     }
487 
488     ierr = PetscOptionsInt("-ksp_max_it","Maximum number of iterations","KSPSetTolerances",ksp->max_it,&ksp->max_it,PETSC_NULL);CHKERRQ(ierr);
489     ierr = PetscOptionsReal("-ksp_rtol","Relative decrease in residual norm","KSPSetTolerances",ksp->rtol,&ksp->rtol,PETSC_NULL);CHKERRQ(ierr);
490     ierr = PetscOptionsReal("-ksp_atol","Absolute value of residual norm","KSPSetTolerances",ksp->atol,&ksp->atol,PETSC_NULL);CHKERRQ(ierr);
491     ierr = PetscOptionsReal("-ksp_divtol","Residual norm increase cause divergence","KSPSetTolerances",ksp->divtol,&ksp->divtol,PETSC_NULL);CHKERRQ(ierr);
492     ierr = PetscOptionsLogical("-ksp_knoll","Use preconditioner applied to b for initial guess","KSPSetInitialGuessKnoll",ksp->guess_knoll,
493                                   &ksp->guess_knoll,PETSC_NULL);CHKERRQ(ierr);
494 
495     ierr = PetscOptionsEList("-ksp_norm_type","KSP Norm type","KSPSetNormType",stype,4,"preconditioned",&indx,&flg);CHKERRQ(ierr);
496     if (flg) {
497       switch (indx) {
498       case 0:
499         ierr = KSPSetNormType(ksp,KSP_NO_NORM);CHKERRQ(ierr);
500         ierr = KSPSetConvergenceTest(ksp,KSPSkipConverged,0);CHKERRQ(ierr);
501         break;
502       case 1:
503         ierr = KSPSetNormType(ksp,KSP_PRECONDITIONED_NORM);CHKERRQ(ierr);
504         break;
505       case 2:
506         ierr = KSPSetNormType(ksp,KSP_UNPRECONDITIONED_NORM);CHKERRQ(ierr);
507         break;
508       case 3:
509         ierr = KSPSetNormType(ksp,KSP_NATURAL_NORM);CHKERRQ(ierr);
510         break;
511       }
512     }
513 
514     ierr = PetscOptionsName("-ksp_diagonal_scale","Diagonal scale matrix before building preconditioner","KSPSetDiagonalScale",&flg);CHKERRQ(ierr);
515     if (flg) {
516       ierr = KSPSetDiagonalScale(ksp,PETSC_TRUE);CHKERRQ(ierr);
517     }
518     ierr = PetscOptionsName("-ksp_diagonal_scale_fix","Fix diagonaled scaled matrix after solve","KSPSetDiagonalScaleFix",&flg);CHKERRQ(ierr);
519     if (flg) {
520       ierr = KSPSetDiagonalScaleFix(ksp,PETSC_TRUE);CHKERRQ(ierr);
521     }
522 
523 
524     ierr = PetscOptionsName("-ksp_constant_null_space","Add constant null space to Krylov solver","KSPSetNullSpace",&flg);CHKERRQ(ierr);
525     if (flg) {
526       MatNullSpace nsp;
527 
528       ierr = MatNullSpaceCreate(ksp->comm,1,0,0,&nsp);CHKERRQ(ierr);
529       ierr = KSPSetNullSpace(ksp,nsp);CHKERRQ(ierr);
530       ierr = MatNullSpaceDestroy(nsp);CHKERRQ(ierr);
531     }
532 
533     /* option is actually checked in KSPSetUp() */
534     if (ksp->nullsp) {
535       ierr = PetscOptionsName("-ksp_test_null_space","Is provided null space correct","None",&flg);CHKERRQ(ierr);
536     }
537 
538     /*
539       Prints reason for convergence or divergence of each linear solve
540     */
541     ierr = PetscOptionsName("-ksp_converged_reason","Print reason for converged or diverged","KSPSolve",&flg);CHKERRQ(ierr);
542     if (flg) {
543       ksp->printreason = PETSC_TRUE;
544     }
545 
546     ierr = PetscOptionsName("-ksp_cancelmonitors","Remove any hardwired monitor routines","KSPClearMonitor",&flg);CHKERRQ(ierr);
547     /* -----------------------------------------------------------------------*/
548     /*
549       Cancels all monitors hardwired into code before call to KSPSetFromOptions()
550     */
551     if (flg) {
552       ierr = KSPClearMonitor(ksp);CHKERRQ(ierr);
553     }
554     /*
555       Prints preconditioned residual norm at each iteration
556     */
557     ierr = PetscOptionsName("-ksp_monitor","Monitor preconditioned residual norm","KSPSetMonitor",&flg);CHKERRQ(ierr);
558     if (flg) {
559       ierr = KSPSetMonitor(ksp,KSPDefaultMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
560     }
561     /*
562       Plots the vector solution
563     */
564     ierr = PetscOptionsName("-ksp_vecmonitor","Monitor solution graphically","KSPSetMonitor",&flg);CHKERRQ(ierr);
565     if (flg) {
566       ierr = KSPSetMonitor(ksp,KSPVecViewMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
567     }
568     /*
569       Prints preconditioned and true residual norm at each iteration
570     */
571     ierr = PetscOptionsName("-ksp_truemonitor","Monitor true (unpreconditioned) residual norm","KSPSetMonitor",&flg);CHKERRQ(ierr);
572     if (flg) {
573       ierr = KSPSetMonitor(ksp,KSPTrueMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
574     }
575     /*
576       Prints extreme eigenvalue estimates at each iteration
577     */
578     ierr = PetscOptionsName("-ksp_singmonitor","Monitor singular values","KSPSetMonitor",&flg);CHKERRQ(ierr);
579     if (flg) {
580       ierr = KSPSetComputeSingularValues(ksp,PETSC_TRUE);CHKERRQ(ierr);
581       ierr = KSPSetMonitor(ksp,KSPSingularValueMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
582     }
583     /*
584       Prints preconditioned residual norm with fewer digits
585     */
586     ierr = PetscOptionsName("-ksp_smonitor","Monitor preconditioned residual norm with fewer digitis","KSPSetMonitor",&flg);CHKERRQ(ierr);
587     if (flg) {
588       ierr = KSPSetMonitor(ksp,KSPDefaultSMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
589     }
590     /*
591       Graphically plots preconditioned residual norm
592     */
593     ierr = PetscOptionsName("-ksp_xmonitor","Monitor graphically preconditioned residual norm","KSPSetMonitor",&flg);CHKERRQ(ierr);
594     if (flg) {
595       ierr = KSPSetMonitor(ksp,KSPLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
596     }
597     /*
598       Graphically plots preconditioned and true residual norm
599     */
600     ierr = PetscOptionsName("-ksp_xtruemonitor","Monitor graphically true residual norm","KSPSetMonitor",&flg);CHKERRQ(ierr);
601     if (flg){
602       ierr = KSPSetMonitor(ksp,KSPLGTrueMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
603     }
604 
605     /* -----------------------------------------------------------------------*/
606 
607     ierr = PetscOptionsLogicalGroupBegin("-ksp_left_pc","Use left preconditioning","KSPSetPreconditionerSide",&flg);CHKERRQ(ierr);
608     if (flg) { ierr = KSPSetPreconditionerSide(ksp,PC_LEFT);CHKERRQ(ierr); }
609     ierr = PetscOptionsLogicalGroup("-ksp_right_pc","Use right preconditioning","KSPSetPreconditionerSide",&flg);CHKERRQ(ierr);
610     if (flg) { ierr = KSPSetPreconditionerSide(ksp,PC_RIGHT);CHKERRQ(ierr);}
611     ierr = PetscOptionsLogicalGroupEnd("-ksp_symmetric_pc","Use symmetric (factorized) preconditioning","KSPSetPreconditionerSide",&flg);CHKERRQ(ierr);
612     if (flg) { ierr = KSPSetPreconditionerSide(ksp,PC_SYMMETRIC);CHKERRQ(ierr);}
613 
614     ierr = PetscOptionsName("-ksp_compute_singularvalues","Compute singular values of preconditioned operator","KSPSetComputeSingularValues",&flg);CHKERRQ(ierr);
615     if (flg) { ierr = KSPSetComputeSingularValues(ksp,PETSC_TRUE);CHKERRQ(ierr); }
616     ierr = PetscOptionsName("-ksp_compute_eigenvalues","Compute eigenvalues of preconditioned operator","KSPSetComputeSingularValues",&flg);CHKERRQ(ierr);
617     if (flg) { ierr = KSPSetComputeSingularValues(ksp,PETSC_TRUE);CHKERRQ(ierr); }
618     ierr = PetscOptionsName("-ksp_plot_eigenvalues","Scatter plot extreme eigenvalues","KSPSetComputeSingularValues",&flg);CHKERRQ(ierr);
619     if (flg) { ierr = KSPSetComputeSingularValues(ksp,PETSC_TRUE);CHKERRQ(ierr); }
620 
621     if (ksp->ops->setfromoptions) {
622       ierr = (*ksp->ops->setfromoptions)(ksp);CHKERRQ(ierr);
623     }
624     ierr = PetscOptionsName("-ksp_view","View linear solver parameters","KSPView",&flg);CHKERRQ(ierr);
625   ierr = PetscOptionsEnd();CHKERRQ(ierr);
626   PetscFunctionReturn(0);
627 }
628 
629 #undef __FUNCT__
630 #define __FUNCT__ "KSPRegister"
631 /*@C
632   KSPRegister - See KSPRegisterDynamic()
633 
634   Level: advanced
635 @*/
636 PetscErrorCode KSPRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(KSP))
637 {
638   PetscErrorCode ierr;
639   char fullname[PETSC_MAX_PATH_LEN];
640 
641   PetscFunctionBegin;
642   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
643   ierr = PetscFListAdd(&KSPList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
644   PetscFunctionReturn(0);
645 }
646 
647 #undef __FUNCT__
648 #define __FUNCT__ "KSPSetNullSpace"
649 /*@C
650   KSPSetNullSpace - Sets the null space of the operator
651 
652   Collective on KSP
653 
654   Input Parameters:
655 +  ksp - the Krylov space object
656 -  nullsp - the null space of the operator
657 
658   Level: advanced
659 
660 .seealso: KSPSetOperators(), MatNullSpaceCreate(), KSPGetNullSpace()
661 @*/
662 PetscErrorCode KSPSetNullSpace(KSP ksp,MatNullSpace nullsp)
663 {
664   PetscErrorCode ierr;
665 
666   PetscFunctionBegin;
667   if (ksp->nullsp) {
668     ierr = MatNullSpaceDestroy(ksp->nullsp);CHKERRQ(ierr);
669   }
670   ksp->nullsp = nullsp;
671   ierr = PetscObjectReference((PetscObject)ksp->nullsp);CHKERRQ(ierr);
672   PetscFunctionReturn(0);
673 }
674 
675 #undef __FUNCT__
676 #define __FUNCT__ "KSPGetNullSpace"
677 /*@C
678   KSPGetNullSpace - Gets the null space of the operator
679 
680   Collective on KSP
681 
682   Input Parameters:
683 +  ksp - the Krylov space object
684 -  nullsp - the null space of the operator
685 
686   Level: advanced
687 
688 .seealso: KSPSetOperators(), MatNullSpaceCreate(), KSPSetNullSpace()
689 @*/
690 PetscErrorCode KSPGetNullSpace(KSP ksp,MatNullSpace *nullsp)
691 {
692   PetscFunctionBegin;
693   *nullsp = ksp->nullsp;
694   PetscFunctionReturn(0);
695 }
696 
697