Skip to content

Commit 6e33f5b

Browse files
authored
Merge pull request #1143 from AayushSabharwal/as/hc-bug
fix: handle bug in `remove_denominators`
2 parents 0890a35 + 3d08d21 commit 6e33f5b

File tree

3 files changed

+89
-21
lines changed

3 files changed

+89
-21
lines changed

ext/CatalystHomotopyContinuationExtension.jl

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -6,7 +6,7 @@ import DynamicPolynomials
66
import ModelingToolkit as MT
77
import HomotopyContinuation as HC
88
import Setfield: @set
9-
import Symbolics: unwrap, wrap, Rewriters, symtype, issym
9+
import Symbolics: unwrap, wrap, Rewriters, symtype, issym, maketerm, BasicSymbolic, metadata
1010
using Symbolics: iscall
1111

1212
# Creates and exports hc_steady_states function.

ext/CatalystHomotopyContinuationExtension/homotopy_continuation_extension.jl

Lines changed: 54 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -10,8 +10,8 @@ Arguments:
1010
- `ps`: The parameter values for which we want to find the steady states.
1111
- `filter_negative=true`: If set to true, solutions with any species concentration <neg_thres is removed from the output.
1212
- `neg_thres=-1e-20`: Determine the minimum values for which a species concentration is to be considered non-negative. Species concentrations ``> neg_thres`` but `< 0.0` are set to `0.0`.
13-
- `u0=nothing`: Initial conditions for which we want to find the steady states. For systems with conservation laws this are required to compute conserved quantities. Initial conditions are not required for all species, only those involved in conserved quantities (if this set is unknown, it is recommended to provide initial conditions for all species).
14-
- `kwargs...`: any additional arguments (like `show_progress= true`) are passed into HomotopyContinuation.jl's `solve` call.
13+
- `u0=nothing`: Initial conditions for which we want to find the steady states. For systems with conservation laws this are required to compute conserved quantities. Initial conditions are not required for all species, only those involved in conserved quantities (if this set is unknown, it is recommended to provide initial conditions for all species).
14+
- `kwargs...`: any additional arguments (like `show_progress= true`) are passed into HomotopyContinuation.jl's `solve` call.
1515
1616
Examples
1717
```@repl
@@ -48,6 +48,8 @@ end
4848

4949
# For a given reaction system, parameter values, and initial conditions, find the polynomial that HC solves to find steady states.
5050
function steady_state_polynomial(rs::ReactionSystem, ps, u0)
51+
# Creates the appropriate nonlinear system, and converts parameters to a form that can
52+
# be substituted in later.
5153
rs = Catalyst.expand_registered_functions(rs)
5254
ns = complete(convert(NonlinearSystem, rs;
5355
remove_conserved = true, remove_conserved_warn = false))
@@ -56,10 +58,15 @@ function steady_state_polynomial(rs::ReactionSystem, ps, u0)
5658
p_vals = ModelingToolkit.varmap_to_vars(pre_varmap, parameters(ns);
5759
defaults = ModelingToolkit.defaults(ns))
5860
p_dict = Dict(parameters(ns) .=> p_vals)
59-
eqs_pars_funcs = vcat(equations(ns), conservedequations(rs))
60-
eqs = map(eq -> substitute(eq.rhs - eq.lhs, p_dict), eqs_pars_funcs)
61-
eqs_intexp = make_int_exps.(eqs)
62-
ss_poly = Catalyst.to_multivariate_poly(remove_denominators.(eqs_intexp))
61+
62+
# Step by step convert the equation to something HC can work on (adds conserved equations,
63+
# inserts parameter values and put everything on a side, converts e.g. 2.0 to 2, remove
64+
# denominators to avoid rational polynomials).
65+
eqs = vcat(equations(ns), conservedequations(rs))
66+
eqs = [substitute(eq.rhs - eq.lhs, p_dict) for eq in eqs]
67+
eqs = make_int_exps.(eqs)
68+
eqs = [remove_denominators(common_denominator(eq)) for eq in eqs]
69+
ss_poly = Catalyst.to_multivariate_poly(eqs)
6370
return poly_type_convert(ss_poly)
6471
end
6572

@@ -78,6 +85,47 @@ function ___make_int_exps(expr)
7885
end
7986
end
8087

88+
# Converts an expression of the form p(x)/q(x) + r(x)/s(x) to t(x)/u(x) (i.e. puts everything
89+
# above a single denominator, which is what `remove_denominators` is able to simplify away).
90+
function common_denominator(expr)
91+
iscall(expr) || return expr
92+
if operation(expr) == /
93+
num, den = arguments(expr)
94+
num = common_denominator(num)
95+
den = common_denominator(den)
96+
return num / den
97+
end
98+
if operation(expr) == +
99+
num = 0
100+
den = 1
101+
for arg in arguments(expr)
102+
arg = common_denominator(arg)
103+
if iscall(arg) && operation(arg) == /
104+
n, d = arguments(arg)
105+
else
106+
n = arg
107+
d = 1
108+
end
109+
num = num * d + den * n
110+
den *= d
111+
end
112+
return num / den
113+
end
114+
if operation(expr) == ^
115+
base, pow = arguments(expr)
116+
base = common_denominator(base)
117+
if iscall(base) && operation(base) == /
118+
num, den = arguments(base)
119+
else
120+
num, den = base, 1
121+
end
122+
num ^= pow
123+
den ^= pow
124+
return num / den
125+
end
126+
return maketerm(BasicSymbolic, operation(expr), map(common_denominator, arguments(expr)), metadata(expr))
127+
end
128+
81129
# If the input is a fraction, removes the denominator.
82130
function remove_denominators(expr)
83131
s_expr = simplify_fractions(expr)

test/extensions/homotopy_continuation.jl

Lines changed: 34 additions & 14 deletions
Original file line numberDiff line numberDiff line change
@@ -14,7 +14,7 @@ include("../test_functions.jl")
1414
# Tests for Symbolics initial condition input.
1515
# Tests for different types (Symbol/Symbolics) for parameters and initial conditions.
1616
# Tests that attempts to find steady states of system with conservation laws, while u0 is not provided, gives an error.
17-
let
17+
let
1818
# Creates the model.
1919
rs = @reaction_network begin
2020
(k1,k2), X1 <--> X2
@@ -37,7 +37,7 @@ end
3737
# Tests for Symbol parameter input.
3838
# Tests that passing kwargs to HC.solve does not error and have an effect (i.e. modifying the seed
3939
# slightly modified the output in some way).
40-
let
40+
let
4141
wilhelm_2009_model = @reaction_network begin
4242
k1, Y --> 2X
4343
k2, 2X --> X + Y
@@ -61,44 +61,44 @@ end
6161
# Tests where input ps/u0 are tuples with mixed types.
6262
let
6363
rs_1 = @reaction_network begin
64-
@parameters kX1=1.0 kX2=2.0 kY1=12345.0
64+
@parameters kX1=1.0 kX2=2.0 kY1=12345.0
6565
@species X1(t)=0.1 X2(t)=0.2 Y1(t)=12345.0
6666
(kX1,kX2), X1 <--> X2
6767
(kY1,kY2), Y1 <--> Y2
6868
(kZ1,kZ2), Z1 <--> Z2
6969
end
7070
ps = (:kY1 => 1.0, :kY2 => 3, :kZ1 => 1.0, :kZ2 => 4.0)
7171
u0_1 = (:Y1 => 1.0, :Y2 => 3, :Z1 => 10, :Z2 =>40.0)
72-
72+
7373
ss_1 = sort(hc_steady_states(rs_1, ps; u0 = u0_1, show_progress = false, seed = 0x000004d1), by = sol->sol[1])
7474
@test ss_1 [[0.2, 0.1, 3.0, 1.0, 40.0, 10.0]]
75-
75+
7676
rs_2 = @reaction_network begin
77-
@parameters kX1=1.0 kX2=2.0 kY1=12345.0
77+
@parameters kX1=1.0 kX2=2.0 kY1=12345.0
7878
@species C2(t)=0.1 C1(t)=0.2 B2(t)=12345.0
7979
(kX1,kX2), C2 <--> C1
8080
(kY1,kY2), B2 <--> B1
8181
(kZ1,kZ2), A2 <--> A1
8282
end
8383
u0_2 = [:B2 => 1.0, :B1 => 3.0, :A2 => 10.0, :A1 =>40.0]
84-
84+
8585
ss_2 = sort(hc_steady_states(rs_2, ps; u0 = u0_2, show_progress = false, seed = 0x000004d1), by = sol->sol[1])
8686
@test ss_1 ss_2
8787
end
8888

8989
# Tests that non-scalar reaction rates work.
90-
# Tests that rational polynomial steady state systems work.
90+
# Tests that rational polynomial steady state systems work.
9191
# Tests that Hill function is correctly expanded even if nested.
9292
# Test filter_negative=false works.
9393
# Tests than non-integer exponents throws an error.
94-
let
94+
let
9595
rs = @reaction_network begin
9696
v*(0.1/v + hill(X,1,K,n)), 0 --> X
9797
d, X --> 0
9898
end
9999
ps = Dict([:v => 5.0, :K => 2.5, :n => 3, :d => 1.0])
100100
sss = hc_steady_states(rs, ps; filter_negative = false, show_progress = false, seed = 0x000004d1)
101-
101+
102102
@test length(sss) == 4
103103
for ss in sss
104104
@test ps[:v]*(0.1/ps[:v] + ss[1]^ps[:n]/(ss[1]^ps[:n] + ps[:K]^ps[:n])) - ps[:d]*ss[1] 0.0 atol = 1e-12
@@ -108,6 +108,26 @@ let
108108
@test_throws Exception hc_steady_states(rs, ps; show_progress = false, seed = 0x000004d1)
109109
end
110110

111+
# Checks that the correct steady states are found for system with multiple, known, steady states.
112+
# Checks where (activating/repressing) Hill function is written out/using `hillar`.
113+
# The system is known to have (exactly five steady states for the given parameter set.
114+
let
115+
# Finds the model steady states.
116+
rs = @reaction_network begin
117+
0.01 + hillar(X,Y,1.0,Kx,3), ∅ --> X
118+
0.01 + (Y^3) / (X^3 + Y^3 + Ky^3), ∅ --> Y
119+
1.0, (X, Y) -->
120+
end
121+
ps = [:Kx => 0.3, :Ky => 0.5]
122+
sss = hc_steady_states(rs, ps)
123+
124+
# Checks that the steady states are correct, all unique, and that 5 (known number) is found.
125+
@test length(sss) == 5
126+
@test allunique(sss)
127+
for ss in sss
128+
@test f_eval(rs, ss, ps, 0.0) [0.0, 0.0] atol = 1e-12 rtol = 1e-12
129+
end
130+
end
111131

112132
### Other Tests ###
113133

@@ -131,22 +151,22 @@ let
131151
end
132152

133153
# Checks that `hc_steady_states` cannot be applied to non-complete `ReactionSystems`s.
134-
let
154+
let
135155
# Create model.
136156
incomplete_network = @network_component begin
137157
(p, d), 0 <--> X
138158
end
139159
p_start = [:p => 1.0, :d => 0.2]
140-
160+
141161
# Computes bifurcation diagram.
142162
@test_throws Exception hc_steady_states(incomplete_network, p_start; show_progress = false, seed = 0x000004d1)
143163
end
144164

145165
# Tests that non-autonomous system throws an error
146-
let
166+
let
147167
rs = @reaction_network begin
148168
(k,t), 0 <--> X
149169
end
150170
ps = [:k => 1.0]
151171
@test_throws Exception hc_steady_states(rs, ps; show_progress = false, seed = 0x000004d1)
152-
end
172+
end

0 commit comments

Comments
 (0)