Skip to content

NonlinearSolve.jl integration

AIBECS ships with its own steady-state solver, CTKAlg(), which has zero extra dependencies and is the default for every tutorial. For users who want to tap into the broader SciML ecosystem, AIBECS also provides an opt-in integration with NonlinearSolve.jl and LinearSolve.jl. Both are weak dependencies: the integration only activates once you load them yourself.

julia
using AIBECS
using NonlinearSolve
using LinearSolve

We use the toy Primeau_2x2x2 circulation here so the docs build stays fast; any AIBECS circulation works the same way.

julia
grd, T = Primeau_2x2x2.load()
(, sparse([1, 2, 3, 2, 5, 1, 3, 1, 4, 2, 4, 5], [1, 1, 1, 2, 2, 3, 3, 4, 4, 5, 5, 5], [4.50923388447017e-9, -5.881609414526308e-10, -3.921072943017539e-9, 9.802682357543846e-10, -5.6015327757393415e-11, -3.921072943017539e-9, 3.921072943017539e-9, -5.881609414526308e-10, 3.360919665443605e-11, -3.9210729430175384e-10, -3.360919665443605e-11, 5.6015327757393415e-11], 5, 5))

Build the standard ideal-age problem (see the Ideal age tutorial for the full derivation).

julia
z = depthvec(grd)
function G(x, p)
    @unpack τ, z₀ = p
    return @. 1 - x / τ * (z  z₀)
end

struct IdealAgeParameters{U} <: AbstractParameters{U}
    τ::U
    z₀::U
end

p = IdealAgeParameters(1.0, z[1])
F = AIBECSFunction(T, G)
nb = sum(iswet(grd))
x_init = zeros(nb)
prob = SteadyStateProblem(F, x_init, p)
SteadyStateProblem with uType Vector{Float64}. In-place: false
u0: 5-element Vector{Float64}:
 0.0
 0.0
 0.0
 0.0
 0.0

Solving with NonlinearSolve

Wrap the SteadyStateProblem with AIBECS.nonlinearproblem before calling solve. The helper performs the SciML conversion to a NonlinearProblem and attaches the sparse Jacobian buffer type (jac_prototype) that LinearSolve's UMFPACK / KLU factorisations expect.

julia
nlprob = AIBECS.nonlinearproblem(prob)
sol = solve(nlprob, NewtonRaphson(linsolve = UMFPACKFactorization()))
retcode: Success
u: 5-element Vector{Float64}:
 28.999999877858574
  8.00000001235138
  1.0000001097900413
  4.76060173862469e10
  1.7852256524842587e10

Or, equivalently, use the AIBECS-recommended default — same Newton step, same UMFPACK factorisation:

julia
sol = solve(nlprob, AIBECS.recommended_nlalg())
retcode: Success
u: 5-element Vector{Float64}:
 28.999999877858574
  8.00000001235138
  1.0000001097900413
  4.76060173862469e10
  1.7852256524842587e10

Sticking with CTKAlg

CTKAlg() continues to work exactly as before — no extra packages required:

julia
sol = solve(prob, CTKAlg())
retcode: Success
u: 5-element Vector{Float64}:
 28.999999877858574
  8.00000001235138
  1.0000001097900413
  4.76060173862469e10
  1.7852256524842587e10

See the solvers explanation page for guidance on when each solver path is appropriate and which algorithm/linear-solver combinations are currently verified to work.


This page was generated using Literate.jl.