SciML / NonlinearSolve.jl

High-performance and differentiation-enabled nonlinear solvers (Newton methods), bracketed rootfinding (bisection, Falsi), with sparsity and Newton-Krylov support.
https://docs.sciml.ai/NonlinearSolve/stable/
MIT License
240 stars 42 forks source link

New Jacobian cache causes copy bug #473

Closed gdalle closed 1 month ago

gdalle commented 1 month ago

Describe the bug 🐞

The new NonlinearSolve v3.15 (with DI in it) causes a BoundaryValueDiffEq solve to error or segfault.

Expected behavior

The solve should run without errors, as it does with NonlinearSolve 3.14.

Minimal Reproducible Example 👇

using BoundaryValueDiffEq

struct EmbeddedTorus
    R::Float64
    r::Float64
end

function affine_connection!(M::EmbeddedTorus, Zc, i, a, Xc, Yc)
    θ = a[1] .+ i[1]
    sinθ, cosθ = sincos(θ)
    Γ¹₂₂ = (M.R + M.r * cosθ) * sinθ / M.r
    Γ²₁₂ = -M.r * sinθ / (M.R + M.r * cosθ)

    Zc[1] = Xc[2] * Γ¹₂₂ * Yc[2]
    Zc[2] = Γ²₁₂ * (Xc[1] * Yc[2] + Xc[2] * Yc[1])
    return Zc
end

function bc1!(residual, u, p, t)
    mid = div(length(u[1]), 2)
    residual[1:mid] = u[1][1:mid] - a1
    return residual[(mid + 1):end] = u[end][1:mid] - a2
end

function chart_log_problem!(du, u, params, t)
    M, i = params
    mid = div(length(u), 2)
    a = u[1:mid]
    dx = u[(mid + 1):end]
    ddx = similar(dx)
    affine_connection!(M, ddx, i, a, dx, dx)
    ddx .*= -1
    du[1:mid] .= dx
    du[(mid + 1):end] .= ddx
    return du
end

M = EmbeddedTorus(3, 2)
a1 = [0.5, -1.2]
a2 = [-0.5, 0.3]
i = (0, 0)
solver = MIRK4()
dt = 0.05
tspan = (0.0, 1.0)
u0 = [vcat(a1, zero(a1)), vcat(a2, zero(a1))]
bvp1 = BVProblem(chart_log_problem!, bc1!, u0, tspan, (M, i))
sol1 = solve(bvp1, solver; dt=dt)

Error & Stacktrace ⚠️

julia> sol1 = solve(bvp1, solver; dt=dt)
ERROR: BoundsError: attempt to access 8×8 Matrix{Float64} at index [1:64]
Stacktrace:
  [1] throw_boundserror(A::Matrix{Float64}, I::Tuple{UnitRange{Int64}})
    @ Base ./abstractarray.jl:737
  [2] checkbounds
    @ ./abstractarray.jl:702 [inlined]
  [3] _copyto_impl!(dest::Matrix{Float64}, doffs::Int64, src::Matrix{Float64}, soffs::Int64, n::Int64)
    @ Base ./array.jl:374
  [4] copyto!
    @ ./array.jl:368 [inlined]
  [5] copyto!
    @ ./array.jl:388 [inlined]
  [6] __set_lincache_A
    @ ~/.julia/packages/NonlinearSolve/EVJf5/src/internal/linear_solve.jl:231 [inlined]
  [7] __update_A!
    @ ~/.julia/packages/NonlinearSolve/EVJf5/src/internal/linear_solve.jl:221 [inlined]
  [8] __update_A!
    @ ~/.julia/packages/NonlinearSolve/EVJf5/src/internal/linear_solve.jl:201 [inlined]
  [9] (::NonlinearSolve.LinearSolverCache{…})(; A::Matrix{…}, b::Vector{…}, linu::Vector{…}, du::Vector{…}, p::Nothing, weight::Nothing, cachedata::Nothing, reuse_A_if_factorization::Bool, verbose::Bool, kwargs::@Kwargs{…})
    @ NonlinearSolve ~/.julia/packages/NonlinearSolve/EVJf5/src/internal/linear_solve.jl:126
 [10] LinearSolverCache
    @ ~/.julia/packages/NonlinearSolve/EVJf5/src/internal/linear_solve.jl:120 [inlined]
 [11] __internal_solve!(cache::NonlinearSolve.NewtonDescentCache{…}, J::Matrix{…}, fu::Vector{…}, u::Vector{…}, idx::Val{…}; skip_solve::Bool, new_jacobian::Bool, kwargs::@Kwargs{…})
    @ NonlinearSolve ~/.julia/packages/NonlinearSolve/EVJf5/src/descent/newton.jl:84
 [12] __internal_solve! (repeats 2 times)
    @ ~/.julia/packages/NonlinearSolve/EVJf5/src/descent/newton.jl:74 [inlined]
 [13] __step!(cache::NonlinearSolve.GeneralizedFirstOrderAlgorithmCache{…}; recompute_jacobian::Nothing, kwargs::@Kwargs{})
    @ NonlinearSolve ~/.julia/packages/NonlinearSolve/EVJf5/src/core/generalized_first_order.jl:245
 [14] __step!
    @ ~/.julia/packages/NonlinearSolve/EVJf5/src/core/generalized_first_order.jl:226 [inlined]
 [15] #step!#131
    @ ~/.julia/packages/NonlinearSolve/EVJf5/src/core/generic.jl:50 [inlined]
 [16] step!
    @ ~/.julia/packages/NonlinearSolve/EVJf5/src/core/generic.jl:45 [inlined]
 [17] solve!(cache::NonlinearSolve.GeneralizedFirstOrderAlgorithmCache{…})
    @ NonlinearSolve ~/.julia/packages/NonlinearSolve/EVJf5/src/core/generic.jl:13
 [18] #__solve#130
    @ ~/.julia/packages/NonlinearSolve/EVJf5/src/core/generic.jl:4 [inlined]
 [19] __solve
    @ ~/.julia/packages/NonlinearSolve/EVJf5/src/core/generic.jl:1 [inlined]
 [20] macro expansion
    @ ~/.julia/packages/NonlinearSolve/EVJf5/src/default.jl:282 [inlined]
 [21] __solve(::NonlinearProblem{…}, ::NonlinearSolvePolyAlgorithm{…}; stats::SciMLBase.NLStats, alias_u0::Bool, verbose::Bool, kwargs::@Kwargs{…})
    @ NonlinearSolve ~/.julia/packages/NonlinearSolve/EVJf5/src/default.jl:255
 [22] __perform_mirk_iteration(cache::BoundaryValueDiffEq.MIRKCache{…}, abstol::Float64, adaptive::Bool; nlsolve_kwargs::@NamedTuple{}, kwargs::@Kwargs{…})
    @ BoundaryValueDiffEq ~/.julia/packages/BoundaryValueDiffEq/YN0of/src/solve/mirk.jl:175
 [23] __perform_mirk_iteration
    @ ~/.julia/packages/BoundaryValueDiffEq/YN0of/src/solve/mirk.jl:171 [inlined]
 [24] solve!(cache::BoundaryValueDiffEq.MIRKCache{…})
    @ BoundaryValueDiffEq ~/.julia/packages/BoundaryValueDiffEq/YN0of/src/solve/mirk.jl:152
 [25] #__solve#296
    @ ~/.julia/packages/BoundaryValueDiffEq/YN0of/src/BoundaryValueDiffEq.jl:51 [inlined]
 [26] __solve
    @ ~/.julia/packages/BoundaryValueDiffEq/YN0of/src/BoundaryValueDiffEq.jl:49 [inlined]
 [27] #solve_call#44
    @ ~/.julia/packages/DiffEqBase/vscDj/src/solve.jl:612 [inlined]
 [28] solve_call
    @ ~/.julia/packages/DiffEqBase/vscDj/src/solve.jl:569 [inlined]
 [29] solve_up(prob::BVProblem{…}, sensealg::Nothing, u0::Vector{…}, p::Tuple{…}, args::MIRK4{…}; kwargs::@Kwargs{…})
    @ DiffEqBase ~/.julia/packages/DiffEqBase/vscDj/src/solve.jl:1084
 [30] solve_up
    @ ~/.julia/packages/DiffEqBase/vscDj/src/solve.jl:1070 [inlined]
 [31] solve(prob::BVProblem{…}, args::MIRK4{…}; sensealg::Nothing, u0::Nothing, p::Nothing, wrap::Val{…}, kwargs::@Kwargs{…})
    @ DiffEqBase ~/.julia/packages/DiffEqBase/vscDj/src/solve.jl:1007
 [32] top-level scope
    @ ~/Work/GitHub/Julia/Sandbox/test.jl:47
