15#ifndef dealii_petsc_snes_h
16#define dealii_petsc_snes_h
20#ifdef DEAL_II_WITH_PETSC
29# include <petscsnes.h>
33# if defined(DEAL_II_HAVE_CXX20)
259 typename AMatrixType = PMatrixType>
262 std::constructible_from<
265 std::constructible_from<
268 std::constructible_from<AMatrixType, Mat>))
283 const MPI_Comm mpi_comm = PETSC_COMM_WORLD);
295 operator SNES()
const;
364 solve(VectorType &x, PMatrixType &P);
377 solve(VectorType &x, AMatrixType &A, PMatrixType &P);
387 std::function<void(
const VectorType &x, VectorType &res)>
residual;
398 std::function<void(
const VectorType &x, AMatrixType &A, PMatrixType &P)>
413 std::function<void(
const VectorType &x,
414 const unsigned int step_number,
444 std::function<void(
const VectorType &src, VectorType &dst)>
std::string options_prefix
int max_n_function_evaluations
void add_parameters(ParameterHandler &prm)
int maximum_non_linear_iterations
real_type absolute_tolerance
std::string snes_linesearch_type
real_type relative_tolerance
std::function< void(const VectorType &src, VectorType &dst)> solve_with_jacobian
void set_matrices(AMatrixType &A, PMatrixType &P)
std::function< void(const VectorType &x, real_type &energy_value)> energy
SmartPointer< PMatrixType, NonlinearSolver > P
unsigned int solve(VectorType &x, AMatrixType &A, PMatrixType &P)
std::function< void(const VectorType &x, VectorType &res)> residual
unsigned int solve(VectorType &x)
std::exception_ptr pending_exception
NonlinearSolver(const NonlinearSolverData &data=NonlinearSolverData(), const MPI_Comm mpi_comm=PETSC_COMM_WORLD)
MPI_Comm get_mpi_communicator() const
void set_matrix(PMatrixType &P)
unsigned int solve(VectorType &x, PMatrixType &P)
std::function< void(const VectorType &x, AMatrixType &A, PMatrixType &P)> jacobian
std::function< void(const VectorType &x)> setup_jacobian
SmartPointer< AMatrixType, NonlinearSolver > A
std::function< void(const VectorType &x, const unsigned int step_number, const real_type f_norm)> monitor
void reinit(const NonlinearSolverData &data)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_CXX20_REQUIRES(condition)
#define DEAL_II_NAMESPACE_CLOSE
NonlinearSolverData(const std::string &options_prefix="", const std::string &snes_type="", const std::string &snes_linesearch_type="", const real_type absolute_tolerance=0, const real_type relative_tolerance=0, const real_type step_tolerance=0, const int maximum_non_linear_iterations=-1, const int max_n_function_evaluations=-1)