Skip to content

Commit 0902c63

Browse files
committed
Use Jacobian Based Algorithms if user supplies custom Jacobian
1 parent c6c875f commit 0902c63

File tree

2 files changed

+34
-14
lines changed

2 files changed

+34
-14
lines changed

Project.toml

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1,7 +1,7 @@
11
name = "NonlinearSolve"
22
uuid = "8913a72c-1f9b-4ce2-8d82-65094dcecaec"
33
authors = ["SciML"]
4-
version = "2.8.0"
4+
version = "2.9.0"
55

66
[deps]
77
ADTypes = "47edcb42-4c32-4615-8424-f2b9edc5f35b"

src/default.jl

Lines changed: 33 additions & 13 deletions
Original file line numberDiff line numberDiff line change
@@ -166,7 +166,7 @@ end
166166

167167
"""
168168
RobustMultiNewton(; concrete_jac = nothing, linsolve = nothing, precs = DEFAULT_PRECS,
169-
adkwargs...)
169+
adkwargs...)
170170
171171
A polyalgorithm focused on robustness. It uses a mixture of Newton methods with different
172172
globalizing techniques (trust region updates, line searches, etc.) in order to find a
@@ -212,7 +212,7 @@ end
212212

213213
"""
214214
FastShortcutNonlinearPolyalg(; concrete_jac = nothing, linsolve = nothing,
215-
precs = DEFAULT_PRECS, adkwargs...)
215+
precs = DEFAULT_PRECS, must_use_jacobian::Val = Val(false), adkwargs...)
216216
217217
A polyalgorithm focused on balancing speed and robustness. It first tries less robust methods
218218
for more performance and then tries more robust techniques if the faster ones fail.
@@ -238,15 +238,25 @@ for more performance and then tries more robust techniques if the faster ones fa
238238
[LinearSolve.jl documentation](https://docs.sciml.ai/LinearSolve/stable/).
239239
"""
240240
function FastShortcutNonlinearPolyalg(; concrete_jac = nothing, linsolve = nothing,
241-
precs = DEFAULT_PRECS, adkwargs...)
242-
algs = (GeneralKlement(; linsolve, precs),
243-
GeneralBroyden(),
244-
NewtonRaphson(; concrete_jac, linsolve, precs, adkwargs...),
245-
NewtonRaphson(; concrete_jac, linsolve, precs, linesearch = BackTracking(),
246-
adkwargs...),
247-
TrustRegion(; concrete_jac, linsolve, precs, adkwargs...),
248-
TrustRegion(; concrete_jac, linsolve, precs,
249-
radius_update_scheme = RadiusUpdateSchemes.Bastin, adkwargs...))
241+
precs = DEFAULT_PRECS, must_use_jacobian::Val{JAC} = Val(false),
242+
adkwargs...) where {JAC}
243+
if JAC
244+
algs = (NewtonRaphson(; concrete_jac, linsolve, precs, adkwargs...),
245+
NewtonRaphson(; concrete_jac, linsolve, precs, linesearch = BackTracking(),
246+
adkwargs...),
247+
TrustRegion(; concrete_jac, linsolve, precs, adkwargs...),
248+
TrustRegion(; concrete_jac, linsolve, precs,
249+
radius_update_scheme = RadiusUpdateSchemes.Bastin, adkwargs...))
250+
else
251+
algs = (GeneralKlement(; linsolve, precs),
252+
GeneralBroyden(),
253+
NewtonRaphson(; concrete_jac, linsolve, precs, adkwargs...),
254+
NewtonRaphson(; concrete_jac, linsolve, precs, linesearch = BackTracking(),
255+
adkwargs...),
256+
TrustRegion(; concrete_jac, linsolve, precs, adkwargs...),
257+
TrustRegion(; concrete_jac, linsolve, precs,
258+
radius_update_scheme = RadiusUpdateSchemes.Bastin, adkwargs...))
259+
end
250260
return NonlinearSolvePolyAlgorithm(algs, Val(:NLS))
251261
end
252262

@@ -288,12 +298,22 @@ end
288298

289299
## Defaults
290300

301+
## TODO: In the long run we want to use an `Assumptions` API like LinearSolve to specify
302+
## the conditioning of the Jacobian and such
303+
## Defaults to a fast and robust poly algorithm in most cases. If the user went through
304+
## the trouble of specifying a custom jacobian function, we should use algorithms that
305+
## can use that!
306+
291307
function SciMLBase.__init(prob::NonlinearProblem, ::Nothing, args...; kwargs...)
292-
return SciMLBase.__init(prob, FastShortcutNonlinearPolyalg(), args...; kwargs...)
308+
must_use_jacobian = Val(prob.f.jac !== nothing)
309+
return SciMLBase.__init(prob, FastShortcutNonlinearPolyalg(; must_use_jacobian),
310+
args...; kwargs...)
293311
end
294312

295313
function SciMLBase.__solve(prob::NonlinearProblem, ::Nothing, args...; kwargs...)
296-
return SciMLBase.__solve(prob, FastShortcutNonlinearPolyalg(), args...; kwargs...)
314+
must_use_jacobian = Val(prob.f.jac !== nothing)
315+
return SciMLBase.__solve(prob, FastShortcutNonlinearPolyalg(; must_use_jacobian),
316+
args...; kwargs...)
297317
end
298318

299319
function SciMLBase.__init(prob::NonlinearLeastSquaresProblem, ::Nothing, args...; kwargs...)

0 commit comments

Comments
 (0)