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