Symbolic ODE Solving
While not all ODEs have an analytical solution, symbolic ODE solving is provided by Symbolics.jl for subsets of cases for which known analytical solutions can be obtained. These expressions can then be merged with other techniques in order to accelerate code or gain a deeper understanding of real-world systems.
Merging Symbolic ODEs with Numerical Methods: ModelingToolkit.jl
If you are looking to merge Symbolics.jl manipulations with numerical solvers such as DifferentialEquations.jl, then we highly recommend checking out the ModelingToolkit.jl system. It represents systems of differentible-algebraic equations (DAEs, and extension to ODEs) for which a sophisticated symbolic analysis pipeline is used to generate highly efficient code. This mtkcompile
pipeline makes use of all tricks from Symbolics.jl and many that are more specific to numerical code generation, such as Pantelides index reduction and tearing of nonlinear systems, and it will analytically solve subsets of the ODE system as finds possible. Thus if you are attempting to use Symbolics to pre-solve some parts of an ODE analytically, we recommend allowing ModelingToolkit.jl to do this optimization.
Note that if ModelingToolkit is able to analytically solve the equation, it will give an ODEProblem
where prob.u0 === nothing
, and then running solve
on the ODEProblem
will give a numerical ODESolution
object that on-demand uses the analytical solution to generate any plots or other artifacts. The analytical solution can be investigated symbolically using observed(sys)
.
Symbolically Solving ODEs
This area is currently under heavy development. More solvers will be available in the near future.
Symbolics.LinearODE
— TypeRepresents a linear ordinary differential equation of the form:
dⁿx/dtⁿ + pₙ(t)(dⁿ⁻¹x/dtⁿ⁻¹) + ... + p₂(t)(dx/dt) + p₁(t)x = q(t)
Fields
x
: dependent variablet
: independent variablep
: coefficient functions oft
ordered in increasing order (p₁, p₂, ...)q
: right hand side function oft
, without anyx
Examples
julia> using Symbolics
julia> @variables x, t
2-element Vector{Num}:
x
t
julia> eq = LinearODE(x, t, [1, 2, 3], 3exp(4t))
(Dt^3)x + (3)(Dt^2)x + (2)(Dt^1)x + (1)(Dt^0)x ~ 3exp(4t)
Symbolics.symbolic_solve_ode
— Functionsymbolic_solve_ode(eq::LinearODE)
Symbolically solve a linear ordinary differential equation
Arguments
- eq: a
LinearODE
to solve
Returns
Symbolic solution to the ODE
Supported Methods
- first-order integrating factor
- constant coefficient homogeneous solutions (can handle repeated and complex characteristic roots)
- exponential and resonant response formula particular solutions (for any linear combination of
exp
,sin
,cos
, orexp
timessin
orcos
(e.g.e^2t * cos(-t) + e^-3t + sin(5t))
) - method of undetermined coefficients particular solutions
- linear combinations of above particular solutions
Examples
julia> using Symbolics; import Nemo, SymPy
julia> @variables x, t
2-element Vector{Num}:
x
t
# Integrating Factor (note that SymPy is required for integration)
julia> symbolic_solve_ode(LinearODE(x, t, [5/t], 7t))
(C₁ + t^7) / (t^5)
# Constant Coefficients and RRF (note that Nemo is required to find characteristic roots)
julia> symbolic_solve_ode(LinearODE(x, t, [9, -6], 4exp(3t)))
C₁*exp(3t) + C₂*t*exp(3t) + (2//1)*(t^2)*exp(3t)
julia> symbolic_solve_ode(LinearODE(x, t, [6, 5], 2exp(-t)*cos(t)))
C₁*exp(-2t) + C₂*exp(-3t) + (1//5)*cos(t)*exp(-t) + (3//5)*exp(-t)*sin(t)
# Method of Undetermined Coefficients
julia> symbolic_solve_ode(LinearODE(x, t, [-3, 2], 2t - 5))
(11//9) - (2//3)*t + C₁*exp(t) + C₂*exp(-3t)
symbolic_solve_ode(expr::Equation, x, t)
Symbolically solve an ODE
Arguments
- expr: a symbolic ODE
- x: dependent variable
- t: independent variable
Supported Methods
- all methods of solving linear ODEs mentioned for
symbolic_solve_ode(eq::LinearODE)
- Clairaut's equation
- Bernoulli equations
Examples
julia> using Symbolics; import Nemo
julia> @variables x, t
2-element Vector{Num}:
x
t
julia> Dt = Differential(t)
Differential(t)
# LinearODE (via constant coefficients and RRF)
julia> symbolic_solve_ode(9t*x - 6*Dt(x) ~ 4exp(3t), x, t)
C₁*exp(3t) + C₂*t*exp(3t) + (2//1)*(t^2)*exp(3t)
# Clairaut's equation
julia> symbolic_solve_ode(x ~ Dt(x)*t - ((Dt(x))^3), x, t)
C₁*t - (C₁^3)
# Bernoulli equations
julia> symbolic_solve_ode(Dt(x) + (4//t)*x ~ t^3 * x^2, x, t)
1 / (C₁*(t^4) - (t^4)*log(t))
Continuous Dynamical Systems
Symbolics.solve_linear_system
— Functionsolve_linear_system(A::Matrix{<:Number}, x0::Vector{<:Number}, t::Num)
Solve linear continuous dynamical system of differential equations of the form Ax = x' with initial condition x0
Arguments
A
: matrix of coefficientsx0
: initial conditions vectort
: independent variable
Returns
vector of symbolic solutions
Examples
!!! note uses method symbolic_solve
, so packages Nemo
and Groebner
are often required
julia> @variables t
1-element Vector{Num}:
t
julia> solve_linear_system([1 0; 0 -1], [1, -1], t) # requires Nemo
2-element Vector{Num}:
exp(t)
-exp(-t)
julia> solve_linear_system([-3 4; -2 3], [7, 2], t) # requires Groebner
2-element Vector{Num}:
(10//1)*exp(-t) - (3//1)*exp(t)
(5//1)*exp(-t) - (3//1)*exp(t)
SymPy
Symbolics.sympy_ode_solve
— Functionsympy_ode_solve(expr, func, var)
Solves ODE expr = 0 for function func w.r.t. var using SymPy.
Arguments
expr
: Symbolics expression representing ODE (set to 0).func
: Symbolics function (e.g., f(x)).var
: Independent Symbolics variable.
Returns
Symbolics solution(s).
Example
@variables x
@syms f(x)
expr = Symbolics.Derivative(f, x) - 2*f
sol = sympy_ode_solve(expr, f, x) # Returns C1*exp(2*x)