petsc-3.13.6 2020-09-29
Report Typos and Errors

SNESSetFromOptions

Sets various SNES and KSP parameters from user options.

Synopsis

#include "petscsnes.h"  
PetscErrorCode  SNESSetFromOptions(SNES snes)
Collective on SNES

Input Parameter

snes -the SNES context

Options Database Keys

-snes_type <type> - newtonls, newtontr, ngmres, ncg, nrichardson, qn, vi, fas, SNESType for complete list
-snes_stol - convergence tolerance in terms of the norm of the change in the solution between steps
-snes_atol <abstol> - absolute tolerance of residual norm
-snes_rtol <rtol> - relative decrease in tolerance norm from initial
-snes_divergence_tolerance <divtol> - if the residual goes above divtol*rnorm0, exit with divergence
-snes_force_iteration <force> - force SNESSolve() to take at least one iteration
-snes_max_it <max_it> - maximum number of iterations
-snes_max_funcs <max_funcs> - maximum number of function evaluations
-snes_max_fail <max_fail> - maximum number of line search failures allowed before stopping, default is none
-snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops
-snes_lag_preconditioner <lag> - how often preconditioner is rebuilt (use -1 to never rebuild)
-snes_lag_jacobian <lag> - how often Jacobian is rebuilt (use -1 to never rebuild)
-snes_trtol <trtol> - trust region tolerance
-snes_no_convergence_test - skip convergence test in nonlinear solver; hence iterations will continue until max_it or some other criterion is reached. Saves expense of convergence test
-snes_monitor [ascii][:filename][:viewer format] - prints residual norm at each iteration. if no filename given prints to stdout
-snes_monitor_solution [ascii binary draw][:filename][:viewer format] - plots solution at each iteration
-snes_monitor_residual [ascii binary draw][:filename][:viewer format] - plots residual (not its norm) at each iteration
-snes_monitor_solution_update [ascii binary draw][:filename][:viewer format] - plots update to solution at each iteration
-snes_monitor_lg_residualnorm - plots residual norm at each iteration
-snes_monitor_lg_range - plots residual norm at each iteration
-snes_fd - use finite differences to compute Jacobian; very slow, only for testing
-snes_fd_color - use finite differences with coloring to compute Jacobian
-snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
-snes_converged_reason - print the reason for convergence/divergence after each solve
-npc_snes_type <type> - the SNES type to use as a nonlinear preconditioner

Options Database for Eisenstat-Walker method

-snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
-snes_ksp_ew_version ver - version of Eisenstat-Walker method
-snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
-snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
-snes_ksp_ew_gamma <gamma> - Sets gamma
-snes_ksp_ew_alpha <alpha> - Sets alpha
-snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
-snes_ksp_ew_threshold <threshold> - Sets threshold

Notes

To see all options, run your program with the -help option or consult the users manual

Notes

SNES supports three approaches for computing (approximate) Jacobians: user provided via SNESSetJacobian(), matrix free, and computing explictly with finite differences and coloring using MatFDColoring. It is also possible to use automatic differentiation and the MatFDColoring object.

See Also

SNESSetOptionsPrefix(), SNESResetFromOptions(), SNES, SNESCreate()

Level

beginner

Location

src/snes/interface/snes.c

Examples

src/snes/tutorials/ex1.c.html
src/snes/tutorials/ex2.c.html
src/snes/tutorials/ex3.c.html
src/snes/tutorials/ex5.c.html
src/snes/tutorials/ex9.c.html
src/snes/tutorials/ex12.c.html
src/snes/tutorials/ex14.c.html
src/snes/tutorials/ex15.c.html
src/snes/tutorials/ex18.c.html
src/snes/tutorials/ex19.c.html
src/snes/tutorials/ex21.c.html

Implementations

SNESSetFromOptions_Composite in src/snes/impls/composite/snescomposite.c
SNESSetFromOptions_FAS in src/snes/impls/fas/fas.c
SNESSetFromOptions_NGS in src/snes/impls/gs/snesgs.c
SNESSetFromOptions_NEWTONLS in src/snes/impls/ls/ls.c
SNESSetFromOptions_MS in src/snes/impls/ms/ms.c
SNESSetFromOptions_Multiblock in src/snes/impls/multiblock/multiblock.c
SNESSetFromOptions_NASM in src/snes/impls/nasm/nasm.c
SNESSetFromOptions_NCG in src/snes/impls/ncg/snesncg.c
SNESSetFromOptions_Anderson in src/snes/impls/ngmres/anderson.c
SNESSetFromOptions_NGMRES in src/snes/impls/ngmres/snesngmres.c
SNESSetFromOptions_Patch in src/snes/impls/patch/snespatch.c
SNESSetFromOptions_QN in src/snes/impls/qn/qn.c
SNESSetFromOptions_NRichardson in src/snes/impls/richardson/snesrichardson.c
SNESSetFromOptions_Shell in src/snes/impls/shell/snesshell.c
SNESSetFromOptions_NEWTONTR in src/snes/impls/tr/tr.c
SNESSetFromOptions_VINEWTONSSLS in src/snes/impls/vi/ss/viss.c
SNESSetFromOptions_VI in src/snes/impls/vi/vi.c

Index of all SNES routines
Table of Contents for all manual pages
Index of all manual pages