Skip to content

Commit f41cdf5

Browse files
Merge pull request #3776 from fchen121/iss3707
Transformation function to turn FDEs into ODEs
2 parents 3c1eac3 + 0411e45 commit f41cdf5

File tree

6 files changed

+331
-3
lines changed

6 files changed

+331
-3
lines changed

Project.toml

Lines changed: 3 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -135,6 +135,7 @@ ModelingToolkitStandardLibrary = "2.20"
135135
Moshi = "0.3"
136136
NaNMath = "0.3, 1"
137137
NonlinearSolve = "4.3"
138+
ODEInterfaceDiffEq = "3.13.4"
138139
OffsetArrays = "1"
139140
OrderedCollections = "1"
140141
OrdinaryDiffEq = "6.82.0"
@@ -184,6 +185,7 @@ LinearSolve = "7ed4a6bd-45f5-4d41-b270-4a48e9bafcae"
184185
Logging = "56ddb016-857b-54e1-b83d-db4d58db5568"
185186
ModelingToolkitStandardLibrary = "16a59e39-deab-5bd0-87e4-056b12336739"
186187
NonlinearSolve = "8913a72c-1f9b-4ce2-8d82-65094dcecaec"
188+
ODEInterfaceDiffEq = "09606e27-ecf5-54fc-bb29-004bd9f985bf"
187189
Optimization = "7f7a1694-90dd-40f0-9382-eb1efda571ba"
188190
OptimizationBase = "bca83a33-5cc9-4baa-983d-23429ab6bcbb"
189191
OptimizationMOI = "fd9f6733-72f4-499f-8506-86b2bdd0dea1"
@@ -206,4 +208,4 @@ Sundials = "c3572dad-4567-51f8-b174-8c6c989267f4"
206208
Test = "8dfed614-e22c-5e08-85e1-65c5234f0b40"
207209

208210
[targets]
209-
test = ["AmplNLWriter", "BenchmarkTools", "BoundaryValueDiffEqMIRK", "BoundaryValueDiffEqAscher", "ControlSystemsBase", "DataInterpolations", "DelayDiffEq", "NonlinearSolve", "ForwardDiff", "Ipopt", "Ipopt_jll", "ModelingToolkitStandardLibrary", "Optimization", "OptimizationOptimJL", "OptimizationMOI", "OrdinaryDiffEq", "OrdinaryDiffEqCore", "OrdinaryDiffEqDefault", "REPL", "Random", "ReferenceTests", "SafeTestsets", "StableRNGs", "Statistics", "SteadyStateDiffEq", "Test", "StochasticDiffEq", "Sundials", "StochasticDelayDiffEq", "Pkg", "JET", "OrdinaryDiffEqNonlinearSolve", "Logging", "OptimizationBase", "LinearSolve"]
211+
test = ["AmplNLWriter", "BenchmarkTools", "BoundaryValueDiffEqMIRK", "BoundaryValueDiffEqAscher", "ControlSystemsBase", "DataInterpolations", "DelayDiffEq", "NonlinearSolve", "ForwardDiff", "Ipopt", "Ipopt_jll", "ModelingToolkitStandardLibrary", "Optimization", "OptimizationOptimJL", "OptimizationMOI", "ODEInterfaceDiffEq", "OrdinaryDiffEq", "OrdinaryDiffEqCore", "OrdinaryDiffEqDefault", "REPL", "Random", "ReferenceTests", "SafeTestsets", "StableRNGs", "Statistics", "SteadyStateDiffEq", "Test", "StochasticDiffEq", "Sundials", "StochasticDelayDiffEq", "Pkg", "JET", "OrdinaryDiffEqNonlinearSolve", "Logging", "OptimizationBase", "LinearSolve"]

docs/src/API/model_building.md

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -219,6 +219,8 @@ symbolic analysis.
219219

