scirs2_optimize/unconstrained/
powell.rs

1//! Powell's method for unconstrained optimization
2
3use crate::error::OptimizeError;
4use crate::unconstrained::result::OptimizeResult;
5use crate::unconstrained::{Bounds, Options};
6use scirs2_core::ndarray::{Array1, ArrayView1};
7
8/// Implements Powell's method for unconstrained optimization with optional bounds support
9#[allow(dead_code)]
10pub fn minimize_powell<F, S>(
11    mut fun: F,
12    x0: Array1<f64>,
13    options: &Options,
14) -> Result<OptimizeResult<S>, OptimizeError>
15where
16    F: FnMut(&ArrayView1<f64>) -> S,
17    S: Into<f64> + Clone,
18{
19    // Get options or use defaults
20    let ftol = options.ftol;
21    let max_iter = options.max_iter;
22    let bounds = options.bounds.as_ref();
23
24    // Initialize variables
25    let n = x0.len();
26    let mut x = x0.to_owned();
27
28    // If bounds are provided, ensure x0 is within bounds
29    if let Some(bounds) = bounds {
30        bounds.project(x.as_slice_mut().unwrap());
31    }
32
33    let mut f = fun(&x.view()).into();
34
35    // Initialize the set of directions as the standard basis
36    let mut directions = Vec::with_capacity(n);
37    for i in 0..n {
38        let mut e_i = Array1::zeros(n);
39        e_i[i] = 1.0;
40        directions.push(e_i);
41    }
42
43    // Counters
44    let mut iter = 0;
45    let mut nfev = 1; // Initial function evaluation
46
47    // Powell's main loop
48    while iter < max_iter {
49        let x_old = x.clone();
50        let f_old = f;
51
52        // Keep track of the greatest function reduction
53        let mut f_reduction_max = 0.0;
54        let mut reduction_idx = 0;
55
56        // Line search along each direction
57        for (i, u) in directions.iter().enumerate().take(n) {
58            let f_before = f;
59
60            // Line search along direction u, respecting bounds
61            let (alpha, f_min) = line_search_powell(&mut fun, &x, u, f, &mut nfev, bounds);
62
63            // Update current position and function value
64            x = &x + &(alpha * u);
65            f = f_min;
66
67            // Update maximum reduction tracker
68            let reduction = f_before - f;
69            if reduction > f_reduction_max {
70                f_reduction_max = reduction;
71                reduction_idx = i;
72            }
73        }
74
75        // Compute the new direction
76        let new_dir = &x - &x_old;
77
78        // Check if the new direction is zero (happens if the point hits a bound and can't move)
79        let new_dir_norm = new_dir.mapv(|x| x.powi(2)).sum().sqrt();
80        if new_dir_norm < 1e-8 {
81            // We're likely at a bound constraint and can't make progress
82            break;
83        }
84
85        // Extra line search along the extrapolated direction
86        let (alpha, f_min) = line_search_powell(&mut fun, &x, &new_dir, f, &mut nfev, bounds);
87        x = &x + &(alpha * &new_dir);
88        f = f_min;
89
90        // Only *now* check for convergence
91        if 2.0 * (f_old - f) <= ftol * (f_old.abs() + f.abs() + 1e-10) {
92            break;
93        }
94
95        // Keep the basis full rank.
96        // If the extrapolated displacement is numerically zero we would
97        // lose a basis direction; just keep the old one instead.
98        if new_dir.iter().any(|v| v.abs() > 1e-12) {
99            // Update the set of directions by replacing the direction of greatest reduction
100            directions[reduction_idx] = new_dir;
101        }
102
103        iter += 1;
104    }
105
106    // Final check for bounds
107    if let Some(bounds) = bounds {
108        bounds.project(x.as_slice_mut().unwrap());
109    }
110
111    // Use original function for final value
112    let final_fun = fun(&x.view());
113
114    // Create and return result
115    Ok(OptimizeResult {
116        x,
117        fun: final_fun,
118        nit: iter,
119        func_evals: nfev,
120        nfev,
121        success: iter < max_iter,
122        message: if iter < max_iter {
123            "Optimization terminated successfully.".to_string()
124        } else {
125            "Maximum iterations reached.".to_string()
126        },
127        jacobian: None,
128        hessian: None,
129    })
130}
131
132/// Calculate the range of the line search parameter to respect bounds.
133///
134/// For a point x and direction p, find a_min and a_max such that:
135/// x + a * p stays within the bounds for all a in [a_min, a_max].
136#[allow(dead_code)]
137fn line_bounds(x: &Array1<f64>, direction: &Array1<f64>, bounds: Option<&Bounds>) -> (f64, f64) {
138    // If no bounds are provided, use unbounded line search
139    if bounds.is_none() {
140        return (f64::NEG_INFINITY, f64::INFINITY);
141    }
142
143    let bounds = bounds.unwrap();
144
145    // Start with unbounded range
146    let mut a_min = f64::NEG_INFINITY;
147    let mut a_max = f64::INFINITY;
148
149    // For each dimension, calculate the range restriction
150    for i in 0..x.len() {
151        let xi = x[i];
152        let pi = direction[i];
153
154        if pi.abs() < 1e-16 {
155            // No movement in this dimension
156            continue;
157        }
158
159        // Lower bound constraint: x_i + a * p_i >= lb_i
160        if let Some(lb) = bounds.lower[i] {
161            let a_lb = (lb - xi) / pi;
162            if pi > 0.0 {
163                a_min = a_min.max(a_lb);
164            } else {
165                a_max = a_max.min(a_lb);
166            }
167        }
168
169        // Upper bound constraint: x_i + a * p_i <= ub_i
170        if let Some(ub) = bounds.upper[i] {
171            let a_ub = (ub - xi) / pi;
172            if pi > 0.0 {
173                a_max = a_max.min(a_ub);
174            } else {
175                a_min = a_min.max(a_ub);
176            }
177        }
178    }
179
180    // Ensure a valid range exists
181    if a_min > a_max {
182        // No feasible movement, set range to zero
183        (0.0, 0.0)
184    } else {
185        (a_min, a_max)
186    }
187}
188
189/// Helper function for line search in Powell's method with bounds support
190/// One-dimensional minimization along `x + α·direction`.
191#[allow(dead_code)]
192fn line_search_powell<F, S>(
193    fun: &mut F,
194    x: &Array1<f64>,
195    direction: &Array1<f64>,
196    f_x: f64,
197    nfev: &mut usize,
198    bounds: Option<&Bounds>,
199) -> (f64, f64)
200where
201    F: FnMut(&ArrayView1<f64>) -> S,
202    S: Into<f64>,
203{
204    // Get bounds on the line search parameter
205    let (a_min, a_max) = line_bounds(x, direction, bounds);
206
207    // Degenerate direction ⇒ no movement
208    if direction.iter().all(|v| v.abs() <= 1e-16) {
209        return (0.0, f_x);
210    }
211
212    // helper ϕ(α)
213    let mut phi = |alpha: f64| -> f64 {
214        let y = x + &(direction * alpha);
215        // Project the point onto bounds if needed
216        let mut y_bounded = y.clone();
217        if let Some(bounds) = bounds {
218            bounds.project(y_bounded.as_slice_mut().unwrap());
219        }
220        *nfev += 1;
221        fun(&y_bounded.view()).into()
222    };
223
224    // --------------------------------------------------------------
225    // 1) Find the downhill direction and an initial bracket
226    // --------------------------------------------------------------
227    let golden = 1.618_033_988_75_f64; // φ
228    let mut step = 1.0; // Δ
229    let mut a = f64::max(0.0, a_min); // Start from 0 or a_min if it's positive
230    let mut fa = f_x; // ϕ(0)
231
232    // probe +Δ and –Δ
233    let mut b = f64::min(step, a_max);
234    let mut fb = phi(b);
235    if fb > fa {
236        // uphill → try –Δ
237        b = f64::max(-step, a_min);
238        fb = phi(b);
239        if fb > fa {
240            // no downhill yet: shrink Δ until we find one
241            for _ in 0..20 {
242                step *= 0.5;
243                if step < 1e-12 {
244                    break;
245                } // give up – extremely flat
246                b = f64::min(step, a_max);
247                fb = phi(b);
248                if fb < fa {
249                    break;
250                }
251                b = f64::max(-step, a_min);
252                fb = phi(b);
253                if fb < fa {
254                    break;
255                }
256            }
257        }
258    }
259
260    // if we *still* have no downhill point, stay put
261    if fb >= fa {
262        return (0.0, f_x);
263    }
264
265    // at this point 'a = 0' is higher, 'b' is lower
266    // grow the interval until we are uphill again
267    let mut c = f64::min(b + golden * (b - a), a_max);
268    let mut fc = phi(c);
269    for _ in 0..50 {
270        if fc > fb {
271            break;
272        } // bracket found
273        a = b;
274        fa = fb;
275        b = c;
276        fb = fc;
277        c = f64::min(b + golden * (b - a), a_max);
278        fc = phi(c);
279    }
280
281    // sort so that a < b < c
282    if a > c {
283        std::mem::swap(&mut a, &mut c);
284        std::mem::swap(&mut fa, &mut fc);
285    }
286
287    // --------------------------------------------------------------
288    // 2) Golden-section search inside the bracket
289    // --------------------------------------------------------------
290    let mut lo = a;
291    let mut hi = c;
292    let mut x1 = hi - (hi - lo) / golden;
293    let mut x2 = lo + (hi - lo) / golden;
294    let mut f1 = phi(x1);
295    let mut f2 = phi(x2);
296
297    const IT_MAX: usize = 100;
298    const TOL: f64 = 1e-8;
299
300    for _ in 0..IT_MAX {
301        if (hi - lo).abs() < TOL {
302            let alpha = 0.5 * (hi + lo);
303            return (alpha, phi(alpha)); // φ counts this eval
304        }
305        if f1 < f2 {
306            hi = x2;
307            x2 = x1;
308            f2 = f1;
309            x1 = hi - (hi - lo) / golden;
310            f1 = phi(x1);
311        } else {
312            lo = x1;
313            x1 = x2;
314            f1 = f2;
315            x2 = lo + (hi - lo) / golden;
316            f2 = phi(x2);
317        }
318    }
319
320    // fall-back: return the best of the two interior points
321    if f1 < f2 {
322        (x1, f1)
323    } else {
324        (x2, f2)
325    }
326}
327
328#[cfg(test)]
329mod tests {
330    use super::*;
331    use approx::assert_abs_diff_eq;
332
333    #[test]
334    fn test_powell_simple() {
335        let quadratic = |x: &ArrayView1<f64>| -> f64 { x[0] * x[0] + x[1] * x[1] };
336
337        let x0 = Array1::from_vec(vec![1.0, 1.0]);
338        let options = Options::default();
339
340        let result = minimize_powell(quadratic, x0, &options).unwrap();
341
342        assert!(result.success);
343        assert_abs_diff_eq!(result.x[0], 0.0, epsilon = 1e-6);
344        assert_abs_diff_eq!(result.x[1], 0.0, epsilon = 1e-6);
345    }
346
347    #[test]
348    fn test_powell_rosenbrock() {
349        let rosenbrock = |x: &ArrayView1<f64>| -> f64 {
350            let a = 1.0;
351            let b = 100.0;
352            (a - x[0]).powi(2) + b * (x[1] - x[0].powi(2)).powi(2)
353        };
354
355        let x0 = Array1::from_vec(vec![0.0, 0.0]);
356        let options = Options::default();
357
358        let result = minimize_powell(rosenbrock, x0, &options).unwrap();
359
360        assert!(result.success);
361        assert_abs_diff_eq!(result.x[0], 1.0, epsilon = 1e-3);
362        assert_abs_diff_eq!(result.x[1], 1.0, epsilon = 1e-3);
363    }
364
365    #[test]
366    fn test_powell_with_bounds() {
367        let quadratic =
368            |x: &ArrayView1<f64>| -> f64 { (x[0] - 2.0).powi(2) + (x[1] - 3.0).powi(2) };
369
370        let x0 = Array1::from_vec(vec![0.0, 0.0]);
371        let mut options = Options::default();
372
373        // Constrain solution to [0, 1] x [0, 1]
374        let bounds = Bounds::new(&[(Some(0.0), Some(1.0)), (Some(0.0), Some(1.0))]);
375        options.bounds = Some(bounds);
376
377        let result = minimize_powell(quadratic, x0, &options).unwrap();
378
379        assert!(result.success);
380        // The optimal point (2, 3) is outside the bounds, so we should get (1, 1)
381        assert_abs_diff_eq!(result.x[0], 1.0, epsilon = 1e-6);
382        assert_abs_diff_eq!(result.x[1], 1.0, epsilon = 1e-6);
383    }
384}