SDESystem
System Constructors
ModelingToolkit.SDESystem
— Typestruct SDESystem <: ModelingToolkit.AbstractODESystem
A system of stochastic differential equations.
Fields
tag
: A tag for the system. If two systems have the same tag, then they are structurally identical.
eqs
: The expressions defining the drift term.noiseeqs
: The expressions defining the diffusion term.iv
: Independent variable.unknowns
: Dependent variables. Must not contain the independent variable.ps
: Parameter variables. Must not contain the independent variable.tspan
: Time span.var_to_name
: Array variables.ctrls
: Control parameters (some subset ofps
).observed
: Observed equations.tgrad
: Time-derivative matrix. Note: this field will not be defined untilcalculate_tgrad
is called on the system.
jac
: Jacobian matrix. Note: this field will not be defined untilcalculate_jacobian
is called on the system.
ctrl_jac
: Control Jacobian matrix. Note: this field will not be defined untilcalculate_control_jacobian
is called on the system.
Wfact
: Note: this field will not be defined untilgenerate_factorized_W
is called on the system.
Wfact_t
: Note: this field will not be defined untilgenerate_factorized_W
is called on the system.
name
: The name of the system.
description
: A description of the system.
systems
: The internal systems. These are required to have unique names.
defaults
: The default values to use when initial conditions and/or parameters are not supplied inODEProblem
.
guesses
: The guesses to use as the initial conditions for the initialization system.
initializesystem
: The system for performing the initialization.
initialization_eqs
: Extra equations to be enforced during the initialization sequence.
connector_type
: Type of the system.
continuous_events
: AVector{SymbolicContinuousCallback}
that model events. The integrator will use root finding to guarantee that it steps at each zero crossing.
discrete_events
: AVector{SymbolicDiscreteCallback}
that models events. Symbolic analog toSciMLBase.DiscreteCallback
that executes an affect when a given condition is true at the end of an integration step.
parameter_dependencies
: Topologically sorted parameter dependency equations, where all symbols are parameters and the LHS is a single parameter.
assertions
: Mapping of conditions which should be true throughout the solution process to corresponding error messages. These will be added to the equations when callingdebug_system
.
metadata
: Metadata for the system, to be used by downstream packages.
gui_metadata
: Metadata for MTK GUI.
namespacing
: If false, thensys.x
no longer performs namespacing.
complete
: If true, denotes the model will not be modified any further.
index_cache
: Cached data for fast symbolic indexing.
parent
: The hierarchical parent system before simplification.
is_scalar_noise
: Signal for whether the noise equations should be treated as a scalar process. This should only betrue
whennoiseeqs isa Vector
.
is_dde
: A boolean indicating if the givenODESystem
represents a system of DDEs.
isscheduled
tearing_state
Example
using ModelingToolkit
using ModelingToolkit: t_nounits as t, D_nounits as D
@parameters σ ρ β
@variables x(t) y(t) z(t)
eqs = [D(x) ~ σ*(y-x),
D(y) ~ x*(ρ-z)-y,
D(z) ~ x*y - β*z]
noiseeqs = [0.1*x,
0.1*y,
0.1*z]
@named de = SDESystem(eqs,noiseeqs,t,[x,y,z],[σ,ρ,β]; tspan = (0, 1000.0))
To convert an ODESystem
to an SDESystem
directly:
ode = ODESystem(eqs,t,[x,y,z],[σ,ρ,β])
sde = SDESystem(ode, noiseeqs)
Composition and Accessor Functions
get_eqs(sys)
orequations(sys)
: The equations that define the SDE.get_unknowns(sys)
orunknowns(sys)
: The set of unknowns in the SDE.get_ps(sys)
orparameters(sys)
: The parameters of the SDE.get_iv(sys)
: The independent variable of the SDE.continuous_events(sys)
: The set of continuous events in the SDE.discrete_events(sys)
: The set of discrete events in the SDE.alg_equations(sys)
: The algebraic equations (i.e. that does not contain a differential) that defines the ODE.get_alg_eqs(sys)
: The algebraic equations (i.e. that does not contain a differential) that defines the ODE. Only returns equations of the current-level system.diff_equations(sys)
: The differential equations (i.e. that contain a differential) that defines the ODE.get_diff_eqs(sys)
: The differential equations (i.e. that contain a differential) that defines the ODE. Only returns equations of the current-level system.has_alg_equations(sys)
: Returnstrue
if the ODE contains any algebraic equations (i.e. that does not contain a differential).has_alg_eqs(sys)
: Returnstrue
if the ODE contains any algebraic equations (i.e. that does not contain a differential). Only considers the current-level system.has_diff_equations(sys)
: Returnstrue
if the ODE contains any differential equations (i.e. that does contain a differential).has_diff_eqs(sys)
: Returnstrue
if the ODE contains any differential equations (i.e. that does contain a differential). Only considers the current-level system.
Transformations
ModelingToolkit.structural_simplify
— Functionstructural_simplify(sys; ...)
structural_simplify(
sys,
io;
additional_passes,
simplify,
split,
allow_symbolic,
allow_parameter,
conservative,
fully_determined,
kwargs...
)
Structurally simplify algebraic equations in a system and compute the topological sort of the observed equations in sys
.
Optional Arguments:
- optional argument
io
may take a tuple(inputs, outputs)
. This will convert allinputs
to parameters and allow them to be unconnected, i.e., simplification will allow models wheren_unknowns = n_equations - n_inputs
.
Optional Keyword Arguments:
- When
simplify=true
, thesimplify
function will be applied during the tearing process. allow_symbolic=false
,allow_parameter=true
, andconservative=false
limit the coefficient types during tearing. In particular,conservative=true
limits tearing to only solve for trivial linear systems where the coefficient has the absolute value of $1$.fully_determined=true
controls whether or not an error will be thrown if the number of equations don't match the number of inputs, outputs, and equations.sort_eqs=true
controls whether equations are sorted lexicographically before simplification or not.
Missing docstring for alias_elimination
. Check Documenter's build log for details.
ModelingToolkit.Girsanov_transform
— FunctionGirsanov_transform(sys::SDESystem, u; θ0) -> SDESystem
Measure transformation method that allows for a reduction in the variance of an estimator Exp(g(X_t))
. Input: Original SDE system and symbolic function u(t,x)
with scalar output that defines the adjustable parameters d
in the Girsanov transformation. Optional: initial condition for θ0
. Output: Modified SDESystem with additional component θ_t
and initial value θ0
, as well as the weight θ_t/θ0
as observed equation, such that the estimator Exp(g(X_t)θ_t/θ0)
has a smaller variance.
Reference: Kloeden, P. E., Platen, E., & Schurz, H. (2012). Numerical solution of SDE through computer experiments. Springer Science & Business Media.
Example
using ModelingToolkit
using ModelingToolkit: t_nounits as t, D_nounits as D
@parameters α β
@variables x(t) y(t) z(t)
eqs = [D(x) ~ α*x]
noiseeqs = [β*x]
@named de = SDESystem(eqs,noiseeqs,t,[x],[α,β])
# define u (user choice)
u = x
θ0 = 0.1
g(x) = x[1]^2
demod = ModelingToolkit.Girsanov_transform(de, u; θ0=0.1)
u0modmap = [
x => x0
]
parammap = [
α => 1.5,
β => 1.0
]
probmod = SDEProblem(complete(demod),u0modmap,(0.0,1.0),parammap)
ensemble_probmod = EnsembleProblem(probmod;
output_func = (sol,i) -> (g(sol[x,end])*sol[demod.weight,end],false),
)
simmod = solve(ensemble_probmod,EM(),dt=dt,trajectories=numtraj)
Analyses
Applicable Calculation and Generation Functions
ModelingToolkit.calculate_jacobian
— Functioncalculate_jacobian(sys::AbstractSystem)
Calculate the Jacobian matrix of a system.
Returns a matrix of Num
instances. The result from the first call will be cached in the system object.
ModelingToolkit.calculate_tgrad
— Functioncalculate_tgrad(sys::AbstractTimeDependentSystem)
Calculate the time gradient of a system.
Returns a vector of Num
instances. The result from the first call will be cached in the system object.
ModelingToolkit.calculate_factorized_W
— Functioncalculate_factorized_W(sys::AbstractSystem)
Calculate the factorized W-matrix of a system.
Returns a matrix of Num
instances. The result from the first call will be cached in the system object.
ModelingToolkit.generate_jacobian
— Functiongenerate_jacobian(sys::AbstractSystem, dvs = unknowns(sys), ps = parameters(sys),
expression = Val{true}; sparse = false, kwargs...)
Generates a function for the Jacobian matrix of a system. Extra arguments control the arguments to the internal build_function
call.
ModelingToolkit.generate_tgrad
— Functiongenerate_tgrad(sys::AbstractTimeDependentSystem, dvs = unknowns(sys), ps = parameters(sys),
expression = Val{true}; kwargs...)
Generates a function for the time gradient of a system. Extra arguments control the arguments to the internal build_function
call.
ModelingToolkit.generate_factorized_W
— Functiongenerate_factorized_W(sys::AbstractSystem, dvs = unknowns(sys), ps = parameters(sys),
expression = Val{true}; sparse = false, kwargs...)
Generates a function for the factorized W matrix of a system. Extra arguments control the arguments to the internal build_function
call.
Missing docstring for jacobian_sparsity
. Check Documenter's build log for details.
Problem Constructors
SciMLBase.SDEFunction
— MethodDiffEqBase.SDEFunction{iip}(sys::SDESystem, dvs = sys.unknowns, ps = sys.ps;
version = nothing, tgrad = false, sparse = false,
jac = false, Wfact = false, kwargs...) where {iip}
Create an SDEFunction
from the SDESystem
. The arguments dvs
and ps
are used to set the order of the dependent variable and parameter vectors, respectively.
SciMLBase.SDEProblem
— MethodDiffEqBase.SDEProblem{iip}(sys::SDESystem, u0map, tspan, p = parammap;
version = nothing, tgrad = false,
jac = false, Wfact = false,
checkbounds = false, sparse = false,
sparsenoise = sparse,
skipzeros = true, fillzeros = true,
linenumbers = true, parallel = SerialForm(),
kwargs...)
Generates an SDEProblem from an SDESystem and allows for automatically symbolically calculating numerical enhancements.
Expression Constructors
ModelingToolkit.SDEFunctionExpr
— TypeDiffEqBase.SDEFunctionExpr{iip}(sys::AbstractODESystem, dvs = unknowns(sys),
ps = parameters(sys);
version = nothing, tgrad = false,
jac = false, Wfact = false,
skipzeros = true, fillzeros = true,
sparse = false,
kwargs...) where {iip}
Create a Julia expression for an SDEFunction
from the SDESystem
. The arguments dvs
and ps
are used to set the order of the dependent variable and parameter vectors, respectively.
ModelingToolkit.SDEProblemExpr
— TypeDiffEqBase.SDEProblemExpr{iip}(sys::AbstractODESystem, u0map, tspan,
parammap = DiffEqBase.NullParameters();
version = nothing, tgrad = false,
jac = false, Wfact = false,
checkbounds = false, sparse = false,
linenumbers = true, parallel = SerialForm(),
kwargs...) where {iip}
Generates a Julia expression for constructing an ODEProblem from an ODESystem and allows for automatically symbolically calculating numerical enhancements.