Skip to content

Commit 7a9c17f

Browse files
Merge pull request #3595 from aml5600/andrew/change-of-vars-units
Change of Variables with units
2 parents 8c6b5ad + 89a375a commit 7a9c17f

File tree

2 files changed

+33
-4
lines changed

2 files changed

+33
-4
lines changed

src/systems/diffeqs/basic_transformations.jl

Lines changed: 16 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -126,10 +126,23 @@ function change_independent_variable(
126126
# Set up intermediate and final variables for the transformation
127127
iv1name = nameof(iv1) # e.g. :t
128128
iv2name = nameof(operation(iv2_of_iv1)) # e.g. :u
129-
iv2, = @independent_variables $iv2name # e.g. u
130-
iv1_of_iv2, = GlobalScope.(@variables $iv1name(iv2)) # inverse, e.g. t(u), global because iv1 has no namespacing in sys
131129
D1 = Differential(iv1) # e.g. d/d(t)
132-
div2_of_iv1 = GlobalScope(default_toterm(D1(iv2_of_iv1))) # e.g. uˍt(t)
130+
131+
# construct new terms, e.g:
132+
# iv2 -> u
133+
# iv1_of_iv2 -> t(u), (inverse, global because iv1 has no namespacing in sys)
134+
# div2_of_iv1 -> uˍt(t)
135+
iv2_unit = getmetadata(iv2_of_iv1, VariableUnit, nothing)
136+
if isnothing(iv2_unit)
137+
iv2, = @independent_variables $iv2name
138+
iv1_of_iv2, = GlobalScope.(@variables $iv1name(iv2))
139+
div2_of_iv1 = GlobalScope(default_toterm(D1(iv2_of_iv1)))
140+
else
141+
iv2, = @independent_variables $iv2name [unit = iv2_unit]
142+
iv1_of_iv2, = GlobalScope.(@variables $iv1name(iv2) [unit = get_unit(iv1)])
143+
div2_of_iv1 = GlobalScope(diff2term_with_unit(D1(iv2_of_iv1), iv1))
144+
end
145+
133146
div2_of_iv2 = substitute(div2_of_iv1, iv1 => iv2) # e.g. uˍt(u)
134147
div2_of_iv2_of_iv1 = substitute(div2_of_iv2, iv2 => iv2_of_iv1) # e.g. uˍt(u(t))
135148

test/basic_transformations.jl

Lines changed: 17 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1,4 +1,4 @@
1-
using ModelingToolkit, OrdinaryDiffEq, DataInterpolations, Test
1+
using ModelingToolkit, OrdinaryDiffEq, DataInterpolations, DynamicQuantities, Test
22

33
@independent_variables t
44
D = Differential(t)
@@ -215,3 +215,19 @@ end
215215
M = ODESystem([D(x(t)) ~ x(t - 1)], t; name = :M)
216216
@test_throws "DDE" change_independent_variable(M, x(t))
217217
end
218+
219+
@testset "Change independent variable w/ units (free fall with 2nd order horizontal equation)" begin
220+
@independent_variables t_units [unit = u"s"]
221+
D_units = Differential(t_units)
222+
@variables x(t_units) [unit = u"m"] y(t_units) [unit = u"m"]
223+
@parameters g = 9.81 [unit = u"m * s^-2"] # gravitational acceleration
224+
Mt = ODESystem([D_units(D_units(y)) ~ -g, D_units(D_units(x)) ~ 0], t_units; name = :M) # gives (x, y) as function of t, ...
225+
Mx = change_independent_variable(Mt, x; add_old_diff = true) # ... but we want y as a function of x
226+
Mx = structural_simplify(Mx; allow_symbolic = true)
227+
Dx = Differential(Mx.x)
228+
u0 = [Mx.y => 0.0, Dx(Mx.y) => 1.0, Mx.t_units => 0.0, Mx.xˍt_units => 10.0]
229+
prob = ODEProblem(Mx, u0, (0.0, 20.0), []) # 1 = dy/dx = (dy/dt)/(dx/dt) means equal initial horizontal and vertical velocities
230+
sol = solve(prob, Tsit5(); reltol = 1e-5)
231+
# compare to analytical solution (x(t) = v*t, y(t) = v*t - g*t^2/2)
232+
@test all(isapprox.(sol[Mx.y], sol[Mx.x - g * (Mx.t_units)^2 / 2]; atol = 1e-10))
233+
end

0 commit comments

Comments
 (0)