|
| 1 | +# [Interactive Simulation and Plotting](@id interactive_brusselator) |
| 2 | + |
| 3 | +Catalyst can utilize the [GLMakie.jl](https://github.com/JuliaPlots/GLMakie.jl) package for creating interactive visualizations of your reaction network dynamics. This tutorial provides a step-by-step guide to creating an interactive visualization of the Brusselator model, building upon the basic [Brusselator](@ref basic_CRN_library_brusselator) example. |
| 4 | + |
| 5 | + |
| 6 | +## [Setting up the Brusselator model](@id setup_brusselator) |
| 7 | + |
| 8 | +Let's again use the oscillating Brusselator model, extending the basic simulation [plotting](@ref simulation_plotting) workflow we saw earlier. |
| 9 | + |
| 10 | +```@example interactive_brusselator; continued = true |
| 11 | +using Catalyst |
| 12 | +using OrdinaryDiffEq |
| 13 | +using GLMakie |
| 14 | +GLMakie.activate!(inline = true, visible = false) # hide |
| 15 | +
|
| 16 | +# Define the Brusselator model |
| 17 | +brusselator = @reaction_network begin |
| 18 | + A, ∅ → X |
| 19 | + 1, 2X + Y → 3X |
| 20 | + B, X → Y |
| 21 | + 1, X → ∅ |
| 22 | +end |
| 23 | +
|
| 24 | +# Initial parameter values and conditions |
| 25 | +p = [:A => 1.0, :B => 4.0] |
| 26 | +u0 = [:X => 1.0, :Y => 0.0] |
| 27 | +tspan = (0.0, 50.0) |
| 28 | +
|
| 29 | +oprob = ODEProblem(brusselator, u0, tspan, p) |
| 30 | +
|
| 31 | +# Function to solve the ODE |
| 32 | +function solve_brusselator(A, B, X0, Y0, prob = oprob) |
| 33 | + p = [:A => A, :B => B] |
| 34 | + u0 = [:X => X0, :Y => Y0] |
| 35 | + newprob = remake(prob, p=p, u0=u0) |
| 36 | + solve(newprob, Tsit5(), saveat = 0.1) |
| 37 | +end |
| 38 | +``` |
| 39 | +This code sets up our Brusselator model using Catalyst.jl's `@reaction_network` macro. We also define initial parameters, initial conditions, create an `ODEProblem`, and define a function to solve the ODE with given parameters. Setting `saveat = 0.1` in the call to `solve` ensures the solution is saved with the desired temporal frequency we want for our later plots. |
| 40 | + |
| 41 | +!!! note |
| 42 | + Be sure to set `saveat` to a value that is appropriate for your system; otherwise, the size of the solution can change during interactivity, which will cause dimension mismatch errors once we add our interactive elements. |
| 43 | + |
| 44 | +## [Basic static plotting](@id basic_static_plotting) |
| 45 | + |
| 46 | +Let's start by creating a basic plot of our Brusselator model: |
| 47 | + |
| 48 | +```@example interactive_brusselator |
| 49 | +# Create the main figure |
| 50 | +fig = Figure(size = (800, 600), fontsize = 18); |
| 51 | +
|
| 52 | +# Create an axis for the plot |
| 53 | +ax = Axis(fig[1, 1], |
| 54 | + title = "Brusselator Model", |
| 55 | + xlabel = "Time", |
| 56 | + ylabel = "Concentration") |
| 57 | +
|
| 58 | +# Solve the ODE |
| 59 | +sol = solve_brusselator(1.0, 4.0, 1.0, 0.0) |
| 60 | +
|
| 61 | +# Plot the solution |
| 62 | +lines!(ax, sol.t, sol[:X], label = "X", color = :blue, linewidth = 3) |
| 63 | +lines!(ax, sol.t, sol[:Y], label = "Y", color = :red, linewidth = 3) |
| 64 | +
|
| 65 | +# Add a legend |
| 66 | +axislegend(ax, position = :rt) |
| 67 | +
|
| 68 | +# Display the figure |
| 69 | +fig |
| 70 | +``` |
| 71 | + |
| 72 | +The plot shows the concentrations of species X and Y over time. Notice the oscillatory behavior characteristic of the Brusselator model. |
| 73 | + |
| 74 | +## [Adding interactivity](@id adding_interactivity) |
| 75 | + |
| 76 | +Now, let's add interactivity to our plot using Observables and sliders. We'll build this up step by step. |
| 77 | + |
| 78 | +### [Creating Observables](@id creating_observables) |
| 79 | + |
| 80 | +Observables are a key concept in reactive programming and are central to how Makie.jl creates interactive visualizations. You can read more about them [here](https://docs.makie.org/stable/explanations/observables). |
| 81 | + |
| 82 | +```@example interactive_brusselator; continued = true |
| 83 | +# Create observables for parameters and initial conditions |
| 84 | +A = Observable(1.0) |
| 85 | +B = Observable(4.0) |
| 86 | +X0 = Observable(1.0) |
| 87 | +Y0 = Observable(0.0) |
| 88 | +``` |
| 89 | + |
| 90 | +An Observable is a container for a value that can change over time. When the value changes, any dependent computations are automatically updated. |
| 91 | + |
| 92 | +### [Adding sliders and connecting to Observables](@id adding_sliders) |
| 93 | + |
| 94 | +Let's add [sliders](https://docs.makie.org/stable/reference/blocks/slider) that will control our Observables: |
| 95 | + |
| 96 | +```@example interactive_brusselator; continued = true |
| 97 | +# Create the main figure |
| 98 | +fig = Figure(size = (800, 600), fontsize = 18); |
| 99 | +
|
| 100 | +# Create layout for plot and sliders |
| 101 | +plot_layout = fig[1, 1] = GridLayout() |
| 102 | +slider_layout = fig[2, 1] = GridLayout() |
| 103 | +
|
| 104 | +# Create sliders |
| 105 | +slider_A = Slider(slider_layout[1, 1], range = 0.0:0.01:5.0, startvalue = to_value(A)) # to_value(A) unwraps the Observable to a value |
| 106 | +slider_B = Slider(slider_layout[2, 1], range = 0.0:0.01:5.0, startvalue = to_value(B)) |
| 107 | +slider_X0 = Slider(slider_layout[3, 1], range = 0.0:0.01:5.0, startvalue = to_value(X0)) |
| 108 | +slider_Y0 = Slider(slider_layout[4, 1], range = 0.0:0.01:5.0, startvalue = to_value(Y0)) |
| 109 | +
|
| 110 | +# Add labels for sliders |
| 111 | +Label(slider_layout[1, 1, Left()], "A") |
| 112 | +Label(slider_layout[2, 1, Left()], "B") |
| 113 | +Label(slider_layout[3, 1, Left()], "X₀") |
| 114 | +Label(slider_layout[4, 1, Left()], "Y₀") |
| 115 | +
|
| 116 | +# Connect the values of the sliders to the observables |
| 117 | +connect!(A, slider_A.value) |
| 118 | +connect!(B, slider_B.value) |
| 119 | +connect!(X0, slider_X0.value) |
| 120 | +connect!(Y0, slider_Y0.value) |
| 121 | +``` |
| 122 | + |
| 123 | +These sliders allow us to interactively change the parameters A and B, as well as the initial conditions X₀ and Y₀. |
| 124 | + |
| 125 | +### [Creating a reactive plot](@id reactive_plot) |
| 126 | + |
| 127 | +Now, let's create a plot that reacts to changes in our sliders: |
| 128 | + |
| 129 | +```@example interactive_brusselator |
| 130 | +# Create an axis for the plot |
| 131 | +ax = Axis(plot_layout[1, 1], |
| 132 | + title = "Brusselator Model", |
| 133 | + xlabel = "Time", |
| 134 | + ylabel = "Concentration") |
| 135 | +
|
| 136 | +# Create an observable for the solution |
| 137 | +# The `@lift` macro is used to create an Observable that depends on the observables `A`, `B`, `X0`, and `Y0`, and automatically updates when any of these observables change |
| 138 | +solution = @lift(solve_brusselator($A, $B, $X0, $Y0)) |
| 139 | +
|
| 140 | +# Plot the solution |
| 141 | +# We don't use the ODESolution plot recipe here, as you've seen in the previous examples where only the solution and an `idxs` argument was passed to the plot method, because we are passing in an Observable wrapping the solution |
| 142 | +lines!(ax, lift(sol -> sol.t, solution), lift(sol -> sol[:X], solution), label = "X", color = :blue, linewidth = 3) # `lift` can either be used as a function or a macro |
| 143 | +lines!(ax, lift(sol -> sol.t, solution), lift(sol -> sol[:Y], solution), label = "Y", color = :red, linewidth = 3) |
| 144 | +
|
| 145 | +# Add a legend |
| 146 | +axislegend(ax, position = :rt) |
| 147 | +
|
| 148 | +# Display the figure |
| 149 | +fig |
| 150 | +``` |
| 151 | + |
| 152 | +This plot will now update in real-time as you move the sliders, allowing for interactive exploration of the Brusselator's behavior under different conditions. |
| 153 | + |
| 154 | +## [Adding a phase plot](@id adding_phase_plot) |
| 155 | + |
| 156 | +To gain more insight into the system's behavior, let's enhance our visualization by adding a phase plot, along with some other improvements: |
| 157 | + |
| 158 | +```@example interactive_brusselator |
| 159 | +# Create the main figure |
| 160 | +fig = Figure(size = (1200, 800), fontsize = 18); |
| 161 | +
|
| 162 | +# Create main layout: plots on top, sliders at bottom |
| 163 | +plot_grid = fig[1, 1] = GridLayout() |
| 164 | +slider_grid = fig[2, 1] = GridLayout() |
| 165 | +
|
| 166 | +# Create sub-grids for plots |
| 167 | +time_plot = plot_grid[1, 1] = GridLayout() |
| 168 | +phase_plot = plot_grid[1, 2] = GridLayout() |
| 169 | +
|
| 170 | +# Create axes for the time series plot and phase plot |
| 171 | +ax_time = Axis(time_plot[1, 1], |
| 172 | + title = "Brusselator Model - Time Series", |
| 173 | + xlabel = "Time", |
| 174 | + ylabel = "Concentration") |
| 175 | +
|
| 176 | +ax_phase = Axis(phase_plot[1, 1], |
| 177 | + title = "Brusselator Model - Phase Plot", |
| 178 | + xlabel = "X", |
| 179 | + ylabel = "Y") |
| 180 | +
|
| 181 | +# Create sub-grids for sliders |
| 182 | +param_grid = slider_grid[1, 1] = GridLayout() |
| 183 | +ic_grid = slider_grid[1, 2] = GridLayout() |
| 184 | +
|
| 185 | +# Create observables for parameters and initial conditions |
| 186 | +A = Observable{Float64}(1.0) # We can specify the type of the Observable value, which can help with type stability and performance |
| 187 | +B = Observable{Float64}(4.0) |
| 188 | +X0 = Observable{Float64}(1.0) |
| 189 | +Y0 = Observable{Float64}(0.0) |
| 190 | +
|
| 191 | +# Create sliders with labels and group titles |
| 192 | +Label(param_grid[1, 1:2], "Parameters", fontsize = 22) |
| 193 | +slider_A = Slider(param_grid[2, 2], range = 0.0:0.01:5.0, startvalue = to_value(A)) |
| 194 | +slider_B = Slider(param_grid[3, 2], range = 0.0:0.01:5.0, startvalue = to_value(B)) |
| 195 | +Label(param_grid[2, 1], "A") |
| 196 | +Label(param_grid[3, 1], "B") |
| 197 | +
|
| 198 | +Label(ic_grid[1, 1:2], "Initial Conditions", fontsize = 22) |
| 199 | +slider_X0 = Slider(ic_grid[2, 2], range = 0.0:0.01:5.0, startvalue = to_value(X0)) |
| 200 | +slider_Y0 = Slider(ic_grid[3, 2], range = 0.0:0.01:5.0, startvalue = to_value(Y0)) |
| 201 | +Label(ic_grid[2, 1], "X₀") |
| 202 | +Label(ic_grid[3, 1], "Y₀") |
| 203 | +
|
| 204 | +# Connect sliders to observables |
| 205 | +connect!(A, slider_A.value) |
| 206 | +connect!(B, slider_B.value) |
| 207 | +connect!(X0, slider_X0.value) |
| 208 | +connect!(Y0, slider_Y0.value) |
| 209 | +
|
| 210 | +# Create an observable for the solution. |
| 211 | +solution = @lift(solve_brusselator($A, $B, $X0, $Y0)) |
| 212 | +
|
| 213 | +# Plot the time series |
| 214 | +lines!(ax_time, lift(sol -> sol.t, solution), lift(sol -> sol[:X], solution), label = "X", color = :blue, linewidth = 3) |
| 215 | +lines!(ax_time, lift(sol -> sol.t, solution), lift(sol -> sol[:Y], solution), label = "Y", color = :red, linewidth = 3) |
| 216 | +
|
| 217 | +# Plot the phase plot |
| 218 | +phase_plot_obj = lines!(ax_phase, lift(sol -> sol[:X], solution), lift(sol -> sol[:Y], solution), |
| 219 | + color = lift(sol -> sol.t, solution), colormap = :viridis) |
| 220 | +
|
| 221 | +# Add a colorbar for the phase plot |
| 222 | +Colorbar(phase_plot[1, 2], phase_plot_obj, label = "Time") |
| 223 | +
|
| 224 | +# Add legends |
| 225 | +axislegend(ax_time, position = :rt) |
| 226 | +
|
| 227 | +# Adjust layout to your liking |
| 228 | +colgap!(plot_grid, 20) |
| 229 | +rowgap!(fig.layout, 20) |
| 230 | +colgap!(param_grid, 10) |
| 231 | +colgap!(ic_grid, 10) |
| 232 | +
|
| 233 | +# Display the figure |
| 234 | +#fig |
| 235 | +``` |
| 236 | + |
| 237 | +This will create a visualization with both time series and phase plots: |
| 238 | + |
| 239 | + |
| 240 | + |
| 241 | +## [Common plotting options](@id common_makie_plotting_options) |
| 242 | + |
| 243 | +Various plotting options can be provided as optional arguments to the `lines!` command. Common options include: |
| 244 | +- `linewidth` or `lw`: Determine plot line widths. |
| 245 | +- `linestyle`: Determines plot line style. |
| 246 | +- `color`: Determines the line colors. |
| 247 | +- `label`: Determines label texts displayed in the legend. |
| 248 | + |
| 249 | +For example: |
| 250 | + |
| 251 | +```julia |
| 252 | +lines!(ax_time, lift(sol -> sol.t, solution), lift(sol -> sol[:X], solution), |
| 253 | + label = "X", color = :green, linewidth = 2, linestyle = :dash) |
| 254 | +``` |
| 255 | + |
| 256 | +## [Extending the interactive visualization](@id extending_interactive_visualization) |
| 257 | + |
| 258 | +You can further extend this visualization by: |
| 259 | +- Adding other interactive elements, such as [buttons](https://docs.makie.org/stable/reference/blocks/button) or [dropdown menus](https://docs.makie.org/stable/reference/blocks/menu) to control different aspects of the simulation or visualization. |
| 260 | +- Adding additonal axes to the plot, such as plotting the derivatives of the species. |
| 261 | +- Color coding the slider and slider labels to match the plot colors. |
| 262 | + |
| 263 | + |
0 commit comments