220220
```@docs
221221
liouville_transform
222+
fractional_to_ordinary
223+
linear_fractional_to_ordinary
222224
change_of_variables
223225
stochastic_integral_transform
224226
Girsanov_transform

src/ModelingToolkit.jl

Lines changed: 2 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -269,7 +269,8 @@ export isinput, isoutput, getbounds, hasbounds, getguess, hasguess, isdisturbanc
269269
subset_tunables
270270
export liouville_transform, change_independent_variable, substitute_component,
271271
add_accumulations, noise_to_brownians, Girsanov_transform, change_of_variables,
272-
respecialize
272+
fractional_to_ordinary, linear_fractional_to_ordinary
273+
export respecialize
273274
export PDESystem
274275
export Differential, expand_derivatives, @derivatives
275276
export Equation, ConstrainedEquation

src/systems/diffeqs/basic_transformations.jl

Lines changed: 236 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -185,6 +185,242 @@ function change_of_variables(
185185
return new_sys
186186
end
187187

188+
"""
189+
Generates the system of ODEs to find solution to FDEs.
190+
191+
Example:
192+
193+
```julia
194+
@independent_variables t
195+
@variables x(t)
196+
D = Differential(t)
197+
tspan = (0., 1.)
198+
199+
α = 0.5
200+
eqs = (9*gamma(1 + α)/4) - (3*t^(4 - α/2)*gamma(5 + α/2)/gamma(5 - α/2))
201+
eqs += (gamma(9)*t^(8 - α)/gamma(9 - α)) + (3/2*t^(α/2)-t^4)^3 - x^(3/2)
202+
sys = fractional_to_ordinary(eqs, x, α, 10^-7, 1)
203+
204+
prob = ODEProblem(sys, [], tspan)
205+
sol = solve(prob, radau5(), abstol = 1e-10, reltol = 1e-10)
206+
```
207+
"""
208+
function fractional_to_ordinary(
209+
eqs, variables, alphas, epsilon, T;
210+
initials = 0, additional_eqs = [], iv = only(@independent_variables t), matrix=false
211+
)
212+
D = Differential(iv)
213+
i = 0
214+
all_eqs = Equation[]
215+
all_def = Pair[]
216+
217+
function fto_helper(sub_eq, sub_var, α; initial=0)
218+
alpha_0 = α
219+
220+
if> 1)
221+
coeff = 1/- 1)
222+
m = 2
223+
while- m > 0)
224+
coeff /= α - m
225+
m += 1
226+
end
227+
alpha_0 = α - m + 1
228+
end
229+
230+
δ = (gamma(alpha_0+1) * epsilon)^(1/alpha_0)
231+
a = pi/2*(1-(1-alpha_0)/((2-alpha_0) * log(epsilon^-1)))
232+
h = 2*pi*a / log(1 + (2/epsilon * (cos(a))^(alpha_0 - 1)))
233+
234+
x_sub = (gamma(2-alpha_0) * epsilon)^(1/(1-alpha_0))
235+
x_sup = -log(gamma(1-alpha_0) * epsilon)
236+
M = floor(Int, log(x_sub / T) / h)
237+
N = ceil(Int, log(x_sup / δ) / h)
238+
239+
function c_i(index)
240+
h * sin(pi * alpha_0) / pi * exp((1-alpha_0)*h*index)
241+
end
242+
243+
function γ_i(index)
244+
exp(h * index)
245+
end
246+
247+
new_eqs = Equation[]
248+
def = Pair[]
249+
250+
if matrix
251+
new_z = Symbol(, :_, i)
252+
i += 1
253+
γs = diagm([γ_i(index) for index in M:N-1])
254+
cs = [c_i(index) for index in M:N-1]
255+
256+
if< 1)
257+
new_z = only(@variables $new_z(iv)[1:N-M])
258+
new_eq = D(new_z) ~ -γs*new_z .+ sub_eq
259+
rhs = dot(cs, new_z) + initial
260+
push!(def, new_z=>zeros(N-M))
261+
else
262+
new_z = only(@variables $new_z(iv)[1:N-M, 1:m])
263+
new_eq = D(new_z) ~ -γs*new_z + hcat(fill(sub_eq, N-M, 1), collect(new_z[:, 1:m-1]*diagm(1:m-1)))
264+
rhs = coeff*sum(cs[i]*new_z[i, m] for i in 1:N-M)
265+
for (index, value) in enumerate(initial)
266+
rhs += value * iv^(index - 1) / gamma(index)
267+
end
268+
push!(def, new_z=>zeros(N-M, m))
269+
end
270+
push!(new_eqs, new_eq)
271+
else
272+
if< 1)
273+
rhs = initial
274+
for index in range(M, N-1; step=1)
275+
new_z = Symbol(, :_, i)
276+
i += 1
277+
new_z = ModelingToolkit.unwrap(only(@variables $new_z(iv)))
278+
new_eq = D(new_z) ~ sub_eq - γ_i(index)*new_z
279+
push!(new_eqs, new_eq)
280+
push!(def, new_z=>0)
281+
rhs += c_i(index)*new_z
282+
end
283+
else
284+
rhs = 0
285+
for (index, value) in enumerate(initial)
286+
rhs += value * iv^(index - 1) / gamma(index)
287+
end
288+
for index in range(M, N-1; step=1)
289+
new_z = Symbol(, :_, i)
290+
i += 1
291+
γ = γ_i(index)
292+
base = sub_eq
293+
for k in range(1, m; step=1)
294+
new_z = Symbol(, :_, index-M, :_, k)
295+
new_z = ModelingToolkit.unwrap(only(@variables $new_z(iv)))
296+
new_eq = D(new_z) ~ base - γ*new_z
297+
base = k * new_z
298+
push!(new_eqs, new_eq)
299+
push!(def, new_z=>0)
300+
end
301+
rhs += coeff*c_i(index)*new_z
302+
end
303+
end
304+
end
305+
push!(new_eqs, sub_var ~ rhs)
306+
return (new_eqs, def)
307+
end
308+
309+
for (eq, cur_var, alpha, init) in zip(eqs, variables, alphas, initials)
310+
(new_eqs, def) = fto_helper(eq, cur_var, alpha; initial=init)
311+
append!(all_eqs, new_eqs)
312+
append!(all_def, def)
313+
end
314+
append!(all_eqs, additional_eqs)
315+
@named sys = System(all_eqs, iv; defaults=all_def)
316+
return mtkcompile(sys)
317+
end
318+
319+
"""
320+
Generates the system of ODEs to find solution to FDEs.
321+
322+
Example:
323+
324+
```julia
325+
@independent_variables t
326+
@variables x_0(t)
327+
D = Differential(t)
328+
tspan = (0., 5000.)
329+
330+
function expect(t)
331+
return sqrt(2) * sin(t + pi/4)
332+
end
333+
334+
sys = linear_fractional_to_ordinary([3, 2.5, 2, 1, .5, 0], [1, 1, 1, 4, 1, 4], 6*cos(t), 10^-5, 5000; initials=[1, 1, -1])
335+
prob = ODEProblem(sys, [], tspan)
336+
sol = solve(prob, radau5(), abstol = 1e-5, reltol = 1e-5)
337+
```
338+
"""
339+
function linear_fractional_to_ordinary(
340+
degrees, coeffs, rhs, epsilon, T;
341+
initials = 0, symbol = :x, iv = only(@independent_variables t), matrix=false
342+
)
343+
previous = Symbol(symbol, :_, 0)
344+
previous = ModelingToolkit.unwrap(only(@variables $previous(iv)))
345+
@variables x_0(iv)
346+
D = Differential(iv)
347+
i = 0
348+
all_eqs = Equation[]
349+
all_def = Pair[]
350+
351+
function fto_helper(sub_eq, α)
352+
δ = (gamma+1) * epsilon)^(1/α)
353+
a = pi/2*(1-(1-α)/((2-α) * log(epsilon^-1)))
354+
h = 2*pi*a / log(1 + (2/epsilon * (cos(a))^- 1)))
355+
356+
x_sub = (gamma(2-α) * epsilon)^(1/(1-α))
357+
x_sup = -log(gamma(1-α) * epsilon)
358+
M = floor(Int, log(x_sub / T) / h)
359+
N = ceil(Int, log(x_sup / δ) / h)
360+
361+
function c_i(index)
362+
h * sin(pi * α) / pi * exp((1-α)*h*index)
363+
end
364+
365+
function γ_i(index)
366+
exp(h * index)
367+
end
368+
369+
new_eqs = Equation[]
370+
def = Pair[]
371+
if matrix
372+
new_z = Symbol(, :_, i)
373+
i += 1
374+
γs = diagm([γ_i(index) for index in M:N-1])
375+
cs = [c_i(index) for index in M:N-1]
376+
377+
new_z = only(@variables $new_z(iv)[1:N-M])
378+
new_eq = D(new_z) ~ -γs*new_z .+ sub_eq
379+
sum = dot(cs, new_z)
380+
push!(def, new_z=>zeros(N-M))
381+
push!(new_eqs, new_eq)
382+
else
383+
sum = 0
384+
for index in range(M, N-1; step=1)
385+
new_z = Symbol(, :_, i)
386+
i += 1
387+
new_z = ModelingToolkit.unwrap(only(@variables $new_z(iv)))
388+
new_eq = D(new_z) ~ sub_eq - γ_i(index)*new_z
389+
push!(new_eqs, new_eq)
390+
push!(def, new_z=>0)
391+
sum += c_i(index)*new_z
392+
end
393+
end
394+
return (new_eqs, def, sum)
395+
end
396+
397+
for i in range(1, ceil(Int, degrees[1]); step=1)
398+
new_x = Symbol(symbol, :_, i)
399+
new_x = ModelingToolkit.unwrap(only(@variables $new_x(iv)))
400+
push!(all_eqs, D(previous) ~ new_x)
401+
push!(all_def, previous => initials[i])
402+
previous = new_x
403+
end
404+
405+
new_rhs = -rhs
406+
for (degree, coeff) in zip(degrees, coeffs)
407+
rounded = ceil(Int, degree)
408+
new_x = Symbol(symbol, :_, rounded)
409+
new_x = ModelingToolkit.unwrap(only(@variables $new_x(iv)))
410+
if isinteger(degree)
411+
new_rhs += coeff * new_x
412+
else
413+
(new_eqs, def, sum) = fto_helper(new_x, rounded - degree)
414+
append!(all_eqs, new_eqs)
415+
append!(all_def, def)
416+
new_rhs += coeff * sum
417+
end
418+
end
419+
push!(all_eqs, 0 ~ new_rhs)
420+
@named sys = System(all_eqs, iv; defaults=all_def)
421+
return mtkcompile(sys)
422+
end
423+
188424
"""
189425
change_independent_variable(
190426
sys::System, iv, eqs = [];

test/fractional_to_ordinary.jl

Lines changed: 86 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,86 @@
1+
using ModelingToolkit, OrdinaryDiffEq, ODEInterfaceDiffEq, SpecialFunctions, LinearAlgebra
2+
using Test
3+
4+
# Testing for α < 1
5+
# Uses example 1 from Section 7 of https://arxiv.org/pdf/2506.04188
6+
@independent_variables t
7+
@variables x(t)
8+
D = Differential(t)
9+
tspan = (0., 1.)
10+
timepoint = [0., 0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1.]
11+
12+
function expect(t, α)
13+
return (3/2*t^/2) - t^4)^2
14+
end
15+
16+
α = 0.5
17+
eqs = (9*gamma(1 + α)/4) - (3*t^(4 - α/2)*gamma(5 + α/2)/gamma(5 - α/2))
18+
eqs += (gamma(9)*t^(8 - α)/gamma(9 - α)) + (3/2*t^/2)-t^4)^3 - x^(3/2)
19+
sys = fractional_to_ordinary(eqs, x, α, 10^-7, 1)
20+
21+
prob = ODEProblem(sys, [], tspan)
22+
sol = solve(prob, radau5(), saveat=timepoint, abstol = 1e-10, reltol = 1e-10)
23+
24+
for time in 0:0.1:1
25+
@test isapprox(expect(time, α), sol(time, idxs=x), atol=1e-7)
26+
time += 0.1
27+
end
28+
29+
α = 0.3
30+
eqs = (9*gamma(1 + α)/4) - (3*t^(4 - α/2)*gamma(5 + α/2)/gamma(5 - α/2))
31+
eqs += (gamma(9)*t^(8 - α)/gamma(9 - α)) + (3/2*t^/2)-t^4)^3 - x^(3/2)
32+
sys = fractional_to_ordinary(eqs, x, α, 10^-7, 1; matrix=true)
33+
34+
prob = ODEProblem(sys, [], tspan)
35+
sol = solve(prob, radau5(), saveat=timepoint, abstol = 1e-10, reltol = 1e-10)
36+
37+
for time in 0:0.1:1
38+
@test isapprox(expect(time, α), sol(time, idxs=x), atol=1e-7)
39+
end
40+
41+
α = 0.9
42+
eqs = (9*gamma(1 + α)/4) - (3*t^(4 - α/2)*gamma(5 + α/2)/gamma(5 - α/2))
43+
eqs += (gamma(9)*t^(8 - α)/gamma(9 - α)) + (3/2*t^/2)-t^4)^3 - x^(3/2)
44+
sys = fractional_to_ordinary(eqs, x, α, 10^-7, 1)
45+
46+
prob = ODEProblem(sys, [], tspan)
47+
sol = solve(prob, radau5(), saveat=timepoint, abstol = 1e-10, reltol = 1e-10)
48+
49+
for time in 0:0.1:1
50+
@test isapprox(expect(time, α), sol(time, idxs=x), atol=1e-7)
51+
end
52+
53+
# Testing for example 2 of Section 7
54+
@independent_variables t
55+
@variables x(t) y(t)
56+
D = Differential(t)
57+
tspan = (0., 220.)
58+
59+
sys = fractional_to_ordinary([1 - 4*x + x^2 * y, 3*x - x^2 * y], [x, y], [1.3, 0.8], 10^-8, 220; initials=[[1.2, 1], 2.8], matrix=true)
60+
prob = ODEProblem(sys, [], tspan)
61+
sol = solve(prob, radau5(), abstol = 1e-8, reltol = 1e-8)
62+
63+
@test isapprox(1.0097684171, sol(220, idxs=x), atol=1e-5)
64+
@test isapprox(2.1581264031, sol(220, idxs=y), atol=1e-5)
65+
66+
#Testing for example 3 of Section 7
67+
@independent_variables t
68+
@variables x_0(t)
69+
D = Differential(t)
70+
tspan = (0., 5000.)
71+
72+
function expect(t)
73+
return sqrt(2) * sin(t + pi/4)
74+
end
75+
76+
sys = linear_fractional_to_ordinary([3, 2.5, 2, 1, .5, 0], [1, 1, 1, 4, 1, 4], 6*cos(t), 10^-5, 5000; initials=[1, 1, -1])
77+
prob = ODEProblem(sys, [], tspan)
78+
sol = solve(prob, radau5(), abstol = 1e-5, reltol = 1e-5)
79+
80+
@test isapprox(expect(5000), sol(5000, idxs=x_0), atol=1e-5)
81+
82+
msys = linear_fractional_to_ordinary([3, 2.5, 2, 1, .5, 0], [1, 1, 1, 4, 1, 4], 6*cos(t), 10^-5, 5000; initials=[1, 1, -1], matrix=true)
83+
mprob = ODEProblem(sys, [], tspan)
84+
msol = solve(prob, radau5(), abstol = 1e-5, reltol = 1e-5)
85+
86+
@test isapprox(expect(5000), msol(5000, idxs=x_0), atol=1e-5)

test/runtests.jl

Lines changed: 2 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -100,9 +100,10 @@ end
100100
@safetestset "Subsystem replacement" include("substitute_component.jl")
101101
@safetestset "Linearization Tests" include("linearize.jl")
102102
@safetestset "LinearProblem Tests" include("linearproblem.jl")
103+
@safetestset "Fractional Differential Equations Tests" include("fractional_to_ordinary.jl")
103104
end
104105
end
105-
106+
106107
if GROUP == "All" || GROUP == "SymbolicIndexingInterface"
107108
@safetestset "SymbolicIndexingInterface test" include("symbolic_indexing_interface.jl")
108109
@safetestset "SciML Problem Input Test" include("sciml_problem_inputs.jl")

0 commit comments

Comments
 (0)