:orphan: # SNESSetFromOptions Sets various `SNES` and `KSP` parameters from user options. ## Synopsis ``` #include "petscsnes.h" PetscErrorCode SNESSetFromOptions(SNES snes) ``` Collective ## Input Parameter - ***snes -*** the `SNES` context ## Options Database Keys - ***-snes_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 -*** absolute tolerance of residual norm - ***-snes_rtol -*** relative decrease in tolerance norm from initial - ***-snes_divergence_tolerance -*** if the residual goes above divtol*rnorm0, exit with divergence - ***-snes_force_iteration -*** force SNESSolve() to take at least one iteration - ***-snes_max_it -*** maximum number of iterations - ***-snes_max_funcs -*** maximum number of function evaluations - ***-snes_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 -*** how often preconditioner is rebuilt (use -1 to never rebuild) - ***-snes_lag_preconditioner_persists -*** retains the -snes_lag_preconditioner information across multiple SNESSolve() - ***-snes_lag_jacobian -*** how often Jacobian is rebuilt (use -1 to never rebuild) - ***-snes_lag_jacobian_persists -*** retains the -snes_lag_jacobian information across multiple SNESSolve() - ***-snes_tr_tol -*** trust region tolerance - ***-snes_convergence_test -*** convergence test in nonlinear solver. default `SNESConvergedDefault()`. skip `SNESConvergedSkip()` means continue iterating until max_it or some other criterion is reached, saving expense of convergence test. correct_pressure S`NESConvergedCorrectPressure()` has special handling of a pressure null space. - ***-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_monitor_pause_final -*** Pauses all monitor drawing after the solver ends - ***-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 -*** the SNES type to use as a nonlinear preconditioner - ***-snes_test_jacobian -*** compare the user provided Jacobian with one computed via finite differences to check for errors. If a threshold is given, display only those entries whose difference is greater than the threshold. - ***-snes_test_jacobian_view -*** display the user provided Jacobian, the finite difference Jacobian and the difference between them to help users detect the location of errors in the user provided Jacobian. ## Options Database Keys 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 -*** Sets rtol0 - ***-snes_ksp_ew_rtolmax -*** Sets rtolmax - ***-snes_ksp_ew_gamma -*** Sets gamma - ***-snes_ksp_ew_alpha -*** Sets alpha - ***-snes_ksp_ew_alpha2 -*** Sets alpha2 - ***-snes_ksp_ew_threshold -*** Sets threshold ## Notes To see all options, run your program with the -help option or consult the users manual `SNES` supports three approaches for computing (approximate) Jacobians: user provided via `SNESSetJacobian()`, matrix free, and computing explicitly with finite differences and coloring using `MatFDColoring`. It is also possible to use automatic differentiation and the `MatFDColoring` object. ## See Also [](ch_snes), `SNESType`, `SNESSetOptionsPrefix()`, `SNESResetFromOptions()`, `SNES`, `SNESCreate()` ## Level beginner ## Location src/snes/interface/snes.c ## Examples src/snes/tutorials/ex1.c
src/snes/tutorials/ex12.c
src/snes/tutorials/ex13.c
src/snes/tutorials/ex14.c
src/snes/tutorials/ex15.c
src/snes/tutorials/ex16.c
src/snes/tutorials/ex17.c
src/snes/tutorials/ex18.c
src/snes/tutorials/ex19.c
src/snes/tutorials/ex1f.F90
src/snes/tutorials/ex2.c
## 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_NEWTONTRDC in src/snes/impls/ntrdc/ntrdc.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
--- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/snes/interface/snes.c) [Index of all SNES routines](index.md) [Table of Contents for all manual pages](/manualpages/index.md) [Index of all manual pages](/manualpages/singleindex.md)