Skip to content

Commit e1c34ec

Browse files
authored
Merge pull request #2218 from SciML/myb/vis
Make introspecting `dummy_derivative` easier
2 parents 4aa29a0 + fd58156 commit e1c34ec

File tree

5 files changed

+15
-10
lines changed

5 files changed

+15
-10
lines changed

src/structural_transformation/bipartite_tearing/modia_tearing.jl

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -82,6 +82,7 @@ function tear_graph_modia(structure::SystemStructure, isder::F = nothing,
8282
var_eq_matching = complete(var_eq_matching,
8383
max(length(var_eq_matching),
8484
maximum(x -> x isa Int ? x : 0, var_eq_matching)))
85+
full_var_eq_matching = copy(var_eq_matching)
8586
var_sccs::Vector{Union{Vector{Int}, Int}} = find_var_sccs(graph, var_eq_matching)
8687
vargraph = DiCMOBiGraph{true}(graph)
8788
ict = IncrementalCycleTracker(vargraph; dir = :in)
@@ -110,6 +111,5 @@ function tear_graph_modia(structure::SystemStructure, isder::F = nothing,
110111
empty!(ieqs)
111112
empty!(filtered_vars)
112113
end
113-
114-
return var_eq_matching
114+
return var_eq_matching, full_var_eq_matching
115115
end

src/structural_transformation/partial_state_selection.jl

Lines changed: 10 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -168,15 +168,15 @@ function partial_state_selection_graph!(structure::SystemStructure, var_eq_match
168168
end
169169

170170
function dummy_derivative_graph!(state::TransformationState, jac = nothing;
171-
state_priority = nothing, kwargs...)
171+
state_priority = nothing, log = Val(false), kwargs...)
172172
state.structure.solvable_graph === nothing && find_solvables!(state; kwargs...)
173173
complete!(state.structure)
174174
var_eq_matching = complete(pantelides!(state))
175-
dummy_derivative_graph!(state.structure, var_eq_matching, jac, state_priority)
175+
dummy_derivative_graph!(state.structure, var_eq_matching, jac, state_priority, log)
176176
end
177177

178178
function dummy_derivative_graph!(structure::SystemStructure, var_eq_matching, jac,
179-
state_priority)
179+
state_priority, ::Val{log} = Val(false)) where {log}
180180
@unpack eq_to_diff, var_to_diff, graph = structure
181181
diff_to_eq = invview(eq_to_diff)
182182
diff_to_var = invview(var_to_diff)
@@ -338,7 +338,7 @@ function dummy_derivative_graph!(structure::SystemStructure, var_eq_matching, ja
338338
end
339339
end
340340

341-
var_eq_matching = tear_graph_modia(structure, isdiffed,
341+
var_eq_matching, full_var_eq_matching = tear_graph_modia(structure, isdiffed,
342342
Union{Unassigned, SelectedState};
343343
varfilter = can_eliminate)
344344
for v in eachindex(var_eq_matching)
@@ -348,5 +348,10 @@ function dummy_derivative_graph!(structure::SystemStructure, var_eq_matching, ja
348348
var_eq_matching[v] = SelectedState()
349349
end
350350

351-
return var_eq_matching
351+
if log
352+
candidates = can_eliminate.(1:ndsts(graph))
353+
return var_eq_matching, full_var_eq_matching, candidates
354+
else
355+
return var_eq_matching
356+
end
352357
end

src/structural_transformation/symbolics_tearing.jl

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -564,7 +564,7 @@ function tearing(state::TearingState; kwargs...)
564564
@unpack graph = state.structure
565565
algvars = BitSet(findall(v -> isalgvar(state.structure, v), 1:ndsts(graph)))
566566
aeqs = algeqs(state.structure)
567-
var_eq_matching′ = tear_graph_modia(state.structure;
567+
var_eq_matching′, = tear_graph_modia(state.structure;
568568
varfilter = var -> var in algvars,
569569
eqfilter = eq -> eq in aeqs)
570570
var_eq_matching = Matching{Union{Unassigned, SelectedState}}(var_eq_matching′)

test/ccompile.jl

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -14,4 +14,4 @@ p = rand(1)
1414
_t = rand()
1515
f(du, u, p, _t)
1616
f2(du2, u, p, _t)
17-
@test du == du2
17+
@test du du2

test/structural_transformation/tearing.jl

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -169,7 +169,7 @@ prob.f(du, u, pr, tt)
169169
# test the initial guess is respected
170170
@named sys = ODESystem(eqs, t, defaults = Dict(z => Inf))
171171
infprob = ODAEProblem(structural_simplify(sys), [x => 1.0], (0, 1.0), [p => 0.2])
172-
@test_throws DomainError infprob.f(du, u, pr, tt)
172+
@test_throws Any infprob.f(du, u, pr, tt)
173173

174174
sol1 = solve(prob, Tsit5())
175175
sol2 = solve(ODEProblem{false}((u, p, t) -> [-asin(u[1] - pr * t)],

0 commit comments

Comments
 (0)