Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Use PolyAlgorithms for Nonlinear Solves #137

Merged
merged 3 commits into from
Nov 7, 2023
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
2 changes: 1 addition & 1 deletion Project.toml
Original file line number Diff line number Diff line change
@@ -1,6 +1,6 @@
name = "BoundaryValueDiffEq"
uuid = "764a87c0-6b3e-53db-9096-fe964310641d"
version = "5.4.0"
version = "5.4.1"

[deps]
ADTypes = "47edcb42-4c32-4615-8424-f2b9edc5f35b"
Expand Down
26 changes: 12 additions & 14 deletions src/algorithms.jl
Original file line number Diff line number Diff line change
Expand Up @@ -3,14 +3,15 @@ abstract type BoundaryValueDiffEqAlgorithm <: SciMLBase.AbstractBVPAlgorithm end
abstract type AbstractMIRK <: BoundaryValueDiffEqAlgorithm end

"""
Shooting(ode_alg; nlsolve = NewtonRaphson(), jac_alg = BVPJacobianAlgorithm())
Shooting(ode_alg = nothing; nlsolve = nothing, jac_alg = BVPJacobianAlgorithm())

Single shooting method, reduces BVP to an initial value problem and solves the IVP.

## Arguments

- `ode_alg`: ODE algorithm to use for solving the IVP. Any solver which conforms to the
SciML `ODEProblem` interface can be used!
SciML `ODEProblem` interface can be used! (Defaults to `nothing` which will use
poly-algorithm if `DifferentialEquations.jl` is loaded else this must be supplied)

## Keyword Arguments

Expand Down Expand Up @@ -39,7 +40,7 @@ function concretize_jacobian_algorithm(alg::Shooting, prob)
return Shooting(alg.ode_alg, alg.nlsolve, BVPJacobianAlgorithm(diffmode))
end

function Shooting(ode_alg; nlsolve = NewtonRaphson(), jac_alg = nothing)
function Shooting(ode_alg = nothing; nlsolve = nothing, jac_alg = nothing)
jac_alg === nothing && (jac_alg = __propagate_nlsolve_ad_to_jac_alg(nlsolve))
return Shooting(ode_alg, nlsolve, jac_alg)
end
Expand All @@ -61,7 +62,7 @@ function __propagate_nlsolve_ad_to_jac_alg(nlsolve::N) where {N}
end

"""
MultipleShooting(nshoots::Int, ode_alg; nlsolve = NewtonRaphson(),
MultipleShooting(nshoots::Int, ode_alg = nothing; nlsolve = nothing,
grid_coarsening = true, jac_alg = BVPJacobianAlgorithm())

Multiple Shooting method, reduces BVP to an initial value problem and solves the IVP.
Expand All @@ -71,7 +72,8 @@ Significantly more stable than Single Shooting.

- `nshoots`: Number of shooting points.
- `ode_alg`: ODE algorithm to use for solving the IVP. Any solver which conforms to the
SciML `ODEProblem` interface can be used!
SciML `ODEProblem` interface can be used! (Defaults to `nothing` which will use
poly-algorithm if `DifferentialEquations.jl` is loaded else this must be supplied)

## Keyword Arguments

Expand Down Expand Up @@ -121,7 +123,7 @@ function update_nshoots(alg::MultipleShooting, nshoots::Int)
alg.grid_coarsening)
end

function MultipleShooting(nshoots::Int, ode_alg; nlsolve = NewtonRaphson(),
function MultipleShooting(nshoots::Int, ode_alg = nothing; nlsolve = nothing,
grid_coarsening = true, jac_alg = BVPJacobianAlgorithm())
@assert grid_coarsening isa Bool || grid_coarsening isa Function ||
grid_coarsening isa AbstractVector{<:Integer} ||
Expand All @@ -141,7 +143,7 @@ for order in (2, 3, 4, 5, 6)
"""
$($alg)(; nlsolve = NewtonRaphson(), jac_alg = BVPJacobianAlgorithm())

$($order)th order Monotonic Implicit Runge Kutta method, with Newton Raphson nonlinear solver as default.
$($order)th order Monotonic Implicit Runge Kutta method.

## Keyword Arguments

Expand Down Expand Up @@ -173,13 +175,9 @@ for order in (2, 3, 4, 5, 6)
pages={479-497}
}
"""
struct $(alg){N, J <: BVPJacobianAlgorithm} <: AbstractMIRK
nlsolve::N
jac_alg::J
end

function $(alg)(; nlsolve = NewtonRaphson(), jac_alg = BVPJacobianAlgorithm())
return $(alg)(nlsolve, jac_alg)
Base.@kwdef struct $(alg){N, J <: BVPJacobianAlgorithm} <: AbstractMIRK
nlsolve::N = nothing
jac_alg::J = BVPJacobianAlgorithm()
end
end
end
Expand Down
2 changes: 1 addition & 1 deletion test/mirk/nonlinear_least_squares.jl
Original file line number Diff line number Diff line change
Expand Up @@ -2,7 +2,7 @@ using BoundaryValueDiffEq, LinearAlgebra, Test

@testset "Overconstrained BVP" begin
SOLVERS = [mirk(; nlsolve) for mirk in (MIRK4, MIRK5, MIRK6),
nlsolve in (LevenbergMarquardt(), GaussNewton())]
nlsolve in (LevenbergMarquardt(), GaussNewton(), nothing)]

# OOP MP-BVP
f1(u, p, t) = [u[2], -u[1]]
Expand Down
2 changes: 2 additions & 0 deletions test/shooting/nonlinear_least_squares.jl
Original file line number Diff line number Diff line change
Expand Up @@ -2,8 +2,10 @@ using BoundaryValueDiffEq, LinearAlgebra, LinearSolve, OrdinaryDiffEq, Test

@testset "Overconstrained BVP" begin
SOLVERS = [
Shooting(Tsit5()),
Shooting(Tsit5(); nlsolve = LevenbergMarquardt()),
Shooting(Tsit5(); nlsolve = GaussNewton()),
MultipleShooting(10, Tsit5()),
MultipleShooting(10, Tsit5(); nlsolve = LevenbergMarquardt()),
MultipleShooting(10, Tsit5(); nlsolve = GaussNewton())]

Expand Down
Loading