Skip to content

adding linsolve option for trust region method solver #278

Open
@kpobrien

Description

Hi, thanks for writing and maintaining this package! I'd like to get input on adding a linsolve option to the trust region method solver, identical to the linsolve option already available for the newton method solver.

A compelling use case is to perform a KLU factorization or other non-default factorizations before solving the linear system. For some problems, this leads to improvements in speed and memory usage. The proposed changes don't change the behavior or speed when used with default options. For example:

default linear solver

@time out=nlsolve(od,method = :trust_region,x)
14.935686 seconds (25.79 k allocations: 20.709 GiB, 1.04% gc time)
12.432642 seconds (25.83 k allocations: 20.709 GiB, 0.39% gc time)
15.237596 seconds (25.93 k allocations: 20.709 GiB, 0.35% gc time)

default behavior with changes below

@time out=nlsolve(od,method = :trust_region,x)
13.542611 seconds (25.98 k allocations: 20.709 GiB, 0.38% gc time)
15.107630 seconds (25.96 k allocations: 20.709 GiB, 0.36% gc time)
13.196935 seconds (25.96 k allocations: 20.709 GiB, 0.40% gc time)

perform KLU before solving the linear system with changes below

@time out=nlsolve(od,method = :trust_region,x,linsolve=(x, A, b) -> copyto!(x,klu(A)\b))
8.797286 seconds (28.68 k allocations: 13.647 GiB, 0.92% gc time)
8.853783 seconds (28.64 k allocations: 13.647 GiB, 0.55% gc time)
8.579454 seconds (28.58 k allocations: 13.647 GiB, 0.67% gc time)

Here are the diffs for the two files I changed, nlsolve.jl and trust_region.jl. I can make a PR if that's easier.

diff solvers/trust_region.jl ~/.julia/packages/NLsolve/gJL1I/src/solvers/trust_region.jl
48c48
<                  r, d, J, delta::Real)
---
>                  r, d, J, delta::Real, linsolve)
51c51
<         copyto!(p_i, J \ vec(r)) # Gauss-Newton step
---
>         linsolve(p_i, J, vec(r)) # Gauss-Newton step
116a117
>                        linsolve,
170c171
<         dogleg!(cache.p, cache.p_c, cache.pi, cache.r, cache.d, jacobian(df), delta)
---
>         dogleg!(cache.p, cache.p_c, cache.pi, cache.r, cache.d, jacobian(df), delta, linsolve)
234,235c235,237
<                       cache = NewtonTrustRegionCache(df)) where T
<     trust_region_(df, initial_x, convert(real(T), xtol), convert(real(T), ftol), iterations, store_trace, show_trace, extended_trace, convert(real(T), factor), autoscale, cache)
---
>                       cache = NewtonTrustRegionCache(df);
>                       linsolve=(x, A, b) -> copyto!(x, A\b)) where T
>     trust_region_(df, initial_x, convert(real(T), xtol), convert(real(T), ftol), iterations, store_trace, show_trace, extended_trace, convert(real(T), factor), autoscale, linsolve, cache)
diff nlsolve/nlsolve.jl ~/.julia/packages/NLsolve/gJL1I/src/nlsolve/nlsolve.jl
28c28
<                      autoscale)
---
>                      autoscale; linsolve=linsolve)

Activity

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment

Metadata

Assignees

No one assigned

    Labels

    No labels
    No labels

    Type

    No type

    Projects

    No projects

    Milestone

    No milestone

    Relationships

    None yet

    Development

    No branches or pull requests

    Issue actions