Skip to content

Commit 222acf1

Browse files
Merge pull request SciML#266 from isaacsas/drop_fastmath
drop fastmaths
2 parents 87d93d0 + 846818b commit 222acf1

File tree

5 files changed

+24
-26
lines changed

5 files changed

+24
-26
lines changed

src/JumpProcesses.jl

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -8,6 +8,7 @@ using DataStructures, PoissonRandom, Random, ArrayInterfaceCore
88
using FunctionWrappers, UnPack
99
using Graphs
1010
using SciMLBase: SciMLBase
11+
using Base.FastMath: add_fast
1112

1213
import DiffEqBase: DiscreteCallback, init, solve, solve!, plot_indices, initialize!
1314
import Base: size, getindex, setindex!, length, similar, show, merge!, merge

src/aggregators/direct.jl

Lines changed: 9 additions & 10 deletions
Original file line numberDiff line numberDiff line change
@@ -58,16 +58,16 @@ end
5858
# calculate the next jump / jump time
5959
function generate_jumps!(p::DirectJumpAggregation, integrator, u, params, t)
6060
p.sum_rate, ttnj = time_to_next_jump(p, u, params, t)
61-
@fastmath p.next_jump_time = t + ttnj
61+
p.next_jump_time = add_fast(t, ttnj)
6262
@inbounds p.next_jump = searchsortedfirst(p.cur_rates, rand(p.rng) * p.sum_rate)
6363
nothing
6464
end
6565

6666
######################## SSA specific helper routines ########################
6767

6868
# tuple-based constant jumps
69-
@fastmath function time_to_next_jump(p::DirectJumpAggregation{T, S, F1, F2, RNG}, u, params,
70-
t) where {T, S, F1 <: Tuple, F2 <: Tuple, RNG}
69+
function time_to_next_jump(p::DirectJumpAggregation{T, S, F1, F2, RNG}, u, params,
70+
t) where {T, S, F1 <: Tuple, F2 <: Tuple, RNG}
7171
prev_rate = zero(t)
7272
new_rate = zero(t)
7373
cur_rates = p.cur_rates
@@ -77,7 +77,7 @@ end
7777
idx = get_num_majumps(majumps)
7878
@inbounds for i in 1:idx
7979
new_rate = evalrxrate(u, i, majumps)
80-
cur_rates[i] = new_rate + prev_rate
80+
cur_rates[i] = add_fast(new_rate, prev_rate)
8181
prev_rate = cur_rates[i]
8282
end
8383

@@ -87,7 +87,7 @@ end
8787
idx += 1
8888
fill_cur_rates(u, params, t, cur_rates, idx, rates...)
8989
@inbounds for i in idx:length(cur_rates)
90-
cur_rates[i] = cur_rates[i] + prev_rate
90+
cur_rates[i] = add_fast(cur_rates[i], prev_rate)
9191
prev_rate = cur_rates[i]
9292
end
9393
end
@@ -108,9 +108,8 @@ end
108108
end
109109

110110
# function wrapper-based constant jumps
111-
@fastmath function time_to_next_jump(p::DirectJumpAggregation{T, S, F1, F2, RNG}, u, params,
112-
t) where {T, S, F1 <: AbstractArray,
113-
F2 <: AbstractArray, RNG}
111+
function time_to_next_jump(p::DirectJumpAggregation{T, S, F1, F2, RNG}, u, params,
112+
t) where {T, S, F1 <: AbstractArray, F2 <: AbstractArray, RNG}
114113
prev_rate = zero(t)
115114
new_rate = zero(t)
116115
cur_rates = p.cur_rates
@@ -120,7 +119,7 @@ end
120119
idx = get_num_majumps(majumps)
121120
@inbounds for i in 1:idx
122121
new_rate = evalrxrate(u, i, majumps)
123-
cur_rates[i] = new_rate + prev_rate
122+
cur_rates[i] = add_fast(new_rate, prev_rate)
124123
prev_rate = cur_rates[i]
125124
end
126125

@@ -129,7 +128,7 @@ end
129128
rates = p.rates
130129
@inbounds for i in 1:length(p.rates)
131130
new_rate = rates[i](u, params, t)
132-
cur_rates[idx] = new_rate + prev_rate
131+
cur_rates[idx] = add_fast(new_rate, prev_rate)
133132
prev_rate = cur_rates[idx]
134133
idx += 1
135134
end

src/aggregators/frm.jl

Lines changed: 6 additions & 8 deletions
Original file line numberDiff line numberDiff line change
@@ -75,7 +75,7 @@ end
7575
######################## SSA specific helper routines ########################
7676