Some type information was truncated. Use `show(err)` to see complete types.

Environment (please complete the following information):

julia> Pkg.status()
Status `~/Work/GitHub/Julia/Sandbox/Project.toml`
  [764a87c0] BoundaryValueDiffEq v5.10.0
  [8913a72c] NonlinearSolve v3.15.0
```julia julia> Pkg.status(; mode = PKGMODE_MANIFEST) Status `~/Work/GitHub/Julia/Sandbox/Manifest.toml` [47edcb42] ADTypes v1.9.0 [7d9f7c33] Accessors v0.1.38 [79e6a3ab] Adapt v4.0.4 [ec485272] ArnoldiMethod v0.4.0 [4fba245c] ArrayInterface v7.16.0 [4c555306] ArrayLayouts v1.10.3 [aae01518] BandedMatrices v1.7.5 [62783981] BitTwiddlingConvenienceFunctions v0.1.6 [764a87c0] BoundaryValueDiffEq v5.10.0 [2a0fbf3d] CPUSummary v0.2.6 [d360d2e6] ChainRulesCore v1.25.0 [fb6a15b2] CloseOpenIntervals v0.1.13 [38540f10] CommonSolve v0.2.4 [bbf7d656] CommonSubexpressions v0.3.1 [f70d9fcc] CommonWorldInvalidations v1.0.0 [34da2185] Compat v4.16.0 [a33af91c] CompositionsBase v0.1.2 [2569d6c7] ConcreteStructs v0.2.3 [187b0558] ConstructionBase v1.5.8 [adafc99b] CpuId v0.3.1 [9a962f9c] DataAPI v1.16.0 [864edb3b] DataStructures v0.18.20 [e2d170a0] DataValueInterfaces v1.0.0 [2b5f629d] DiffEqBase v6.156.1 [163ba53b] DiffResults v1.1.0 [b552c78f] DiffRules v1.15.1 ⌃ [a0c0ee7d] DifferentiationInterface v0.6.5 [ffbed154] DocStringExtensions v0.9.3 [4e289a0a] EnumX v1.0.4 [f151be2c] EnzymeCore v0.8.4 [d4d017d3] ExponentialUtilities v1.26.1 [e2ba6199] ExprTools v0.1.10 ⌅ [6b7a57c9] Expronicon v0.8.5 [9d29842c] FastAlmostBandedMatrices v0.1.3 [7034ab61] FastBroadcast v0.3.5 [9aa1b823] FastClosures v0.3.2 [29a986be] FastLapackInterface v2.0.4 [1a297f60] FillArrays v1.13.0 [6a86dc24] FiniteDiff v2.24.0 [f6369f11] ForwardDiff v0.10.36 [069b7b12] FunctionWrappers v1.1.3 [77dc65aa] FunctionWrappersWrappers v0.1.3 [46192b85] GPUArraysCore v0.1.6 [c145ed77] GenericSchur v0.5.4 [86223c79] Graphs v1.12.0 [3e5b6fbb] HostCPUFeatures v0.1.17 [615f187c] IfElse v0.1.1 [d25df0c9] Inflate v0.1.5 [3587e190] InverseFunctions v0.1.17 [92d709cd] IrrationalConstants v0.2.2 [82899510] IteratorInterfaceExtensions v1.0.0 [692b3bcd] JLLWrappers v1.6.0 [ef3ab10e] KLU v0.6.0 [ba0b0d4f] Krylov v0.9.6 [10f19ff3] LayoutPointers v0.1.17 [5078a376] LazyArrays v2.2.1 [87fe0de2] LineSearch v0.1.2 [d3d80556] LineSearches v7.3.0 [7ed4a6bd] LinearSolve v2.35.0 `dev/LinearSolve` [2ab3a3ac] LogExpFunctions v0.3.28 [bdcacae8] LoopVectorization v0.12.171 [d8e11817] MLStyle v0.4.17 [1914dd2f] MacroTools v0.5.13 [d125e4d3] ManualMemory v0.1.8 [a3b82374] MatrixFactorizations v3.0.1 [bb5d69b7] MaybeInplace v0.1.4 [46d2c3a1] MuladdMacro v0.2.4 [d41bc354] NLSolversBase v7.8.3 [77ba4419] NaNMath v1.0.2 [8913a72c] NonlinearSolve v3.15.0 [6fe1bfb0] OffsetArrays v1.14.1 [bac558e1] OrderedCollections v1.6.3 [1dea7af3] OrdinaryDiffEq v6.89.0 [89bda076] OrdinaryDiffEqAdamsBashforthMoulton v1.1.0 [6ad6398a] OrdinaryDiffEqBDF v1.1.2 [bbf590c4] OrdinaryDiffEqCore v1.6.0 [50262376] OrdinaryDiffEqDefault v1.1.0 [4302a76b] OrdinaryDiffEqDifferentiation v1.1.0 [9286f039] OrdinaryDiffEqExplicitRK v1.1.0 [e0540318] OrdinaryDiffEqExponentialRK v1.1.0 [becaefa8] OrdinaryDiffEqExtrapolation v1.1.0 [5960d6e9] OrdinaryDiffEqFIRK v1.1.1 [101fe9f7] OrdinaryDiffEqFeagin v1.1.0 [d3585ca7] OrdinaryDiffEqFunctionMap v1.1.1 [d28bc4f8] OrdinaryDiffEqHighOrderRK v1.1.0 [9f002381] OrdinaryDiffEqIMEXMultistep v1.1.0 [521117fe] OrdinaryDiffEqLinear v1.1.0 [1344f307] OrdinaryDiffEqLowOrderRK v1.2.0 [b0944070] OrdinaryDiffEqLowStorageRK v1.2.1 [127b3ac7] OrdinaryDiffEqNonlinearSolve v1.2.1 [c9986a66] OrdinaryDiffEqNordsieck v1.1.0 [5dd0a6cf] OrdinaryDiffEqPDIRK v1.1.0 [5b33eab2] OrdinaryDiffEqPRK v1.1.0 [04162be5] OrdinaryDiffEqQPRK v1.1.0 [af6ede74] OrdinaryDiffEqRKN v1.1.0 [43230ef6] OrdinaryDiffEqRosenbrock v1.2.0 [2d112036] OrdinaryDiffEqSDIRK v1.1.0 [669c94d9] OrdinaryDiffEqSSPRK v1.2.0 [e3e12d00] OrdinaryDiffEqStabilizedIRK v1.1.0 [358294b1] OrdinaryDiffEqStabilizedRK v1.1.0 [fa646aed] OrdinaryDiffEqSymplecticRK v1.1.0 [b1df2697] OrdinaryDiffEqTsit5 v1.1.0 [79d7bb75] OrdinaryDiffEqVerner v1.1.1 [65ce6f38] PackageExtensionCompat v1.0.2 [d96e819e] Parameters v0.12.3 [f517fe37] Polyester v0.7.16 [1d0040c9] PolyesterWeave v0.2.2 [d236fae5] PreallocationTools v0.4.24 [aea7be01] PrecompileTools v1.2.1 [21216c6a] Preferences v1.4.3 [3cdcf5f2] RecipesBase v1.3.4 [731186ca] RecursiveArrayTools v3.27.0 [f2c3362d] RecursiveFactorization v0.2.23 [189a3867] Reexport v1.2.2 [ae029012] Requires v1.3.0 [7e49a35a] RuntimeGeneratedFunctions v0.5.13 [94e857df] SIMDTypes v0.1.0 [476501e8] SLEEFPirates v0.6.43 [0bca4576] SciMLBase v2.55.0 [19f34311] SciMLJacobianOperators v0.1.0 [c0aeaf25] SciMLOperators v0.3.10 [53ae85a6] SciMLStructures v1.5.0 [efcf1570] Setfield v1.1.1 [727e6d20] SimpleNonlinearSolve v1.12.3 [699a6c99] SimpleTraits v0.9.4 [ce78b400] SimpleUnPack v1.1.0 [9f842d2f] SparseConnectivityTracer v0.6.6 [47a9eef4] SparseDiffTools v2.22.0 [0a514795] SparseMatrixColorings v0.4.6 [e56a9233] Sparspak v0.3.9 [276daf66] SpecialFunctions v2.4.0 [aedffcd0] Static v1.1.1 [0d7ed370] StaticArrayInterface v1.8.0 [90137ffa] StaticArrays v1.9.7 [1e83bf80] StaticArraysCore v1.4.3 [7792a7ef] StrideArraysCore v0.5.7 [2efcf032] SymbolicIndexingInterface v0.3.31 [3783bdb8] TableTraits v1.0.1 [bd369af6] Tables v1.12.0 [8290d209] ThreadingUtilities v0.5.2 [a759f4b9] TimerOutputs v0.5.24 [d5829a12] TriangularSolve v0.2.1 [410a4b4d] Tricks v0.1.9 [781d530d] TruncatedStacktraces v1.4.0 [3a884ed6] UnPack v1.0.2 [3d5dd08c] VectorizationBase v0.21.70 [19fa3120] VertexSafeGraphs v0.2.0 [1d5cc7b8] IntelOpenMP_jll v2024.2.1+0 [856f044c] MKL_jll v2024.2.0+0 [efe28fd5] OpenSpecFun_jll v0.5.5+0 [1317d2d5] oneTBB_jll v2021.12.0+0 [0dad84c5] ArgTools v1.1.1 [56f22d72] Artifacts [2a0f44e3] Base64 [ade2ca70] Dates [8ba89e20] Distributed [f43a241f] Downloads v1.6.0 [7b1f6079] FileWatching [9fa8497b] Future [b77e0a4c] InteractiveUtils [4af54fe1] LazyArtifacts [b27032c2] LibCURL v0.6.4 [76f85450] LibGit2 [8f399da3] Libdl [37e2e46d] LinearAlgebra [56ddb016] Logging [d6f4376e] Markdown [a63ad114] Mmap [ca575930] NetworkOptions v1.2.0 [44cfe95a] Pkg v1.10.0 [de0858da] Printf [3fa0cd96] REPL [9a3f8284] Random [ea8e919c] SHA v0.7.0 [9e88b42a] Serialization [1a1011a3] SharedArrays [6462fe0b] Sockets [2f01184e] SparseArrays v1.10.0 [10745b16] Statistics v1.10.0 [fa267f1f] TOML v1.0.3 [a4e569a6] Tar v1.10.0 [8dfed614] Test [cf7118a7] UUIDs [4ec0a83e] Unicode [e66e0078] CompilerSupportLibraries_jll v1.1.1+0 [deac9b47] LibCURL_jll v8.4.0+0 [e37daf67] LibGit2_jll v1.6.4+0 [29816b5a] LibSSH2_jll v1.11.0+1 [c8ffd9c3] MbedTLS_jll v2.28.2+1 [14a3606d] MozillaCACerts_jll v2023.1.10 [4536629a] OpenBLAS_jll v0.3.23+4 [05823500] OpenLibm_jll v0.8.1+2 [bea87d4a] SuiteSparse_jll v7.2.1+1 [83775a58] Zlib_jll v1.2.13+1 [8e850b90] libblastrampoline_jll v5.11.0+0 [8e850ede] nghttp2_jll v1.52.0+1 [3f19e933] p7zip_jll v17.4.0+2 Info Packages marked with ⌃ and ⌅ have new versions available. Those with ⌃ may be upgradable, but those with ⌅ are restricted by compatibility constraints from upgrading. To see why use `status --outdated -m` ```
julia> versioninfo()
Julia Version 1.10.5
Commit 6f3fdf7b362 (2024-08-27 14:19 UTC)
Build Info:
  Official https://julialang.org/ release
