Extension to the SciML ecosystem
The SciML ecosystem provides a rich set of tools to solve (non)-linear equations, differential equations and inverse problems. We provide an interface (in the form of Package extensions) to export the the derived harmonic equations computed with harmonic balance method or krylov-bogoliubov method to the SciML ecosystem.
ModeligToolkit.jl
The ModelingToolkit.jl
(MTK) package provides a symbolic framework for defining and simplifying mathematical models. Through, MTK SciML provides a symbolic interface for their ecosystem
SciMLBase.ODEProblem Method
ODEProblem(
eom::Union{DifferentialEquation, HarmonicEquation},
u0,
tspan::Tuple,
p::AbstractDict;
in_place,
kwargs...
) -> Any
Creates and ModelingToolkit.ODEProblem from a DifferentialEquation or HarmonicEquation.
Example
using ModelingToolkit, StaticArrays
@variables α ω ω0 F γ t x(t)
diff_eq = DifferentialEquation(
d(x, t, 2) + ω0^2 * x + α * x^3 + γ * d(x, t) ~ F * cos(ω * t), x
)
add_harmonic!(diff_eq, x, ω) #
harmonic_eq = get_harmonic_equations(diff_eq)
# in place (most performant for large systems)
ODEProblem(harmonic_eq, [1.0, 0.0], (0, 100), param)
# out of place (most performant for small systems with StaticArrays)
ODEProblem(
harmonic_eq, [1.0, 0.0], (0, 100), param;
in_place=false, u0_constructor=x -> SVector(x...)
)
ModelingToolkit.ODESystem Type
ODESystem(eom::HarmonicEquation) -> Any
Creates and ModelingToolkit.ODESystem from a HarmonicEquation.
Example
using ModelingToolkit
@variables α ω ω0 F γ t x(t)
diff_eq = DifferentialEquation(
d(x, t, 2) + ω0^2 * x + α * x^3 + γ * d(x, t) ~ F * cos(ω * t), x
)
add_harmonic!(diff_eq, x, ω) #
harmonic_eq = get_harmonic_equations(diff_eq)
sys = ODESystem(harmonic_eq)
param = (α => 1.0, ω0 => 1.1, F => 0.01, γ => 0.01, ω => 1.1)
ODEProblem(sys, [1.0, 0.0], (0, 100), param)
ODESystem(diff_eq::DifferentialEquation) -> Any
Creates and ModelingToolkit.ODESystem from a DifferentialEquation.
Example
using ModelingToolkit
@variables α ω ω0 F γ t x(t)
diff_eq = DifferentialEquation(
d(x, t, 2) + ω0^2 * x + α * x^3 + γ * d(x, t) ~ F * cos(ω * t), x
)
sys = ODESystem(diff_eq)
param = (α => 1.0, ω0 => 1.1, F => 0.01, γ => 0.01, ω => 1.1)
ODEProblem(sys, [1.0, 0.0], (0, 100), param)
SciMLBase.SteadyStateProblem Type
SteadyStateProblem(
eom::HarmonicEquation,
u0,
p::AbstractDict;
in_place,
kwargs...
) -> Any
Creates and ModelingToolkit.SteadyStateProblem from a HarmonicEquation.
Example
using ModelingToolkit, StaticArrays
@variables α ω ω0 F γ t x(t)
diff_eq = DifferentialEquation(
d(x, t, 2) + ω0^2 * x + α * x^3 + γ * d(x, t) ~ F * cos(ω * t), x
)
add_harmonic!(diff_eq, x, ω) #
harmonic_eq = get_harmonic_equations(diff_eq)
SteadyStateProblem(harmonic_eq, [1.0, 0.0], param)
SciMLBase.NonlinearProblem Type
NonlinearProblem(
eom::HarmonicEquation,
u0,
p::AbstractDict;
in_place,
kwargs...
) -> Any
Creates and ModelingToolkit.NonlinearProblem from a HarmonicEquation.
Example
using ModelingToolkit, StaticArrays
@variables α ω ω0 F γ t x(t)
diff_eq = DifferentialEquation(
d(x, t, 2) + ω0^2 * x + α * x^3 + γ * d(x, t) ~ F * cos(ω * t), x
)
add_harmonic!(diff_eq, x, ω) #
harmonic_eq = get_harmonic_equations(diff_eq)
NonlinearProblem(harmonic_eq, [1.0, 0.0], param)