7777
# mass action jumps
78-
@fastmath function next_ma_jump(p::FRMJumpAggregation, u, params, t)
78+
function next_ma_jump(p::FRMJumpAggregation, u, params, t)
7979
ttnj = typemax(typeof(t))
8080
nextrx = zero(Int)
8181
majumps = p.ma_jumps
@@ -91,9 +91,8 @@ end
9191
end
9292

9393
# tuple-based constant jumps
94-
@fastmath function next_constant_rate_jump(p::FRMJumpAggregation{T, S, F1, F2, RNG}, u,
95-
params,
96-
t) where {T, S, F1 <: Tuple, F2 <: Tuple, RNG}
94+
function next_constant_rate_jump(p::FRMJumpAggregation{T, S, F1, F2, RNG}, u, params,
95+
t) where {T, S, F1 <: Tuple, F2 <: Tuple, RNG}
9796
ttnj = typemax(typeof(t))
9897
nextrx = zero(Int)
9998
if !isempty(p.rates)
@@ -111,10 +110,9 @@ end
111110
end
112111

113112
# function wrapper-based constant jumps
114-
@fastmath function next_constant_rate_jump(p::FRMJumpAggregation{T, S, F1, F2, RNG}, u,
115-
params,
116-
t) where {T, S, F1 <: AbstractArray,
117-
F2 <: AbstractArray, RNG}
113+
function next_constant_rate_jump(p::FRMJumpAggregation{T, S, F1, F2, RNG}, u, params,
114+
t) where {T, S, F1 <: AbstractArray, F2 <: AbstractArray,
115+
RNG}
118116
ttnj = typemax(typeof(t))
119117
nextrx = zero(Int)
120118
if !isempty(p.rates)

src/massaction_rates.jl

Lines changed: 6 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -3,8 +3,8 @@
33
# stochiometric coefficient.
44
###############################################################################
55

6-
@inline @fastmath function evalrxrate(speciesvec::AbstractVector{T}, rxidx::S,
7-
majump::MassActionJump{U, V, W, X})::R where
6+
@inline function evalrxrate(speciesvec::AbstractVector{T}, rxidx::S,
7+
majump::MassActionJump{U, V, W, X})::R where
88
{T, S, R, U <: AbstractVector{R}, V, W, X}
99
val = one(T)
1010
@inbounds for specstoch in majump.reactant_stoch[rxidx]
@@ -19,17 +19,17 @@
1919
@inbounds return val * majump.scaled_rates[rxidx]
2020
end
2121

22-
@inline @fastmath function executerx!(speciesvec::AbstractVector{T}, rxidx::S,
23-
majump::M) where {T, S, M <: AbstractMassActionJump}
22+
@inline function executerx!(speciesvec::AbstractVector{T}, rxidx::S,
23+
majump::M) where {T, S, M <: AbstractMassActionJump}
2424
@inbounds net_stoch = majump.net_stoch[rxidx]
2525
@inbounds for specstoch in net_stoch
2626
speciesvec[specstoch[1]] += specstoch[2]
2727
end
2828
nothing
2929
end
3030

31-
@inline @fastmath function executerx(speciesvec::SVector{T}, rxidx::S,
32-
majump::M) where {T, S, M <: AbstractMassActionJump}
31+
@inline function executerx(speciesvec::SVector{T}, rxidx::S,
32+
majump::M) where {T, S, M <: AbstractMassActionJump}
3333
@inbounds net_stoch = majump.net_stoch[rxidx]
3434
@inbounds for specstoch in net_stoch
3535
speciesvec = setindex(speciesvec, speciesvec[specstoch[1]] + specstoch[2],

test/ssa_callback_test.jl

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -25,13 +25,13 @@ function fuel_affect!(integrator)
2525
end
2626
cb = DiscreteCallback(condition, fuel_affect!, save_positions = (false, true))
2727

28-
sol = solve(jump_prob, SSAStepper(), callback = cb, tstops = [5])
28+
sol = solve(jump_prob, SSAStepper(); callback = cb, tstops = [5])
2929
@test sol.t[1:2] == [0.0, 5.0] # no jumps between t=0 and t=5
3030
@test sol(5 + 1e-10) == [100, 0] # state just after fueling before any decays can happen
3131

3232
# test can pass callbacks via JumpProblem
3333
jump_prob2 = JumpProblem(prob, Direct(), jump; rng = rng, callback = cb)
34-
sol2 = solve(jump_prob2, SSAStepper(), tstops = [5])
34+
sol2 = solve(jump_prob2, SSAStepper(); tstops = [5])
3535
@test sol2.t[1:2] == [0.0, 5.0] # no jumps between t=0 and t=5
3636
@test sol2(5 + 1e-10) == [100, 0] # state just after fueling before any decays can happen
3737

0 commit comments

Comments
 (0)