Platform Info:
  OS: Linux (x86_64-linux-gnu)
  CPU: 12 × Intel(R) Core(TM) i7-8850H CPU @ 2.60GHz
  WORD_SIZE: 64
  LIBM: libopenlibm
  LLVM: libLLVM-15.0.7 (ORCJIT, skylake)
Threads: 1 default, 0 interactive, 1 GC (on 12 virtual cores)
Environment:
  LD_LIBRARY_PATH = :/home/guillaume/Software/gurobi1002/linux64/lib
  JULIA_EDITOR = code
  JULIA_NUM_THREADS = 1
  JULIA_CONDAPKG_BACKEND = System
  JULIA_CONDAPKG_EXE = /home/guillaume/miniforge3/bin/mamba

Additional context

The bug is probably not deterministic. Sometimes it doesn't even happen, and sometimes it downright causes a segfault.

By pkg> dev-ing NonlinearSolve and adding some logging before line src/internal/linear_solve.jl:231, I obtained some more information about the matrices lincache.A and new_A:

┌ Warning: Trying to perform `copyto!(lincache.A::Matrix{Float64}, new_A::Matrix{Float64})
│   size(lincache.A) = (8, 8)
│   size(new_A) = (8, 8)
│   length(lincache.A) = 4
│   length(new_A) = 64
└ @ NonlinearSolve ~/.julia/dev/NonlinearSolve/src/internal/linear_solve.jl:231

In other words, the copyto! fails even though both matrices have the same size, because the one from lincache doesn't have the right length. Note that length(lincache.A) is not always 1.

gdalle commented 1 month ago

Further investigation: when I run the tests with

julia --check-bounds=yes

I get a different error:

julia> sol1 = solve(bvp1, solver; dt=dt)
ERROR: BoundsError: attempt to access 4-element UnitRange{Int64} at index [5]
Stacktrace:
  [1] throw_boundserror(A::UnitRange{Int64}, I::Int64)
    @ Base ./abstractarray.jl:737
  [2] getindex
    @ ./range.jl:930 [inlined]
  [3] _colorediteration!
    @ ~/.julia/packages/FiniteDiff/wm6iC/src/iteration_utils.jl:4 [inlined]
  [4] forwarddiff_color_jacobian!(J::SubArray{…}, f::BoundaryValueDiffEq.__Fix3{…}, x::Vector{…}, jac_cache::ForwardColorJacCache{…})
    @ SparseDiffTools ~/.julia/packages/SparseDiffTools/RGOJw/src/differentiation/compute_jacobian_ad.jl:392
  [5] sparse_jacobian!
    @ ~/.julia/packages/SparseDiffTools/RGOJw/src/highlevel/forward_mode.jl:63 [inlined]
  [6] __mirk_mpoint_jacobian!(J::Matrix{…}, ::BandedMatrices.BandedMatrix{…}, x::Vector{…}, bc_diffmode::AutoForwardDiff{…}, nonbc_diffmode::AutoSparse{…}, bc_diffcache::SparseDiffTools.ForwardDiffJacobianCache{…}, nonbc_diffcache::SparseDiffTools.ForwardDiffJacobianCache{…}, loss_bc::BoundaryValueDiffEq.__Fix3{…}, loss_collocation::BoundaryValueDiffEq.__Fix3{…}, resid_bc::Vector{…}, resid_collocation::Vector{…}, L::Int64)
    @ BoundaryValueDiffEq ~/.julia/packages/BoundaryValueDiffEq/YN0of/src/solve/mirk.jl:390
  [7] #268
    @ ~/.julia/packages/BoundaryValueDiffEq/YN0of/src/solve/mirk.jl:371 [inlined]
  [8] JacobianCache
    @ ~/.julia/dev/NonlinearSolve/src/internal/jacobian.jl:149 [inlined]
  [9] JacobianCache
    @ ~/.julia/dev/NonlinearSolve/src/internal/jacobian.jl:120 [inlined]
 [10] __step!(cache::NonlinearSolve.GeneralizedFirstOrderAlgorithmCache{…}; recompute_jacobian::Nothing, kwargs::@Kwargs{})
    @ NonlinearSolve ~/.julia/dev/NonlinearSolve/src/core/generalized_first_order.jl:230
 [11] __step!
    @ ~/.julia/dev/NonlinearSolve/src/core/generalized_first_order.jl:226 [inlined]
 [12] #step!#131
    @ ~/.julia/dev/NonlinearSolve/src/core/generic.jl:50 [inlined]
 [13] step!
    @ ~/.julia/dev/NonlinearSolve/src/core/generic.jl:45 [inlined]
 [14] solve!(cache::NonlinearSolve.GeneralizedFirstOrderAlgorithmCache{…})
    @ NonlinearSolve ~/.julia/dev/NonlinearSolve/src/core/generic.jl:13
 [15] #__solve#130
    @ ~/.julia/dev/NonlinearSolve/src/core/generic.jl:4 [inlined]
 [16] __solve
    @ ~/.julia/dev/NonlinearSolve/src/core/generic.jl:1 [inlined]
 [17] macro expansion
    @ ~/.julia/dev/NonlinearSolve/src/default.jl:282 [inlined]
 [18] __solve(::NonlinearProblem{…}, ::NonlinearSolvePolyAlgorithm{…}; stats::SciMLBase.NLStats, alias_u0::Bool, verbose::Bool, kwargs::@Kwargs{…})
    @ NonlinearSolve ~/.julia/dev/NonlinearSolve/src/default.jl:255
 [19] __perform_mirk_iteration(cache::BoundaryValueDiffEq.MIRKCache{…}, abstol::Float64, adaptive::Bool; nlsolve_kwargs::@NamedTuple{}, kwargs::@Kwargs{…})
    @ BoundaryValueDiffEq ~/.julia/packages/BoundaryValueDiffEq/YN0of/src/solve/mirk.jl:175
 [20] __perform_mirk_iteration
    @ ~/.julia/packages/BoundaryValueDiffEq/YN0of/src/solve/mirk.jl:171 [inlined]
 [21] solve!(cache::BoundaryValueDiffEq.MIRKCache{…})
    @ BoundaryValueDiffEq ~/.julia/packages/BoundaryValueDiffEq/YN0of/src/solve/mirk.jl:152
 [22] #__solve#296
    @ ~/.julia/packages/BoundaryValueDiffEq/YN0of/src/BoundaryValueDiffEq.jl:51 [inlined]
 [23] __solve
    @ ~/.julia/packages/BoundaryValueDiffEq/YN0of/src/BoundaryValueDiffEq.jl:49 [inlined]
 [24] #solve_call#44
    @ ~/.julia/packages/DiffEqBase/vscDj/src/solve.jl:612 [inlined]
 [25] solve_call
    @ ~/.julia/packages/DiffEqBase/vscDj/src/solve.jl:569 [inlined]
 [26] solve_up(prob::BVProblem{…}, sensealg::Nothing, u0::Vector{…}, p::Tuple{…}, args::MIRK4{…}; kwargs::@Kwargs{…})
    @ DiffEqBase ~/.julia/packages/DiffEqBase/vscDj/src/solve.jl:1084
 [27] solve_up
    @ ~/.julia/packages/DiffEqBase/vscDj/src/solve.jl:1070 [inlined]
 [28] solve(prob::BVProblem{…}, args::MIRK4{…}; sensealg::Nothing, u0::Nothing, p::Nothing, wrap::Val{…}, kwargs::@Kwargs{…})
    @ DiffEqBase ~/.julia/packages/DiffEqBase/vscDj/src/solve.jl:1007
 [29] top-level scope
    @ ~/Work/GitHub/Julia/Sandbox/test.jl:47
Some type information was truncated. Use `show(err)` to see complete types.

So the one I was observing before might have been due to incorrect use of @inbounds in FiniteDiff:

https://github.com/JuliaDiff/FiniteDiff.jl/blob/b574440c5d54708aa23a90c48f4a000dd1c6d915/src/iteration_utils.jl#L1-L7

gdalle commented 1 month ago

And in this same part of FiniteDiff, the array we're incorrectly accessing is rows_index because the matrix has more columns than rows. This is right before the error occurs:

i = 5
rows_index = 1:4
cols_index = 1:8
vfx = [1.0, 0.0, 0.0, 0.0]
ERROR: BoundsError: attempt to access 4-element UnitRange{Int64} at index [5]
ChrisRackauckas commented 1 month ago

[18] __solve(::NonlinearProblem{…}, ::NonlinearSolvePolyAlgorithm{…}; stats::SciMLBase.NLStats, alias_u0::Bool, verbose::Bool, kwargs::@Kwargs{…}) @ NonlinearSolve ~/.julia/dev/NonlinearSolve/src/default.jl:255

it should be square? What are the dimensions?

gdalle commented 1 month ago
size(J) = (4, 8)
size(sparsity) = (4, 8)
typeof(J) = SubArray{Float64, 2, Matrix{Float64}, Tuple{UnitRange{Int64}, Base.Slice{Base.OneTo{Int64}}}, false}
typeof(sparsity) = BandedMatrices.BandedMatrix{Float64, Matrix{Float64}, Base.OneTo{Int64}}
ChrisRackauckas commented 1 month ago

oh it's differencing the banded subproblem, I see.

ChrisRackauckas commented 1 month ago

Shouldn't this end up in the bandedmatrix iteration https://github.com/JuliaArrays/ArrayInterface.jl/blob/master/ext/ArrayInterfaceBandedMatricesExt.jl#L56-L62?

gdalle commented 1 month ago

It probably doesn't because J is not itself banded. But even when we miss that dispatch we should still get the correct answer and not an indexing error.

gdalle commented 1 month ago

I think there is something fundamentally wrong with the default FiniteDiff._colorediteration!, so if someone can tell me what each argument is for we can try to look into it.

ChrisRackauckas commented 1 month ago

Agreed, but it looks like there's two errors here. One is that NonlinearSolve changes are now missing the BandedMatrix dispatches which isn't great, but then the second is that the dense fallback should be corrected.

gdalle commented 1 month ago

Only the second one is an error, the first one is a suboptimality. I guess we need to ask @avik-pal where the jacobian prototype is built for this one, and how it ends up being a SubArray

ChrisRackauckas commented 1 month ago

https://github.com/JuliaDiff/FiniteDiff.jl/blob/b574440c5d54708aa23a90c48f4a000dd1c6d915/ext/FiniteDiffBandedMatricesExt.jl#L13 if you relax Jac::BandedMatrices.BandedMatrix, does it hit it and the error goes away?

avik-pal commented 1 month ago

Partial bisection done

git bisect start
# status: waiting for both good and bad commits
# good: [f667683bb6c5bf78a8cef5defa52f887c526fc71] feat: integrate SciMLJacobianOperators into NonlinearSolve
git bisect good f667683bb6c5bf78a8cef5defa52f887c526fc71
# status: waiting for bad commit, 1 good commit known
# bad: [d56476c69c749a4b575b47d58c03cf366ab7134c] fix: use DI preparation result when initializing Jacobian (#472)
git bisect bad d56476c69c749a4b575b47d58c03cf366ab7134c
# skip: [fe1bd3dd6ec4923afad60453600316038b3cd9f1] chore: bump minimum versions
git bisect skip fe1bd3dd6ec4923afad60453600316038b3cd9f1
# skip: [59b7d022bc282a04f858b973954b26fc52792771] fix: update bruss testing
git bisect skip 59b7d022bc282a04f858b973954b26fc52792771
# skip: [92c59dbb2d360a4a2ed968668395473be57fc81e] chore: run formatter
git bisect skip 92c59dbb2d360a4a2ed968668395473be57fc81e
# skip: [058b87d5af338ad7e2049f2f27c95b476798d32b] feat: remove Setfield dep
git bisect skip 058b87d5af338ad7e2049f2f27c95b476798d32b
# skip: [3f59f4949829b74d78b9eea8de834dd4f9239152] fix: autodiff cannot be nothing
git bisect skip 3f59f4949829b74d78b9eea8de834dd4f9239152
# skip: [40eb6c867fbed99368e7005ae1ba4246083269eb] chore: mark master as DEV (DON'T TAG)
git bisect skip 40eb6c867fbed99368e7005ae1ba4246083269eb
# skip: [adac560c3f0f737b5c854d61642a182078bb50a1] fix: __value_derivative removal from line searches
git bisect skip adac560c3f0f737b5c854d61642a182078bb50a1
# skip: [b2bb25252e3526b753b3852dfbfff9cc386e3f7a] fix: missing qualifier
git bisect skip b2bb25252e3526b753b3852dfbfff9cc386e3f7a
# skip: [3b12e73e945941ac8d0d4eebe038f1c8ae15d5cf] fix: unwrap sparse AD for derivative call
git bisect skip 3b12e73e945941ac8d0d4eebe038f1c8ae15d5cf
# skip: [43df907989a67b707e433753b09ef3bcd3d08a4c] test: structured jacobians
git bisect skip 43df907989a67b707e433753b09ef3bcd3d08a4c
# skip: [0f2a62f00334d680096e55c9a8ecdc8012819aef] fix: minor test fixes
git bisect skip 0f2a62f00334d680096e55c9a8ecdc8012819aef
# good: [bdacb9ca9327b9d5d2fe6a8f0c08f44ac206f8ec] docs: disable ILU for now
git bisect good bdacb9ca9327b9d5d2fe6a8f0c08f44ac206f8ec
# bad: [b6221b7d61027339ec5fc50bb1f3f90398ce607e] refactor: reorder imports
git bisect bad b6221b7d61027339ec5fc50bb1f3f90398ce607e
# skip: [ba45097e841f9875a7e53d94a1d86f1ffcfd1425] fix: sparse jacobian
git bisect skip ba45097e841f9875a7e53d94a1d86f1ffcfd1425
# good: [6201fcc2f49e2ec2019eb61491f30bc9664e0dc5] fix: DI now works with ReverseDiff
git bisect good 6201fcc2f49e2ec2019eb61491f30bc9664e0dc5
gdalle commented 1 month ago

Yes, if I relax that I get a lot of these warnings but at least the solve works:

┌ Warning: Potential Rank Deficient Matrix Detected. Attempting to solve using Pivoted QR Factorization.
└ @ NonlinearSolve ~/.julia/dev/NonlinearSolve/src/internal/linear_solve.jl:156
ChrisRackauckas commented 1 month ago

Okay cool, one mystery solved. And that looks optimal too, since it's using the 2D indexing. That's worth a PR.

The warnings, we should probably pass verbose false there but they make sense.

Now for the fallback, vfx is the vector form of f(x). colorvec should make sense. rows_index, cols_index = ArrayInterface.findstructralnz(sparsity)

https://github.com/JuliaArrays/ArrayInterface.jl/blob/2e2b81bfb53e9e9fc41f26b927fb980afcbfa9f5/src/ArrayInterface.jl#L330-L335

i.e. they are defined the same as findnz(::SparseMatrixCSC), which is that, in column order, (rows_index[i], cols_index[i]) is the (i,j) index of the ith non-zero element.

For dense then we should specialize on Matrix and do a fast path if square. If non-square, these should be:

julia> n = 4
4

julia> m = 5
5

julia> A = rand(4,5)
4×5 Matrix{Float64}:
 0.0732819  0.236601  0.971541  0.0239939  0.834559
 0.297407   0.472245  0.963563  0.841125   0.296606
 0.467831   0.297999  0.760833  0.87123    0.309619
 0.781913   0.879559  0.871359  0.120851   0.0749186

julia> col_idxs = vec(stack([i*ones(Int,n) for i in 1:m]))
20-element Vector{Int64}:
 1
 1
 1
 1
 2
 2
 2
 ⋮
 4
 4
 5
 5
 5
 5

julia> row_idxs = vec(stack([1:m for i in 1:n]))
20-element Vector{Int64}:
 1
 2
 3
 4
 5
 1
 2
 ⋮
 5
 1
 2
 3
 4
 5
avik-pal commented 1 month ago
git bisect start
# status: waiting for both good and bad commits
# good: [f667683bb6c5bf78a8cef5defa52f887c526fc71] feat: integrate SciMLJacobianOperators into NonlinearSolve
git bisect good f667683bb6c5bf78a8cef5defa52f887c526fc71
# status: waiting for bad commit, 1 good commit known
# bad: [d56476c69c749a4b575b47d58c03cf366ab7134c] fix: use DI preparation result when initializing Jacobian (#472)
git bisect bad d56476c69c749a4b575b47d58c03cf366ab7134c
# skip: [fe1bd3dd6ec4923afad60453600316038b3cd9f1] chore: bump minimum versions
git bisect skip fe1bd3dd6ec4923afad60453600316038b3cd9f1
# skip: [59b7d022bc282a04f858b973954b26fc52792771] fix: update bruss testing
git bisect skip 59b7d022bc282a04f858b973954b26fc52792771
# skip: [92c59dbb2d360a4a2ed968668395473be57fc81e] chore: run formatter
git bisect skip 92c59dbb2d360a4a2ed968668395473be57fc81e
# skip: [058b87d5af338ad7e2049f2f27c95b476798d32b] feat: remove Setfield dep
git bisect skip 058b87d5af338ad7e2049f2f27c95b476798d32b
# skip: [3f59f4949829b74d78b9eea8de834dd4f9239152] fix: autodiff cannot be nothing
git bisect skip 3f59f4949829b74d78b9eea8de834dd4f9239152
# skip: [40eb6c867fbed99368e7005ae1ba4246083269eb] chore: mark master as DEV (DON'T TAG)
git bisect skip 40eb6c867fbed99368e7005ae1ba4246083269eb
# skip: [adac560c3f0f737b5c854d61642a182078bb50a1] fix: __value_derivative removal from line searches
git bisect skip adac560c3f0f737b5c854d61642a182078bb50a1
# skip: [b2bb25252e3526b753b3852dfbfff9cc386e3f7a] fix: missing qualifier
git bisect skip b2bb25252e3526b753b3852dfbfff9cc386e3f7a
# skip: [3b12e73e945941ac8d0d4eebe038f1c8ae15d5cf] fix: unwrap sparse AD for derivative call
git bisect skip 3b12e73e945941ac8d0d4eebe038f1c8ae15d5cf
# skip: [43df907989a67b707e433753b09ef3bcd3d08a4c] test: structured jacobians
git bisect skip 43df907989a67b707e433753b09ef3bcd3d08a4c
# skip: [0f2a62f00334d680096e55c9a8ecdc8012819aef] fix: minor test fixes
git bisect skip 0f2a62f00334d680096e55c9a8ecdc8012819aef
# good: [bdacb9ca9327b9d5d2fe6a8f0c08f44ac206f8ec] docs: disable ILU for now
git bisect good bdacb9ca9327b9d5d2fe6a8f0c08f44ac206f8ec
# bad: [b6221b7d61027339ec5fc50bb1f3f90398ce607e] refactor: reorder imports
git bisect bad b6221b7d61027339ec5fc50bb1f3f90398ce607e
# skip: [ba45097e841f9875a7e53d94a1d86f1ffcfd1425] fix: sparse jacobian
git bisect skip ba45097e841f9875a7e53d94a1d86f1ffcfd1425
# good: [6201fcc2f49e2ec2019eb61491f30bc9664e0dc5] fix: DI now works with ReverseDiff
git bisect good 6201fcc2f49e2ec2019eb61491f30bc9664e0dc5
# bad: [33ac4501df94150288f0e7c15a5aff9ab6c20792] fix: update minimum compats
git bisect bad 33ac4501df94150288f0e7c15a5aff9ab6c20792
# good: [05aa3db8351e2da44729e99795bd2b30f7c7a601] fix: remove stale dep
git bisect good 05aa3db8351e2da44729e99795bd2b30f7c7a601
# bad: [a8e60d5f379d89d73d5d4509b3512768e55365a4] feat: use DI for dense jacobians
git bisect bad a8e60d5f379d89d73d5d4509b3512768e55365a4
# first bad commit: [a8e60d5f379d89d73d5d4509b3512768e55365a4] feat: use DI for dense jacobians

found it [a8e60d5f379d89d73d5d4509b3512768e55365a4] feat: use DI for dense jacobians

gdalle commented 1 month ago

Trying to figure out why the prototype isn't banded. We get a SubArray because of the views here https://github.com/SciML/BoundaryValueDiffEq.jl/blob/0a2b4ac73413fcf3ab1514607a4e816fe4d9ccc3/src/solve/mirk.jl#L386-L393 I guess the question to @avik-pal is where does the creation of cache.J happen for this example?

ChrisRackauckas commented 1 month ago

found it [https://github.com/SciML/NonlinearSolve.jl/commit/a8e60d5f379d89d73d5d4509b3512768e55365a4] feat: use DI for dense jacobians

That makes a lot of sense. It must've changed the jac or sparsity type here. We can:

  1. Figure out why the type of Jac changed
  2. Expand the banded matrix iteration dispatch, it doesn't actually need to be writing into a banded matrix, a Matrix that is actually banded is fine.
  3. Isolate, fix and test color iteration for dense non-square matrices. This isn't strictly required to fix our tests, but it's a case that has been identified from this that clearly needs work.
gdalle commented 1 month ago

For 1 I have a hunch: if you're using DI with a dense backend to compute your prototype, you're always gonna get something that is the result of stack, hence a Matrix. If you wanted to preserve the sparsity pattern (in this case banded) you would need to compute the prototype with an AutoSparse. @avik-pal does that resonate with you?

avik-pal commented 1 month ago
infil> f.jac_prototype
8×8 AlmostBandedMatrix{Float64} with bandwidths (5, 4) and fill rank 4 with data 8×8 BandedMatrix{Float64} with bandwidths (5, 4) and fill 4×8 Matrix{Float64}:
 1.0  1.0  1.0  1.0  1.0  0.0  0.0  0.0
 1.0  1.0  1.0  1.0  1.0  1.0  0.0  0.0
 1.0  1.0  1.0  1.0  1.0  1.0  1.0  0.0
 1.0  1.0  1.0  1.0  1.0  1.0  1.0  1.0
 1.0  1.0  1.0  1.0  1.0  1.0  1.0  1.0
 1.0  1.0  1.0  1.0  1.0  1.0  1.0  1.0
  ⋅   1.0  1.0  1.0  1.0  1.0  1.0  1.0
  ⋅    ⋅   1.0  1.0  1.0  1.0  1.0  1.0

infil> similar(f.jac_prototype)
8×8 Matrix{Float64}:
 2.122e-314    1.2732e-313  1.4854e-313  1.9098e-313  2.3342e-313   2.54639e-313  2.75859e-313  7.0e-323
 4.24399e-314  3.0e-323     1.4854e-313  1.9098e-313  2.3342e-313   2.54639e-313  2.75859e-313  3.18299e-313
 8.48798e-314  1.4854e-313  1.6976e-313  2.122e-313   2.54639e-313  2.54639e-313  2.75859e-313  3.18299e-313
 8.48798e-314  1.4854e-313  1.6976e-313  2.122e-313   2.54639e-313  2.54639e-313  2.97079e-313  3.18299e-313
 8.48798e-314  1.4854e-313  1.6976e-313  2.122e-313   2.54639e-313  2.75859e-313  2.97079e-313  3.18299e-313
 8.48798e-314  1.4854e-313  1.6976e-313  2.122e-313   2.54639e-313  2.75859e-313  2.97079e-313  3.18299e-313
 1.061e-313    1.4854e-313  1.6976e-313  2.122e-313   2.54639e-313  2.75859e-313  2.97079e-313  7.4e-323
 1.2732e-313   1.4854e-313  1.9098e-313  2.3342e-313  2.54639e-313  2.75859e-313  2.97079e-313  7.4e-323
gdalle commented 1 month ago

Okay then that's not on me :rofl: good find Avik!

avik-pal commented 1 month ago

I have always hated Julia's fallback semantics for similar :angry:

ChrisRackauckas commented 1 month ago

Oh I see what's going on here.

https://github.com/JuliaDiff/SparseDiffTools.jl/blob/a984510bf1611985b833f227a08b536d72627419/src/differentiation/compute_jacobian_ad.jl#L348-L353

https://github.com/JuliaDiff/FiniteDiff.jl/blob/b574440c5d54708aa23a90c48f4a000dd1c6d915/ext/FiniteDiffBandedMatricesExt.jl#L11

Is an internal setup to say "are you going to use structural_nz anyways? Don't make the arrays if you don't have to". It checks the sparsity pattern and says, "banded matrix? We have one of those at home https://github.com/JuliaDiff/FiniteDiff.jl/blob/b574440c5d54708aa23a90c48f4a000dd1c6d915/ext/FiniteDiffBandedMatricesExt.jl#L13-L27 which effectively just uses a trivial rows/cols, so don't spend the time building that array". Then Jac isn't a banded matrix so it doesn't use that dispatch 😅.

gdalle commented 1 month ago

So what's the plan of action here? If we define a generic findstructralnz for AbstractMatrix then we no longer need these two wrong lines, right?

https://github.com/JuliaDiff/SparseDiffTools.jl/blob/a984510bf1611985b833f227a08b536d72627419/src/differentiation/compute_jacobian_ad.jl#L351-L352

But it would maybe send a signal that dense matrices also have a structure

ChrisRackauckas commented 1 month ago

Those should probably just be nothing since they are free. Basically, if you need the structuralnz, compute it once. If you don't need it, then they are nothing. If you need something like that inside of the specialization because it's "free", so banded matrix and Matrix can just

        rows_index = 1:size(J, 1)
        cols_index = 1:size(J, 2)

on the inside which is non-allocating so it's fine to do it there. nothing is just a clearer sentinel.

ChrisRackauckas commented 1 month ago

Or we can just make sure it returns a correct non-allocating iterator. What would be bad is to allocate some big arrays (O(nonzeros)) that we don't end up using because the iteration scheme uses an analytical solution though. That would be a bit wasteful because for a matrix it's as big as the matrix itself.

gdalle commented 1 month ago

Here's an easy solution:

julia> using Base.Iterators

julia> using BenchmarkTools

julia> function findstructralnz(x::AbstractMatrix)
           cartind = vec(CartesianIndices(x))
           return Iterators.map(first ∘ Tuple, cartind), Iterators.map(last ∘ Tuple, cartind)
       end
findstructralnz (generic function with 1 method)

julia> collect.(findstructralnz(rand(2, 3)))
([1, 2, 1, 2, 1, 2], [1, 1, 2, 2, 3, 3])

julia> @btime findstructralnz($(rand(2, 3)));
  28.162 ns (0 allocations: 0 bytes)
ChrisRackauckas commented 1 month ago

👍that looks pretty decent.

gdalle commented 1 month ago

If that's okay I would like to let you insert it in ArrayInterface and modify SparseDiffTools though, I don't know these code bases and I don't want to accidentally break something

gdalle commented 1 month ago

Have a nice weekend!

ChrisRackauckas commented 1 month ago

One day I'll remember what weekends are 😅

gdalle commented 1 month ago

I think the least disruptive fix is this one:

gdalle commented 1 month ago

And to avoid the diagnosis issues we've ran into here:

ChrisRackauckas commented 1 month ago

And https://github.com/JuliaDiff/FiniteDiff.jl/pull/194 is the other part, and https://github.com/SciML/FastAlmostBandedMatrices.jl/pull/16 is the better similar. So I think all facets of this are getting handled.

ChrisRackauckas commented 1 month ago

Merged, yay.