40 Real & current_residual,
46 if ((current_residual < last_residual) ||
61 (current_residual > last_residual &&
68 substepdivision =
std::min(0.5, last_residual/current_residual);
69 substepdivision =
std::max(substepdivision, tol*2.);
72 substepdivision = 0.5;
74 newton_iterate.
add (bx * (1.-substepdivision),
76 newton_iterate.
close();
77 bx *= substepdivision;
89 current_residual = rhs.
l2_norm();
92 << current_residual << std::endl;
96 (current_residual > last_residual)))
103 libmesh_convergence_failure();
123 Real x = bx, w = bx, v = bx;
126 Real fx = current_residual,
127 fw = current_residual,
128 fv = current_residual;
131 const unsigned int max_i = 20;
134 const Real golden_ratio = 1.-(std::sqrt(5.)-1.)/2.;
136 for (
unsigned int i=1; i <= max_i; i++)
138 Real xm = (ax+cx)*0.5;
140 Real tol2 = 2.0 * tol1;
143 if (
std::abs(x-xm) <= (tol2 - 0.5 * (cx - ax)))
151 Real r = (x-w)*(fx-fv);
152 Real q = (x-v)*(fx-fw);
153 Real p = (x-v)*q-(x-w)*r;
164 e = x >= xm ? ax-x : cx-x;
165 d = golden_ratio * e;
171 if (x+d-ax < tol2 || cx-(x+d) < tol2)
172 d =
SIGN(tol1, xm - x);
178 e = x >= xm ? ax-x : cx-x;
179 d = golden_ratio * e;
185 newton_iterate.
add (bx - u, linear_solution);
186 newton_iterate.
close();
209 fv = fw; fw = fx; fx = fu;
217 if (fu <= fw || w == x)
222 else if (fu <= fv || v == x || v == w)
231 libMesh::out <<
"Warning! Too many iterations used in Brent line search!" 239 require_residual_reduction(true),
240 require_finite_residual(true),
241 brent_line_search(true),
242 track_linear_convergence(false),
244 linear_tolerance_multiplier(1e-3),
284 LOG_SCOPE(
"solve()",
"NewtonSolver");
291 std::unique_ptr<NumericVector<Number>> linear_solution_ptr = newton_iterate.
zero_clone();
295 newton_iterate.
close();
296 linear_solution.
close();
299 #ifdef LIBMESH_ENABLE_CONSTRAINTS 329 <<
" with residual Not-a-Number" 331 libmesh_convergence_failure();
340 <<
" meets absolute tolerance " 349 if (current_residual == 0)
363 Real last_residual = current_residual;
375 << current_residual << std::endl;
378 current_linear_tolerance =
std::min (current_linear_tolerance,
395 linear_solution.
zero();
399 << current_linear_tolerance << std::endl;
402 const std::pair<unsigned int, Real> rval =
404 linear_solution, rhs, current_linear_tolerance,
412 if (linear_c_reason < 0)
417 libMesh::out <<
"Linear solver failed during Newton step, dropping out." 426 #ifdef LIBMESH_ENABLE_CONSTRAINTS 428 (
_system, &linear_solution,
true);
431 const unsigned int linear_steps = rval.first;
435 const bool linear_solve_finished =
439 libMesh::out <<
"Linear solve finished, step " << linear_steps
440 <<
", residual " << rval.second
449 newton_iterate.
add (-1., linear_solution);
450 newton_iterate.
close();
455 norm_total = newton_iterate.
l2_norm();
458 newton_iterate, norm_total,
474 current_residual = rhs.l2_norm();
477 << current_residual << std::endl;
481 linear_solve_finished &&
482 current_residual <= last_residual))
486 norm_delta, linear_solve_finished &&
487 current_residual <= last_residual);
496 last_residual, current_residual,
497 newton_iterate, linear_solution);
498 norm_delta *= steplength;
510 libMesh::out <<
" Nonlinear solver reached maximum step " 512 << current_residual << std::endl;
520 libmesh_convergence_failure();
526 norm_total = newton_iterate.
l2_norm();
534 << norm_delta / norm_total
535 <<
", |du| = " << norm_delta
541 linear_solve_finished))
545 norm_delta / steplength,
546 linear_solve_finished);
553 #ifdef LIBMESH_ENABLE_CONSTRAINTS 573 bool linear_solve_finished)
576 bool has_converged =
false;
582 has_converged =
true;
590 has_converged =
true;
594 if (!linear_solve_finished)
596 return has_converged;
603 has_converged =
true;
611 has_converged =
true;
614 return has_converged;
619 Real current_residual,
621 bool linear_solve_finished)
626 libMesh::out <<
" Nonlinear solver converged, step " << step_num
627 <<
", residual " << current_residual
634 << current_residual <<
" > " 642 libMesh::out <<
" Nonlinear solver converged, step " << step_num
643 <<
", residual reduction " 658 if (!linear_solve_finished)
664 libMesh::out <<
" Nonlinear solver converged, step " << step_num
665 <<
", absolute step size " 683 libMesh::out <<
" Nonlinear solver converged, step " << step_num
684 <<
", relative step size " bool continue_after_max_iterations
virtual void init() override
NewtonSolver(sys_type &system)
Real minimum_linear_tolerance
std::unique_ptr< LinearSolver< Number > > _linear_solver
Real absolute_step_tolerance
virtual std::unique_ptr< NumericVector< T > > zero_clone() const =0
unsigned int max_nonlinear_iterations
Real absolute_residual_tolerance
virtual unsigned int solve() override
NumericVector< Number > * rhs
static const Real TOLERANCE
Real initial_linear_tolerance
long double max(long double a, double b)
unsigned int _solve_result
unsigned int _inner_iterations
virtual void assembly(bool, bool, bool=false, bool=false)
virtual Real l2_norm() const =0
unsigned int max_linear_iterations
bool require_finite_residual
std::unique_ptr< NumericVector< Number > > solution
bool track_linear_convergence
Real linear_tolerance_multiplier
const SparseMatrix< Number > * request_matrix(const std::string &mat_name) const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
bool continue_after_backtrack_failure
bool libmesh_isnan(float a)
virtual void reinit() override
unsigned int _outer_iterations
SparseMatrix< Number > * matrix
bool test_convergence(Real current_residual, Real step_norm, bool linear_solve_finished)
bool require_residual_reduction
bool on_command_line(std::string arg)
const std::string & name() const
virtual void add(const numeric_index_type i, const T value)=0
Real line_search(Real tol, Real last_residual, Real ¤t_residual, NumericVector< Number > &newton_iterate, const NumericVector< Number > &linear_solution)
void print_convergence(unsigned int step_num, Real current_residual, Real step_norm, bool linear_solve_finished)
OStreamProxy out(std::cout)
const DofMap & get_dof_map() const
long double min(long double a, double b)
Real relative_residual_tolerance
Real relative_step_tolerance
std::unique_ptr< LinearSolutionMonitor > linear_solution_monitor
void enforce_constraints_exactly(const System &system, NumericVector< Number > *v=nullptr, bool homogeneous=false) const
Manages consistently variables, degrees of freedom, coefficient vectors, and matrices for implicit sy...