Skip to content
This repository was archived by the owner on May 15, 2025. It is now read-only.

Commit bc57490

Browse files
committed
Add custom jac in SimpleTrustRegion
1 parent 411ad53 commit bc57490

File tree

2 files changed

+10
-5
lines changed

2 files changed

+10
-5
lines changed

src/trustRegion.jl

Lines changed: 5 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -112,7 +112,11 @@ function SciMLBase.__solve(prob::NonlinearProblem,
112112
real(oneunit(eltype(T))) * (eps(real(one(eltype(T)))))^(4 // 5)
113113
rtol = reltol !== nothing ? reltol : eps(real(one(eltype(T))))^(4 // 5)
114114

115-
if alg_autodiff(alg)
115+
116+
if DiffEqBase.has_jac(prob.f)
117+
∇f = prob.f.jac(x, prob.p)
118+
F = f(x)
119+
elseif alg_autodiff(alg)
116120
F, ∇f = value_derivative(f, x)
117121
elseif x isa AbstractArray
118122
F = f(x)

test/basictests.jl

Lines changed: 5 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -491,7 +491,7 @@ for alg in (BATCHED_BROYDEN_SOLVERS..., BATCHED_LBROYDEN_SOLVERS...)
491491
@test abs.(sol.u) sqrt.(p)
492492
end
493493

494-
## Specified Jacobian
494+
## User specified Jacobian
495495

496496
f, u0 = (u, p) -> u .* u .- p, randn(3)
497497

@@ -501,6 +501,7 @@ p = [2.0, 1.0, 5.0];
501501

502502
probN = NonlinearProblem(NonlinearFunction(f, jac = f_jac), u0, p)
503503

504-
sol = solve(probN, SimpleNewtonRaphson())
505-
506-
@test abs.(sol.u) sqrt.(p)
504+
for alg in (SimpleNewtonRaphson(), SimpleTrustRegion())
505+
sol = solve(probN, alg)
506+
@test abs.(sol.u) sqrt.(p)
507+
end

0 commit comments

Comments
 (0)