Skip to content

Commit 8445499

Browse files
Merge pull request #921 from SciML/noisesize_error
Add a better error message for noise size incompatability
2 parents 76dcda3 + f90ed89 commit 8445499

File tree

2 files changed

+49
-1
lines changed

2 files changed

+49
-1
lines changed

src/solve.jl

Lines changed: 33 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -187,6 +187,34 @@ function Base.showerror(io::IO, e::NonSolverError)
187187
println(io, TruncatedStacktraces.VERBOSE_MSG)
188188
end
189189

190+
const NOISE_SIZE_MESSAGE = """
191+
Noise sizes are incompatible. The expected number of noise terms in the defined
192+
`noise_rate_prototype` does not match the number of noise terms in the defined
193+
`AbstractNoiseProcess`. Please ensure that
194+
size(prob.noise_rate_prototype,2) == length(prob.noise.W[1]).
195+
196+
Note: Noise process definitions require that users specify `u0`, and this value is
197+
directly used in the definition. For example, if `noise = WienerProcess(0.0,0.0)`,
198+
then the noise process is a scalar with `u0=0.0`. If `noise = WienerProcess(0.0,[0.0])`,
199+
then the noise process is a vector with `u0=0.0`. If `noise_rate_prototype = zeros(2,4)`,
200+
then the noise process must be a 4-dimensional process, for example
201+
`noise = WienerProcess(0.0,zeros(4))`. This error is a sign that the user definition
202+
of `noise_rate_prototype` and `noise` are not aligned in this manner and the definitions should
203+
be double checked.
204+
"""
205+
206+
struct NoiseSizeIncompatabilityError <: Exception
207+
prototypesize::Int
208+
noisesize::Int
209+
end
210+
211+
function Base.showerror(io::IO, e::NoiseSizeIncompatabilityError)
212+
println(io, NOISE_SIZE_MESSAGE)
213+
println(io, "size(prob.noise_rate_prototype,2) = $(e.prototypesize)")
214+
println(io, "length(prob.noise.W[1]) = $(e.noisesize)")
215+
println(io, TruncatedStacktraces.VERBOSE_MSG)
216+
end
217+
190218
const PROBSOLVER_PAIRING_MESSAGE = """
191219
Incompatible problem+solver pairing.
192220
For example, this can occur if an ODE solver is passed with an SDEProblem.
@@ -1278,6 +1306,11 @@ function check_prob_alg_pairing(prob, alg)
12781306
throw(DirectAutodiffError())
12791307
end
12801308

1309+
if prob isa SDEProblem && prob.noise_rate_prototype !== nothing &&
1310+
prob.noise !== nothing && size(prob.noise_rate_prototype,2) != length(prob.noise.W[1])
1311+
throw(NoiseSizeIncompatabilityError(size(prob.noise_rate_prototype,2), length(prob.noise.W[1])))
1312+
end
1313+
12811314
# Complex number support comes before arbitrary number support for a more direct
12821315
# error message.
12831316
if !SciMLBase.allowscomplex(alg)

test/downstream/solve_error_handling.jl

Lines changed: 16 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1,4 +1,4 @@
1-
using OrdinaryDiffEq, Test, Sundials
1+
using OrdinaryDiffEq, StochasticDiffEq, Test, Sundials
22

33
f(u, p, t) = 2u
44
u0 = 0.5
@@ -48,3 +48,18 @@ fmm = ODEFunction((du, u, t) -> nothing, mass_matrix = zeros(0, 0))
4848
prob = ODEProblem(fmm, nothing, (0.0, 1.0))
4949
sol = solve(prob, Tsit5())
5050
@test isa(sol, DiffEqBase.ODESolution)
51+
52+
f(du, u, p, t) = du .= 1.01u
53+
function g(du, u, p, t)
54+
du[1, 1] = 0.3u[1]
55+
du[1, 2] = 0.6u[1]
56+
du[1, 3] = 0.9u[1]
57+
du[1, 4] = 0.12u[1]
58+
du[2, 1] = 1.2u[2]
59+
du[2, 2] = 0.2u[2]
60+
du[2, 3] = 0.3u[2]
61+
du[2, 4] = 1.8u[2]
62+
end
63+
64+
prob = SDEProblem(f, g, randn(ComplexF64,2), (0.0, 1.0), noise_rate_prototype =complex(zeros(2, 4)),noise=StochasticDiffEq.RealWienerProcess(0.0,zeros(3)))
65+
@test_throws DiffEqBase.NoiseSizeIncompatabilityError solve(prob, LambaEM())

0 commit comments

Comments
 (0)