From 821b62146579c7bcfe8f9676329e593070d87490 Mon Sep 17 00:00:00 2001 From: Chris Rackauckas Date: Mon, 22 Aug 2022 05:09:20 -0400 Subject: [PATCH] add precompile snooping --- Project.toml | 2 ++ src/Sundials.jl | 46 ++++++++++++++++++++++++++++++++++++++++++++++ 2 files changed, 48 insertions(+) diff --git a/Project.toml b/Project.toml index 94e1e9b..7e008f8 100644 --- a/Project.toml +++ b/Project.toml @@ -13,6 +13,7 @@ Logging = "56ddb016-857b-54e1-b83d-db4d58db5568" Reexport = "189a3867-3050-52da-a836-e630ba90ab69" SparseArrays = "2f01184e-e22b-5df5-ae63-d93ebab69eaf" Sundials_jll = "fb77eaff-e24c-56d4-86b1-d163f2edb164" +SnoopPrecompile = "66db9d55-30c0-4569-8b51-7e840670fc0c" [compat] CEnum = "0.2, 0.3, 0.4" @@ -20,6 +21,7 @@ DataStructures = "0.18" DiffEqBase = "6.44" Reexport = "0.2, 1.0" Sundials_jll = "5.2" +SnoopPrecompile = "1" julia = "1.6" [extras] diff --git a/src/Sundials.jl b/src/Sundials.jl index c10ae61..c7f072d 100644 --- a/src/Sundials.jl +++ b/src/Sundials.jl @@ -83,6 +83,52 @@ include("common_interface/integrator_types.jl") include("common_interface/integrator_utils.jl") include("common_interface/solve.jl") +import SnoopPrecompile + +SnoopPrecompile.@precompile_all_calls begin + function lorenz(du, u, p, t) + du[1] = 10.0(u[2] - u[1]) + du[2] = u[1] * (28.0 - u[3]) - u[2] + du[3] = u[1] * u[2] - (8 / 3) * u[3] + end + + function lorenz_oop(u, p, t) + [10.0(u[2] - u[1]), u[1] * (28.0 - u[3]) - u[2], u[1] * u[2] - (8 / 3) * u[3]] + end + + solver_list = [ + ARKODE(), CVODE_Adams(), CVODE_BDF(), + ] + + prob_list = [ODEProblem(lorenz, [1.0; 0.0; 0.0], (0.0, 1.0)); + ODEProblem{true, false}(lorenz, [1.0; 0.0; 0.0], (0.0, 1.0)); + ODEProblem{true, false}(lorenz, [1.0; 0.0; 0.0], (0.0, 1.0), Float64[]); + ODEProblem(lorenz_oop, [1.0; 0.0; 0.0], (0.0, 1.0)); + ODEProblem{false, false}(lorenz_oop, [1.0; 0.0; 0.0], (0.0, 1.0)); + ODEProblem{false, false}(lorenz_oop, [1.0; 0.0; 0.0], (0.0, 1.0), + Float64[])] + + for prob in prob_list, solver in solver_list + solve(prob, solver)(5.0) + end + + prob_list = nothing + + function f(out, du, u, p, t) + out[1] = -0.04u[1] + 1e4 * u[2] * u[3] - du[1] + out[2] = +0.04u[1] - 3e7 * u[2]^2 - 1e4 * u[2] * u[3] - du[2] + out[3] = u[1] + u[2] + u[3] - 1.0 + end + u₀ = [1.0, 0, 0] + du₀ = [-0.04, 0.04, 0.0] + tspan = (0.0, 100000.0) + differential_vars = [true, true, false] + prob = DAEProblem(f, du₀, u₀, tspan, differential_vars = differential_vars) + sol = solve(prob, IDA()) + + prob = nothing +end + ################################################################## # Deprecations ##################################################################