API
SimulationLogs.SimulationLogs — ModuleSimulationLogsLog internal variables in a DifferentialEquations.jl simulation with the @log macro and use the get_log function to retrieve the logged variables.
SimulationLogs.Logged — Typestruct Logged{T,N,A,S<:AbstractTimeseriesSolution{T,N,A}} <: AbstractTimeseriesSolution{T,N,A}
Logged(sol::AbstractTimeseriesSolution, log::SimulationLog)Logged differential equation solution. All properties of the underlying solution can be accessed with the getproperty or 'dot' accessing. The inner SimulationLog can be accessed with .log. A Logged solution can be used the same way as its underlying solution (i.e. indexed with sol[i] and interpolated with sol(t))
SimulationLogs.SimulationLog — TypeSimulationLog(; value_dict::Dict{Symbol,Vector}=Dict{Symbol,Vector}(), is_active::Bool=false)A log of saved variables from a simulation. Variables can be accessed with dot notation or the getproperty function.
SimulationLogs.activate! — Functionactivate!(log::SimulationLog=GLOBAL_LOG)Activate a SimulationLog. If no log is given, activate the GLOBAL_LOG.
SimulationLogs.deactivate! — Functiondeactivate!(log::SimulationLog=GLOBAL_LOG)Deactivate a SimulationLog. If no log is given, deactivate the GLOBAL_LOG.
SimulationLogs.get_log — Methodget_log(sol; callback=nothing)
get_log(sol, t; callback=nothing)Get the variables logged by @log from an ODESolution. If a callback or callback set was used to change parameters during the simulation, it must be passed in through the keyword callback to obtain correct results. When doing this, be sure to reset any parameters that changed over the course of the simulation to their starting values.
Alternatively, replace your solve call with logged_solve to handle this all automatically.
See also: logged_solve
SimulationLogs.is_active — Methodis_active(log::SimulationLog)See if simulation is currently active (able to be logged to).
SimulationLogs.logged_solve — Methodlogged_solve(prob, args...; kwargs...)Create a Logged ODE solution whose logged variables can be accessed through the .log property.
See also: get_log
SimulationLogs.reset! — Functionreset!(log::SimulationLog=GLOBAL_LOG)Reset the value dictionary of a SimulationLog. If no log is given, the GLOBAL_LOG will be reset
SimulationLogs.scope — Functionscope(sol, varnames)
scope!(sol, varnames)Plot the variables in varnames from an ODESolution sol. If varnames is a vector, the variables will be plotted against time on the same axis. If varnames is a tuple, the variables will be plotted against each other.
SimulationLogs.scope! — Functionscope(sol, varnames)
scope!(sol, varnames)Plot the variables in varnames from an ODESolution sol. If varnames is a vector, the variables will be plotted against time on the same axis. If varnames is a tuple, the variables will be plotted against each other.
SimulationLogs.value_dict — Methodvalue_dict(log::SimulationLog)Get the value dictionary of a SimulationLog.
SimulationLogs.@log — Macro@log var_name = val
@log var_name exprLog var_name to the GLOBAL_LOG. If @log is placed on an assignment for var_name, the variable will also be created in the local scope. If @log is placed with a variable name before an expression, the expression will run in the local scope without creating the variable in that scope.
Example
function lorenz!(du, u, p, t)
@log a = u[2]-u[1]
@log b u[3]+a
du[1] = p[1]*a
du[2] = u[1]*(p[2]-u[3]) - u[2]
du[3] = u[1]*u[2] - p[3]*u[3]
endIn this example, a is evaluated into the scope of the lorenz! function and is able to be used within that scope. b, however is not evaluated into the lorenz! scope, so no b variable is created in that scope. Both a and b will be logged.