Skip to content
New issue

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

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

Already on GitHub? Sign in to your account

Change the plot recipes to use the SymbolicIndexingInterface #572

Merged
merged 2 commits into from
Dec 27, 2023

Conversation

ChrisRackauckas
Copy link
Member

Fixes a whole lot of small issues.

@ChrisRackauckas
Copy link
Member Author

using ModelingToolkit
using LinearAlgebra
using OrdinaryDiffEq
using Plots

function pendulum!(du, u, p, t)
    x, dx, y, dy, T = u
    g, L = p
    du[1] = dx
    du[2] = T * x
    du[3] = dy
    du[4] = T * y - g
    du[5] = x^2 + y^2 - L^2
    return nothing
end
pendulum_fun! = ODEFunction(pendulum!, mass_matrix = Diagonal([1, 1, 1, 1, 0]))
u0 = [1.0, 0, 0, 0, 0]
p = [9.8, 1]
tspan = (0, 10.0)
pendulum_prob = ODEProblem(pendulum_fun!, u0, tspan, p)
traced_sys = modelingtoolkitize(pendulum_prob)
pendulum_sys = structural_simplify(dae_index_lowering(traced_sys))
prob = ODAEProblem(pendulum_sys, [], tspan)
sol = solve(prob, Tsit5(), abstol = 1e-8, reltol = 1e-8)
t = independent_variables(pendulum_sys)[1]

plot(sol); savefig("plot1.png")
plot(sol, idxs = states(traced_sys)); savefig("plot2.png")
plot(sol, idxs = states(traced_sys)[4]); savefig("plot3.png")
plot(sol, idxs = states(traced_sys)[5]); savefig("plot4.png")
plot(sol, idxs = states(traced_sys)[4:5]); savefig("plot5.png")
plot(sol, idxs = tuple(states(traced_sys)[4:5]...), plotdensity=10); savefig("plot6.png")
plot(sol, idxs = tuple(states(traced_sys)[4:5]...)); savefig("plot7.png")
plot(sol, idxs = (0,states(traced_sys)[4])); savefig("plot8.png")
plot(sol, idxs = (0,states(traced_sys)[5])); savefig("plot9.png")
plot(sol, idxs = (t,states(traced_sys)[4])); savefig("plot10.png")
plot(sol, idxs = (t,states(traced_sys)[5])); savefig("plot11.png")
plot(sol, idxs = (0,states(traced_sys)[4],states(traced_sys)[5])); savefig("plot12.png")
plot(sol, idxs = (states(traced_sys)[3],states(traced_sys)[4],states(traced_sys)[5])); savefig("plot13.png")
plot(sol, idxs = (states(traced_sys)[3],0,states(traced_sys)[5])); savefig("plot14.png")
plot(sol, idxs = (0,states(traced_sys)[4])); savefig("plot15.png")
plot(sol, idxs = (t,states(traced_sys)[4],states(traced_sys)[5])); savefig("plot16.png")
plot(sol, idxs = (states(traced_sys)[3],t,states(traced_sys)[5])); savefig("plot17.png")

function lorenz(du, u, p, t)
    du[1] = p[1] * (u[2] - u[1])
    du[2] = u[1] * (p[2] - u[3]) - u[2]
    du[3] = u[1] * u[2] - p[3] * u[3]
end

u0 = [1.0, 5.0, 10.0]
tspan = (0.0, 100.0)
p = (10.0, 28.0, 8 / 3)
prob = ODEProblem(lorenz, u0, tspan, p)
sol = solve(prob, Tsit5()); savefig("plot18.png")
xyzt = plot(sol, plotdensity = 10000, lw = 1.5); savefig("plot19.png")
xy = plot(sol, plotdensity = 10000, idxs = (1, 2)); savefig("plot20.png")
xz = plot(sol, plotdensity = 10000, idxs = (1, 3)); savefig("plot21.png")
yz = plot(sol, plotdensity = 10000, idxs = (2, 3)); savefig("plot22.png")
xyz = plot(sol, plotdensity = 10000, idxs = (1, 2, 3)); savefig("plot23.png")
plot(plot(xyzt, xyz), plot(xy, xz, yz, layout = (1, 3), w = 1), layout = (2, 1)); savefig("plot23.png")

f(t, x, y, z) = (t, sqrt(x^2 + y^2 + z^2))
plot(sol, idxs = (f, 0, 1, 2, 3)); savefig("plot25.png")
plot!(sol); savefig("plot26.png")

f(x, y, z) = (sqrt(x^2 + y^2 + z^2), x)
plot(sol, idxs = (f, 1, 2, 3)); savefig("plot27.png")
plot!(sol, idxs = (1,2,3)); savefig("plot28.png")

# Linear ODE which starts at 0.5 and solves from t=0.0 to t=1.0
prob = ODEProblem((u, p, t) -> 1.01u, 0.5, (0.0, 1.0))
integrator = init(prob, Tsit5(); dt = 1 // 2^(4), tstops = [0.5])
plot(integrator)
for i in integrator
    display(plot!(integrator, idxs = (0, 1), legend = false))
end
step!(integrator);
plot!(integrator, idxs = (0, 1), legend = false);
savefig("iteratorplot.png")

plot1
plot2
plot3
plot4
plot5
plot6
plot7
plot8
plot9
plot10
plot11
plot12
plot13
plot14
plot15
plot16
plot17
plot18
plot19
plot20
plot21
plot22
plot23
plot25
plot26
plot27
plot28
iteratorplot

Copy link

codecov bot commented Dec 27, 2023

Codecov Report

Attention: 66 lines in your changes are missing coverage. Please review.

Comparison is base (1ccfa8d) 5.26% compared to head (db52700) 30.49%.

Files Patch % Lines
src/solutions/solution_interface.jl 40.90% 52 Missing ⚠️
src/integrator_interface.jl 0.00% 14 Missing ⚠️
Additional details and impacted files
@@             Coverage Diff             @@
##           master     #572       +/-   ##
===========================================
+ Coverage    5.26%   30.49%   +25.22%     
===========================================
  Files          53       53               
  Lines        4061     4020       -41     
===========================================
+ Hits          214     1226     +1012     
+ Misses       3847     2794     -1053     

☔ View full report in Codecov by Sentry.
📢 Have feedback on the report? Share it here.

@ChrisRackauckas ChrisRackauckas merged commit d7bc515 into master Dec 27, 2023
32 of 41 checks passed
@ChrisRackauckas ChrisRackauckas deleted the plotrecipe_symbolic branch December 27, 2023 20:52
Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Labels
None yet
Projects
None yet
Development

Successfully merging this pull request may close these issues.

1 participant