diff --git a/README.md b/README.md index 7e41d4b..32be5fd 100644 --- a/README.md +++ b/README.md @@ -307,6 +307,17 @@ By default, no linesearch is performed. **Note:** it is assumed that a passed linesearch function will at least update the solution vector and evaluate the function at the new point. +If `linesearch=LineSearches.Static()` (the default), an additional parameter +`apply_step!` (default value `(x, x_old, newton_step)->(x .= x_old .+ newton_step)`) can be used to +define a problem-specific function to update the `x` value after the Newton iteration. +For example, `apply_step! = (x, x_old, newton_step)->(x .= x_old .+ newton_step; x .= max.(x, 1e-80))` +will enforce a constraint `x[i] >= 1e-80`. + +If the optional parameter `always_step = false` (the default), no Newton iteration is performed +if the solution is already converged for the supplied initial `x` values, and the supplied initial `x` is returned +in the `zero` field of `SolverResults`. If `always_step = true`, a Newton step is performed +before testing for convergence. + ## Anderson acceleration This method is selected with `method = :anderson`. diff --git a/src/nlsolve/nlsolve.jl b/src/nlsolve/nlsolve.jl index 9bc06e8..2d1ca43 100644 --- a/src/nlsolve/nlsolve.jl +++ b/src/nlsolve/nlsolve.jl @@ -9,6 +9,8 @@ function nlsolve(df::Union{NonDifferentiable, OnceDifferentiable}, extended_trace::Bool = false, linesearch = LineSearches.Static(), linsolve=(x, A, b) -> copyto!(x, A\b), + apply_step! = (x, x_old, newton_step)->(x .= x_old .+ newton_step), + always_step = false, factor::Real = one(real(eltype(initial_x))), autoscale::Bool = true, m::Integer = 10, @@ -21,7 +23,7 @@ function nlsolve(df::Union{NonDifferentiable, OnceDifferentiable}, end if method == :newton newton(df, initial_x, xtol, ftol, iterations, - store_trace, show_trace, extended_trace, linesearch; linsolve=linsolve) + store_trace, show_trace, extended_trace, linesearch; linsolve=linsolve, apply_step! =apply_step!, always_step=always_step) elseif method == :trust_region trust_region(df, initial_x, xtol, ftol, iterations, store_trace, show_trace, extended_trace, factor, diff --git a/src/solvers/newton.jl b/src/solvers/newton.jl index d9a424b..4f16a78 100644 --- a/src/solvers/newton.jl +++ b/src/solvers/newton.jl @@ -42,6 +42,8 @@ function newton_(df::OnceDifferentiable, extended_trace::Bool, linesearch, linsolve, + apply_step!, + always_step, cache = NewtonCache(df)) where T n = length(initial_x) copyto!(cache.x, initial_x) @@ -52,7 +54,12 @@ function newton_(df::OnceDifferentiable, x_converged, f_converged = assess_convergence(initial_x, cache.xold, value(df), NaN, ftol) stopped = any(isnan, cache.x) || any(isnan, value(df)) ? true : false - converged = x_converged || f_converged + if always_step + converged = false + else + converged = x_converged || f_converged + end + x_ls = copy(cache.x) tr = SolverTrace() tracing = store_trace || show_trace || extended_trace @@ -110,7 +117,8 @@ function newton_(df::OnceDifferentiable, if linesearch isa Static - x_ls .= cache.x .+ cache.p + # default is x_ls .= cache.x .+ cache.p + apply_step!(x_ls, cache.x, cache.p) value_jacobian!(df, x_ls) alpha, ϕalpha = one(real(T)), value(dfo) else @@ -142,6 +150,8 @@ function newton(df::OnceDifferentiable, extended_trace::Bool, linesearch, cache = NewtonCache(df); - linsolve=(x, A, b) -> copyto!(x, A\b)) where T - newton_(df, initial_x, convert(real(T), xtol), convert(real(T), ftol), iterations, store_trace, show_trace, extended_trace, linesearch, linsolve, cache) + linsolve=(x, A, b) -> copyto!(x, A\b), + apply_step! = (x, x_old, newton_step)->(x .= x_old .+ newton_step), + always_step::Bool) where T + newton_(df, initial_x, convert(real(T), xtol), convert(real(T), ftol), iterations, store_trace, show_trace, extended_trace, linesearch, linsolve, apply_step!, always_step, cache) end