PETScSNESSolver.h¶
Note
The documentation on this page was automatically extracted from the DOLFIN C++ code and may need to be edited or expanded.
-
class
PETScSNESSolver¶ Parent class(es)
This class implements methods for solving nonlinear systems via PETSc’s SNES interface. It includes line search and trust region techniques for globalising the convergence of the nonlinear iteration.
-
PETScSNESSolver(std::string nls_type = "default")¶ Create SNES solver for a particular method
-
std::pair<std::size_t, bool>
solve(NonlinearProblem &nonlinear_problem, GenericVector &x, const GenericVector &lb, const GenericVector &ub)¶ Solve a nonlinear variational inequality with bound constraints
- Arguments
- nonlinear_function (
NonlinearProblem) - The nonlinear problem.
- x (
GenericVector) - The vector.
- lb (
GenericVector) - The lower bound.
- ub (
GenericVector) - The upper bound.
- nonlinear_function (
- Returns
- std::pair<std::size_t, bool>
- Pair of number of Newton iterations, and whether iteration converged)
-
std::pair<std::size_t, bool>
solve(NonlinearProblem &nonlinear_function, GenericVector &x)¶ Solve abstract nonlinear problem \(F(x) = 0\) for given \(F\) and Jacobian \(\dfrac{\partial F}{\partial x}\).
- Arguments
- nonlinear_function (
NonlinearProblem) - The nonlinear problem.
- x (
GenericVector) - The vector.
- nonlinear_function (
- Returns
- std::pair<std::size_t, bool>
- Pair of number of Newton iterations, and whether iteration converged)
-
static std::vector<std::pair<std::string, std::string>>
methods()¶ Return a list of available solver methods
-
static Parameters
default_parameters()¶ Default parameter values
-
boost::shared_ptr<SNES>
snes() const¶ Return PETSc SNES pointer
-
void
init(const std::string &method)¶ Initialize SNES solver
-
void
set_linear_solver_parameters()¶ Update the linear solver parameters
-
static PetscErrorCode
FormFunction(SNES snes, Vec x, Vec f, void *ctx)¶ The callback for PETSc to compute F, the nonlinear residual
-
static PetscErrorCode
FormJacobian(SNES snes, Vec x, Mat *A, Mat *B, MatStructure *flag, void *ctx)¶ The callback for PETSc to compute A, the Jacobian
-
void
set_bounds(GenericVector &x)¶ Set the bounds on the problem from the parameters, if desired Here, x is passed in as a model vector from which we make our Vecs that tell PETSc the bounds if the “sign” parameter is used.
-