Actual source code: itcl.c
petsc-3.4.5 2014-06-29
2: /*
3: Code for setting KSP options from the options database.
4: */
6: #include <petsc-private/kspimpl.h> /*I "petscksp.h" I*/
8: extern PetscBool KSPRegisterAllCalled;
12: /*@C
13: KSPSetOptionsPrefix - Sets the prefix used for searching for all
14: KSP options in the database.
16: Logically Collective on KSP
18: Input Parameters:
19: + ksp - the Krylov context
20: - prefix - the prefix string to prepend to all KSP option requests
22: Notes:
23: A hyphen (-) must NOT be given at the beginning of the prefix name.
24: The first character of all runtime options is AUTOMATICALLY the
25: hyphen.
27: For example, to distinguish between the runtime options for two
28: different KSP contexts, one could call
29: .vb
30: KSPSetOptionsPrefix(ksp1,"sys1_")
31: KSPSetOptionsPrefix(ksp2,"sys2_")
32: .ve
34: This would enable use of different options for each system, such as
35: .vb
36: -sys1_ksp_type gmres -sys1_ksp_rtol 1.e-3
37: -sys2_ksp_type bcgs -sys2_ksp_rtol 1.e-4
38: .ve
40: Level: advanced
42: .keywords: KSP, set, options, prefix, database
44: .seealso: KSPAppendOptionsPrefix(), KSPGetOptionsPrefix()
45: @*/
46: PetscErrorCode KSPSetOptionsPrefix(KSP ksp,const char prefix[])
47: {
52: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
53: PCSetOptionsPrefix(ksp->pc,prefix);
54: PetscObjectSetOptionsPrefix((PetscObject)ksp,prefix);
55: return(0);
56: }
60: /*@C
61: KSPAppendOptionsPrefix - Appends to the prefix used for searching for all
62: KSP options in the database.
64: Logically Collective on KSP
66: Input Parameters:
67: + ksp - the Krylov context
68: - prefix - the prefix string to prepend to all KSP option requests
70: Notes:
71: A hyphen (-) must NOT be given at the beginning of the prefix name.
72: The first character of all runtime options is AUTOMATICALLY the hyphen.
74: Level: advanced
76: .keywords: KSP, append, options, prefix, database
78: .seealso: KSPSetOptionsPrefix(), KSPGetOptionsPrefix()
79: @*/
80: PetscErrorCode KSPAppendOptionsPrefix(KSP ksp,const char prefix[])
81: {
86: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
87: PCAppendOptionsPrefix(ksp->pc,prefix);
88: PetscObjectAppendOptionsPrefix((PetscObject)ksp,prefix);
89: return(0);
90: }
94: /*@
95: KSPGetTabLevel - Gets the number of tabs that ASCII output used by ksp.
97: Not Collective
99: Input Parameter:
100: . ksp - a KSP object.
102: Output Parameter:
103: . tab - the number of tabs
105: Level: developer
107: Notes: this is used in conjunction with KSPSetTabLevel() to manage the output from the KSP and its PC coherently.
110: .seealso: KSPSetTabLevel()
112: @*/
113: PetscErrorCode KSPGetTabLevel(KSP ksp,PetscInt *tab)
114: {
119: PetscObjectGetTabLevel((PetscObject)ksp, tab);
120: return(0);
121: }
125: /*@
126: KSPSetTabLevel - Sets the number of tabs that ASCII output for the ksp andn its pc will use.
128: Not Collective
130: Input Parameters:
131: + ksp - a KSP object
132: - tab - the number of tabs
134: Level: developer
136: Notes: this is used to manage the output from KSP and PC objects that are imbedded in other objects,
137: for example, the KSP object inside a SNES object. By indenting each lower level further the heirarchy
138: of objects is very clear. By setting the KSP object's tab level with KSPSetTabLevel() its PC object
139: automatically receives the same tab level, so that whatever objects the pc might create are tabbed
140: appropriately, too.
142: .seealso: KSPGetTabLevel()
143: @*/
144: PetscErrorCode KSPSetTabLevel(KSP ksp, PetscInt tab)
145: {
150: PetscObjectSetTabLevel((PetscObject)ksp, tab);
151: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
152: /* Do we need a PCSetTabLevel()? */
153: PetscObjectSetTabLevel((PetscObject)ksp->pc, tab);
154: return(0);
155: }
159: /*@C
160: KSPSetUseFischerGuess - Use the Paul Fischer algorithm, see KSPFischerGuessCreate()
162: Logically Collective on KSP
164: Input Parameters:
165: + ksp - the Krylov context
166: . model - use model 1, model 2 or 0 to turn it off
167: - size - size of subspace used to generate initial guess
169: Options Database:
170: . -ksp_fischer_guess <model,size> - uses the Fischer initial guess generator for repeated linear solves
172: Level: advanced
174: .keywords: KSP, set, options, prefix, database
176: .seealso: KSPSetOptionsPrefix(), KSPAppendOptionsPrefix(), KSPSetUseFischerGuess(), KSPSetFischerGuess(), KSPGetFischerInitialGuess()
177: @*/
178: PetscErrorCode KSPSetUseFischerGuess(KSP ksp,PetscInt model,PetscInt size)
179: {
186: KSPFischerGuessDestroy(&ksp->guess);
187: if (model == 1 || model == 2) {
188: KSPFischerGuessCreate(ksp,model,size,&ksp->guess);
189: KSPFischerGuessSetFromOptions(ksp->guess);
190: } else if (model != 0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Model must be 1 or 2 (or 0 to turn off guess generation)");
191: return(0);
192: }
196: /*@C
197: KSPSetFischerGuess - Use the Paul Fischer algorithm created by KSPFischerGuessCreate()
199: Logically Collective on KSP
201: Input Parameters:
202: + ksp - the Krylov context
203: - guess - the object created with KSPFischerGuessCreate()
205: Level: advanced
207: Notes: this allows a single KSP to be used with several different initial guess generators (likely for different linear
208: solvers, see KSPSetPC()).
210: This increases the reference count of the guess object, you must destroy the object with KSPFischerGuessDestroy()
211: before the end of the program.
213: .keywords: KSP, set, options, prefix, database
215: .seealso: KSPSetOptionsPrefix(), KSPAppendOptionsPrefix(), KSPSetUseFischerGuess(), KSPSetFischerGuess(), KSPGetFischerGuess()
216: @*/
217: PetscErrorCode KSPSetFischerGuess(KSP ksp,KSPFischerGuess guess)
218: {
223: KSPFischerGuessDestroy(&ksp->guess);
224: ksp->guess = guess;
225: if (guess) guess->refcnt++;
226: return(0);
227: }
231: /*@C
232: KSPGetFischerGuess - Gets the initial guess generator set with either KSPSetFischerGuess() or KSPCreateFischerGuess()/KSPSetFischerGuess()
234: Not Collective
236: Input Parameters:
237: . ksp - the Krylov context
239: Output Parameters:
240: . guess - the object
242: Level: developer
244: .keywords: KSP, set, options, prefix, database
246: .seealso: KSPSetOptionsPrefix(), KSPAppendOptionsPrefix(), KSPSetUseFischerGuess(), KSPSetFischerGuess()
247: @*/
248: PetscErrorCode KSPGetFischerGuess(KSP ksp,KSPFischerGuess *guess)
249: {
251: *guess = ksp->guess;
252: return(0);
253: }
257: /*@C
258: KSPGetOptionsPrefix - Gets the prefix used for searching for all
259: KSP options in the database.
261: Not Collective
263: Input Parameters:
264: . ksp - the Krylov context
266: Output Parameters:
267: . prefix - pointer to the prefix string used is returned
269: Notes: On the fortran side, the user should pass in a string 'prifix' of
270: sufficient length to hold the prefix.
272: Level: advanced
274: .keywords: KSP, set, options, prefix, database
276: .seealso: KSPSetOptionsPrefix(), KSPAppendOptionsPrefix()
277: @*/
278: PetscErrorCode KSPGetOptionsPrefix(KSP ksp,const char *prefix[])
279: {
284: PetscObjectGetOptionsPrefix((PetscObject)ksp,prefix);
285: return(0);
286: }
290: /*@
291: KSPSetFromOptions - Sets KSP options from the options database.
292: This routine must be called before KSPSetUp() if the user is to be
293: allowed to set the Krylov type.
295: Collective on KSP
297: Input Parameters:
298: . ksp - the Krylov space context
300: Options Database Keys:
301: + -ksp_max_it - maximum number of linear iterations
302: . -ksp_rtol rtol - relative tolerance used in default determination of convergence, i.e.
303: if residual norm decreases by this factor than convergence is declared
304: . -ksp_atol abstol - absolute tolerance used in default convergence test, i.e. if residual
305: norm is less than this then convergence is declared
306: . -ksp_divtol tol - if residual norm increases by this factor than divergence is declared
307: . -ksp_converged_use_initial_residual_norm - see KSPDefaultConvergedSetUIRNorm()
308: . -ksp_converged_use_min_initial_residual_norm - see KSPDefaultConvergedSetUMIRNorm()
309: . -ksp_norm_type - none - skip norms used in convergence tests (useful only when not using
310: $ convergence test (say you always want to run with 5 iterations) to
311: $ save on communication overhead
312: $ preconditioned - default for left preconditioning
313: $ unpreconditioned - see KSPSetNormType()
314: $ natural - see KSPSetNormType()
315: . -ksp_check_norm_iteration it - do not compute residual norm until iteration number it (does compute at 0th iteration)
316: $ works only for PCBCGS, PCIBCGS and and PCCG
317: -ksp_lag_norm - compute the norm of the residual for the ith iteration on the i+1 iteration; this means that one can use
318: $ the norm of the residual for convergence test WITHOUT an extra MPI_Allreduce() limiting global synchronizations.
319: $ This will require 1 more iteration of the solver than usual.
320: . -ksp_fischer_guess <model,size> - uses the Fischer initial guess generator for repeated linear solves
321: . -ksp_constant_null_space - assume the operator (matrix) has the constant vector in its null space
322: . -ksp_test_null_space - tests the null space set with KSPSetNullSpace() to see if it truly is a null space
323: . -ksp_knoll - compute initial guess by applying the preconditioner to the right hand side
324: . -ksp_monitor_cancel - cancel all previous convergene monitor routines set
325: . -ksp_monitor <optional filename> - print residual norm at each iteration
326: . -ksp_monitor_lg_residualnorm - plot residual norm at each iteration
327: . -ksp_monitor_solution - plot solution at each iteration
328: - -ksp_monitor_singular_value - monitor extremem singular values at each iteration
330: Notes:
331: To see all options, run your program with the -help option
332: or consult <A href="../../docs/manual.pdf#nameddest=Chapter 4 KSP: Linear Equations Solvers">KSP chapter of the users manual</A>.
334: Level: beginner
336: .keywords: KSP, set, from, options, database
338: .seealso: KSPSetUseFischerInitialGuess()
340: @*/
341: PetscErrorCode KSPSetFromOptions(KSP ksp)
342: {
344: PetscInt indx;
345: const char *convtests[] = {"default","skip"};
346: char type[256], monfilename[PETSC_MAX_PATH_LEN];
347: PetscViewer monviewer;
348: PetscBool flg,flag;
349: PetscInt model[2]={0,0},nmax;
350: KSPNormType normtype;
351: PCSide pcside;
352: void *ctx;
356: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
357: PCSetFromOptions(ksp->pc);
359: if (!KSPRegisterAllCalled) {KSPRegisterAll();}
360: PetscObjectOptionsBegin((PetscObject)ksp);
361: PetscOptionsList("-ksp_type","Krylov method","KSPSetType",KSPList,(char*)(((PetscObject)ksp)->type_name ? ((PetscObject)ksp)->type_name : KSPGMRES),type,256,&flg);
362: if (flg) {
363: KSPSetType(ksp,type);
364: }
365: /*
366: Set the type if it was never set.
367: */
368: if (!((PetscObject)ksp)->type_name) {
369: KSPSetType(ksp,KSPGMRES);
370: }
372: PetscOptionsInt("-ksp_max_it","Maximum number of iterations","KSPSetTolerances",ksp->max_it,&ksp->max_it,NULL);
373: PetscOptionsReal("-ksp_rtol","Relative decrease in residual norm","KSPSetTolerances",ksp->rtol,&ksp->rtol,NULL);
374: PetscOptionsReal("-ksp_atol","Absolute value of residual norm","KSPSetTolerances",ksp->abstol,&ksp->abstol,NULL);
375: PetscOptionsReal("-ksp_divtol","Residual norm increase cause divergence","KSPSetTolerances",ksp->divtol,&ksp->divtol,NULL);
377: flag = PETSC_FALSE;
378: PetscOptionsBool("-ksp_converged_use_initial_residual_norm","Use initial residual residual norm for computing relative convergence","KSPDefaultConvergedSetUIRNorm",flag,&flag,NULL);
379: if (flag) {KSPDefaultConvergedSetUIRNorm(ksp);}
380: flag = PETSC_FALSE;
381: PetscOptionsBool("-ksp_converged_use_min_initial_residual_norm","Use minimum of initial residual norm and b for computing relative convergence","KSPDefaultConvergedSetUMIRNorm",flag,&flag,NULL);
382: if (flag) {KSPDefaultConvergedSetUMIRNorm(ksp);}
383: KSPGetInitialGuessNonzero(ksp,&flag);
384: PetscOptionsBool("-ksp_initial_guess_nonzero","Use the contents of the solution vector for initial guess","KSPSetInitialNonzero",flag,&flag,&flg);
385: if (flg) {
386: KSPSetInitialGuessNonzero(ksp,flag);
387: }
389: PetscOptionsBool("-ksp_knoll","Use preconditioner applied to b for initial guess","KSPSetInitialGuessKnoll",ksp->guess_knoll,&ksp->guess_knoll,NULL);
390: PetscOptionsBool("-ksp_error_if_not_converged","Generate error if solver does not converge","KSPSetErrorIfNotConverged",ksp->errorifnotconverged,&ksp->errorifnotconverged,NULL);
391: nmax = 2;
392: PetscOptionsIntArray("-ksp_fischer_guess","Use Paul Fischer's algorithm for initial guess","KSPSetUseFischerGuess",model,&nmax,&flag);
393: if (flag) {
394: if (nmax != 2) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Must pass in model,size as arguments");
395: KSPSetUseFischerGuess(ksp,model[0],model[1]);
396: }
398: PetscOptionsEList("-ksp_convergence_test","Convergence test","KSPSetConvergenceTest",convtests,2,"default",&indx,&flg);
399: if (flg) {
400: switch (indx) {
401: case 0:
402: KSPDefaultConvergedCreate(&ctx);
403: KSPSetConvergenceTest(ksp,KSPDefaultConverged,ctx,KSPDefaultConvergedDestroy);
404: break;
405: case 1: KSPSetConvergenceTest(ksp,KSPSkipConverged,NULL,NULL); break;
406: }
407: }
409: KSPSetUpNorms_Private(ksp,&normtype,&pcside);
410: PetscOptionsEnum("-ksp_norm_type","KSP Norm type","KSPSetNormType",KSPNormTypes,(PetscEnum)normtype,(PetscEnum*)&normtype,&flg);
411: if (flg) { KSPSetNormType(ksp,normtype); }
413: PetscOptionsInt("-ksp_check_norm_iteration","First iteration to compute residual norm","KSPSetCheckNormIteration",ksp->chknorm,&ksp->chknorm,NULL);
415: flag = ksp->lagnorm;
416: PetscOptionsBool("-ksp_lag_norm","Lag the calculation of the residual norm","KSPSetLagNorm",flag,&flag,&flg);
417: if (flg) {
418: KSPSetLagNorm(ksp,flag);
419: }
421: KSPGetDiagonalScale(ksp,&flag);
422: PetscOptionsBool("-ksp_diagonal_scale","Diagonal scale matrix before building preconditioner","KSPSetDiagonalScale",flag,&flag,&flg);
423: if (flg) {
424: KSPSetDiagonalScale(ksp,flag);
425: }
426: KSPGetDiagonalScaleFix(ksp,&flag);
427: PetscOptionsBool("-ksp_diagonal_scale_fix","Fix diagonally scaled matrix after solve","KSPSetDiagonalScaleFix",flag,&flag,&flg);
428: if (flg) {
429: KSPSetDiagonalScaleFix(ksp,flag);
430: }
432: flg = PETSC_FALSE;
433: PetscOptionsBool("-ksp_constant_null_space","Add constant null space to Krylov solver","KSPSetNullSpace",flg,&flg,NULL);
434: if (flg) {
435: MatNullSpace nsp;
437: MatNullSpaceCreate(PetscObjectComm((PetscObject)ksp),PETSC_TRUE,0,0,&nsp);
438: KSPSetNullSpace(ksp,nsp);
439: MatNullSpaceDestroy(&nsp);
440: }
442: /* option is actually checked in KSPSetUp(), just here so goes into help message */
443: if (ksp->nullsp) {
444: PetscOptionsName("-ksp_test_null_space","Is provided null space correct","None",&flg);
445: }
447: /*
448: Prints reason for convergence or divergence of each linear solve
449: */
450: flg = PETSC_FALSE;
451: PetscOptionsBool("-ksp_converged_reason","Print reason for converged or diverged","KSPSolve",flg,&flg,NULL);
452: if (flg) ksp->printreason = PETSC_TRUE;
454: flg = PETSC_FALSE;
455: PetscOptionsBool("-ksp_monitor_cancel","Remove any hardwired monitor routines","KSPMonitorCancel",flg,&flg,NULL);
456: /* -----------------------------------------------------------------------*/
457: /*
458: Cancels all monitors hardwired into code before call to KSPSetFromOptions()
459: */
460: if (flg) {
461: KSPMonitorCancel(ksp);
462: }
463: /*
464: Prints preconditioned residual norm at each iteration
465: */
466: PetscOptionsString("-ksp_monitor","Monitor preconditioned residual norm","KSPMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
467: if (flg) {
468: PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ksp),monfilename,&monviewer);
469: KSPMonitorSet(ksp,KSPMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
470: }
471: /*
472: Prints preconditioned residual norm at each iteration
473: */
474: PetscOptionsString("-ksp_monitor_range","Monitor percent of residual entries more than 10 percent of max","KSPMonitorRange","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
475: if (flg) {
476: PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ksp),monfilename,&monviewer);
477: KSPMonitorSet(ksp,KSPMonitorRange,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
478: }
479: PetscObjectTypeCompare((PetscObject)ksp->pc,PCKSP,&flg);
480: PetscObjectTypeCompare((PetscObject)ksp->pc,PCBJACOBI,&flag);
481: if (flg || flag) {
482: /* A hack for using dynamic tolerance in preconditioner */
483: PetscOptionsString("-sub_ksp_dynamic_tolerance","Use dynamic tolerance for PC if PC is a KSP","KSPMonitorDynamicTolerance","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
484: if (flg) {
485: KSPDynTolCtx *scale = NULL;
486: PetscReal defaultv = 1.0;
487: PetscMalloc(1*sizeof(KSPDynTolCtx),&scale);
488: scale->bnrm = -1.0;
489: scale->coef = defaultv;
490: PetscOptionsReal("-sub_ksp_dynamic_tolerance_param","Parameter of dynamic tolerance for inner PCKSP","KSPMonitorDynamicToleranceParam",defaultv,&(scale->coef),&flg);
491: KSPMonitorSet(ksp,KSPMonitorDynamicTolerance,scale,KSPMonitorDynamicToleranceDestroy);
492: }
493: }
494: /*
495: Plots the vector solution
496: */
497: flg = PETSC_FALSE;
498: PetscOptionsBool("-ksp_monitor_solution","Monitor solution graphically","KSPMonitorSet",flg,&flg,NULL);
499: if (flg) {
500: KSPMonitorSet(ksp,KSPMonitorSolution,NULL,NULL);
501: }
502: /*
503: Prints preconditioned and true residual norm at each iteration
504: */
505: PetscOptionsString("-ksp_monitor_true_residual","Monitor true residual norm","KSPMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
506: if (flg) {
507: PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ksp),monfilename,&monviewer);
508: KSPMonitorSet(ksp,KSPMonitorTrueResidualNorm,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
509: }
510: /*
511: Prints with max norm at each iteration
512: */
513: PetscOptionsString("-ksp_monitor_max","Monitor true residual max norm","KSPMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
514: if (flg) {
515: PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ksp),monfilename,&monviewer);
516: KSPMonitorSet(ksp,KSPMonitorTrueResidualMaxNorm,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
517: }
518: /*
519: Prints extreme eigenvalue estimates at each iteration
520: */
521: PetscOptionsString("-ksp_monitor_singular_value","Monitor singular values","KSPMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
522: if (flg) {
523: KSPSetComputeSingularValues(ksp,PETSC_TRUE);
524: PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ksp),monfilename,&monviewer);
525: KSPMonitorSet(ksp,KSPMonitorSingularValue,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
526: }
527: /*
528: Prints preconditioned residual norm with fewer digits
529: */
530: PetscOptionsString("-ksp_monitor_short","Monitor preconditioned residual norm with fewer digits","KSPMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
531: if (flg) {
532: PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ksp),monfilename,&monviewer);
533: KSPMonitorSet(ksp,KSPMonitorDefaultShort,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
534: }
535: /*
536: Calls Python function
537: */
538: PetscOptionsString("-ksp_monitor_python","Use Python function","KSPMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);
539: if (flg) {PetscPythonMonitorSet((PetscObject)ksp,monfilename);}
540: /*
541: Graphically plots preconditioned residual norm
542: */
543: flg = PETSC_FALSE;
544: PetscOptionsBool("-ksp_monitor_lg_residualnorm","Monitor graphically preconditioned residual norm","KSPMonitorSet",flg,&flg,NULL);
545: if (flg) {
546: PetscDrawLG ctx;
548: KSPMonitorLGResidualNormCreate(0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,&ctx);
549: KSPMonitorSet(ksp,KSPMonitorLGResidualNorm,ctx,(PetscErrorCode (*)(void**))KSPMonitorLGResidualNormDestroy);
550: }
551: /*
552: Graphically plots preconditioned and true residual norm
553: */
554: flg = PETSC_FALSE;
555: PetscOptionsBool("-ksp_monitor_lg_true_residualnorm","Monitor graphically true residual norm","KSPMonitorSet",flg,&flg,NULL);
556: if (flg) {
557: PetscDrawLG ctx;
559: KSPMonitorLGTrueResidualNormCreate(PetscObjectComm((PetscObject)ksp),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,&ctx);
560: KSPMonitorSet(ksp,KSPMonitorLGTrueResidualNorm,ctx,(PetscErrorCode (*)(void**))KSPMonitorLGTrueResidualNormDestroy);
561: }
562: /*
563: Graphically plots preconditioned residual norm and range of residual element values
564: */
565: flg = PETSC_FALSE;
566: PetscOptionsBool("-ksp_monitor_lg_range","Monitor graphically range of preconditioned residual norm","KSPMonitorSet",flg,&flg,NULL);
567: if (flg) {
568: PetscViewer ctx;
570: PetscViewerDrawOpen(PetscObjectComm((PetscObject)ksp),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,&ctx);
571: KSPMonitorSet(ksp,KSPMonitorLGRange,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);
572: }
574: #if defined(PETSC_HAVE_AMS)
575: /*
576: Publish convergence information using AMS
577: */
578: flg = PETSC_FALSE;
579: PetscOptionsBool("-ksp_monitor_ams","Publish KSP progress using AMS","KSPMonitorSet",flg,&flg,NULL);
580: if (flg) {
581: void *ctx;
582: KSPMonitorAMSCreate(ksp,NULL,&ctx);
583: KSPMonitorSet(ksp,KSPMonitorAMS,ctx,KSPMonitorAMSDestroy);
584: KSPSetComputeSingularValues(ksp,PETSC_TRUE);
585: }
586: #endif
588: /* -----------------------------------------------------------------------*/
589: KSPSetUpNorms_Private(ksp,&normtype,&pcside);
590: PetscOptionsEnum("-ksp_pc_side","KSP preconditioner side","KSPSetPCSide",PCSides,(PetscEnum)pcside,(PetscEnum*)&pcside,&flg);
591: if (flg) {KSPSetPCSide(ksp,pcside);}
593: flg = PETSC_FALSE;
594: PetscOptionsBool("-ksp_compute_singularvalues","Compute singular values of preconditioned operator","KSPSetComputeSingularValues",flg,&flg,NULL);
595: if (flg) { KSPSetComputeSingularValues(ksp,PETSC_TRUE); }
596: flg = PETSC_FALSE;
597: PetscOptionsBool("-ksp_compute_eigenvalues","Compute eigenvalues of preconditioned operator","KSPSetComputeSingularValues",flg,&flg,NULL);
598: if (flg) { KSPSetComputeSingularValues(ksp,PETSC_TRUE); }
599: flg = PETSC_FALSE;
600: PetscOptionsBool("-ksp_plot_eigenvalues","Scatter plot extreme eigenvalues","KSPSetComputeSingularValues",flg,&flg,NULL);
601: if (flg) { KSPSetComputeSingularValues(ksp,PETSC_TRUE); }
603: #if defined(PETSC_HAVE_AMS)
604: {
605: PetscBool set;
606: flg = PETSC_FALSE;
607: PetscOptionsBool("-ksp_ams_block","Block for AMS memory snooper at end of KSPSolve","PetscObjectAMSBlock",((PetscObject)ksp)->amspublishblock,&flg,&set);
608: if (set) {
609: PetscObjectAMSSetBlock((PetscObject)ksp,flg);
610: }
611: }
612: #endif
614: if (ksp->ops->setfromoptions) {
615: (*ksp->ops->setfromoptions)(ksp);
616: }
617: /* actually check in setup this is just here so goes into help message */
618: PetscOptionsName("-ksp_view","View linear solver parameters","KSPView",&flg);
620: /* process any options handlers added with PetscObjectAddOptionsHandler() */
621: PetscObjectProcessOptionsHandlers((PetscObject)ksp);
622: PetscOptionsEnd();
623: return(0);
624: }