-
-
Notifications
You must be signed in to change notification settings - Fork 259
Expand file tree
/
Copy pathDelayDiffEq.jl
More file actions
134 lines (115 loc) · 4.7 KB
/
DelayDiffEq.jl
File metadata and controls
134 lines (115 loc) · 4.7 KB
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
module DelayDiffEq
import Reexport: @reexport, Reexport
@reexport using SciMLBase
import OrdinaryDiffEqCore, OrdinaryDiffEqNonlinearSolve, OrdinaryDiffEqDifferentiation,
OrdinaryDiffEqRosenbrock
import OrdinaryDiffEqDefault: OrdinaryDiffEqDefault
import OrdinaryDiffEqFunctionMap: OrdinaryDiffEqFunctionMap
using BinaryHeaps: BinaryMinHeap
using LinearAlgebra: opnorm, I
using Logging: @logmsg
using RecursiveArrayTools: copyat_or_push!, recursivecopy, recursivecopy!,
recursive_bottom_eltype, recursive_unitless_bottom_eltype,
recursive_unitless_eltype
using ForwardDiff: ForwardDiff
import ArrayInterface
import SimpleNonlinearSolve
import SymbolicIndexingInterface as SII
using SciMLBase: AbstractDDEAlgorithm, AbstractDDEIntegrator, AbstractODEIntegrator,
DEIntegrator
using SciMLBase: AbstractSDDEProblem, SDDEProblem
using Base: deleteat!
import FastBroadcast: @..
using OrdinaryDiffEqNonlinearSolve: NLAnderson, NLFunctional
using OrdinaryDiffEqCore: AbstractNLSolverCache, SlowConvergence,
alg_extrapolates, alg_maximum_order, initialize!, DEVerbosity
using OrdinaryDiffEqCore: StochasticDiffEqAlgorithm,
StochasticDiffEqRODEAlgorithm,
StochasticDiffEqConstantCache, StochasticDiffEqMutableCache
using OrdinaryDiffEqRosenbrock: RosenbrockMutableCache
using OrdinaryDiffEqFunctionMap: FunctionMap
# using OrdinaryDiffEqDifferentiation: resize_grad_config!, resize_jac_config!
using DiffEqBase: is_diagonal_noise
# Explicit imports for functions
using OrdinaryDiffEqCore: AutoSwitch, CompositeAlgorithm
using OrdinaryDiffEqDefault: DefaultODEAlgorithm
using SciMLBase: CallbackSet, DAEProblem, DDEProblem, DESolution, ODEProblem, ReturnCode,
VectorContinuousCallback, addat!, addat_non_user_cache!,
deleteat_non_user_cache!,
du_cache, full_cache, get_tmp_cache, isinplace,
reeval_internals_due_to_modification!,
reinit!, resize_non_user_cache!, savevalues!, u_cache, user_cache,
step!, terminate!, u_modified!, get_proposed_dt, set_proposed_dt!,
auto_dt_reset!,
add_tstop!, add_saveat!, get_du, get_du!, addsteps!,
change_t_via_interpolation!, isadaptive
using DiffEqBase: initialize!
import DiffEqBase
using SciMLLogging: AbstractVerbosityPreset, None, @SciMLMessage
import SciMLBase
const SDEAlgUnion = Union{StochasticDiffEqAlgorithm, StochasticDiffEqRODEAlgorithm}
# Internal hook functions for SDDE support. These are overloaded by the
# StochasticDiffEqCore extension to provide actual implementations.
# Calling them without the extension loaded gives a clear error.
function _sde_alg_cache end
function _create_sdde_noise end
export Discontinuity, MethodOfSteps
include("discontinuity_type.jl")
include("functionwrapper.jl")
include("integrators/type.jl")
include("integrators/utils.jl")
include("integrators/interface.jl")
include("fpsolve/type.jl")
include("fpsolve/fpsolve.jl")
include("fpsolve/utils.jl")
include("fpsolve/functional.jl")
include("cache_utils.jl")
include("interpolants.jl")
include("history_function.jl")
include("algorithms.jl")
include("track.jl")
include("alg_utils.jl")
include("solve.jl")
include("utils.jl")
# Default solver for DDEProblems (no algorithm specified)
function SciMLBase.__solve(prob::DDEProblem; kwargs...)
return SciMLBase.solve(prob, MethodOfSteps(DefaultODEAlgorithm()); kwargs...)
end
function SciMLBase.__init(prob::DDEProblem; kwargs...)
return DiffEqBase.init(prob, MethodOfSteps(DefaultODEAlgorithm()); kwargs...)
end
# Auto-wrap bare ODE/SDE algorithms in MethodOfSteps for DDE/SDDE problems.
# Allows: solve(DDEProblem(...), Tsit5()) instead of MethodOfSteps(Tsit5())
# and: solve(SDDEProblem(...), EM()) instead of MethodOfSteps(EM())
const OrdinaryDiffEqAlgorithm = OrdinaryDiffEqCore.OrdinaryDiffEqAlgorithm
function SciMLBase.__solve(
prob::SciMLBase.AbstractDDEProblem,
alg::OrdinaryDiffEqAlgorithm, args...; kwargs...
)
return SciMLBase.__solve(prob, MethodOfSteps(alg), args...; kwargs...)
end
function SciMLBase.__init(
prob::SciMLBase.AbstractDDEProblem,
alg::OrdinaryDiffEqAlgorithm, args...; kwargs...
)
return SciMLBase.__init(prob, MethodOfSteps(alg), args...; kwargs...)
end
function DiffEqBase.check_prob_alg_pairing(prob::DDEProblem, alg::OrdinaryDiffEqAlgorithm)
return nothing
end
function SciMLBase.__solve(
prob::AbstractSDDEProblem,
alg::SDEAlgUnion, args...; kwargs...
)
return SciMLBase.__solve(prob, MethodOfSteps(alg), args...; kwargs...)
end
function SciMLBase.__init(
prob::AbstractSDDEProblem,
alg::SDEAlgUnion, args...; kwargs...
)
return SciMLBase.__init(prob, MethodOfSteps(alg), args...; kwargs...)
end
function DiffEqBase.check_prob_alg_pairing(prob::SDDEProblem, alg::SDEAlgUnion)
return nothing
end
end # module