SNES
The SNES (Scalable Nonlinear Equations Solvers) module provides methods for solving nonlinear systems of the form F(x) = 0. It builds on KSP for the linear solves within Newton-like methods.
Overview
SNES provides:
- Newton methods: Newton line search, Newton trust region
- Quasi-Newton: L-BFGS, Broyden
- Nonlinear Richardson: With various line search strategies
- FAS multigrid: Full Approximation Scheme for nonlinear problems
- Composite solvers: Combine multiple nonlinear solvers
Creating a SNES Solver
# Basic creation
snes = SNES(petsclib, MPI.COMM_WORLD)
# With options
snes = SNES(petsclib, MPI.COMM_WORLD;
snes_type = "newtonls",
snes_rtol = 1e-8,
snes_max_it = 50
)Setting the Nonlinear Function
Define the residual function F(x):
function residual!(fx, snes, x)
# Compute F(x) and store in fx
fx[1] = x[1]^2 + x[2] - 1
fx[2] = x[1] + x[2]^2 - 1
return 0
end
setfunction!(snes, residual!, f_vec)Setting the Jacobian
Define the Jacobian J = dF/dx:
function jacobian!(J, snes, x)
# Fill Jacobian matrix
J[1, 1] = 2*x[1]
J[1, 2] = 1.0
J[2, 1] = 1.0
J[2, 2] = 2*x[2]
assemble!(J)
return 0
end
setjacobian!(snes, jacobian!, J, J) # (J, P) where P is preconditioner matrixUsing a DM
For PDE problems, associate the SNES with a DM:
setDM!(snes, dm)
# Get the DM from SNES
dm = getDMDA(snes)Solving
# Solve with initial guess x
solve!(x, snes)
# Get solution vector
sol = get_solution(snes)Common Solver Options
Nonlinear Solver Types (snes_type)
newtonls- Newton with line search (default)newtontr- Newton with trust regionnrichardson- Nonlinear Richardsonqn- Quasi-Newton (L-BFGS)fas- Full Approximation Scheme multigrid
Convergence Options
snes_rtol- Relative tolerancesnes_atol- Absolute tolerancesnes_stol- Step tolerancesnes_max_it- Maximum iterationssnes_monitor- Print residual each iteration
Line Search Options
snes_linesearch_type-bt(backtracking),basic,l2,cp
Example: Full Setup
petsclib = PETSc.petsclibs[1]
PETSc.initialize(petsclib)
snes = SNES(petsclib, MPI.COMM_WORLD;
snes_monitor = true,
ksp_type = "gmres",
pc_type = "ilu"
)
setfunction!(snes, residual!, f)
setjacobian!(snes, jacobian!, J, J)
setfromoptions!(snes)
solve!(x, snes)Functions
PETSc.SNES — Method
SNES(petsclib, comm::MPI.Comm; prefix="", options...)Create a PETSc nonlinear solver (SNES) context on the communicator comm.
Arguments
petsclib: The PETSc library instancecomm::MPI.Comm: MPI communicatorprefix::String: Optional prefix for command-line optionsoptions...: Additional PETSc options as keyword arguments
If comm has size 1, the garbage collector will handle cleanup automatically. Otherwise, the user is responsible for calling destroy.
External Links
- PETSc Manual:
SNES/SNESCreate
- PETSc Manual:
SNES/SNESSetFromOptions
PETSc.destroy — Method
destroy(snes::AbstractPetscSNES)Destroy a SNES (nonlinear solver) object and release associated resources.
This function is typically called automatically via finalizers when the object is garbage collected, but can be called explicitly to free resources immediately.
External Links
- PETSc Manual:
SNES/SNESDestroy
PETSc.getDMDA — Method
dm = getDMDA(snes::AbstractPetscSNES)Get dmda for snes
The returned dmda is owned by the snes
External Links
- PETSc Manual:
SNES/SNESGetDM
PETSc.setDM! — Method
setDM!(snes::AbstractPetscSNES, dm::AbstractDM)Set dm for snes
External Links
- PETSc Manual:
SNES/SNESSetDM
PETSc.setfunction! — Method
setfunction!(snes, f!, vec)
setfunction!(f!, snes, vec)Set the residual function f! for the nonlinear solver snes.
The function f! will be called as f!(fx, snes, x) where:
fx: Output vector to store the residual F(x)snes: The SNES contextx: Input vector with current solution
The vec argument is a template vector used for the residual.
External Links
- PETSc Manual:
SNES/SNESSetFunction
PETSc.setjacobian! — Function
setjacobian!(
snes::AbstractSNES,
updateJ!::Function,
J::AbstractMat,
P::AbstractMat = J
)
setjacobian!(
updateJ!::Function,
snes::AbstractSNES,
J::AbstractMat,
P::AbstractMat = J
)Define updateJ! to be the function that updates the Jacobian of the snes.
If J == P then a call to updateJ!(J, snes, x) should set the elements of the PETSc Jacobian (approximation).
If J ≠ P then a call to updateJ!(J, P, snes, x) should set the elements of the PETSc Jacobian (approximation) and preconditioning matrix P.
If you set snes.user_ctx, then updateJ! may optionally accept that as an additional last argument:
updateJ!(J, snes, x, user_ctx)whenJ == PupdateJ!(J, P, snes, x, user_ctx)whenJ ≠ P
External Links
- PETSc Manual:
SNES/SNESSetJacobian