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
use super::{RootFinder, Stats};
use crate::StrError;
impl RootFinder {
/// Employs Brent's method to find a single root of an equation
///
/// See: <https://mathworld.wolfram.com/BrentsMethod.html>
///
/// See also: <https://en.wikipedia.org/wiki/Brent%27s_method>
///
/// # Input
///
/// * `xa` -- lower bound such that `f(xa) × f(xb) < 0`
/// * `xb` -- upper bound such that `f(xa) × f(xb) < 0`
/// * `args` -- extra arguments for the callback function
/// * `f` -- callback function implementing `f(x)` as `f(x, args)`; it returns `f @ x` or it may return an error.
///
/// **Note:** `xa` must be different from `xb`
///
/// # Output
///
/// Returns `(xo, stats)` where:
///
/// * `xo` -- is the root such that `f(xo) = 0`
/// * `stats` -- some statistics regarding the computations
///
/// # Examples
///
/// ```
/// use russell_lab::*;
///
/// fn main() -> Result<(), StrError> {
/// let args = &mut 0;
/// let solver = RootFinder::new();
/// let (xa, xb) = (-4.0, 0.0);
/// let (xo, stats) = solver.brent(xa, xb, args, |x, _| Ok(4.0 - x * x))?;
/// println!("\nroot = {:?}", xo);
/// println!("\n{}", stats);
/// approx_eq(xo, -2.0, 1e-14);
/// Ok(())
/// }
/// ```
pub fn brent<F, A>(&self, xa: f64, xb: f64, args: &mut A, mut f: F) -> Result<(f64, Stats), StrError>
where
F: FnMut(f64, &mut A) -> Result<f64, StrError>,
{
// Based on ZEROIN C math library: <http://www.netlib.org/c/>
// By: Oleg Keselyov <oleg@ponder.csci.unt.edu, oleg@unt.edu> May 23, 1991
//
// G.Forsythe, M.Malcolm, C.Moler, Computer methods for mathematical
// computations. M., Mir, 1980, p.180 of the Russian edition
//
// The function makes use of the bisection procedure combined with
// the linear or quadric inverse interpolation.
// At every step program operates on three abscissae - a, b, and c.
//
// * b - the last and the best approximation to the root
// * a - the last but one approximation
// * c - the last but one or even earlier approximation than a that
// 1. |f(b)| <= |f(c)|
// 2. f(b) and f(c) have opposite signs, i.e. b and c confine the root
//
// At every step Zeroin selects one of the two new approximations, the
// former being obtained by the bisection procedure and the latter
// resulting in the interpolation (if a,b, and c are all different
// the quadric interpolation is utilized, otherwise the linear one).
// If the latter (i.e. obtained by the interpolation) point is
// reasonable (i.e. lies within the current interval `[b,c]` not being
// too close to the boundaries) it is accepted. The bisection result
// is used in the other case. Therefore, the range of uncertainty is
// ensured to be reduced at least by the factor 1.6
// check
if f64::abs(xa - xb) < 10.0 * f64::EPSILON {
return Err("xa must be different from xb");
}
// allocate stats struct
let mut stats = Stats::new();
// initialization
let (mut a, mut b) = (xa, xb);
let (mut fa, mut fb) = (f(a, args)?, f(b, args)?);
let mut c = a;
let mut fc = fa;
stats.n_function += 2;
// check
if fa * fb >= -f64::EPSILON {
return Err("xa and xb must bracket the root and f(xa) × f(xb) < 0");
}
// solve
let mut converged = false;
for _ in 0..self.brent_max_iterations {
stats.n_iterations += 1;
// old step
let step_old = b - a;
// swap data
if f64::abs(fc) < f64::abs(fb) {
a = b;
b = c;
c = a;
fa = fb;
fb = fc;
fc = fa;
}
// tolerance
let tol = 2.0 * f64::EPSILON * f64::abs(b) + self.brent_tolerance / 2.0;
// new step
let mut step_new = (c - b) / 2.0;
// converged?
if f64::abs(step_new) <= tol || fb == 0.0 {
converged = true;
break;
}
// perform interpolation
if f64::abs(step_old) >= tol && f64::abs(fa) > f64::abs(fb) {
// delta
let del = c - b;
// linear interpolation
let (mut p, mut q) = if a == c {
let t0 = fb / fa;
(del * t0, 1.0 - t0)
// quadratic inverse interpolation
} else {
let t0 = fa / fc;
let t1 = fb / fc;
let t2 = fb / fa;
(
t2 * (del * t0 * (t0 - t1) - (b - a) * (t1 - 1.0)),
(t0 - 1.0) * (t1 - 1.0) * (t2 - 1.0),
)
};
// fix the sign of p and q
if p > 0.0 {
q = -q;
} else {
p = -p;
}
// update step
if p < (0.75 * del * q - f64::abs(tol * q) / 2.0) && p < f64::abs(step_old * q / 2.0) {
step_new = p / q;
}
}
// limit the step
if f64::abs(step_new) < tol {
if step_new > 0.0 {
step_new = tol;
} else {
step_new = -tol;
}
}
// update a
a = b;
fa = fb;
// update b
b += step_new;
fb = f(b, args)?;
stats.n_function += 1;
// update c
if (fb > 0.0 && fc > 0.0) || (fb < 0.0 && fc < 0.0) {
c = a;
fc = fa;
}
}
// check
if !converged {
return Err("brent solver failed to converge");
}
// done
stats.stop_sw_total();
Ok((b, stats))
}
}
////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
#[cfg(test)]
mod tests {
use crate::algo::testing::get_test_functions;
use crate::approx_eq;
use crate::{NoArgs, RootFinder};
#[test]
fn brent_find_captures_errors_1() {
let f = |x, _: &mut NoArgs| Ok(x * x - 1.0);
let args = &mut 0;
let mut solver = RootFinder::new();
assert_eq!(f(1.0, args).unwrap(), 0.0);
assert_eq!(
solver.brent(-0.5, -0.5, args, f).err(),
Some("xa must be different from xb")
);
assert_eq!(
solver.brent(-0.5, -0.5 - 10.0 * f64::EPSILON, args, f).err(),
Some("xa and xb must bracket the root and f(xa) × f(xb) < 0")
);
solver.brent_max_iterations = 0;
assert_eq!(
solver.brent(0.0, 2.0, args, f).err(),
Some("brent solver failed to converge")
);
}
#[test]
fn brent_find_captures_errors_2() {
struct Args {
count: usize,
target: usize,
}
let f = |x, args: &mut Args| {
let res = if args.count == args.target {
Err("stop")
} else {
Ok(x * x - 1.0)
};
args.count += 1;
res
};
let args = &mut Args { count: 0, target: 0 };
let solver = RootFinder::new();
// first function call
assert_eq!(solver.brent(-0.5, 2.0, args, f).err(), Some("stop"));
// second function call
args.count = 0;
args.target = 1;
assert_eq!(solver.brent(-0.5, 2.0, args, f).err(), Some("stop"));
// third function call
args.count = 0;
args.target = 2;
assert_eq!(solver.brent(-0.5, 2.0, args, f).err(), Some("stop"));
}
#[test]
fn brent_find_works() {
let args = &mut 0;
let solver = RootFinder::new();
for test in &get_test_functions() {
println!("\n===================================================================");
println!("\n{}", test.name);
if let Some(bracket) = test.root1 {
let (xo, stats) = solver.brent(bracket.a, bracket.b, args, test.f).unwrap();
println!("\nxo = {:?}", xo);
println!("\n{}", stats);
approx_eq(xo, bracket.xo, 1e-11);
approx_eq((test.f)(xo, args).unwrap(), 0.0, test.tol_root);
}
if let Some(bracket) = test.root2 {
let (xo, stats) = solver.brent(bracket.a, bracket.b, args, test.f).unwrap();
println!("\nxo = {:?}", xo);
println!("\n{}", stats);
approx_eq(xo, bracket.xo, 1e-11);
approx_eq((test.f)(xo, args).unwrap(), 0.0, test.tol_root);
}
if let Some(bracket) = test.root3 {
let (xo, stats) = solver.brent(bracket.a, bracket.b, args, test.f).unwrap();
println!("\nxo = {:?}", xo);
println!("\n{}", stats);
approx_eq(xo, bracket.xo, 1e-13);
approx_eq((test.f)(xo, args).unwrap(), 0.0, test.tol_root);
}
}
println!("\n===================================================================\n");
}
}