petsc-3.7.5 2017-01-01
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_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

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 Users-Manual: ch_snes

Keywords

SNES, nonlinear, set, options, database

See Also

SNESSetOptionsPrefix()

Level:beginner
Location:
src/snes/interface/snes.c
Index of all SNES routines
Table of Contents for all manual pages
Index of all manual pages

Examples

src/snes/examples/tutorials/ex1.c.html
src/snes/examples/tutorials/ex2.c.html
src/snes/examples/tutorials/ex3.c.html
src/snes/examples/tutorials/ex5.c.html
src/snes/examples/tutorials/ex5s.c.html
src/snes/examples/tutorials/ex7.c.html
src/snes/examples/tutorials/ex12.c.html
src/snes/examples/tutorials/ex14.c.html
src/snes/examples/tutorials/ex15.c.html
src/snes/examples/tutorials/ex18.c.html
src/snes/examples/tutorials/ex19.c.html