Skip to content

fix: use Diagonal for diagonal mass matrices #3615

New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Merged
merged 1 commit into from
May 9, 2025
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
5 changes: 5 additions & 0 deletions src/systems/diffeqs/abstractodesystem.jl
Original file line number Diff line number Diff line change
Expand Up @@ -284,6 +284,9 @@ function calculate_massmatrix(sys::AbstractODESystem; simplify = false)
end
M = simplify ? ModelingToolkit.simplify.(M) : M
# M should only contain concrete numbers
if isdiag(M)
M = Diagonal(M)
end
Comment on lines +287 to +289
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

It should always be diagonal?

Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Technically users can manually create an ODESystem with non-diagonal mass matrix, complete and codegen it.

Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I see

M == I ? I : M
end

Expand Down Expand Up @@ -410,6 +413,8 @@ function DiffEqBase.ODEFunction{iip, specialize}(sys::AbstractODESystem,
SparseArrays.sparse(M)
elseif u0 === nothing || M === I
M
elseif M isa Diagonal
Diagonal(ArrayInterface.restructure(u0, diag(M)))
else
ArrayInterface.restructure(u0 .* u0', M)
end
Expand Down
10 changes: 9 additions & 1 deletion src/systems/diffeqs/sdesystem.jl
Original file line number Diff line number Diff line change
Expand Up @@ -652,7 +652,15 @@ function DiffEqBase.SDEFunction{iip, specialize}(sys::SDESystem, dvs = unknowns(
W_prototype = nothing
end

_M = (u0 === nothing || M == I) ? M : ArrayInterface.restructure(u0 .* u0', M)
_M = if sparse && !(u0 === nothing || M === I)
SparseArrays.sparse(M)
elseif u0 === nothing || M === I
M
elseif M isa Diagonal
Diagonal(ArrayInterface.restructure(u0, diag(M)))
else
ArrayInterface.restructure(u0 .* u0', M)
end

observedfun = ObservedFunctionCache(
sys; eval_expression, eval_module, checkbounds = get(kwargs, :checkbounds, false), cse)
Expand Down
21 changes: 20 additions & 1 deletion test/mass_matrix.jl
Original file line number Diff line number Diff line change
@@ -1,4 +1,4 @@
using OrdinaryDiffEq, ModelingToolkit, Test, LinearAlgebra
using OrdinaryDiffEq, ModelingToolkit, Test, LinearAlgebra, StaticArrays
using ModelingToolkit: t_nounits as t, D_nounits as D, MTKParameters

@variables y(t)[1:3]
Expand All @@ -12,13 +12,18 @@ eqs = [D(y[1]) ~ -k[1] * y[1] + k[3] * y[2] * y[3],
sys = complete(sys)
@test_throws ArgumentError ODESystem(eqs, y[1])
M = calculate_massmatrix(sys)
@test M isa Diagonal
@test M == [1 0 0
0 1 0
0 0 0]

prob_mm = ODEProblem(sys, [y => [1.0, 0.0, 0.0]], (0.0, 1e5),
[k => [0.04, 3e7, 1e4]])
@test prob_mm.f.mass_matrix isa Diagonal{Float64, Vector{Float64}}
sol = solve(prob_mm, Rodas5(), reltol = 1e-8, abstol = 1e-8)
prob_mm = ODEProblem(sys, SA[y => [1.0, 0.0, 0.0]], (0.0, 1e5),
[k => [0.04, 3e7, 1e4]])
@test prob_mm.f.mass_matrix isa Diagonal{Float64, SVector{3, Float64}}

function rober(du, u, p, t)
y₁, y₂, y₃ = u
Expand All @@ -43,3 +48,17 @@ eqs = [D(y[1]) ~ y[1], D(y[2]) ~ y[2], D(y[3]) ~ y[3]]
@named sys = ODESystem(eqs, t, collect(y), [k])

@test calculate_massmatrix(sys) === I

@testset "Mass matrix `isa Diagonal` for `SDEProblem`" begin
eqs = [D(y[1]) ~ -k[1] * y[1] + k[3] * y[2] * y[3],
D(y[2]) ~ k[1] * y[1] - k[3] * y[2] * y[3] - k[2] * y[2]^2,
0 ~ y[1] + y[2] + y[3] - 1]

@named sys = ODESystem(eqs, t, collect(y), [k])
@named sys = SDESystem(sys, [1, 1, 0])
sys = complete(sys)
prob = SDEProblem(sys, [y => [1.0, 0.0, 0.0]], (0.0, 1e5), [k => [0.04, 3e7, 1e4]])
@test prob.f.mass_matrix isa Diagonal{Float64, Vector{Float64}}
prob = SDEProblem(sys, SA[y => [1.0, 0.0, 0.0]], (0.0, 1e5), [k => [0.04, 3e7, 1e4]])
@test prob.f.mass_matrix isa Diagonal{Float64, SVector{3, Float64}}
end
Loading