forked from josdejong/mathjs
-
Notifications
You must be signed in to change notification settings - Fork 0
Expand file tree
/
Copy pathsolveODE.js
More file actions
322 lines (290 loc) · 10.8 KB
/
solveODE.js
File metadata and controls
322 lines (290 loc) · 10.8 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
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
import { isUnit, isNumber, isBigNumber } from '../../utils/is.js'
import { factory } from '../../utils/factory.js'
const name = 'solveODE'
const dependencies = [
'typed',
'add',
'subtract',
'multiply',
'divide',
'max',
'map',
'abs',
'isPositive',
'isNegative',
'larger',
'smaller',
'matrix',
'bignumber',
'unaryMinus'
]
export const createSolveODE = /* #__PURE__ */ factory(name, dependencies, (
{
typed,
add,
subtract,
multiply,
divide,
max,
map,
abs,
isPositive,
isNegative,
larger,
smaller,
matrix,
bignumber,
unaryMinus
}
) => {
/**
* Numerical Integration of Ordinary Differential Equations
*
* Two variable step methods are provided:
* - "RK23": Bogacki–Shampine method
* - "RK45": Dormand-Prince method RK5(4)7M (default)
*
* The arguments are expected as follows.
*
* - `func` should be the forcing function `f(t, y)`
* - `tspan` should be a vector of two numbers or units `[tStart, tEnd]`
* - `y0` the initial state values, should be a scalar or a flat array
* - `options` should be an object with the following information:
* - `method` ('RK45'): ['RK23', 'RK45']
* - `tol` (1e-3): Numeric tolerance of the method, the solver keeps the error estimates less than this value
* - `firstStep`: Initial step size
* - `minStep`: minimum step size of the method
* - `maxStep`: maximum step size of the method
* - `minDelta` (0.2): minimum ratio of change for the step
* - `maxDelta` (5): maximum ratio of change for the step
* - `maxIter` (1e4): maximum number of iterations
*
* The returned value is an object with `{t, y}` please note that even though `t` means time, it can represent any other independant variable like `x`:
* - `t` an array of size `[n]`
* - `y` the states array can be in two ways
* - **if `y0` is a scalar:** returns an array-like of size `[n]`
* - **if `y0` is a flat array-like of size [m]:** returns an array like of size `[n, m]`
*
* Syntax:
*
* math.solveODE(func, tspan, y0)
* math.solveODE(func, tspan, y0, options)
*
* Examples:
*
* function func(t, y) {return y}
* const tspan = [0, 4]
* const y0 = 1
* math.solveODE(func, tspan, y0)
* math.solveODE(func, tspan, [1, 2])
* math.solveODE(func, tspan, y0, { method:"RK23", maxStep:0.1 })
*
* See also:
*
* derivative, simplifyCore
*
* @param {function} func The forcing function f(t,y)
* @param {Array | Matrix} tspan The time span
* @param {number | BigNumber | Unit | Array | Matrix} y0 The initial value
* @param {Object} [options] Optional configuration options
* @return {Object} Return an object with t and y values as arrays
*/
function _rk (butcherTableau) {
// generates an adaptive runge kutta method from it's butcher tableau
return function (f, tspan, y0, options) {
// adaptive runge kutta methods
const wrongTSpan = !((tspan.length === 2) && (tspan.every(isNumOrBig) || tspan.every(isUnit)))
if (wrongTSpan) {
throw new Error('"tspan" must be an Array of two numeric values or two units [tStart, tEnd]')
}
const t0 = tspan[0] // initial time
const tf = tspan[1] // final time
const isForwards = larger(tf, t0)
const firstStep = options.firstStep
if (firstStep !== undefined && !isPositive(firstStep)) {
throw new Error('"firstStep" must be positive')
}
const maxStep = options.maxStep
if (maxStep !== undefined && !isPositive(maxStep)) {
throw new Error('"maxStep" must be positive')
}
const minStep = options.minStep
if (minStep && isNegative(minStep)) {
throw new Error('"minStep" must be positive or zero')
}
const timeVars = [t0, tf, firstStep, minStep, maxStep].filter(x => x !== undefined)
if (!(timeVars.every(isNumOrBig) || timeVars.every(isUnit))) {
throw new Error('Inconsistent type of "t" dependant variables')
}
const steps = 1 // divide time in this number of steps
const tol = options.tol ? options.tol : 1e-4 // define a tolerance (must be an option)
const minDelta = options.minDelta ? options.minDelta : 0.2
const maxDelta = options.maxDelta ? options.maxDelta : 5
const maxIter = options.maxIter ? options.maxIter : 10_000 // stop inifite evaluation if something goes wrong
const hasBigNumbers = [t0, tf, ...y0, maxStep, minStep].some(isBigNumber)
const [a, c, b, bp] = hasBigNumbers
? [
bignumber(butcherTableau.a),
bignumber(butcherTableau.c),
bignumber(butcherTableau.b),
bignumber(butcherTableau.bp)
]
: [butcherTableau.a, butcherTableau.c, butcherTableau.b, butcherTableau.bp]
let h = firstStep
? isForwards ? firstStep : unaryMinus(firstStep)
: divide(subtract(tf, t0), steps) // define the first step size
const t = [t0] // start the time array
const y = [y0] // start the solution array
const deltaB = subtract(b, bp) // b - bp
let n = 0
let iter = 0
const ongoing = _createOngoing(isForwards)
const trimStep = _createTrimStep(isForwards)
// iterate unitil it reaches either the final time or maximum iterations
while (ongoing(t[n], tf)) {
const k = []
// trim the time step so that it doesn't overshoot
h = trimStep(t[n], tf, h)
// calculate the first value of k
k.push(f(t[n], y[n]))
// calculate the rest of the values of k
for (let i = 1; i < c.length; ++i) {
k.push(
f(
add(t[n], multiply(c[i], h)),
add(y[n], multiply(h, a[i], k))
)
)
}
// estimate the error by comparing solutions of different orders
const TE = max(
abs(
map(multiply(deltaB, k), (X) =>
isUnit(X) ? X.value : X
)
)
)
if (TE < tol && tol / TE > 1 / 4) {
// push solution if within tol
t.push(add(t[n], h))
y.push(add(y[n], multiply(h, b, k)))
n++
}
// estimate the delta value that will affect the step size
let delta = 0.84 * (tol / TE) ** (1 / 5)
if (smaller(delta, minDelta)) {
delta = minDelta
} else if (larger(delta, maxDelta)) {
delta = maxDelta
}
delta = hasBigNumbers ? bignumber(delta) : delta
h = multiply(h, delta)
if (maxStep && larger(abs(h), maxStep)) {
h = isForwards ? maxStep : unaryMinus(maxStep)
} else if (minStep && smaller(abs(h), minStep)) {
h = isForwards ? minStep : unaryMinus(minStep)
}
iter++
if (iter > maxIter) {
throw new Error('Maximum number of iterations reached, try changing options')
}
}
return { t, y }
}
}
function _rk23 (f, tspan, y0, options) {
// Bogacki–Shampine method
// Define the butcher table
const a = [
[],
[1 / 2],
[0, 3 / 4],
[2 / 9, 1 / 3, 4 / 9]
]
const c = [null, 1 / 2, 3 / 4, 1]
const b = [2 / 9, 1 / 3, 4 / 9, 0]
const bp = [7 / 24, 1 / 4, 1 / 3, 1 / 8]
const butcherTableau = { a, c, b, bp }
// Solve an adaptive step size rk method
return _rk(butcherTableau)(f, tspan, y0, options)
}
function _rk45 (f, tspan, y0, options) {
// Dormand Prince method
// Define the butcher tableau
const a = [
[],
[1 / 5],
[3 / 40, 9 / 40],
[44 / 45, -56 / 15, 32 / 9],
[19372 / 6561, -25360 / 2187, 64448 / 6561, -212 / 729],
[9017 / 3168, -355 / 33, 46732 / 5247, 49 / 176, -5103 / 18656],
[35 / 384, 0, 500 / 1113, 125 / 192, -2187 / 6784, 11 / 84]
]
const c = [null, 1 / 5, 3 / 10, 4 / 5, 8 / 9, 1, 1]
const b = [35 / 384, 0, 500 / 1113, 125 / 192, -2187 / 6784, 11 / 84, 0]
const bp = [5179 / 57600, 0, 7571 / 16695, 393 / 640, -92097 / 339200, 187 / 2100, 1 / 40]
const butcherTableau = { a, c, b, bp }
// Solve an adaptive step size rk method
return _rk(butcherTableau)(f, tspan, y0, options)
}
function _solveODE (f, tspan, y0, opt) {
const method = opt.method ? opt.method : 'RK45'
const methods = {
RK23: _rk23,
RK45: _rk45
}
if (method.toUpperCase() in methods) {
const methodOptions = { ...opt } // clone the options object
delete methodOptions.method // delete the method as it won't be needed
return methods[method.toUpperCase()](f, tspan, y0, methodOptions)
} else {
// throw an error indicating there is no such method
const methodsWithQuotes = Object.keys(methods).map(x => `"${x}"`)
// generates a string of methods like: "BDF", "RK23" and "RK45"
const availableMethodsString = `${methodsWithQuotes.slice(0, -1).join(', ')} and ${methodsWithQuotes.slice(-1)}`
throw new Error(`Unavailable method "${method}". Available methods are ${availableMethodsString}`)
}
}
function _createOngoing (isForwards) {
// returns the correct function to test if it's still iterating
return isForwards ? smaller : larger
}
function _createTrimStep (isForwards) {
const outOfBounds = isForwards ? larger : smaller
return function (t, tf, h) {
const next = add(t, h)
return outOfBounds(next, tf) ? subtract(tf, t) : h
}
}
function isNumOrBig (x) {
// checks if it's a number or bignumber
return isBigNumber(x) || isNumber(x)
}
function _matrixSolveODE (f, T, y0, options) {
// receives matrices and returns matrices
const sol = _solveODE(f, T.toArray(), y0.toArray(), options)
return { t: matrix(sol.t), y: matrix(sol.y) }
}
return typed('solveODE', {
'function, Array, Array, Object': _solveODE,
'function, Matrix, Matrix, Object': _matrixSolveODE,
'function, Array, Array': (f, T, y0) => _solveODE(f, T, y0, {}),
'function, Matrix, Matrix': (f, T, y0) => _matrixSolveODE(f, T, y0, {}),
'function, Array, number | BigNumber | Unit': (f, T, y0) => {
const sol = _solveODE(f, T, [y0], {})
return { t: sol.t, y: sol.y.map((Y) => Y[0]) }
},
'function, Matrix, number | BigNumber | Unit': (f, T, y0) => {
const sol = _solveODE(f, T.toArray(), [y0], {})
return { t: matrix(sol.t), y: matrix(sol.y.map((Y) => Y[0])) }
},
'function, Array, number | BigNumber | Unit, Object': (f, T, y0, options) => {
const sol = _solveODE(f, T, [y0], options)
return { t: sol.t, y: sol.y.map((Y) => Y[0]) }
},
'function, Matrix, number | BigNumber | Unit, Object': (f, T, y0, options) => {
const sol = _solveODE(f, T.toArray(), [y0], options)
return { t: matrix(sol.t), y: matrix(sol.y.map((Y) => Y[0])) }
}
})
})