Go to the documentation of this file.
20 #ifndef LIBMESH_NEWTON_SOLVER_H
21 #define LIBMESH_NEWTON_SOLVER_H
24 #include "libmesh/libmesh_common.h"
25 #include "libmesh/linear_solver.h"
26 #include "libmesh/reference_counted_object.h"
27 #include "libmesh/diff_solver.h"
67 virtual void init ()
override;
73 virtual void reinit ()
override;
79 virtual unsigned int solve ()
override;
164 Real & current_residual,
173 Real current_residual,
175 bool linear_solve_finished);
183 bool linear_solve_finished);
190 #endif // LIBMESH_NEWTON_SOLVER_H
Manages consistently variables, degrees of freedom, coefficient vectors, and matrices for implicit sy...
bool track_linear_convergence
If set to true, check for convergence of the linear solve.
bool require_residual_reduction
If this is set to true, the solver is forced to test the residual after each Newton step,...
virtual void init() override
The initialization function.
The libMesh namespace provides an interface to certain functionality in the library.
std::unique_ptr< LinearSolver< Number > > _linear_solver
The LinearSolver defines the interface used to solve the linear_implicit system.
This is a generic class that defines a solver to handle ImplicitSystem classes, including NonlinearIm...
const sys_type & system() const
double linear_tolerance_multiplier
The tolerance for linear solves is kept below this multiplier (which defaults to 1e-3) times the norm...
const LinearSolver< Number > & get_linear_solver() const
Real line_search(Real tol, Real last_residual, Real ¤t_residual, NumericVector< Number > &newton_iterate, const NumericVector< Number > &linear_solution)
This does a line search in the direction opposite linear_solution to try and minimize the residual of...
void print_convergence(unsigned int step_num, Real current_residual, Real step_norm, bool linear_solve_finished)
This prints output for the convergence criteria based on by the given residual and step size.
This class defines a solver which uses the default libMesh linear solver in a quasiNewton method to h...
virtual void reinit() override
The reinitialization function.
bool test_convergence(Real current_residual, Real step_norm, bool linear_solve_finished)
LinearSolver< Number > & get_linear_solver()
virtual unsigned int solve() override
This method performs a solve, using an inexact Newton-Krylov method with line search.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
NewtonSolver(sys_type &system)
Constructor.
virtual ~NewtonSolver()
Destructor.
bool require_finite_residual
If this is set to true, the solver is forced to test the residual after each Newton step,...
Real minsteplength
If the quasi-Newton step length must be reduced to below this factor to give a residual reduction,...
bool brent_line_search
If require_residual_reduction is true, the solver may reduce step lengths when required.