Parameter optimisation
In this how-to we fit two parameters of the coupled PO₄–POP model (τ_DIP, τ_POP) to PO₄ observations from WorldOceanAtlasTools.jl on the OCCA circulation. Restricting the fit to two parameters lets us also sweep the cost function in a 2D grid around the optimiser trajectory and visualise the loss landscape directly.
We use the SciML Optimization.jl frontend, with F1Method.jl supplying the cached steady-state gradient and Hessian, and share an UMFPACK factorisation between F1Method's adjoint solve and AIBECS's inner Newton step (warm-start trick) so the sweep stays cheap.
The API surface we touch:
@flattenable @bounds @logscaled @prior— declare which parameters are optimised and supply their priorsp2λ/λ2p— transform a parameter struct to/from an unconstrained flat vectorf_and_∇ₓf— build the volume-weighted mismatch objective and its gradient w.r.t. the state vector xF1Method— gradient and Hessian (Pasquier et al., 2019)Optimization.jl— the SciML optimisation frontend, wrapping many optimizers (we will useNewtonTrustRegion).LinearSolve.jl— the SciML linear-solver frontend (we will use the default sparse UMFPACK), and share the factorisation between F1Method and AIBECS.
Re-build the P-model on OCCA
We reuse the source/sink terms from the P-model tutorial verbatim, just swapping the circulation to OCCA so the resolution stays moderate (~100k wet boxes × 2 tracers; OCIM2 is ~200k wet boxes for comparison).
using AIBECS
using JLD2 # required by `OCCA.load`
grd, T_OCIM = OCCA.load()
T_DIP(p) = T_OCIM
T_POP(p) = transportoperator(grd, z -> w(z, p))
function w(z, p)
@unpack w₀, w′ = p
return @. w₀ + w′ * z
end
z = depthvec(grd)
function U(x, p)
@unpack τ_DIP, k, z₀ = p
return @. x / τ_DIP * x / (x + k) * (z ≤ z₀) * (x ≥ 0)
end
function R(x, p)
@unpack τ_POP = p
return x / τ_POP
end
function G_DIP(DIP, POP, p)
@unpack DIP_geo, τ_geo = p
return @. -$U(DIP, p) + $R(POP, p) + (DIP_geo - DIP) / τ_geo
end
function G_POP(DIP, POP, p)
@unpack τ_geo = p
return @. $U(DIP, p) - $R(POP, p) - POP / τ_geo
endG_POP (generic function with 1 method)Optimisable parameter struct
We layer the full optimisation metadata stack on top of the tutorial's PmodelParameters. The leftmost macro consumes the leftmost column (canonical example at src/Parameters.jl). Only τ_DIP and τ_POP are optimised and marked as flattenable. So we only pass a 2-element vector λ to the optimiser.
import AIBECS: @initial_value, initial_value
import AIBECS: @units, units
import AIBECS: @flattenable, flattenable
import AIBECS: @bounds, bounds
import AIBECS: @logscaled, logscaled
import AIBECS: prior
using Unitful: m, d, s, yr, kyr, Myr, mol, mmol, μmol, μM
using Distributions
using Bijectors # activates the `AIBECSBijectorsExt` extension (p2λ, λ2p)
using DataFrames # activates the `AIBECSDataFramesExt` extension (AIBECS.table)
const ∞ = Inf
@initial_value @units @flattenable @bounds @logscaled struct OptParams{U} <: AbstractParameters{U}
w₀::U | 5.0 | m/d | false | (0, ∞) | true
w′::U | 0.035 | m/d/m | false | (0, ∞) | true
τ_DIP::U | 100.0 | d | true | (0, ∞) | true
k::U | 10.0 | μmol/m^3 | false | (0, ∞) | true
z₀::U | 80.0 | m | false | (0, ∞) | false
τ_POP::U | 30.0 | d | true | (0, ∞) | true
τ_geo::U | 1.0 | kyr | false | (0, ∞) | true
DIP_geo::U | 2.12 | mmol/m^3 | false | (0, ∞) | false
endinitial_value (generic function with 17 methods)We apply a penalty to the parameters. For this we derive priors from (lb, ub) = bounds(T, s) and the initial value, matching the convention in test/parameters.jl and GNOM: positive half-line gets a LogNormal, the real line a wide Normal, and finite-interval parameters get a LogitNormal stretched onto (lb, ub). Here we will only use the LogNormal priors but this gives you an idea of how to set up different priors if you want to.
Note
Under the hood AIBECS will perform a change of variables based on Turing.jl's Bijectors.jl extension to transform the constrained parameter space (0, ∞) into an unconstrained one, which is what the optimiser operates in. This is why we can use LogNormal and LogitNormal priors even though the optimiser only sees unconstrained vectors.
function prior(::Type{T}, s::Symbol) where {T <: OptParams}
flattenable(T, s) || return nothing
lb, ub = bounds(T, s)
if (lb, ub) == (0, ∞)
return LogNormal(log(initial_value(T, s)), 1.0)
elseif (lb, ub) == (-∞, ∞)
return Normal(initial_value(T, s), 10.0)
else
m = initial_value(T, s)
f = (m - lb) / (ub - lb)
return LocationScale(lb, ub - lb, LogitNormal(log(f / (1 - f)), 1.0))
end
end
p₀ = OptParams()Main.OptParams{Float64}
w₀ = 5.0 (m d⁻¹)
w′ = 0.035 (d⁻¹)
τ_DIP = 100.0 (d)
k = 10.0 (μmol m⁻³)
z₀ = 80.0 (m)
τ_POP = 30.0 (d)
τ_geo = 1.0 (kyr)
DIP_geo = 2.12 (mmol m⁻³)Two fields (τ_DIP, τ_POP) carry flattenable = true; the rest are held fixed at their initial value. The metadata table:
AIBECS.table(p₀)| Row | Symbol | Value | Initial value | Unit | Bounds | Logscaled | Optimizable |
|---|---|---|---|---|---|---|---|
| Symbol | Float64 | Float64 | FreeUnit… | Tuple… | Bool | Bool | |
| 1 | w₀ | 5.0 | 5.0 | m d⁻¹ | (0, Inf) | true | false |
| 2 | w′ | 0.035 | 0.035 | d⁻¹ | (0, Inf) | true | false |
| 3 | τ_DIP | 100.0 | 100.0 | d | (0, Inf) | true | true |
| 4 | k | 10.0 | 10.0 | μmol m⁻³ | (0, Inf) | true | false |
| 5 | z₀ | 80.0 | 80.0 | m | (0, Inf) | false | false |
| 6 | τ_POP | 30.0 | 30.0 | d | (0, Inf) | true | true |
| 7 | τ_geo | 1.0 | 1.0 | kyr | (0, Inf) | true | false |
| 8 | DIP_geo | 2.12 | 2.12 | mmol m⁻³ | (0, Inf) | false | false |
State function in parameter-type form
Using the AIBECSFunction(Ts, Gs, nb, ::Type{P}) form lets the SciML solvers feed in a flat parameter vector λ; AIBECS will call λ2p(OptParams, λ) internally to reconstruct the struct.
nb = sum(iswet(grd))
F = AIBECSFunction((T_DIP, T_POP), (G_DIP, G_POP), nb, OptParams)
x0 = p₀.DIP_geo * ones(2nb)169322-element Vector{Float64}:
2.12
2.12
2.12
2.12
2.12
2.12
2.12
2.12
2.12
2.12
⋮
2.12
2.12
2.12
2.12
2.12
2.12
2.12
2.12
2.12Load WOA PO₄ observations
WorldOceanAtlasTools.fit_to_grid returns the mean and variance of WOA PO₄ regridded onto our grd as 3D fields. There is a subtle unit gotcha: the WOA 2018 PO₄ file's units attribute is "micromoles_per_kilogram", and WorldOceanAtlasTools.convert_to_SI_unit! does χ_3D .*= ustrip(upreferred(1.0u"μmol/kg")), which gives 1.0e-6 mol/kg — kg is already SI, so upreferred does not convert to mol/m³. The returned array is therefore in mol/kg (mean ≈ 1.6e-6), a factor of ~1000 smaller than the AIBECS model state. We close the gap with a nominal seawater density before anything downstream consumes these arrays.
using WorldOceanAtlasTools
μDIPobs3D, σ²DIPobs3D = WorldOceanAtlasTools.fit_to_grid(grd, "PO₄")
iwet = iswet(grd)
ρ_factor = ustrip(upreferred(1.035u"kg/L")) # = 1035.0 kg/m³
μDIPobs = μDIPobs3D[iwet] .* ρ_factor # mol/kg → mol/m³
σ²DIPobs = σ²DIPobs3D[iwet] .* ρ_factor^2 # (mol/kg)² → (mol/m³)²
v = volumevec(grd)84661-element Vector{Float64}:
4.1259441424340015e11
8.858714556466658e11
9.658719592542885e11
1.0446956966607107e12
1.122246633283809e12
1.1984302852530718e12
1.2731538345237942e12
1.3463262419614636e12
1.417858358258977e12
1.487663032551257e12
⋮
6.642031628003376e13
5.1445109315594125e13
4.172486815218424e13
1.3409911172962852e13
1.824228219687533e13
2.7565991089214062e13
3.717262516144412e13
4.7318477973128164e13
4.2196399175523055e13Build the objective (f, ∇ₓf)
We weight the DIP misfit at 1, ignore POP (no direct observations), and add a small parameter-prior penalty.
ωs = (1.0, 0.0)
ωp = 1.0e-4
μx = (μDIPobs, missing)
σ²x = (σ²DIPobs, missing)
f, ∇ₓf = f_and_∇ₓf(ωs, μx, σ²x, v, ωp, OptParams)(AIBECS.var"#f#generate_f##1"{Main.OptParams, Tuple{Float64, Float64}, Tuple{Vector{Float64}, Missing}, Tuple{Vector{Float64}, Missing}, Vector{Float64}, Float64, AIBECS.var"#tracers#generate_f##0"{Int64, Int64}}((1.0, 0.0), ([0.001993383847718889, 0.0019914651177959006, 0.0019652102196352048, 0.0019188907593353226, 0.0018741692451726304, 0.0018593180609968575, 0.0018762683281708847, 0.0018866130577569658, 0.0018618487300791524, 0.0017936109269207174 … 0.001597349353171885, 0.0016139239530996844, 0.0016146148645877836, 0.0016014571666717527, 0.0015958896704982308, 0.0015948419428572936, 0.0015868619350592294, 0.0015825331862406298, 0.0015800279642719972, 0.001589979759638126], missing), ([4.7499243508165185e-9, 2.3806476583116648e-9, 1.644850324443933e-9, 2.474518014121221e-9, 2.8420192007964507e-9, 1.9923744404030257e-9, 8.150346589335414e-10, 2.974775216581886e-10, 4.021567996605766e-10, 8.840013268111294e-10 … 4.92697369868074e-10, 2.974775216581886e-10, 2.974775216581886e-10, 2.3496760609714243e-9, 2.446816223482753e-9, 2.974775216581886e-10, 8.053041266157682e-10, 8.39703761723315e-10, 1.2798255287465522e-9, 1.7406554311458743e-9], missing), [4.1259441424340015e11, 8.858714556466658e11, 9.658719592542885e11, 1.0446956966607107e12, 1.122246633283809e12, 1.1984302852530718e12, 1.2731538345237942e12, 1.3463262419614636e12, 1.417858358258977e12, 1.487663032551257e12 … 5.551473635870713e13, 6.642031628003376e13, 5.1445109315594125e13, 4.172486815218424e13, 1.3409911172962852e13, 1.824228219687533e13, 2.7565991089214062e13, 3.717262516144412e13, 4.7318477973128164e13, 4.2196399175523055e13], 0.0001, AIBECS.var"#tracers#generate_f##0"{Int64, Int64}(84661, 2)), AIBECS.var"#∇ₓf#generate_∇ₓf##1"{Tuple{Float64, Float64}, Tuple{Vector{Float64}, Missing}, Tuple{Vector{Float64}, Missing}, Vector{Float64}, AIBECS.var"#tracers#generate_∇ₓf##0"{Int64, Int64}}((1.0, 0.0), ([0.001993383847718889, 0.0019914651177959006, 0.0019652102196352048, 0.0019188907593353226, 0.0018741692451726304, 0.0018593180609968575, 0.0018762683281708847, 0.0018866130577569658, 0.0018618487300791524, 0.0017936109269207174 … 0.001597349353171885, 0.0016139239530996844, 0.0016146148645877836, 0.0016014571666717527, 0.0015958896704982308, 0.0015948419428572936, 0.0015868619350592294, 0.0015825331862406298, 0.0015800279642719972, 0.001589979759638126], missing), ([4.7499243508165185e-9, 2.3806476583116648e-9, 1.644850324443933e-9, 2.474518014121221e-9, 2.8420192007964507e-9, 1.9923744404030257e-9, 8.150346589335414e-10, 2.974775216581886e-10, 4.021567996605766e-10, 8.840013268111294e-10 … 4.92697369868074e-10, 2.974775216581886e-10, 2.974775216581886e-10, 2.3496760609714243e-9, 2.446816223482753e-9, 2.974775216581886e-10, 8.053041266157682e-10, 8.39703761723315e-10, 1.2798255287465522e-9, 1.7406554311458743e-9], missing), [4.1259441424340015e11, 8.858714556466658e11, 9.658719592542885e11, 1.0446956966607107e12, 1.122246633283809e12, 1.1984302852530718e12, 1.2731538345237942e12, 1.3463262419614636e12, 1.417858358258977e12, 1.487663032551257e12 … 5.551473635870713e13, 6.642031628003376e13, 5.1445109315594125e13, 4.172486815218424e13, 1.3409911172962852e13, 1.824228219687533e13, 2.7565991089214062e13, 3.717262516144412e13, 4.7318477973128164e13, 4.2196399175523055e13], AIBECS.var"#tracers#generate_∇ₓf##0"{Int64, Int64}(84661, 2), Core.Box(AIBECS.var"#∇ₓf#generate_∇ₓf##1"{Tuple{Float64, Float64}, Tuple{Vector{Float64}, Missing}, Tuple{Vector{Float64}, Missing}, Vector{Float64}, AIBECS.var"#tracers#generate_∇ₓf##0"{Int64, Int64}}(#= circular reference @-2 =#))))F-1 memory cache + warm-started linear solve
F1Method.initialize_mem solves the steady state once and caches the factorised Jacobian; later objective / gradient / Hessian calls reuse it whenever the parameter vector has not moved. We route F1Method's linear algebra through LinearSolve.jl and share its UMFPACK factors with a sibling cache that AIBECS's inner Newton step uses (AIBECScache.cacheval = mem.linear_cache.cacheval) — this makes the cost-surface sweep below ~10× faster. CTKAlg's Shamanskii + Armijo logic refreshes the factors automatically when they go stale.
using F1Method, Optimization, OptimizationOptimJL, ADTypes
using LinearSolve
solver_kwargs = (; maxItNewton = 15)
λ₀ = p2λ(p₀)
mem = F1Method.initialize_mem(F, ∇ₓf, x0, λ₀, CTKAlg();
ad = AutoForwardDiff(),
linsolve = UMFPACKFactorization(),
solver_kwargs...)
AIBECScache = init(
LinearProblem(F.jac(x0, λ2p(OptParams, λ₀)), similar(x0)),
UMFPACKFactorization(),
)
AIBECScache.cacheval = mem.linear_cache.cacheval
AIBECScache.isfresh = false
solver_kwargs = (; solver_kwargs..., linear_cache = AIBECScache)
optfn = F1Method.optimization_function(f, F, ∇ₓf, mem, CTKAlg(); solver_kwargs...)SciMLBase.OptimizationFunction{true, SciMLBase.NoAD, F1MethodOptimizationExt.var"#5#6"{Base.Pairs{Symbol, Any, Nothing, @NamedTuple{maxItNewton::Int64, linear_cache::LinearSolve.LinearCache{SparseMatrixCSC{Float64, Int64}, Vector{Float64}, Vector{Float64}, SciMLBase.NullParameters, LinearSolve.UMFPACKFactorization, SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}, SciMLOperators.IdentityOperator, SciMLOperators.IdentityOperator, Float64, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}, Bool, LinearSolve.LinearSolveAdjoint{Missing}}}}, AIBECS.var"#f#generate_f##1"{Main.OptParams, Tuple{Float64, Float64}, Tuple{Vector{Float64}, Missing}, Tuple{Vector{Float64}, Missing}, Vector{Float64}, Float64, AIBECS.var"#tracers#generate_f##0"{Int64, Int64}}, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#AIBECSFunction##1"{Main.OptParams, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#AIBECSFunction##0"{Main.OptParams, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}, F1Method.F1Cache{LinearSolve.LinearCache{SparseMatrixCSC{Float64, Int64}, Matrix{Float64}, Matrix{Float64}, SciMLBase.NullParameters, LinearSolve.UMFPACKFactorization, SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}, SciMLOperators.IdentityOperator, SciMLOperators.IdentityOperator, Float64, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}, Bool, LinearSolve.LinearSolveAdjoint{Missing}}, Vector{Float64}, Matrix{Float64}, LinearAlgebra.Transpose{Float64, Vector{Float64}}, Vector{Float64}, ADTypes.AutoForwardDiff{nothing, Nothing}}, CTKAlg{LinearSolve.UMFPACKFactorization}}, F1MethodOptimizationExt.var"#7#8"{Base.Pairs{Symbol, Any, Nothing, @NamedTuple{maxItNewton::Int64, linear_cache::LinearSolve.LinearCache{SparseMatrixCSC{Float64, Int64}, Vector{Float64}, Vector{Float64}, SciMLBase.NullParameters, LinearSolve.UMFPACKFactorization, SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}, SciMLOperators.IdentityOperator, SciMLOperators.IdentityOperator, Float64, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}, Bool, LinearSolve.LinearSolveAdjoint{Missing}}}}, AIBECS.var"#f#generate_f##1"{Main.OptParams, Tuple{Float64, Float64}, Tuple{Vector{Float64}, Missing}, Tuple{Vector{Float64}, Missing}, Vector{Float64}, Float64, AIBECS.var"#tracers#generate_f##0"{Int64, Int64}}, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#AIBECSFunction##1"{Main.OptParams, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#AIBECSFunction##0"{Main.OptParams, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}, AIBECS.var"#∇ₓf#generate_∇ₓf##1"{Tuple{Float64, Float64}, Tuple{Vector{Float64}, Missing}, Tuple{Vector{Float64}, Missing}, Vector{Float64}, AIBECS.var"#tracers#generate_∇ₓf##0"{Int64, Int64}}, F1Method.F1Cache{LinearSolve.LinearCache{SparseMatrixCSC{Float64, Int64}, Matrix{Float64}, Matrix{Float64}, SciMLBase.NullParameters, LinearSolve.UMFPACKFactorization, SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}, SciMLOperators.IdentityOperator, SciMLOperators.IdentityOperator, Float64, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}, Bool, LinearSolve.LinearSolveAdjoint{Missing}}, Vector{Float64}, Matrix{Float64}, LinearAlgebra.Transpose{Float64, Vector{Float64}}, Vector{Float64}, ADTypes.AutoForwardDiff{nothing, Nothing}}, CTKAlg{LinearSolve.UMFPACKFactorization}}, Nothing, F1MethodOptimizationExt.var"#9#10"{Base.Pairs{Symbol, Any, Nothing, @NamedTuple{maxItNewton::Int64, linear_cache::LinearSolve.LinearCache{SparseMatrixCSC{Float64, Int64}, Vector{Float64}, Vector{Float64}, SciMLBase.NullParameters, LinearSolve.UMFPACKFactorization, SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}, SciMLOperators.IdentityOperator, SciMLOperators.IdentityOperator, Float64, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}, Bool, LinearSolve.LinearSolveAdjoint{Missing}}}}, AIBECS.var"#f#generate_f##1"{Main.OptParams, Tuple{Float64, Float64}, Tuple{Vector{Float64}, Missing}, Tuple{Vector{Float64}, Missing}, Vector{Float64}, Float64, AIBECS.var"#tracers#generate_f##0"{Int64, Int64}}, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#AIBECSFunction##1"{Main.OptParams, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#AIBECSFunction##0"{Main.OptParams, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}, AIBECS.var"#∇ₓf#generate_∇ₓf##1"{Tuple{Float64, Float64}, Tuple{Vector{Float64}, Missing}, Tuple{Vector{Float64}, Missing}, Vector{Float64}, AIBECS.var"#tracers#generate_∇ₓf##0"{Int64, Int64}}, F1Method.F1Cache{LinearSolve.LinearCache{SparseMatrixCSC{Float64, Int64}, Matrix{Float64}, Matrix{Float64}, SciMLBase.NullParameters, LinearSolve.UMFPACKFactorization, SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}, SciMLOperators.IdentityOperator, SciMLOperators.IdentityOperator, Float64, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}, Bool, LinearSolve.LinearSolveAdjoint{Missing}}, Vector{Float64}, Matrix{Float64}, LinearAlgebra.Transpose{Float64, Vector{Float64}}, Vector{Float64}, ADTypes.AutoForwardDiff{nothing, Nothing}}, CTKAlg{LinearSolve.UMFPACKFactorization}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED_NO_TIME), Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing}(F1MethodOptimizationExt.var"#5#6"{Base.Pairs{Symbol, Any, Nothing, @NamedTuple{maxItNewton::Int64, linear_cache::LinearSolve.LinearCache{SparseMatrixCSC{Float64, Int64}, Vector{Float64}, Vector{Float64}, SciMLBase.NullParameters, LinearSolve.UMFPACKFactorization, SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}, SciMLOperators.IdentityOperator, SciMLOperators.IdentityOperator, Float64, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}, Bool, LinearSolve.LinearSolveAdjoint{Missing}}}}, AIBECS.var"#f#generate_f##1"{Main.OptParams, Tuple{Float64, Float64}, Tuple{Vector{Float64}, Missing}, Tuple{Vector{Float64}, Missing}, Vector{Float64}, Float64, AIBECS.var"#tracers#generate_f##0"{Int64, Int64}}, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#AIBECSFunction##1"{Main.OptParams, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#AIBECSFunction##0"{Main.OptParams, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}, F1Method.F1Cache{LinearSolve.LinearCache{SparseMatrixCSC{Float64, Int64}, Matrix{Float64}, Matrix{Float64}, SciMLBase.NullParameters, LinearSolve.UMFPACKFactorization, SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}, SciMLOperators.IdentityOperator, SciMLOperators.IdentityOperator, Float64, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}, Bool, LinearSolve.LinearSolveAdjoint{Missing}}, Vector{Float64}, Matrix{Float64}, LinearAlgebra.Transpose{Float64, Vector{Float64}}, Vector{Float64}, ADTypes.AutoForwardDiff{nothing, Nothing}}, CTKAlg{LinearSolve.UMFPACKFactorization}}(Base.Pairs{Symbol, Any, Nothing, @NamedTuple{maxItNewton::Int64, linear_cache::LinearSolve.LinearCache{SparseMatrixCSC{Float64, Int64}, Vector{Float64}, Vector{Float64}, SciMLBase.NullParameters, LinearSolve.UMFPACKFactorization, SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}, SciMLOperators.IdentityOperator, SciMLOperators.IdentityOperator, Float64, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}, Bool, LinearSolve.LinearSolveAdjoint{Missing}}}}(:maxItNewton => 15, :linear_cache => LinearSolve.LinearCache{SparseMatrixCSC{Float64, Int64}, Vector{Float64}, Vector{Float64}, SciMLBase.NullParameters, LinearSolve.UMFPACKFactorization, SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}, SciMLOperators.IdentityOperator, SciMLOperators.IdentityOperator, Float64, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}, Bool, LinearSolve.LinearSolveAdjoint{Missing}}(sparse([1, 2, 9551, 9606, 84662, 1, 2, 3, 57, 9552 … 84657, 169318, 84658, 169319, 84659, 169320, 84660, 169321, 84661, 169322], [1, 1, 1, 1, 1, 2, 2, 2, 2, 2 … 169318, 169318, 169319, 169319, 169320, 169320, 169321, 169321, 169322, 169322], [-3.440352272546861e-7, -1.8030621712982388e-10, 2.202759122121374e-7, 1.669794081839543e-8, 1.157407407381655e-7, 8.32821319857397e-8, -4.803223973230515e-7, 2.4345873581518986e-8, -2.346272528650749e-8, 3.294980961099172e-7 … 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7], 169322, 169322), [0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 … 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0], [0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 … 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0], SciMLBase.NullParameters(), LinearSolve.UMFPACKFactorization(true, true), SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}(SparseArrays.UMFPACK.Symbolic{Float64, Int64}(Ptr{Nothing}(0x0000000125718b20)), SparseArrays.UMFPACK.Numeric{Float64, Int64}(Ptr{Nothing}(0x000000030121df80)), 169322, 169322, [0, 5, 12, 19, 26, 33, 40, 47, 54, 61 … 818900, 818902, 818904, 818906, 818908, 818910, 818912, 818914, 818916, 818918], [0, 1, 9550, 9605, 84661, 0, 1, 2, 56, 9551 … 84656, 169317, 84657, 169318, 84658, 169319, 84659, 169320, 84660, 169321], [-3.4050287839078227e-7, -1.8030621712982388e-10, 2.202759122121374e-7, 1.669794081839543e-8, 1.1220839187426168e-7, 8.32821319857397e-8, -4.783272470058407e-7, 2.4345873581518986e-8, -2.346272528650749e-8, 3.294980961099172e-7 … 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7], 0, SparseArrays.UMFPACK.UmfpackWS{Int64}([6167, 646, 2414, 2416, 6147, 6157, 6158, 675, 2402, 6137 … 100589, 98900, 100764, 97371, 98568, 102479, 2971, 98602, 3061, 98610], [0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 … -3.634301428314729e-9, 6.08605475490754e-9, 3.126003909338102e-8, 4.993000031340455e-8, 4.535492077596722e-9, 4.902645702569492e-9, 1.897352028675765e-7, 3.2195576264689603e-9, 2.6818076806826084e-7, -1.3166840863933614e-6]), [1.0, 0.2, 0.2, 0.1, 32.0, 0.0, 0.7, 0.0, 1.0, 0.3, 1.0, 1.0, 0.9, 0.0, 10.0, 0.001, 1.0, 0.5, 0.0, 1.0], [0.0, 169322.0, 818918.0, 16.0, 4.0, 8.0, 8.0, 8.0, 0.0, 0.0 … 0.0, 0.0, -1.0, -1.0, 4.1972978e7, 0.01968699999997625, 0.01968699999997625, -1.0, -1.0, -1.0], ReentrantLock()), false, false, IdentityOperator(169322), IdentityOperator(169322), 1.4901161193847656e-8, 1.4901161193847656e-8, 169322, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}(SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.DebugLevel(), SciMLLogging.Silent(), SciMLLogging.InfoLevel(), SciMLLogging.Silent(), SciMLLogging.WarnLevel(), SciMLLogging.WarnLevel(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.WarnLevel(), SciMLLogging.WarnLevel(), SciMLLogging.WarnLevel()), LinearSolve.OperatorAssumptions{Bool}(true, LinearSolve.OperatorCondition.IllConditioned), LinearSolve.LinearSolveAdjoint{Missing}(missing))), AIBECS.var"#f#generate_f##1"{Main.OptParams, Tuple{Float64, Float64}, Tuple{Vector{Float64}, Missing}, Tuple{Vector{Float64}, Missing}, Vector{Float64}, Float64, AIBECS.var"#tracers#generate_f##0"{Int64, Int64}}((1.0, 0.0), ([0.001993383847718889, 0.0019914651177959006, 0.0019652102196352048, 0.0019188907593353226, 0.0018741692451726304, 0.0018593180609968575, 0.0018762683281708847, 0.0018866130577569658, 0.0018618487300791524, 0.0017936109269207174 … 0.001597349353171885, 0.0016139239530996844, 0.0016146148645877836, 0.0016014571666717527, 0.0015958896704982308, 0.0015948419428572936, 0.0015868619350592294, 0.0015825331862406298, 0.0015800279642719972, 0.001589979759638126], missing), ([4.7499243508165185e-9, 2.3806476583116648e-9, 1.644850324443933e-9, 2.474518014121221e-9, 2.8420192007964507e-9, 1.9923744404030257e-9, 8.150346589335414e-10, 2.974775216581886e-10, 4.021567996605766e-10, 8.840013268111294e-10 … 4.92697369868074e-10, 2.974775216581886e-10, 2.974775216581886e-10, 2.3496760609714243e-9, 2.446816223482753e-9, 2.974775216581886e-10, 8.053041266157682e-10, 8.39703761723315e-10, 1.2798255287465522e-9, 1.7406554311458743e-9], missing), [4.1259441424340015e11, 8.858714556466658e11, 9.658719592542885e11, 1.0446956966607107e12, 1.122246633283809e12, 1.1984302852530718e12, 1.2731538345237942e12, 1.3463262419614636e12, 1.417858358258977e12, 1.487663032551257e12 … 5.551473635870713e13, 6.642031628003376e13, 5.1445109315594125e13, 4.172486815218424e13, 1.3409911172962852e13, 1.824228219687533e13, 2.7565991089214062e13, 3.717262516144412e13, 4.7318477973128164e13, 4.2196399175523055e13], 0.0001, AIBECS.var"#tracers#generate_f##0"{Int64, Int64}(84661, 2)), SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#AIBECSFunction##1"{Main.OptParams, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#AIBECSFunction##0"{Main.OptParams, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}(AIBECS.var"#f#AIBECSFunction##1"{Main.OptParams, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}}(SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}(AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}((Main.T_DIP, Main.T_POP), 1:2, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}((Main.G_DIP, Main.G_POP), AIBECS.var"#tracers#66"{Int64, Int64}(84661, 2)), AIBECS.var"#tracer#67"{Int64, Int64}(84661, 2)), LinearAlgebra.UniformScaling{Bool}(true), nothing, nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}(AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}((Main.T_DIP, Main.T_POP), 2, 1:2), AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}((Main.G_DIP, Main.G_POP), 84661, 2)), nothing, nothing, nothing, nothing, nothing, nothing, nothing, nothing, nothing, SciMLBase.DEFAULT_OBSERVED, nothing, nothing, nothing, nothing)), LinearAlgebra.UniformScaling{Bool}(true), nothing, nothing, AIBECS.var"#jac#AIBECSFunction##0"{Main.OptParams, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}}(SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}(AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}((Main.T_DIP, Main.T_POP), 1:2, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}((Main.G_DIP, Main.G_POP), AIBECS.var"#tracers#66"{Int64, Int64}(84661, 2)), AIBECS.var"#tracer#67"{Int64, Int64}(84661, 2)), LinearAlgebra.UniformScaling{Bool}(true), nothing, nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}(AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}((Main.T_DIP, Main.T_POP), 2, 1:2), AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}((Main.G_DIP, Main.G_POP), 84661, 2)), nothing, nothing, nothing, nothing, nothing, nothing, nothing, nothing, nothing, SciMLBase.DEFAULT_OBSERVED, nothing, nothing, nothing, nothing)), nothing, nothing, nothing, nothing, nothing, nothing, nothing, nothing, nothing, SciMLBase.DEFAULT_OBSERVED, nothing, nothing, nothing, nothing), F1Method.F1Cache{LinearSolve.LinearCache{SparseMatrixCSC{Float64, Int64}, Matrix{Float64}, Matrix{Float64}, SciMLBase.NullParameters, LinearSolve.UMFPACKFactorization, SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}, SciMLOperators.IdentityOperator, SciMLOperators.IdentityOperator, Float64, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}, Bool, LinearSolve.LinearSolveAdjoint{Missing}}, Vector{Float64}, Matrix{Float64}, LinearAlgebra.Transpose{Float64, Vector{Float64}}, Vector{Float64}, ADTypes.AutoForwardDiff{nothing, Nothing}}(LinearSolve.LinearCache{SparseMatrixCSC{Float64, Int64}, Matrix{Float64}, Matrix{Float64}, SciMLBase.NullParameters, LinearSolve.UMFPACKFactorization, SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}, SciMLOperators.IdentityOperator, SciMLOperators.IdentityOperator, Float64, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}, Bool, LinearSolve.LinearSolveAdjoint{Missing}}(sparse([1, 2, 9551, 9606, 84662, 1, 2, 3, 57, 9552 … 84657, 169318, 84658, 169319, 84659, 169320, 84660, 169321, 84661, 169322], [1, 1, 1, 1, 1, 2, 2, 2, 2, 2 … 169318, 169318, 169319, 169319, 169320, 169320, 169321, 169321, 169322, 169322], [-3.4050287839078227e-7, -1.8030621712982388e-10, 2.202759122121374e-7, 1.669794081839543e-8, 1.1220839187426168e-7, 8.32821319857397e-8, -4.783272470058407e-7, 2.4345873581518986e-8, -2.346272528650749e-8, 3.294980961099172e-7 … 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7], 169322, 169322), [-4.512560904168677e-12 8.935618833654465e-13; -6.652536669931426e-12 1.317312568668086e-12; … ; -0.0 -2.135436433303535e-13; -0.0 -1.6871098776291334e-13], [5.7489882607834106e-5 -2.308373390511117e-5; 5.646653132302884e-5 -3.1428196353351436e-5; … ; -4.603312782669057e-8 1.1570159339507784e-6; 5.471653491870314e-9 8.932325243468619e-7], SciMLBase.NullParameters(), LinearSolve.UMFPACKFactorization(true, true), SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}(SparseArrays.UMFPACK.Symbolic{Float64, Int64}(Ptr{Nothing}(0x0000000125718b20)), SparseArrays.UMFPACK.Numeric{Float64, Int64}(Ptr{Nothing}(0x000000030121df80)), 169322, 169322, [0, 5, 12, 19, 26, 33, 40, 47, 54, 61 … 818900, 818902, 818904, 818906, 818908, 818910, 818912, 818914, 818916, 818918], [0, 1, 9550, 9605, 84661, 0, 1, 2, 56, 9551 … 84656, 169317, 84657, 169318, 84658, 169319, 84659, 169320, 84660, 169321], [-3.4050287839078227e-7, -1.8030621712982388e-10, 2.202759122121374e-7, 1.669794081839543e-8, 1.1220839187426168e-7, 8.32821319857397e-8, -4.783272470058407e-7, 2.4345873581518986e-8, -2.346272528650749e-8, 3.294980961099172e-7 … 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7], 0, SparseArrays.UMFPACK.UmfpackWS{Int64}([6167, 646, 2414, 2416, 6147, 6157, 6158, 675, 2402, 6137 … 100589, 98900, 100764, 97371, 98568, 102479, 2971, 98602, 3061, 98610], [0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 … -3.634301428314729e-9, 6.08605475490754e-9, 3.126003909338102e-8, 4.993000031340455e-8, 4.535492077596722e-9, 4.902645702569492e-9, 1.897352028675765e-7, 3.2195576264689603e-9, 2.6818076806826084e-7, -1.3166840863933614e-6]), [1.0, 0.2, 0.2, 0.1, 32.0, 0.0, 0.7, 0.0, 1.0, 0.3, 1.0, 1.0, 0.9, 0.0, 10.0, 0.001, 1.0, 0.5, 0.0, 1.0], [0.0, 169322.0, 818918.0, 16.0, 4.0, 8.0, 8.0, 8.0, 0.0, 0.0 … 0.0, 0.0, -1.0, -1.0, 4.1972978e7, 0.01968699999997625, 0.01968699999997625, -1.0, -1.0, -1.0], ReentrantLock()), false, false, IdentityOperator(169322), IdentityOperator(169322), 1.4901161193847656e-8, 1.4901161193847656e-8, 338644, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}(SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.DebugLevel(), SciMLLogging.Silent(), SciMLLogging.InfoLevel(), SciMLLogging.Silent(), SciMLLogging.WarnLevel(), SciMLLogging.WarnLevel(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.WarnLevel(), SciMLLogging.WarnLevel(), SciMLLogging.WarnLevel()), LinearSolve.OperatorAssumptions{Bool}(true, LinearSolve.OperatorCondition.IllConditioned), LinearSolve.LinearSolveAdjoint{Missing}(missing)), [4.7241543264007345e-5, 6.616497742507608e-5, 0.0007607097140667825, 0.0002001818247681397, 0.00015228345495271658, 0.0001377887715644329, 0.00014231181003870498, 0.00013249304462571178, 0.00012386276033894944, 8.471960438304386e-5 … 1.6870928837112562e-7, 1.107674177977734e-7, 1.6619353326168306e-7, 2.6533564688704694e-7, 1.0705570309942015e-6, 9.402540764874025e-7, 6.90519304703392e-7, 6.642473304457541e-7, 5.535051235122763e-7, 4.372988802814714e-7], [5.7489882607834106e-5 -2.308373390511117e-5; 5.646653132302884e-5 -3.1428196353351436e-5; … ; -4.603312782669057e-8 1.1570159339507784e-6; 5.471653491870314e-9 8.932325243468619e-7], transpose([-2.7297252656153418e-5, -0.00011568624430491445, -0.00011421097914680908, -0.00011716805364084193, -0.00010979266846249825, -0.00016721064603409375, -0.00043737172418684654, -0.0012819265870092442, -0.0009894443952221357, -0.0004643803065352418 … 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0]), [4.605170185988092, 3.4011973816621555], [4.605170185988092, 3.4011973816621555], ADTypes.AutoForwardDiff()), CTKAlg{LinearSolve.UMFPACKFactorization}(LinearSolve.UMFPACKFactorization(true, true))), SciMLBase.NoAD(), F1MethodOptimizationExt.var"#7#8"{Base.Pairs{Symbol, Any, Nothing, @NamedTuple{maxItNewton::Int64, linear_cache::LinearSolve.LinearCache{SparseMatrixCSC{Float64, Int64}, Vector{Float64}, Vector{Float64}, SciMLBase.NullParameters, LinearSolve.UMFPACKFactorization, SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}, SciMLOperators.IdentityOperator, SciMLOperators.IdentityOperator, Float64, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}, Bool, LinearSolve.LinearSolveAdjoint{Missing}}}}, AIBECS.var"#f#generate_f##1"{Main.OptParams, Tuple{Float64, Float64}, Tuple{Vector{Float64}, Missing}, Tuple{Vector{Float64}, Missing}, Vector{Float64}, Float64, AIBECS.var"#tracers#generate_f##0"{Int64, Int64}}, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#AIBECSFunction##1"{Main.OptParams, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#AIBECSFunction##0"{Main.OptParams, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}, AIBECS.var"#∇ₓf#generate_∇ₓf##1"{Tuple{Float64, Float64}, Tuple{Vector{Float64}, Missing}, Tuple{Vector{Float64}, Missing}, Vector{Float64}, AIBECS.var"#tracers#generate_∇ₓf##0"{Int64, Int64}}, F1Method.F1Cache{LinearSolve.LinearCache{SparseMatrixCSC{Float64, Int64}, Matrix{Float64}, Matrix{Float64}, SciMLBase.NullParameters, LinearSolve.UMFPACKFactorization, SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}, SciMLOperators.IdentityOperator, SciMLOperators.IdentityOperator, Float64, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}, Bool, LinearSolve.LinearSolveAdjoint{Missing}}, Vector{Float64}, Matrix{Float64}, LinearAlgebra.Transpose{Float64, Vector{Float64}}, Vector{Float64}, ADTypes.AutoForwardDiff{nothing, Nothing}}, CTKAlg{LinearSolve.UMFPACKFactorization}}(Base.Pairs{Symbol, Any, Nothing, @NamedTuple{maxItNewton::Int64, linear_cache::LinearSolve.LinearCache{SparseMatrixCSC{Float64, Int64}, Vector{Float64}, Vector{Float64}, SciMLBase.NullParameters, LinearSolve.UMFPACKFactorization, SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}, SciMLOperators.IdentityOperator, SciMLOperators.IdentityOperator, Float64, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}, Bool, LinearSolve.LinearSolveAdjoint{Missing}}}}(:maxItNewton => 15, :linear_cache => LinearSolve.LinearCache{SparseMatrixCSC{Float64, Int64}, Vector{Float64}, Vector{Float64}, SciMLBase.NullParameters, LinearSolve.UMFPACKFactorization, SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}, SciMLOperators.IdentityOperator, SciMLOperators.IdentityOperator, Float64, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}, Bool, LinearSolve.LinearSolveAdjoint{Missing}}(sparse([1, 2, 9551, 9606, 84662, 1, 2, 3, 57, 9552 … 84657, 169318, 84658, 169319, 84659, 169320, 84660, 169321, 84661, 169322], [1, 1, 1, 1, 1, 2, 2, 2, 2, 2 … 169318, 169318, 169319, 169319, 169320, 169320, 169321, 169321, 169322, 169322], [-3.440352272546861e-7, -1.8030621712982388e-10, 2.202759122121374e-7, 1.669794081839543e-8, 1.157407407381655e-7, 8.32821319857397e-8, -4.803223973230515e-7, 2.4345873581518986e-8, -2.346272528650749e-8, 3.294980961099172e-7 … 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7], 169322, 169322), [0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 … 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0], [0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 … 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0], SciMLBase.NullParameters(), LinearSolve.UMFPACKFactorization(true, true), SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}(SparseArrays.UMFPACK.Symbolic{Float64, Int64}(Ptr{Nothing}(0x0000000125718b20)), SparseArrays.UMFPACK.Numeric{Float64, Int64}(Ptr{Nothing}(0x000000030121df80)), 169322, 169322, [0, 5, 12, 19, 26, 33, 40, 47, 54, 61 … 818900, 818902, 818904, 818906, 818908, 818910, 818912, 818914, 818916, 818918], [0, 1, 9550, 9605, 84661, 0, 1, 2, 56, 9551 … 84656, 169317, 84657, 169318, 84658, 169319, 84659, 169320, 84660, 169321], [-3.4050287839078227e-7, -1.8030621712982388e-10, 2.202759122121374e-7, 1.669794081839543e-8, 1.1220839187426168e-7, 8.32821319857397e-8, -4.783272470058407e-7, 2.4345873581518986e-8, -2.346272528650749e-8, 3.294980961099172e-7 … 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7], 0, SparseArrays.UMFPACK.UmfpackWS{Int64}([6167, 646, 2414, 2416, 6147, 6157, 6158, 675, 2402, 6137 … 100589, 98900, 100764, 97371, 98568, 102479, 2971, 98602, 3061, 98610], [0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 … -3.634301428314729e-9, 6.08605475490754e-9, 3.126003909338102e-8, 4.993000031340455e-8, 4.535492077596722e-9, 4.902645702569492e-9, 1.897352028675765e-7, 3.2195576264689603e-9, 2.6818076806826084e-7, -1.3166840863933614e-6]), [1.0, 0.2, 0.2, 0.1, 32.0, 0.0, 0.7, 0.0, 1.0, 0.3, 1.0, 1.0, 0.9, 0.0, 10.0, 0.001, 1.0, 0.5, 0.0, 1.0], [0.0, 169322.0, 818918.0, 16.0, 4.0, 8.0, 8.0, 8.0, 0.0, 0.0 … 0.0, 0.0, -1.0, -1.0, 4.1972978e7, 0.01968699999997625, 0.01968699999997625, -1.0, -1.0, -1.0], ReentrantLock()), false, false, IdentityOperator(169322), IdentityOperator(169322), 1.4901161193847656e-8, 1.4901161193847656e-8, 169322, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}(SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.DebugLevel(), SciMLLogging.Silent(), SciMLLogging.InfoLevel(), SciMLLogging.Silent(), SciMLLogging.WarnLevel(), SciMLLogging.WarnLevel(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.WarnLevel(), SciMLLogging.WarnLevel(), SciMLLogging.WarnLevel()), LinearSolve.OperatorAssumptions{Bool}(true, LinearSolve.OperatorCondition.IllConditioned), LinearSolve.LinearSolveAdjoint{Missing}(missing))), AIBECS.var"#f#generate_f##1"{Main.OptParams, Tuple{Float64, Float64}, Tuple{Vector{Float64}, Missing}, Tuple{Vector{Float64}, Missing}, Vector{Float64}, Float64, AIBECS.var"#tracers#generate_f##0"{Int64, Int64}}((1.0, 0.0), ([0.001993383847718889, 0.0019914651177959006, 0.0019652102196352048, 0.0019188907593353226, 0.0018741692451726304, 0.0018593180609968575, 0.0018762683281708847, 0.0018866130577569658, 0.0018618487300791524, 0.0017936109269207174 … 0.001597349353171885, 0.0016139239530996844, 0.0016146148645877836, 0.0016014571666717527, 0.0015958896704982308, 0.0015948419428572936, 0.0015868619350592294, 0.0015825331862406298, 0.0015800279642719972, 0.001589979759638126], missing), ([4.7499243508165185e-9, 2.3806476583116648e-9, 1.644850324443933e-9, 2.474518014121221e-9, 2.8420192007964507e-9, 1.9923744404030257e-9, 8.150346589335414e-10, 2.974775216581886e-10, 4.021567996605766e-10, 8.840013268111294e-10 … 4.92697369868074e-10, 2.974775216581886e-10, 2.974775216581886e-10, 2.3496760609714243e-9, 2.446816223482753e-9, 2.974775216581886e-10, 8.053041266157682e-10, 8.39703761723315e-10, 1.2798255287465522e-9, 1.7406554311458743e-9], missing), [4.1259441424340015e11, 8.858714556466658e11, 9.658719592542885e11, 1.0446956966607107e12, 1.122246633283809e12, 1.1984302852530718e12, 1.2731538345237942e12, 1.3463262419614636e12, 1.417858358258977e12, 1.487663032551257e12 … 5.551473635870713e13, 6.642031628003376e13, 5.1445109315594125e13, 4.172486815218424e13, 1.3409911172962852e13, 1.824228219687533e13, 2.7565991089214062e13, 3.717262516144412e13, 4.7318477973128164e13, 4.2196399175523055e13], 0.0001, AIBECS.var"#tracers#generate_f##0"{Int64, Int64}(84661, 2)), SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#AIBECSFunction##1"{Main.OptParams, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#AIBECSFunction##0"{Main.OptParams, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}(AIBECS.var"#f#AIBECSFunction##1"{Main.OptParams, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}}(SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}(AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}((Main.T_DIP, Main.T_POP), 1:2, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}((Main.G_DIP, Main.G_POP), AIBECS.var"#tracers#66"{Int64, Int64}(84661, 2)), AIBECS.var"#tracer#67"{Int64, Int64}(84661, 2)), LinearAlgebra.UniformScaling{Bool}(true), nothing, nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}(AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}((Main.T_DIP, Main.T_POP), 2, 1:2), AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}((Main.G_DIP, Main.G_POP), 84661, 2)), nothing, nothing, nothing, nothing, nothing, nothing, nothing, nothing, nothing, SciMLBase.DEFAULT_OBSERVED, nothing, nothing, nothing, nothing)), LinearAlgebra.UniformScaling{Bool}(true), nothing, nothing, AIBECS.var"#jac#AIBECSFunction##0"{Main.OptParams, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}}(SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}(AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}((Main.T_DIP, Main.T_POP), 1:2, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}((Main.G_DIP, Main.G_POP), AIBECS.var"#tracers#66"{Int64, Int64}(84661, 2)), AIBECS.var"#tracer#67"{Int64, Int64}(84661, 2)), LinearAlgebra.UniformScaling{Bool}(true), nothing, nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}(AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}((Main.T_DIP, Main.T_POP), 2, 1:2), AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}((Main.G_DIP, Main.G_POP), 84661, 2)), nothing, nothing, nothing, nothing, nothing, nothing, nothing, nothing, nothing, SciMLBase.DEFAULT_OBSERVED, nothing, nothing, nothing, nothing)), nothing, nothing, nothing, nothing, nothing, nothing, nothing, nothing, nothing, SciMLBase.DEFAULT_OBSERVED, nothing, nothing, nothing, nothing), AIBECS.var"#∇ₓf#generate_∇ₓf##1"{Tuple{Float64, Float64}, Tuple{Vector{Float64}, Missing}, Tuple{Vector{Float64}, Missing}, Vector{Float64}, AIBECS.var"#tracers#generate_∇ₓf##0"{Int64, Int64}}((1.0, 0.0), ([0.001993383847718889, 0.0019914651177959006, 0.0019652102196352048, 0.0019188907593353226, 0.0018741692451726304, 0.0018593180609968575, 0.0018762683281708847, 0.0018866130577569658, 0.0018618487300791524, 0.0017936109269207174 … 0.001597349353171885, 0.0016139239530996844, 0.0016146148645877836, 0.0016014571666717527, 0.0015958896704982308, 0.0015948419428572936, 0.0015868619350592294, 0.0015825331862406298, 0.0015800279642719972, 0.001589979759638126], missing), ([4.7499243508165185e-9, 2.3806476583116648e-9, 1.644850324443933e-9, 2.474518014121221e-9, 2.8420192007964507e-9, 1.9923744404030257e-9, 8.150346589335414e-10, 2.974775216581886e-10, 4.021567996605766e-10, 8.840013268111294e-10 … 4.92697369868074e-10, 2.974775216581886e-10, 2.974775216581886e-10, 2.3496760609714243e-9, 2.446816223482753e-9, 2.974775216581886e-10, 8.053041266157682e-10, 8.39703761723315e-10, 1.2798255287465522e-9, 1.7406554311458743e-9], missing), [4.1259441424340015e11, 8.858714556466658e11, 9.658719592542885e11, 1.0446956966607107e12, 1.122246633283809e12, 1.1984302852530718e12, 1.2731538345237942e12, 1.3463262419614636e12, 1.417858358258977e12, 1.487663032551257e12 … 5.551473635870713e13, 6.642031628003376e13, 5.1445109315594125e13, 4.172486815218424e13, 1.3409911172962852e13, 1.824228219687533e13, 2.7565991089214062e13, 3.717262516144412e13, 4.7318477973128164e13, 4.2196399175523055e13], AIBECS.var"#tracers#generate_∇ₓf##0"{Int64, Int64}(84661, 2), Core.Box(AIBECS.var"#∇ₓf#generate_∇ₓf##1"{Tuple{Float64, Float64}, Tuple{Vector{Float64}, Missing}, Tuple{Vector{Float64}, Missing}, Vector{Float64}, AIBECS.var"#tracers#generate_∇ₓf##0"{Int64, Int64}}(#= circular reference @-2 =#))), F1Method.F1Cache{LinearSolve.LinearCache{SparseMatrixCSC{Float64, Int64}, Matrix{Float64}, Matrix{Float64}, SciMLBase.NullParameters, LinearSolve.UMFPACKFactorization, SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}, SciMLOperators.IdentityOperator, SciMLOperators.IdentityOperator, Float64, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}, Bool, LinearSolve.LinearSolveAdjoint{Missing}}, Vector{Float64}, Matrix{Float64}, LinearAlgebra.Transpose{Float64, Vector{Float64}}, Vector{Float64}, ADTypes.AutoForwardDiff{nothing, Nothing}}(LinearSolve.LinearCache{SparseMatrixCSC{Float64, Int64}, Matrix{Float64}, Matrix{Float64}, SciMLBase.NullParameters, LinearSolve.UMFPACKFactorization, SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}, SciMLOperators.IdentityOperator, SciMLOperators.IdentityOperator, Float64, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}, Bool, LinearSolve.LinearSolveAdjoint{Missing}}(sparse([1, 2, 9551, 9606, 84662, 1, 2, 3, 57, 9552 … 84657, 169318, 84658, 169319, 84659, 169320, 84660, 169321, 84661, 169322], [1, 1, 1, 1, 1, 2, 2, 2, 2, 2 … 169318, 169318, 169319, 169319, 169320, 169320, 169321, 169321, 169322, 169322], [-3.4050287839078227e-7, -1.8030621712982388e-10, 2.202759122121374e-7, 1.669794081839543e-8, 1.1220839187426168e-7, 8.32821319857397e-8, -4.783272470058407e-7, 2.4345873581518986e-8, -2.346272528650749e-8, 3.294980961099172e-7 … 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7], 169322, 169322), [-4.512560904168677e-12 8.935618833654465e-13; -6.652536669931426e-12 1.317312568668086e-12; … ; -0.0 -2.135436433303535e-13; -0.0 -1.6871098776291334e-13], [5.7489882607834106e-5 -2.308373390511117e-5; 5.646653132302884e-5 -3.1428196353351436e-5; … ; -4.603312782669057e-8 1.1570159339507784e-6; 5.471653491870314e-9 8.932325243468619e-7], SciMLBase.NullParameters(), LinearSolve.UMFPACKFactorization(true, true), SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}(SparseArrays.UMFPACK.Symbolic{Float64, Int64}(Ptr{Nothing}(0x0000000125718b20)), SparseArrays.UMFPACK.Numeric{Float64, Int64}(Ptr{Nothing}(0x000000030121df80)), 169322, 169322, [0, 5, 12, 19, 26, 33, 40, 47, 54, 61 … 818900, 818902, 818904, 818906, 818908, 818910, 818912, 818914, 818916, 818918], [0, 1, 9550, 9605, 84661, 0, 1, 2, 56, 9551 … 84656, 169317, 84657, 169318, 84658, 169319, 84659, 169320, 84660, 169321], [-3.4050287839078227e-7, -1.8030621712982388e-10, 2.202759122121374e-7, 1.669794081839543e-8, 1.1220839187426168e-7, 8.32821319857397e-8, -4.783272470058407e-7, 2.4345873581518986e-8, -2.346272528650749e-8, 3.294980961099172e-7 … 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7], 0, SparseArrays.UMFPACK.UmfpackWS{Int64}([6167, 646, 2414, 2416, 6147, 6157, 6158, 675, 2402, 6137 … 100589, 98900, 100764, 97371, 98568, 102479, 2971, 98602, 3061, 98610], [0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 … -3.634301428314729e-9, 6.08605475490754e-9, 3.126003909338102e-8, 4.993000031340455e-8, 4.535492077596722e-9, 4.902645702569492e-9, 1.897352028675765e-7, 3.2195576264689603e-9, 2.6818076806826084e-7, -1.3166840863933614e-6]), [1.0, 0.2, 0.2, 0.1, 32.0, 0.0, 0.7, 0.0, 1.0, 0.3, 1.0, 1.0, 0.9, 0.0, 10.0, 0.001, 1.0, 0.5, 0.0, 1.0], [0.0, 169322.0, 818918.0, 16.0, 4.0, 8.0, 8.0, 8.0, 0.0, 0.0 … 0.0, 0.0, -1.0, -1.0, 4.1972978e7, 0.01968699999997625, 0.01968699999997625, -1.0, -1.0, -1.0], ReentrantLock()), false, false, IdentityOperator(169322), IdentityOperator(169322), 1.4901161193847656e-8, 1.4901161193847656e-8, 338644, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}(SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.DebugLevel(), SciMLLogging.Silent(), SciMLLogging.InfoLevel(), SciMLLogging.Silent(), SciMLLogging.WarnLevel(), SciMLLogging.WarnLevel(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.WarnLevel(), SciMLLogging.WarnLevel(), SciMLLogging.WarnLevel()), LinearSolve.OperatorAssumptions{Bool}(true, LinearSolve.OperatorCondition.IllConditioned), LinearSolve.LinearSolveAdjoint{Missing}(missing)), [4.7241543264007345e-5, 6.616497742507608e-5, 0.0007607097140667825, 0.0002001818247681397, 0.00015228345495271658, 0.0001377887715644329, 0.00014231181003870498, 0.00013249304462571178, 0.00012386276033894944, 8.471960438304386e-5 … 1.6870928837112562e-7, 1.107674177977734e-7, 1.6619353326168306e-7, 2.6533564688704694e-7, 1.0705570309942015e-6, 9.402540764874025e-7, 6.90519304703392e-7, 6.642473304457541e-7, 5.535051235122763e-7, 4.372988802814714e-7], [5.7489882607834106e-5 -2.308373390511117e-5; 5.646653132302884e-5 -3.1428196353351436e-5; … ; -4.603312782669057e-8 1.1570159339507784e-6; 5.471653491870314e-9 8.932325243468619e-7], transpose([-2.7297252656153418e-5, -0.00011568624430491445, -0.00011421097914680908, -0.00011716805364084193, -0.00010979266846249825, -0.00016721064603409375, -0.00043737172418684654, -0.0012819265870092442, -0.0009894443952221357, -0.0004643803065352418 … 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0]), [4.605170185988092, 3.4011973816621555], [4.605170185988092, 3.4011973816621555], ADTypes.AutoForwardDiff()), CTKAlg{LinearSolve.UMFPACKFactorization}(LinearSolve.UMFPACKFactorization(true, true))), nothing, F1MethodOptimizationExt.var"#9#10"{Base.Pairs{Symbol, Any, Nothing, @NamedTuple{maxItNewton::Int64, linear_cache::LinearSolve.LinearCache{SparseMatrixCSC{Float64, Int64}, Vector{Float64}, Vector{Float64}, SciMLBase.NullParameters, LinearSolve.UMFPACKFactorization, SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}, SciMLOperators.IdentityOperator, SciMLOperators.IdentityOperator, Float64, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}, Bool, LinearSolve.LinearSolveAdjoint{Missing}}}}, AIBECS.var"#f#generate_f##1"{Main.OptParams, Tuple{Float64, Float64}, Tuple{Vector{Float64}, Missing}, Tuple{Vector{Float64}, Missing}, Vector{Float64}, Float64, AIBECS.var"#tracers#generate_f##0"{Int64, Int64}}, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#AIBECSFunction##1"{Main.OptParams, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#AIBECSFunction##0"{Main.OptParams, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}, AIBECS.var"#∇ₓf#generate_∇ₓf##1"{Tuple{Float64, Float64}, Tuple{Vector{Float64}, Missing}, Tuple{Vector{Float64}, Missing}, Vector{Float64}, AIBECS.var"#tracers#generate_∇ₓf##0"{Int64, Int64}}, F1Method.F1Cache{LinearSolve.LinearCache{SparseMatrixCSC{Float64, Int64}, Matrix{Float64}, Matrix{Float64}, SciMLBase.NullParameters, LinearSolve.UMFPACKFactorization, SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}, SciMLOperators.IdentityOperator, SciMLOperators.IdentityOperator, Float64, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}, Bool, LinearSolve.LinearSolveAdjoint{Missing}}, Vector{Float64}, Matrix{Float64}, LinearAlgebra.Transpose{Float64, Vector{Float64}}, Vector{Float64}, ADTypes.AutoForwardDiff{nothing, Nothing}}, CTKAlg{LinearSolve.UMFPACKFactorization}}(Base.Pairs{Symbol, Any, Nothing, @NamedTuple{maxItNewton::Int64, linear_cache::LinearSolve.LinearCache{SparseMatrixCSC{Float64, Int64}, Vector{Float64}, Vector{Float64}, SciMLBase.NullParameters, LinearSolve.UMFPACKFactorization, SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}, SciMLOperators.IdentityOperator, SciMLOperators.IdentityOperator, Float64, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}, Bool, LinearSolve.LinearSolveAdjoint{Missing}}}}(:maxItNewton => 15, :linear_cache => LinearSolve.LinearCache{SparseMatrixCSC{Float64, Int64}, Vector{Float64}, Vector{Float64}, SciMLBase.NullParameters, LinearSolve.UMFPACKFactorization, SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}, SciMLOperators.IdentityOperator, SciMLOperators.IdentityOperator, Float64, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}, Bool, LinearSolve.LinearSolveAdjoint{Missing}}(sparse([1, 2, 9551, 9606, 84662, 1, 2, 3, 57, 9552 … 84657, 169318, 84658, 169319, 84659, 169320, 84660, 169321, 84661, 169322], [1, 1, 1, 1, 1, 2, 2, 2, 2, 2 … 169318, 169318, 169319, 169319, 169320, 169320, 169321, 169321, 169322, 169322], [-3.440352272546861e-7, -1.8030621712982388e-10, 2.202759122121374e-7, 1.669794081839543e-8, 1.157407407381655e-7, 8.32821319857397e-8, -4.803223973230515e-7, 2.4345873581518986e-8, -2.346272528650749e-8, 3.294980961099172e-7 … 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7], 169322, 169322), [0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 … 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0], [0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 … 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0], SciMLBase.NullParameters(), LinearSolve.UMFPACKFactorization(true, true), SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}(SparseArrays.UMFPACK.Symbolic{Float64, Int64}(Ptr{Nothing}(0x0000000125718b20)), SparseArrays.UMFPACK.Numeric{Float64, Int64}(Ptr{Nothing}(0x000000030121df80)), 169322, 169322, [0, 5, 12, 19, 26, 33, 40, 47, 54, 61 … 818900, 818902, 818904, 818906, 818908, 818910, 818912, 818914, 818916, 818918], [0, 1, 9550, 9605, 84661, 0, 1, 2, 56, 9551 … 84656, 169317, 84657, 169318, 84658, 169319, 84659, 169320, 84660, 169321], [-3.4050287839078227e-7, -1.8030621712982388e-10, 2.202759122121374e-7, 1.669794081839543e-8, 1.1220839187426168e-7, 8.32821319857397e-8, -4.783272470058407e-7, 2.4345873581518986e-8, -2.346272528650749e-8, 3.294980961099172e-7 … 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7], 0, SparseArrays.UMFPACK.UmfpackWS{Int64}([6167, 646, 2414, 2416, 6147, 6157, 6158, 675, 2402, 6137 … 100589, 98900, 100764, 97371, 98568, 102479, 2971, 98602, 3061, 98610], [0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 … -3.634301428314729e-9, 6.08605475490754e-9, 3.126003909338102e-8, 4.993000031340455e-8, 4.535492077596722e-9, 4.902645702569492e-9, 1.897352028675765e-7, 3.2195576264689603e-9, 2.6818076806826084e-7, -1.3166840863933614e-6]), [1.0, 0.2, 0.2, 0.1, 32.0, 0.0, 0.7, 0.0, 1.0, 0.3, 1.0, 1.0, 0.9, 0.0, 10.0, 0.001, 1.0, 0.5, 0.0, 1.0], [0.0, 169322.0, 818918.0, 16.0, 4.0, 8.0, 8.0, 8.0, 0.0, 0.0 … 0.0, 0.0, -1.0, -1.0, 4.1972978e7, 0.01968699999997625, 0.01968699999997625, -1.0, -1.0, -1.0], ReentrantLock()), false, false, IdentityOperator(169322), IdentityOperator(169322), 1.4901161193847656e-8, 1.4901161193847656e-8, 169322, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}(SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.DebugLevel(), SciMLLogging.Silent(), SciMLLogging.InfoLevel(), SciMLLogging.Silent(), SciMLLogging.WarnLevel(), SciMLLogging.WarnLevel(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.WarnLevel(), SciMLLogging.WarnLevel(), SciMLLogging.WarnLevel()), LinearSolve.OperatorAssumptions{Bool}(true, LinearSolve.OperatorCondition.IllConditioned), LinearSolve.LinearSolveAdjoint{Missing}(missing))), AIBECS.var"#f#generate_f##1"{Main.OptParams, Tuple{Float64, Float64}, Tuple{Vector{Float64}, Missing}, Tuple{Vector{Float64}, Missing}, Vector{Float64}, Float64, AIBECS.var"#tracers#generate_f##0"{Int64, Int64}}((1.0, 0.0), ([0.001993383847718889, 0.0019914651177959006, 0.0019652102196352048, 0.0019188907593353226, 0.0018741692451726304, 0.0018593180609968575, 0.0018762683281708847, 0.0018866130577569658, 0.0018618487300791524, 0.0017936109269207174 … 0.001597349353171885, 0.0016139239530996844, 0.0016146148645877836, 0.0016014571666717527, 0.0015958896704982308, 0.0015948419428572936, 0.0015868619350592294, 0.0015825331862406298, 0.0015800279642719972, 0.001589979759638126], missing), ([4.7499243508165185e-9, 2.3806476583116648e-9, 1.644850324443933e-9, 2.474518014121221e-9, 2.8420192007964507e-9, 1.9923744404030257e-9, 8.150346589335414e-10, 2.974775216581886e-10, 4.021567996605766e-10, 8.840013268111294e-10 … 4.92697369868074e-10, 2.974775216581886e-10, 2.974775216581886e-10, 2.3496760609714243e-9, 2.446816223482753e-9, 2.974775216581886e-10, 8.053041266157682e-10, 8.39703761723315e-10, 1.2798255287465522e-9, 1.7406554311458743e-9], missing), [4.1259441424340015e11, 8.858714556466658e11, 9.658719592542885e11, 1.0446956966607107e12, 1.122246633283809e12, 1.1984302852530718e12, 1.2731538345237942e12, 1.3463262419614636e12, 1.417858358258977e12, 1.487663032551257e12 … 5.551473635870713e13, 6.642031628003376e13, 5.1445109315594125e13, 4.172486815218424e13, 1.3409911172962852e13, 1.824228219687533e13, 2.7565991089214062e13, 3.717262516144412e13, 4.7318477973128164e13, 4.2196399175523055e13], 0.0001, AIBECS.var"#tracers#generate_f##0"{Int64, Int64}(84661, 2)), SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#AIBECSFunction##1"{Main.OptParams, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#AIBECSFunction##0"{Main.OptParams, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}(AIBECS.var"#f#AIBECSFunction##1"{Main.OptParams, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}}(SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}(AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}((Main.T_DIP, Main.T_POP), 1:2, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}((Main.G_DIP, Main.G_POP), AIBECS.var"#tracers#66"{Int64, Int64}(84661, 2)), AIBECS.var"#tracer#67"{Int64, Int64}(84661, 2)), LinearAlgebra.UniformScaling{Bool}(true), nothing, nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}(AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}((Main.T_DIP, Main.T_POP), 2, 1:2), AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}((Main.G_DIP, Main.G_POP), 84661, 2)), nothing, nothing, nothing, nothing, nothing, nothing, nothing, nothing, nothing, SciMLBase.DEFAULT_OBSERVED, nothing, nothing, nothing, nothing)), LinearAlgebra.UniformScaling{Bool}(true), nothing, nothing, AIBECS.var"#jac#AIBECSFunction##0"{Main.OptParams, SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}}(SciMLBase.ODEFunction{false, SciMLBase.FullSpecialize, AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}, LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing, Nothing, Nothing, Nothing}(AIBECS.var"#f#71"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, UnitRange{Int64}, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}, AIBECS.var"#tracer#67"{Int64, Int64}}((Main.T_DIP, Main.T_POP), 1:2, AIBECS.var"#G#68"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, AIBECS.var"#tracers#66"{Int64, Int64}}((Main.G_DIP, Main.G_POP), AIBECS.var"#tracers#66"{Int64, Int64}(84661, 2)), AIBECS.var"#tracer#67"{Int64, Int64}(84661, 2)), LinearAlgebra.UniformScaling{Bool}(true), nothing, nothing, AIBECS.var"#jac#81"{AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}, AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}}(AIBECS.var"#T#76"{Tuple{typeof(Main.T_DIP), typeof(Main.T_POP)}, Int64, UnitRange{Int64}}((Main.T_DIP, Main.T_POP), 2, 1:2), AIBECS.var"#72#73"{Tuple{typeof(Main.G_DIP), typeof(Main.G_POP)}, Int64, Int64}((Main.G_DIP, Main.G_POP), 84661, 2)), nothing, nothing, nothing, nothing, nothing, nothing, nothing, nothing, nothing, SciMLBase.DEFAULT_OBSERVED, nothing, nothing, nothing, nothing)), nothing, nothing, nothing, nothing, nothing, nothing, nothing, nothing, nothing, SciMLBase.DEFAULT_OBSERVED, nothing, nothing, nothing, nothing), AIBECS.var"#∇ₓf#generate_∇ₓf##1"{Tuple{Float64, Float64}, Tuple{Vector{Float64}, Missing}, Tuple{Vector{Float64}, Missing}, Vector{Float64}, AIBECS.var"#tracers#generate_∇ₓf##0"{Int64, Int64}}((1.0, 0.0), ([0.001993383847718889, 0.0019914651177959006, 0.0019652102196352048, 0.0019188907593353226, 0.0018741692451726304, 0.0018593180609968575, 0.0018762683281708847, 0.0018866130577569658, 0.0018618487300791524, 0.0017936109269207174 … 0.001597349353171885, 0.0016139239530996844, 0.0016146148645877836, 0.0016014571666717527, 0.0015958896704982308, 0.0015948419428572936, 0.0015868619350592294, 0.0015825331862406298, 0.0015800279642719972, 0.001589979759638126], missing), ([4.7499243508165185e-9, 2.3806476583116648e-9, 1.644850324443933e-9, 2.474518014121221e-9, 2.8420192007964507e-9, 1.9923744404030257e-9, 8.150346589335414e-10, 2.974775216581886e-10, 4.021567996605766e-10, 8.840013268111294e-10 … 4.92697369868074e-10, 2.974775216581886e-10, 2.974775216581886e-10, 2.3496760609714243e-9, 2.446816223482753e-9, 2.974775216581886e-10, 8.053041266157682e-10, 8.39703761723315e-10, 1.2798255287465522e-9, 1.7406554311458743e-9], missing), [4.1259441424340015e11, 8.858714556466658e11, 9.658719592542885e11, 1.0446956966607107e12, 1.122246633283809e12, 1.1984302852530718e12, 1.2731538345237942e12, 1.3463262419614636e12, 1.417858358258977e12, 1.487663032551257e12 … 5.551473635870713e13, 6.642031628003376e13, 5.1445109315594125e13, 4.172486815218424e13, 1.3409911172962852e13, 1.824228219687533e13, 2.7565991089214062e13, 3.717262516144412e13, 4.7318477973128164e13, 4.2196399175523055e13], AIBECS.var"#tracers#generate_∇ₓf##0"{Int64, Int64}(84661, 2), Core.Box(AIBECS.var"#∇ₓf#generate_∇ₓf##1"{Tuple{Float64, Float64}, Tuple{Vector{Float64}, Missing}, Tuple{Vector{Float64}, Missing}, Vector{Float64}, AIBECS.var"#tracers#generate_∇ₓf##0"{Int64, Int64}}(#= circular reference @-2 =#))), F1Method.F1Cache{LinearSolve.LinearCache{SparseMatrixCSC{Float64, Int64}, Matrix{Float64}, Matrix{Float64}, SciMLBase.NullParameters, LinearSolve.UMFPACKFactorization, SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}, SciMLOperators.IdentityOperator, SciMLOperators.IdentityOperator, Float64, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}, Bool, LinearSolve.LinearSolveAdjoint{Missing}}, Vector{Float64}, Matrix{Float64}, LinearAlgebra.Transpose{Float64, Vector{Float64}}, Vector{Float64}, ADTypes.AutoForwardDiff{nothing, Nothing}}(LinearSolve.LinearCache{SparseMatrixCSC{Float64, Int64}, Matrix{Float64}, Matrix{Float64}, SciMLBase.NullParameters, LinearSolve.UMFPACKFactorization, SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}, SciMLOperators.IdentityOperator, SciMLOperators.IdentityOperator, Float64, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}, Bool, LinearSolve.LinearSolveAdjoint{Missing}}(sparse([1, 2, 9551, 9606, 84662, 1, 2, 3, 57, 9552 … 84657, 169318, 84658, 169319, 84659, 169320, 84660, 169321, 84661, 169322], [1, 1, 1, 1, 1, 2, 2, 2, 2, 2 … 169318, 169318, 169319, 169319, 169320, 169320, 169321, 169321, 169322, 169322], [-3.4050287839078227e-7, -1.8030621712982388e-10, 2.202759122121374e-7, 1.669794081839543e-8, 1.1220839187426168e-7, 8.32821319857397e-8, -4.783272470058407e-7, 2.4345873581518986e-8, -2.346272528650749e-8, 3.294980961099172e-7 … 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7], 169322, 169322), [-4.512560904168677e-12 8.935618833654465e-13; -6.652536669931426e-12 1.317312568668086e-12; … ; -0.0 -2.135436433303535e-13; -0.0 -1.6871098776291334e-13], [5.7489882607834106e-5 -2.308373390511117e-5; 5.646653132302884e-5 -3.1428196353351436e-5; … ; -4.603312782669057e-8 1.1570159339507784e-6; 5.471653491870314e-9 8.932325243468619e-7], SciMLBase.NullParameters(), LinearSolve.UMFPACKFactorization(true, true), SparseArrays.UMFPACK.UmfpackLU{Float64, Int64}(SparseArrays.UMFPACK.Symbolic{Float64, Int64}(Ptr{Nothing}(0x0000000125718b20)), SparseArrays.UMFPACK.Numeric{Float64, Int64}(Ptr{Nothing}(0x000000030121df80)), 169322, 169322, [0, 5, 12, 19, 26, 33, 40, 47, 54, 61 … 818900, 818902, 818904, 818906, 818908, 818910, 818912, 818914, 818916, 818918], [0, 1, 9550, 9605, 84661, 0, 1, 2, 56, 9551 … 84656, 169317, 84657, 169318, 84658, 169319, 84659, 169320, 84660, 169321], [-3.4050287839078227e-7, -1.8030621712982388e-10, 2.202759122121374e-7, 1.669794081839543e-8, 1.1220839187426168e-7, 8.32821319857397e-8, -4.783272470058407e-7, 2.4345873581518986e-8, -2.346272528650749e-8, 3.294980961099172e-7 … 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7, 3.858024691358024e-7, -3.8583415722361643e-7], 0, SparseArrays.UMFPACK.UmfpackWS{Int64}([6167, 646, 2414, 2416, 6147, 6157, 6158, 675, 2402, 6137 … 100589, 98900, 100764, 97371, 98568, 102479, 2971, 98602, 3061, 98610], [0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 … -3.634301428314729e-9, 6.08605475490754e-9, 3.126003909338102e-8, 4.993000031340455e-8, 4.535492077596722e-9, 4.902645702569492e-9, 1.897352028675765e-7, 3.2195576264689603e-9, 2.6818076806826084e-7, -1.3166840863933614e-6]), [1.0, 0.2, 0.2, 0.1, 32.0, 0.0, 0.7, 0.0, 1.0, 0.3, 1.0, 1.0, 0.9, 0.0, 10.0, 0.001, 1.0, 0.5, 0.0, 1.0], [0.0, 169322.0, 818918.0, 16.0, 4.0, 8.0, 8.0, 8.0, 0.0, 0.0 … 0.0, 0.0, -1.0, -1.0, 4.1972978e7, 0.01968699999997625, 0.01968699999997625, -1.0, -1.0, -1.0], ReentrantLock()), false, false, IdentityOperator(169322), IdentityOperator(169322), 1.4901161193847656e-8, 1.4901161193847656e-8, 338644, LinearSolve.LinearVerbosity{SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.DebugLevel, SciMLLogging.Silent, SciMLLogging.InfoLevel, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.Silent, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel, SciMLLogging.WarnLevel}(SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.DebugLevel(), SciMLLogging.Silent(), SciMLLogging.InfoLevel(), SciMLLogging.Silent(), SciMLLogging.WarnLevel(), SciMLLogging.WarnLevel(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.Silent(), SciMLLogging.WarnLevel(), SciMLLogging.WarnLevel(), SciMLLogging.WarnLevel()), LinearSolve.OperatorAssumptions{Bool}(true, LinearSolve.OperatorCondition.IllConditioned), LinearSolve.LinearSolveAdjoint{Missing}(missing)), [4.7241543264007345e-5, 6.616497742507608e-5, 0.0007607097140667825, 0.0002001818247681397, 0.00015228345495271658, 0.0001377887715644329, 0.00014231181003870498, 0.00013249304462571178, 0.00012386276033894944, 8.471960438304386e-5 … 1.6870928837112562e-7, 1.107674177977734e-7, 1.6619353326168306e-7, 2.6533564688704694e-7, 1.0705570309942015e-6, 9.402540764874025e-7, 6.90519304703392e-7, 6.642473304457541e-7, 5.535051235122763e-7, 4.372988802814714e-7], [5.7489882607834106e-5 -2.308373390511117e-5; 5.646653132302884e-5 -3.1428196353351436e-5; … ; -4.603312782669057e-8 1.1570159339507784e-6; 5.471653491870314e-9 8.932325243468619e-7], transpose([-2.7297252656153418e-5, -0.00011568624430491445, -0.00011421097914680908, -0.00011716805364084193, -0.00010979266846249825, -0.00016721064603409375, -0.00043737172418684654, -0.0012819265870092442, -0.0009894443952221357, -0.0004643803065352418 … 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0]), [4.605170185988092, 3.4011973816621555], [4.605170185988092, 3.4011973816621555], ADTypes.AutoForwardDiff()), CTKAlg{LinearSolve.UMFPACKFactorization}(LinearSolve.UMFPACKFactorization(true, true))), nothing, nothing, nothing, nothing, nothing, nothing, nothing, nothing, nothing, nothing, SciMLBase.DEFAULT_OBSERVED_NO_TIME, nothing, nothing, nothing, nothing, nothing, nothing, nothing, nothing, nothing, nothing)Run the optimiser
maxiters = 20 is plenty to watch the outer Newton trust-region trace drive the gradient norm down by several decades. We attach a callback that records the parameter vector at each iteration so the diagnostics figure below can overlay the Newton trajectory on the cost-surface contourf.
λ_trajectory = Vector{Vector{Float64}}()
push!(λ_trajectory, copy(λ₀))
trajectory_cb = function (state, _)
push!(λ_trajectory, copy(state.u))
return false
end
prob = OptimizationProblem(optfn, copy(λ₀))
results = solve(prob, NewtonTrustRegion();
maxiters = 20, callback = trajectory_cb, show_trace = true)
p_opt = λ2p(OptParams, results.u)Main.OptParams{Float64}
w₀ = 5.0 (m d⁻¹)
w′ = 0.035 (d⁻¹)
τ_DIP = 398.61539377922907 (d)
k = 10.0 (μmol m⁻³)
z₀ = 80.0 (m)
τ_POP = 12.169557505218696 (d)
τ_geo = 1.0 (kyr)
DIP_geo = 2.12 (mmol m⁻³)Cost-function sweep around the trajectory
Evaluate f(x_ss, λ) on a 10×10 grid in λ-space centred on the trajectory (padded by ±0.2). Each cell solves the steady state from the same cold-start x0, warm-started by the shared UMFPACK cache from the previous section. Takes a couple of minutes.
n_grid = 10
pad = 0.2
λ1_lo, λ1_hi = extrema(λ[1] for λ in λ_trajectory) .+ (-pad, pad)
λ2_lo, λ2_hi = extrema(λ[2] for λ in λ_trajectory) .+ (-pad, pad)
λ1_vec = range(λ1_lo, λ1_hi; length = n_grid)
λ2_vec = range(λ2_lo, λ2_hi; length = n_grid)
F_grid = Matrix{Float64}(undef, n_grid, n_grid)
for j in 1:n_grid, i in 1:n_grid
λ = [λ1_vec[i], λ2_vec[j]]
sol = solve(SteadyStateProblem(F, x0, λ), CTKAlg(); solver_kwargs...)
F_grid[i, j] = f(sol.u, λ)
end┌ Warning: CTKAlg did not converge to a successful retcode.
│ retcode = ReturnCode.MaxIters = 4
│ residual_inf_norm = 2.1483971561731135e-14
│ abstol = 0.0
│ reltol = 3.168808781402895e-14
│ maxItNewton = 15
└ @ AIBECS ~/work/AIBECS.jl/AIBECS.jl/src/overload_solve.jl:164
┌ Warning: CTKAlg did not converge to a successful retcode.
│ retcode = ReturnCode.MaxIters = 4
│ residual_inf_norm = 2.4991252418133468e-15
│ abstol = 0.0
│ reltol = 3.168808781402895e-14
│ maxItNewton = 15
└ @ AIBECS ~/work/AIBECS.jl/AIBECS.jl/src/overload_solve.jl:164
┌ Warning: CTKAlg did not converge to a successful retcode.
│ retcode = ReturnCode.MaxIters = 4
│ residual_inf_norm = 3.784922900580025e-15
│ abstol = 0.0
│ reltol = 3.168808781402895e-14
│ maxItNewton = 15
└ @ AIBECS ~/work/AIBECS.jl/AIBECS.jl/src/overload_solve.jl:164
┌ Warning: CTKAlg did not converge to a successful retcode.
│ retcode = ReturnCode.MaxIters = 4
│ residual_inf_norm = 7.486948351552839e-11
│ abstol = 0.0
│ reltol = 3.168808781402895e-14
│ maxItNewton = 15
└ @ AIBECS ~/work/AIBECS.jl/AIBECS.jl/src/overload_solve.jl:164
┌ Warning: CTKAlg did not converge to a successful retcode.
│ retcode = ReturnCode.MaxIters = 4
│ residual_inf_norm = 1.3309906925042203e-14
│ abstol = 0.0
│ reltol = 3.168808781402895e-14
│ maxItNewton = 15
└ @ AIBECS ~/work/AIBECS.jl/AIBECS.jl/src/overload_solve.jl:164
┌ Warning: CTKAlg did not converge to a successful retcode.
│ retcode = ReturnCode.MaxIters = 4
│ residual_inf_norm = 2.8351611483166536e-12
│ abstol = 0.0
│ reltol = 3.168808781402895e-14
│ maxItNewton = 15
└ @ AIBECS ~/work/AIBECS.jl/AIBECS.jl/src/overload_solve.jl:164
┌ Warning: CTKAlg did not converge to a successful retcode.
│ retcode = ReturnCode.MaxIters = 4
│ residual_inf_norm = 5.532024368145776e-15
│ abstol = 0.0
│ reltol = 3.168808781402895e-14
│ maxItNewton = 15
└ @ AIBECS ~/work/AIBECS.jl/AIBECS.jl/src/overload_solve.jl:164
┌ Warning: CTKAlg did not converge to a successful retcode.
│ retcode = ReturnCode.MaxIters = 4
│ residual_inf_norm = 1.9248858459254306e-14
│ abstol = 0.0
│ reltol = 3.168808781402895e-14
│ maxItNewton = 15
└ @ AIBECS ~/work/AIBECS.jl/AIBECS.jl/src/overload_solve.jl:164
┌ Warning: CTKAlg did not converge to a successful retcode.
│ retcode = ReturnCode.MaxIters = 4
│ residual_inf_norm = 3.1690160794423613e-12
│ abstol = 0.0
│ reltol = 3.168808781402895e-14
│ maxItNewton = 15
└ @ AIBECS ~/work/AIBECS.jl/AIBECS.jl/src/overload_solve.jl:164
┌ Warning: CTKAlg did not converge to a successful retcode.
│ retcode = ReturnCode.MaxIters = 4
│ residual_inf_norm = 2.5973607520649588e-14
│ abstol = 0.0
│ reltol = 3.168808781402895e-14
│ maxItNewton = 15
└ @ AIBECS ~/work/AIBECS.jl/AIBECS.jl/src/overload_solve.jl:164
┌ Warning: CTKAlg did not converge to a successful retcode.
│ retcode = ReturnCode.MaxIters = 4
│ residual_inf_norm = 3.219302917198701e-12
│ abstol = 0.0
│ reltol = 3.168808781402895e-14
│ maxItNewton = 15
└ @ AIBECS ~/work/AIBECS.jl/AIBECS.jl/src/overload_solve.jl:164
┌ Warning: CTKAlg did not converge to a successful retcode.
│ retcode = ReturnCode.MaxIters = 4
│ residual_inf_norm = 3.312083398850626e-14
│ abstol = 0.0
│ reltol = 3.168808781402895e-14
│ maxItNewton = 15
└ @ AIBECS ~/work/AIBECS.jl/AIBECS.jl/src/overload_solve.jl:164
┌ Warning: CTKAlg did not converge to a successful retcode.
│ retcode = ReturnCode.MaxIters = 4
│ residual_inf_norm = 3.1786108611166313e-12
│ abstol = 0.0
│ reltol = 3.168808781402895e-14
│ maxItNewton = 15
└ @ AIBECS ~/work/AIBECS.jl/AIBECS.jl/src/overload_solve.jl:164
┌ Warning: CTKAlg did not converge to a successful retcode.
│ retcode = ReturnCode.MaxIters = 4
│ residual_inf_norm = 4.38177660736921e-14
│ abstol = 0.0
│ reltol = 3.168808781402895e-14
│ maxItNewton = 15
└ @ AIBECS ~/work/AIBECS.jl/AIBECS.jl/src/overload_solve.jl:164
┌ Warning: CTKAlg did not converge to a successful retcode.
│ retcode = ReturnCode.MaxIters = 4
│ residual_inf_norm = 1.0736993472968038e-13
│ abstol = 0.0
│ reltol = 3.168808781402895e-14
│ maxItNewton = 15
└ @ AIBECS ~/work/AIBECS.jl/AIBECS.jl/src/overload_solve.jl:164
┌ Warning: CTKAlg did not converge to a successful retcode.
│ retcode = ReturnCode.MaxIters = 4
│ residual_inf_norm = 1.8001670317891646e-11
│ abstol = 0.0
│ reltol = 3.168808781402895e-14
│ maxItNewton = 15
└ @ AIBECS ~/work/AIBECS.jl/AIBECS.jl/src/overload_solve.jl:164
┌ Warning: CTKAlg did not converge to a successful retcode.
│ retcode = ReturnCode.MaxIters = 4
│ residual_inf_norm = 1.6546650051179805e-12
│ abstol = 0.0
│ reltol = 3.168808781402895e-14
│ maxItNewton = 15
└ @ AIBECS ~/work/AIBECS.jl/AIBECS.jl/src/overload_solve.jl:164
┌ Warning: CTKAlg did not converge to a successful retcode.
│ retcode = ReturnCode.MaxIters = 4
│ residual_inf_norm = 1.1063176324372192e-13
│ abstol = 0.0
│ reltol = 3.168808781402895e-14
│ maxItNewton = 15
└ @ AIBECS ~/work/AIBECS.jl/AIBECS.jl/src/overload_solve.jl:164Diagnostics — setup
We use CairoMakie (see also the Makie how-to) for plotting.
using CairoMakie
using OceanBasins
using KernelDensity
CairoMakie.activate!(type = "png")
prob_opt = SteadyStateProblem(F, x0, p_opt)
DIPopt, _ = state_to_tracers(solve(prob_opt, CTKAlg(); solver_kwargs...).u, grd)
DIPopt_μM = ustrip.(μM, DIPopt .* (mol/m^3))
μWOA_μM = ustrip.(μM, μDIPobs .* (mol/m^3))
diff_μM = DIPopt_μM .- μWOA_μM
depth = ustrip.(grd.depth)
lat = ustrip.(grd.lat)
lon = ustrip.(grd.lon)
maxdepth_round = ceil(maximum(depth) / 1000) * 1000
dip_levels = 0:0.25:3.0
dip_cmap = cgrad(:viridis, length(dip_levels) - 1; categorical = true)
diff_levels = -1.0:0.2:1.0
diff_cmap = cgrad(:RdBu, length(diff_levels) - 1; categorical = true, rev = true)
latticks = (-90:30:90, ["90°S", "60°S", "30°S", "0°", "30°N", "60°N", "90°N"])
lonticks = (0:60:360, ["0°", "60°E", "120°E", "180°", "120°W", "60°W", "0°"])
gcdepthticks = (0:1000:maxdepth_round, string.(Int.(0:1000:maxdepth_round)))
plus_minus_ticks = (-1:0.5:1, ["−1", "−0.5", "0", "+0.5", "+1"])(-1.0:0.5:1.0, ["−1", "−0.5", "0", "+0.5", "+1"])Basin masks, each extended to the pole via its Southern Ocean sector.
OCEANS = oceanpolygons()
lat_grd, lon_grd = latvec(grd), lonvec(grd)
mPAC = ispacific2(lat_grd, lon_grd, OCEANS) .| isSOpacific(lat_grd, lon_grd, OCEANS)
mATL = isatlantic2(lat_grd, lon_grd, OCEANS) .| isSOatlantic(lat_grd, lon_grd, OCEANS)
mIND = isindian2(lat_grd, lon_grd, OCEANS) .| isSOindian(lat_grd, lon_grd, OCEANS)
basin_specs = [
(name = "Pacific", mask = mPAC, color = :dodgerblue),
(name = "Atlantic", mask = mATL, color = :firebrick),
(name = "Indian", mask = mIND, color = :darkorange),
]3-element Vector{@NamedTuple{name::String, mask::BitVector, color::Symbol}}:
(name = "Pacific", mask = [0, 0, 0, 0, 0, 0, 0, 0, 0, 0 … 0, 0, 0, 0, 0, 0, 0, 0, 0, 0], color = :dodgerblue)
(name = "Atlantic", mask = [1, 1, 1, 1, 1, 1, 1, 1, 1, 1 … 1, 1, 1, 1, 1, 1, 1, 1, 1, 1], color = :firebrick)
(name = "Indian", mask = [0, 0, 0, 0, 0, 0, 0, 0, 0, 0 … 0, 0, 0, 0, 0, 0, 0, 0, 0, 0], color = :darkorange)Optimised vs WOA scatter
Each wet box at (WOA, optimised), coloured by its mass-weighted density quantile under a volume-weighted bivariate KDE. Walking the sorted cumulative density once gives both the iso-mass contour levels and the per-point colours.
bandwidth = 0.25 .* KernelDensity.default_bandwidth((μWOA_μM, DIPopt_μM))
d_kde = kde((μWOA_μM, DIPopt_μM); bandwidth, weights = v ./ sum(v))
ik = InterpKDE(d_kde)
dx, dy = step(d_kde.x), step(d_kde.y)
ρ_asc = sort(vec(d_kde.density))
cum_asc = cumsum(ρ_asc) .* (dx * dy)
cum_asc ./= cum_asc[end]
mass_qs = 0.1:0.1:0.9
levels = [ρ_asc[searchsortedfirst(cum_asc, q)] for q in mass_qs]
point_q = [(i = searchsortedlast(ρ_asc, d);
i == 0 ? 0.0 : cum_asc[i]) for d in ik.itp.(μWOA_μM, DIPopt_μM)]
fig = Figure()
ax = Axis(fig[1, 1];
xlabel = "WOA PO₄ (μM)", ylabel = "optimised PO₄ (μM)",
xgridvisible = false, ygridvisible = false,
aspect = DataAspect(), limits = (0, 3, 0, 3))
sc = scatter!(ax, μWOA_μM, DIPopt_μM;
markersize = 4, color = point_q, colormap = :plasma, colorrange = (0, 1))
contour!(ax, d_kde.x, d_kde.y, d_kde.density;
levels = levels, color = (:black, 0.5), linewidth = 1.2)
lines!(ax, [0, 3], [0, 3]; color = :red, linestyle = :dash)
Colorbar(fig[1, 2], sc; label = "density quantile", ticks = mass_qs)
fig
Per-basin depth profiles — model vs WOA
fig = Figure(size = (900, 350))
for (c, basin) in enumerate(basin_specs)
ax = Axis(fig[1, c];
xgridvisible = false, ygridvisible = false,
xlabel = "PO₄ (μM)",
ylabel = c == 1 ? "depth (m)" : "",
yreversed = true, yticks = gcdepthticks,
title = basin.name,
limits = (nothing, (0, maxdepth_round)))
prof_model = AIBECS.horizontalmean(DIPopt_μM, grd, basin.mask)
prof_obs = AIBECS.horizontalmean(μWOA_μM, grd, basin.mask)
lines!(ax, prof_model, depth;
color = basin.color, linewidth = 2, label = "model")
lines!(ax, prof_obs, depth;
color = basin.color, linewidth = 2, linestyle = :dash, label = "WOA")
c == 1 && axislegend(ax; position = :rb, framevisible = false)
c > 1 && hideydecorations!(ax;
ticklabels = true, label = true, ticks = false, grid = false)
end
fig
Per-basin zonal means — model / WOA / model − WOA
row_specs = [
(label = "model", field = DIPopt_μM, levels = dip_levels, cmap = dip_cmap),
(label = "WOA", field = μWOA_μM, levels = dip_levels, cmap = dip_cmap),
(label = "model − WOA", field = diff_μM, levels = diff_levels, cmap = diff_cmap),
]
zm_grid = [
[AIBECS.zonalmean(spec.field, grd, basin.mask) for basin in basin_specs]
for spec in row_specs
]
function wetlatrange(zm, pad = 5)
haswet = [any(!isnan, view(zm, i, :)) for i in eachindex(lat)]
wetidx = findall(haswet)
isempty(wetidx) && return (-90.0, 90.0)
return (max(-90, lat[first(wetidx)] - pad),
min( 90, lat[last(wetidx)] + pad))
end
xlim_per_basin = [wetlatrange(zm_grid[1][c]) for c in eachindex(basin_specs)]
# Documenter `@example` blocks run at module scope, where bare
# `co_v = co` inside a `for` would shadow the outer binding; `Ref` sidesteps it.
co_v = Ref{Any}()
co_d = Ref{Any}()
fig = Figure(size = (1400, 900))
for (r, spec) in enumerate(row_specs), (c, basin) in enumerate(basin_specs)
xlim = xlim_per_basin[c]
ax = Axis(fig[r, c];
backgroundcolor = :lightgray,
xgridvisible = false, ygridvisible = false,
xticksmirrored = true, yticksmirrored = true,
xticks = latticks, yticks = gcdepthticks,
yreversed = true,
ylabel = c == 1 ? "depth (m)" : "",
title = r == 1 ? basin.name : "",
limits = (xlim[1], xlim[2], 0, maxdepth_round))
co = contourf!(ax, lat, depth, zm_grid[r][c];
levels = spec.levels, colormap = spec.cmap, nan_color = :lightgray,
extendlow = spec.cmap[1], extendhigh = spec.cmap[end])
r == 1 && c == 1 && (co_v[] = co)
r == 3 && c == 1 && (co_d[] = co)
c > 1 && hideydecorations!(ax; ticklabels = true, label = true, ticks = false, grid = false)
r < 3 && hidexdecorations!(ax; ticklabels = true, label = true, ticks = false, grid = false)
end
Label(fig[1, 0], "model"; rotation = π/2, font = :bold, fontsize = 16, tellheight = false)
Label(fig[2, 0], "WOA"; rotation = π/2, font = :bold, fontsize = 16, tellheight = false)
Label(fig[3, 0], "model − WOA"; rotation = π/2, font = :bold, fontsize = 16, tellheight = false)
Colorbar(fig[1:2, length(basin_specs) + 1], co_v[];
label = "PO₄ (μM)", ticks = 0:1:3)
Colorbar(fig[3, length(basin_specs) + 1], co_d[];
label = "model − WOA (μM)", ticks = plus_minus_ticks)
colsize!(fig.layout, 0, Auto(0.15))
for (c, xlim) in enumerate(xlim_per_basin)
colsize!(fig.layout, c, Auto(xlim[2] - xlim[1]))
end
colsize!(fig.layout, length(basin_specs) + 1, Auto(0.25))
fig
Horizontal slices at 1000 m — model / WOA / diff
sl_opt = AIBECS.horizontalslice(DIPopt_μM, grd; depth = 1000)
sl_obs = AIBECS.horizontalslice(μWOA_μM, grd; depth = 1000)
sl_diff = sl_opt .- sl_obs
function map_axis(grid_pos; title = "")
Axis(grid_pos;
backgroundcolor = :lightgray,
xgridvisible = false, ygridvisible = false,
xticks = lonticks, yticks = latticks,
title = title, limits = (0, 360, -90, 90))
end
fig = Figure(size = (1400, 500))
ax_m = map_axis(fig[1, 1]; title = "optimised")
co_m = contourf!(ax_m, lon, lat, sl_opt';
levels = dip_levels, colormap = dip_cmap, nan_color = :lightgray,
extendlow = dip_cmap[1], extendhigh = dip_cmap[end])
ax_o = map_axis(fig[1, 2]; title = "WOA")
contourf!(ax_o, lon, lat, sl_obs';
levels = dip_levels, colormap = dip_cmap, nan_color = :lightgray,
extendlow = dip_cmap[1], extendhigh = dip_cmap[end])
ax_d = map_axis(fig[1, 3]; title = "optimised − WOA")
co_dd = contourf!(ax_d, lon, lat, sl_diff';
levels = diff_levels, colormap = diff_cmap, nan_color = :lightgray,
extendlow = diff_cmap[1], extendhigh = diff_cmap[end])
linkyaxes!(ax_m, ax_o, ax_d)
hideydecorations!(ax_o; ticklabels = true, label = false, ticks = false, grid = false)
hideydecorations!(ax_d; ticklabels = true, label = false, ticks = false, grid = false)
cb_v = Colorbar(fig[2, 1:2], co_m; vertical = false, flipaxis = false,
label = "PO₄ @ 1000 m (μM)", ticks = 0:1:3)
cb_v.width = Relative(0.7)
cb_diff = Colorbar(fig[2, 3], co_dd; vertical = false, flipaxis = false,
label = "optimised − WOA (μM)", ticks = plus_minus_ticks)
cb_diff.width = Relative(0.85)
fig
Side-by-side parameter tables, initial vs optimised
AIBECS.table(p₀)| Row | Symbol | Value | Initial value | Unit | Bounds | Logscaled | Optimizable |
|---|---|---|---|---|---|---|---|
| Symbol | Float64 | Float64 | FreeUnit… | Tuple… | Bool | Bool | |
| 1 | w₀ | 5.0 | 5.0 | m d⁻¹ | (0, Inf) | true | false |
| 2 | w′ | 0.035 | 0.035 | d⁻¹ | (0, Inf) | true | false |
| 3 | τ_DIP | 100.0 | 100.0 | d | (0, Inf) | true | true |
| 4 | k | 10.0 | 10.0 | μmol m⁻³ | (0, Inf) | true | false |
| 5 | z₀ | 80.0 | 80.0 | m | (0, Inf) | false | false |
| 6 | τ_POP | 30.0 | 30.0 | d | (0, Inf) | true | true |
| 7 | τ_geo | 1.0 | 1.0 | kyr | (0, Inf) | true | false |
| 8 | DIP_geo | 2.12 | 2.12 | mmol m⁻³ | (0, Inf) | false | false |
AIBECS.table(p_opt)| Row | Symbol | Value | Initial value | Unit | Bounds | Logscaled | Optimizable |
|---|---|---|---|---|---|---|---|
| Symbol | Float64 | Float64 | FreeUnit… | Tuple… | Bool | Bool | |
| 1 | w₀ | 5.0 | 5.0 | m d⁻¹ | (0, Inf) | true | false |
| 2 | w′ | 0.035 | 0.035 | d⁻¹ | (0, Inf) | true | false |
| 3 | τ_DIP | 398.615 | 100.0 | d | (0, Inf) | true | true |
| 4 | k | 10.0 | 10.0 | μmol m⁻³ | (0, Inf) | true | false |
| 5 | z₀ | 80.0 | 80.0 | m | (0, Inf) | false | false |
| 6 | τ_POP | 12.1696 | 30.0 | d | (0, Inf) | true | true |
| 7 | τ_geo | 1.0 | 1.0 | kyr | (0, Inf) | true | false |
| 8 | DIP_geo | 2.12 | 2.12 | mmol m⁻³ | (0, Inf) | false | false |
Extra: cost-surface sweep + Newton trajectory
Mostly diagnostic. Once you've converged, the optimal parameters are the point — but the sweep + trajectory is useful when something goes wrong. Convert the λ-space sweep grid and trajectory to physical (τ_DIP, τ_POP) via λ2p, then a single contourf! + scatterlines! does it.
τDIP_grid = [λ2p(OptParams, [λ1, λ₀[2]]).τ_DIP for λ1 in λ1_vec]
τPOP_grid = [λ2p(OptParams, [λ₀[1], λ2]).τ_POP for λ2 in λ2_vec]
τDIP_traj = [λ2p(OptParams, λ).τ_DIP for λ in λ_trajectory]
τPOP_traj = [λ2p(OptParams, λ).τ_POP for λ in λ_trajectory]
fig = Figure()
ax = Axis(fig[1, 1];
xlabel = "τ_DIP (d)", ylabel = "τ_POP (d)",
xscale = log10, yscale = log10,
title = "√f(p) cost surface + Newton trajectory")
co = contourf!(ax, τDIP_grid, τPOP_grid, sqrt.(F_grid);
levels = 20, colormap = cgrad(:BrBG; rev = true))
scatterlines!(ax, τDIP_traj, τPOP_traj;
color = :red, linewidth = 2, linestyle = :dash, marker = :circle, markersize = 4)
scatter!(ax, [τDIP_traj[end]], [τPOP_traj[end]];
marker = :star5, markersize = 22, color = :yellow,
strokecolor = :black, strokewidth = 1.5)
annotation!(ax, 70, -40, τDIP_traj[1], τPOP_traj[1];
text = "first iterate", path = Ann.Paths.Corner())
annotation!(ax, -70, +40, τDIP_traj[end], τPOP_traj[end];
text = "optimal solution", path = Ann.Paths.Corner(), color = :white)
Colorbar(fig[1, 2], co; label = "√f(p)")
fig
Where to next
This example only fits two parameters of a single-tracer mismatch objective on a moderate-resolution circulation; the multi-tracer, multi-observation fit lives in GNOM. See the roadmap for alternative parameter-estimation backends (Turing, ImplicitDifferentiation, SciMLSensitivity), the table-form f_and_∇ₓf(ωs, ωp, grd, obs, ::Type{P}) for fitting against WorldOceanAtlasTools.observations(...) tables, and other directions where contributions are welcome.
This page was generated using Literate.jl.