Actual source code: itcl.c
petsc-3.3-p1 2012-06-15
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: {
51: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
52: PCSetOptionsPrefix(ksp->pc,prefix);
53: PetscObjectSetOptionsPrefix((PetscObject)ksp,prefix);
54: return(0);
55: }
56:
59: /*@C
60: KSPAppendOptionsPrefix - Appends to the prefix used for searching for all
61: KSP options in the database.
63: Logically Collective on KSP
65: Input Parameters:
66: + ksp - the Krylov context
67: - prefix - the prefix string to prepend to all KSP option requests
69: Notes:
70: A hyphen (-) must NOT be given at the beginning of the prefix name.
71: The first character of all runtime options is AUTOMATICALLY the hyphen.
73: Level: advanced
75: .keywords: KSP, append, options, prefix, database
77: .seealso: KSPSetOptionsPrefix(), KSPGetOptionsPrefix()
78: @*/
79: PetscErrorCode KSPAppendOptionsPrefix(KSP ksp,const char prefix[])
80: {
84: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
85: PCAppendOptionsPrefix(ksp->pc,prefix);
86: PetscObjectAppendOptionsPrefix((PetscObject)ksp,prefix);
87: return(0);
88: }
92: /*@C
93: KSPSetUseFischerGuess - Use the Paul Fischer algorithm, see KSPFischerGuessCreate()
95: Logically Collective on KSP
97: Input Parameters:
98: + ksp - the Krylov context
99: . model - use model 1, model 2 or 0 to turn it off
100: - size - size of subspace used to generate initial guess
102: Options Database:
103: . -ksp_fischer_guess <model,size> - uses the Fischer initial guess generator for repeated linear solves
105: Level: advanced
107: .keywords: KSP, set, options, prefix, database
109: .seealso: KSPSetOptionsPrefix(), KSPAppendOptionsPrefix(), KSPSetUseFischerGuess(), KSPSetFischerGuess(), KSPGetFischerInitialGuess()
110: @*/
111: PetscErrorCode KSPSetUseFischerGuess(KSP ksp,PetscInt model,PetscInt size)
112: {
118: KSPFischerGuessDestroy(&ksp->guess);
119: if (model == 1 || model == 2) {
120: KSPFischerGuessCreate(ksp,model,size,&ksp->guess);
121: KSPFischerGuessSetFromOptions(ksp->guess);
122: } else if (model != 0) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Model must be 1 or 2 (or 0 to turn off guess generation)");
123: return(0);
124: }
128: /*@C
129: KSPSetFischerGuess - Use the Paul Fischer algorithm created by KSPFischerGuessCreate()
131: Logically Collective on KSP
133: Input Parameters:
134: + ksp - the Krylov context
135: - guess - the object created with KSPFischerGuessCreate()
137: Level: advanced
139: Notes: this allows a single KSP to be used with several different initial guess generators (likely for different linear
140: solvers, see KSPSetPC()).
142: This increases the reference count of the guess object, you must destroy the object with KSPFischerGuessDestroy()
143: before the end of the program.
145: .keywords: KSP, set, options, prefix, database
147: .seealso: KSPSetOptionsPrefix(), KSPAppendOptionsPrefix(), KSPSetUseFischerGuess(), KSPSetFischerGuess(), KSPGetFischerGuess()
148: @*/
149: PetscErrorCode KSPSetFischerGuess(KSP ksp,KSPFischerGuess guess)
150: {
154: KSPFischerGuessDestroy(&ksp->guess);
155: ksp->guess = guess;
156: if (guess) guess->refcnt++;
157: return(0);
158: }
162: /*@C
163: KSPGetFischerGuess - Gets the initial guess generator set with either KSPSetFischerGuess() or KSPCreateFischerGuess()/KSPSetFischerGuess()
165: Not Collective
167: Input Parameters:
168: . ksp - the Krylov context
170: Output Parameters:
171: . guess - the object
173: Level: developer
175: .keywords: KSP, set, options, prefix, database
177: .seealso: KSPSetOptionsPrefix(), KSPAppendOptionsPrefix(), KSPSetUseFischerGuess(), KSPSetFischerGuess()
178: @*/
179: PetscErrorCode KSPGetFischerGuess(KSP ksp,KSPFischerGuess *guess)
180: {
182: *guess = ksp->guess;
183: return(0);
184: }
188: /*@C
189: KSPGetOptionsPrefix - Gets the prefix used for searching for all
190: KSP options in the database.
192: Not Collective
194: Input Parameters:
195: . ksp - the Krylov context
197: Output Parameters:
198: . prefix - pointer to the prefix string used is returned
200: Notes: On the fortran side, the user should pass in a string 'prifix' of
201: sufficient length to hold the prefix.
203: Level: advanced
205: .keywords: KSP, set, options, prefix, database
207: .seealso: KSPSetOptionsPrefix(), KSPAppendOptionsPrefix()
208: @*/
209: PetscErrorCode KSPGetOptionsPrefix(KSP ksp,const char *prefix[])
210: {
214: PetscObjectGetOptionsPrefix((PetscObject)ksp,prefix);
215: return(0);
216: }
220: /*@
221: KSPSetFromOptions - Sets KSP options from the options database.
222: This routine must be called before KSPSetUp() if the user is to be
223: allowed to set the Krylov type.
225: Collective on KSP
227: Input Parameters:
228: . ksp - the Krylov space context
230: Options Database Keys:
231: + -ksp_max_it - maximum number of linear iterations
232: . -ksp_rtol rtol - relative tolerance used in default determination of convergence, i.e.
233: if residual norm decreases by this factor than convergence is declared
234: . -ksp_atol abstol - absolute tolerance used in default convergence test, i.e. if residual
235: norm is less than this then convergence is declared
236: . -ksp_divtol tol - if residual norm increases by this factor than divergence is declared
237: . -ksp_converged_use_initial_residual_norm - see KSPDefaultConvergedSetUIRNorm()
238: . -ksp_converged_use_min_initial_residual_norm - see KSPDefaultConvergedSetUMIRNorm()
239: . -ksp_norm_type - none - skip norms used in convergence tests (useful only when not using
240: $ convergence test (say you always want to run with 5 iterations) to
241: $ save on communication overhead
242: $ preconditioned - default for left preconditioning
243: $ unpreconditioned - see KSPSetNormType()
244: $ natural - see KSPSetNormType()
245: . -ksp_check_norm_iteration it - do not compute residual norm until iteration number it (does compute at 0th iteration)
246: $ works only for PCBCGS, PCIBCGS and and PCCG
247: -ksp_lag_norm - compute the norm of the residual for the ith iteration on the i+1 iteration; this means that one can use
248: $ the norm of the residual for convergence test WITHOUT an extra MPI_Allreduce() limiting global synchronizations.
249: $ This will require 1 more iteration of the solver than usual.
250: . -ksp_fischer_guess <model,size> - uses the Fischer initial guess generator for repeated linear solves
251: . -ksp_constant_null_space - assume the operator (matrix) has the constant vector in its null space
252: . -ksp_test_null_space - tests the null space set with KSPSetNullSpace() to see if it truly is a null space
253: . -ksp_knoll - compute initial guess by applying the preconditioner to the right hand side
254: . -ksp_monitor_cancel - cancel all previous convergene monitor routines set
255: . -ksp_monitor <optional filename> - print residual norm at each iteration
256: . -ksp_monitor_draw - plot residual norm at each iteration
257: . -ksp_monitor_solution - plot solution at each iteration
258: - -ksp_monitor_singular_value - monitor extremem singular values at each iteration
260: Notes:
261: To see all options, run your program with the -help option
262: or consult <A href="../../docs/manual.pdf#nameddest=Chapter 4 KSP: Linear Equations Solvers">KSP chapter of the users manual</A>.
264: Level: beginner
266: .keywords: KSP, set, from, options, database
268: .seealso: KSPSetUseFischerInitialGuess()
270: @*/
271: PetscErrorCode KSPSetFromOptions(KSP ksp)
272: {
273: PetscErrorCode ierr;
274: PetscInt indx;
275: const char *convtests[] = {"default","skip"};
276: char type[256], monfilename[PETSC_MAX_PATH_LEN];
277: PetscViewer monviewer;
278: PetscBool flg,flag;
279: PetscInt model[2],nmax;
280: void *ctx;
284: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
285: PCSetFromOptions(ksp->pc);
287: if (!KSPRegisterAllCalled) {KSPRegisterAll(PETSC_NULL);}
288: PetscObjectOptionsBegin((PetscObject)ksp);
289: PetscOptionsList("-ksp_type","Krylov method","KSPSetType",KSPList,(char*)(((PetscObject)ksp)->type_name?((PetscObject)ksp)->type_name:KSPGMRES),type,256,&flg);
290: if (flg) {
291: KSPSetType(ksp,type);
292: }
293: /*
294: Set the type if it was never set.
295: */
296: if (!((PetscObject)ksp)->type_name) {
297: KSPSetType(ksp,KSPGMRES);
298: }
300: PetscOptionsInt("-ksp_max_it","Maximum number of iterations","KSPSetTolerances",ksp->max_it,&ksp->max_it,PETSC_NULL);
301: PetscOptionsReal("-ksp_rtol","Relative decrease in residual norm","KSPSetTolerances",ksp->rtol,&ksp->rtol,PETSC_NULL);
302: PetscOptionsReal("-ksp_atol","Absolute value of residual norm","KSPSetTolerances",ksp->abstol,&ksp->abstol,PETSC_NULL);
303: PetscOptionsReal("-ksp_divtol","Residual norm increase cause divergence","KSPSetTolerances",ksp->divtol,&ksp->divtol,PETSC_NULL);
305: flag = PETSC_FALSE;
306: PetscOptionsBool("-ksp_converged_use_initial_residual_norm","Use initial residual residual norm for computing relative convergence","KSPDefaultConvergedSetUIRNorm",flag,&flag,PETSC_NULL);
307: if (flag) {KSPDefaultConvergedSetUIRNorm(ksp);}
308: flag = PETSC_FALSE;
309: PetscOptionsBool("-ksp_converged_use_min_initial_residual_norm","Use minimum of initial residual norm and b for computing relative convergence","KSPDefaultConvergedSetUMIRNorm",flag,&flag,PETSC_NULL);
310: if (flag) {KSPDefaultConvergedSetUMIRNorm(ksp);}
311: KSPGetInitialGuessNonzero(ksp,&flag);
312: PetscOptionsBool("-ksp_initial_guess_nonzero","Use the contents of the solution vector for initial guess","KSPSetInitialNonzero",flag,&flag,&flg);
313: if (flg) {
314: KSPSetInitialGuessNonzero(ksp,flag);
315: }
317: PetscOptionsBool("-ksp_knoll","Use preconditioner applied to b for initial guess","KSPSetInitialGuessKnoll",ksp->guess_knoll,&ksp->guess_knoll,PETSC_NULL);
318: PetscOptionsBool("-ksp_error_if_not_converged","Generate error if solver does not converge","KSPSetErrorIfNotConverged",ksp->errorifnotconverged,&ksp->errorifnotconverged,PETSC_NULL);
319: nmax = 2;
320: PetscOptionsIntArray("-ksp_fischer_guess","Use Paul Fischer's algorithm for initial guess","KSPSetUseFischerGuess",model,&nmax,&flag);
321: if (flag) {
322: if (nmax != 2) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Must pass in model,size as arguments");
323: KSPSetUseFischerGuess(ksp,model[0],model[1]);
324: }
326: PetscOptionsEList("-ksp_convergence_test","Convergence test","KSPSetConvergenceTest",convtests,2,"default",&indx,&flg);
327: if (flg) {
328: switch (indx) {
329: case 0:
330: KSPDefaultConvergedCreate(&ctx);
331: KSPSetConvergenceTest(ksp,KSPDefaultConverged,ctx,KSPDefaultConvergedDestroy);
332: break;
333: case 1: KSPSetConvergenceTest(ksp,KSPSkipConverged,PETSC_NULL,PETSC_NULL); break;
334: }
335: }
337: PetscOptionsEList("-ksp_norm_type","KSP Norm type","KSPSetNormType",KSPNormTypes,4,"preconditioned",&indx,&flg);
338: if (flg) { KSPSetNormType(ksp,(KSPNormType)indx); }
340: PetscOptionsInt("-ksp_check_norm_iteration","First iteration to compute residual norm","KSPSetCheckNormIteration",ksp->chknorm,&ksp->chknorm,PETSC_NULL);
342: flag = ksp->lagnorm;
343: PetscOptionsBool("-ksp_lag_norm","Lag the calculation of the residual norm","KSPSetLagNorm",flag,&flag,&flg);
344: if (flg) {
345: KSPSetLagNorm(ksp,flag);
346: }
348: KSPGetDiagonalScale(ksp,&flag);
349: PetscOptionsBool("-ksp_diagonal_scale","Diagonal scale matrix before building preconditioner","KSPSetDiagonalScale",flag,&flag,&flg);
350: if (flg) {
351: KSPSetDiagonalScale(ksp,flag);
352: }
353: KSPGetDiagonalScaleFix(ksp,&flag);
354: PetscOptionsBool("-ksp_diagonal_scale_fix","Fix diagonally scaled matrix after solve","KSPSetDiagonalScaleFix",flag,&flag,&flg);
355: if (flg) {
356: KSPSetDiagonalScaleFix(ksp,flag);
357: }
359: flg = PETSC_FALSE;
360: PetscOptionsBool("-ksp_constant_null_space","Add constant null space to Krylov solver","KSPSetNullSpace",flg,&flg,PETSC_NULL);
361: if (flg) {
362: MatNullSpace nsp;
364: MatNullSpaceCreate(((PetscObject)ksp)->comm,PETSC_TRUE,0,0,&nsp);
365: KSPSetNullSpace(ksp,nsp);
366: MatNullSpaceDestroy(&nsp);
367: }
369: /* option is actually checked in KSPSetUp(), just here so goes into help message */
370: if (ksp->nullsp) {
371: PetscOptionsName("-ksp_test_null_space","Is provided null space correct","None",&flg);
372: }
374: /*
375: Prints reason for convergence or divergence of each linear solve
376: */
377: flg = PETSC_FALSE;
378: PetscOptionsBool("-ksp_converged_reason","Print reason for converged or diverged","KSPSolve",flg,&flg,PETSC_NULL);
379: if (flg) {
380: ksp->printreason = PETSC_TRUE;
381: }
383: flg = PETSC_FALSE;
384: PetscOptionsBool("-ksp_monitor_cancel","Remove any hardwired monitor routines","KSPMonitorCancel",flg,&flg,PETSC_NULL);
385: /* -----------------------------------------------------------------------*/
386: /*
387: Cancels all monitors hardwired into code before call to KSPSetFromOptions()
388: */
389: if (flg) {
390: KSPMonitorCancel(ksp);
391: }
392: /*
393: Prints preconditioned residual norm at each iteration
394: */
395: PetscOptionsString("-ksp_monitor","Monitor preconditioned residual norm","KSPMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
396: if (flg) {
397: PetscViewerASCIIOpen(((PetscObject)ksp)->comm,monfilename,&monviewer);
398: KSPMonitorSet(ksp,KSPMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
399: }
400: /*
401: Prints preconditioned residual norm at each iteration
402: */
403: PetscOptionsString("-ksp_monitor_range","Monitor percent of residual entries more than 10 percent of max","KSPMonitorRange","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
404: if (flg) {
405: PetscViewerASCIIOpen(((PetscObject)ksp)->comm,monfilename,&monviewer);
406: KSPMonitorSet(ksp,KSPMonitorRange,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
407: }
408: /*
409: Plots the vector solution
410: */
411: flg = PETSC_FALSE;
412: PetscOptionsBool("-ksp_monitor_solution","Monitor solution graphically","KSPMonitorSet",flg,&flg,PETSC_NULL);
413: if (flg) {
414: KSPMonitorSet(ksp,KSPMonitorSolution,PETSC_NULL,PETSC_NULL);
415: }
416: /*
417: Prints preconditioned and true residual norm at each iteration
418: */
419: PetscOptionsString("-ksp_monitor_true_residual","Monitor true residual norm","KSPMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
420: if (flg) {
421: PetscViewerASCIIOpen(((PetscObject)ksp)->comm,monfilename,&monviewer);
422: KSPMonitorSet(ksp,KSPMonitorTrueResidualNorm,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
423: }
424: /*
425: Prints extreme eigenvalue estimates at each iteration
426: */
427: PetscOptionsString("-ksp_monitor_singular_value","Monitor singular values","KSPMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
428: if (flg) {
429: KSPSetComputeSingularValues(ksp,PETSC_TRUE);
430: PetscViewerASCIIOpen(((PetscObject)ksp)->comm,monfilename,&monviewer);
431: KSPMonitorSet(ksp,KSPMonitorSingularValue,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
432: }
433: /*
434: Prints preconditioned residual norm with fewer digits
435: */
436: PetscOptionsString("-ksp_monitor_short","Monitor preconditioned residual norm with fewer digits","KSPMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
437: if (flg) {
438: PetscViewerASCIIOpen(((PetscObject)ksp)->comm,monfilename,&monviewer);
439: KSPMonitorSet(ksp,KSPMonitorDefaultShort,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
440: }
441: /*
442: Calls Python function
443: */
444: PetscOptionsString("-ksp_monitor_python","Use Python function","KSPMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);
445: if (flg) {PetscPythonMonitorSet((PetscObject)ksp,monfilename);}
446: /*
447: Graphically plots preconditioned residual norm
448: */
449: flg = PETSC_FALSE;
450: PetscOptionsBool("-ksp_monitor_draw","Monitor graphically preconditioned residual norm","KSPMonitorSet",flg,&flg,PETSC_NULL);
451: if (flg) {
452: KSPMonitorSet(ksp,KSPMonitorLG,PETSC_NULL,PETSC_NULL);
453: }
454: /*
455: Graphically plots preconditioned and true residual norm
456: */
457: flg = PETSC_FALSE;
458: PetscOptionsBool("-ksp_monitor_draw_true_residual","Monitor graphically true residual norm","KSPMonitorSet",flg,&flg,PETSC_NULL);
459: if (flg){
460: KSPMonitorSet(ksp,KSPMonitorLGTrueResidualNorm,PETSC_NULL,PETSC_NULL);
461: }
462: /*
463: Graphically plots preconditioned residual norm and range of residual element values
464: */
465: flg = PETSC_FALSE;
466: PetscOptionsBool("-ksp_monitor_range_draw","Monitor graphically range of preconditioned residual norm","KSPMonitorSet",flg,&flg,PETSC_NULL);
467: if (flg) {
468: KSPMonitorSet(ksp,KSPMonitorLGRange,PETSC_NULL,PETSC_NULL);
469: }
471: /*
472: Publishes convergence information using AMS
473: */
474: flg = PETSC_FALSE;
475: PetscOptionsBool("-ksp_monitor_ams","Publish KSP progress using AMS","KSPMonitorSet",flg,&flg,PETSC_NULL);
476: if (flg) {
477: char amscommname[256];
478: void *ctx;
479: PetscSNPrintf(amscommname,sizeof amscommname,"%sksp_monitor_ams",((PetscObject)ksp)->prefix?((PetscObject)ksp)->prefix:"");
480: KSPMonitorAMSCreate(ksp,amscommname,&ctx);
481: KSPMonitorSet(ksp,KSPMonitorAMS,ctx,KSPMonitorAMSDestroy);
482: KSPSetComputeSingularValues(ksp,PETSC_TRUE);
483: }
485: /* -----------------------------------------------------------------------*/
486: PetscOptionsEList("-ksp_pc_side","KSP preconditioner side","KSPSetPCSide",PCSides,3,PCSides[ksp->pc_side],&indx,&flg);
487: if (flg) {KSPSetPCSide(ksp,(PCSide)indx);}
489: flg = PETSC_FALSE;
490: PetscOptionsBool("-ksp_compute_singularvalues","Compute singular values of preconditioned operator","KSPSetComputeSingularValues",flg,&flg,PETSC_NULL);
491: if (flg) { KSPSetComputeSingularValues(ksp,PETSC_TRUE); }
492: flg = PETSC_FALSE;
493: PetscOptionsBool("-ksp_compute_eigenvalues","Compute eigenvalues of preconditioned operator","KSPSetComputeSingularValues",flg,&flg,PETSC_NULL);
494: if (flg) { KSPSetComputeSingularValues(ksp,PETSC_TRUE); }
495: flg = PETSC_FALSE;
496: PetscOptionsBool("-ksp_plot_eigenvalues","Scatter plot extreme eigenvalues","KSPSetComputeSingularValues",flg,&flg,PETSC_NULL);
497: if (flg) { KSPSetComputeSingularValues(ksp,PETSC_TRUE); }
500: if (ksp->ops->setfromoptions) {
501: (*ksp->ops->setfromoptions)(ksp);
502: }
503: /* actually check in setup this is just here so goes into help message */
504: PetscOptionsName("-ksp_view","View linear solver parameters","KSPView",&flg);
506: /* process any options handlers added with PetscObjectAddOptionsHandler() */
507: PetscObjectProcessOptionsHandlers((PetscObject)ksp);
508: PetscOptionsEnd();
509: return(0);
510: }