Skip to main content

scirs2_optimize/unconstrained/
mod.rs

1//! Unconstrained optimization algorithms
2//!
3//! This module provides various algorithms for unconstrained minimization problems.
4
5use crate::error::OptimizeError;
6use scirs2_core::ndarray::{Array1, Array2, ArrayView1};
7use std::fmt;
8
9// Sub-modules
10pub mod adaptive_convergence;
11pub mod advanced_line_search;
12pub mod bfgs;
13pub mod callback_diagnostics;
14pub mod conjugate_gradient;
15pub mod convergence_diagnostics;
16pub mod efficient_sparse;
17pub mod lbfgs;
18pub mod line_search;
19pub mod memory_efficient;
20pub mod memory_efficient_sparse;
21pub mod nelder_mead;
22pub mod newton;
23pub mod powell;
24pub mod quasi_newton;
25pub mod result;
26pub mod robust_convergence;
27pub mod simd_bfgs;
28pub mod sparse_optimization;
29pub mod strong_wolfe;
30pub mod subspace_methods;
31pub mod truncated_newton;
32pub mod trust_region;
33pub mod utils;
34
35// Import result type
36pub use result::OptimizeResult;
37
38// Re-export commonly used items
39pub use adaptive_convergence::{
40    check_convergence_adaptive, create_adaptive_options_for_problem, AdaptationStats,
41    AdaptiveToleranceOptions, AdaptiveToleranceState, ConvergenceStatus,
42};
43pub use advanced_line_search::{
44    advanced_line_search, create_non_monotone_state, AdvancedLineSearchOptions,
45    InterpolationStrategy, LineSearchMethod, LineSearchResult, LineSearchStats,
46};
47pub use bfgs::{minimize_bfgs, minimize_bfgs_no_grad};
48pub use callback_diagnostics::{
49    minimize_with_diagnostics, optimize_with_diagnostics, CallbackInfo, CallbackResult,
50    DiagnosticOptimizer, OptimizationCallback,
51};
52pub use conjugate_gradient::minimize_conjugate_gradient;
53pub use convergence_diagnostics::{
54    ConvergenceDiagnostics, DiagnosticCollector, DiagnosticOptions, DiagnosticWarning,
55    ExportFormat, IterationDiagnostic, LineSearchDiagnostic, PerformanceMetrics, ProblemAnalysis,
56    ProblemDifficulty, WarningSeverity,
57};
58pub use efficient_sparse::{
59    minimize_efficient_sparse_newton, EfficientSparseOptions, SparsityInfo,
60};
61pub use lbfgs::{minimize_lbfgs, minimize_lbfgsb};
62pub use memory_efficient::{
63    create_memory_efficient_optimizer, minimize_memory_efficient_lbfgs, MemoryOptions,
64};
65pub use memory_efficient_sparse::{
66    create_advanced_scale_optimizer, minimize_advanced_scale, AdvancedScaleOptions,
67};
68pub use nelder_mead::minimize_nelder_mead;
69pub use newton::minimize_newton_cg;
70pub use powell::minimize_powell;
71pub use quasi_newton::{minimize_dfp, minimize_quasi_newton, minimize_sr1, UpdateFormula};
72pub use robust_convergence::{
73    create_robust_options_for_problem, RobustConvergenceOptions, RobustConvergenceResult,
74    RobustConvergenceState,
75};
76pub use simd_bfgs::{minimize_simd_bfgs, minimize_simd_bfgs_default, SimdBfgsOptions};
77pub use sparse_optimization::{
78    auto_detect_sparsity, minimize_sparse_bfgs, SparseOptimizationOptions,
79};
80pub use strong_wolfe::{
81    create_strong_wolfe_options_for_method, strong_wolfe_line_search, StrongWolfeOptions,
82    StrongWolfeResult,
83};
84pub use subspace_methods::{
85    minimize_adaptive_subspace, minimize_block_coordinate_descent,
86    minimize_cyclical_coordinate_descent, minimize_random_coordinate_descent,
87    minimize_random_subspace, minimize_subspace, SubspaceMethod, SubspaceOptions,
88};
89pub use truncated_newton::{
90    minimize_truncated_newton, minimize_trust_region_newton, Preconditioner, TruncatedNewtonOptions,
91};
92pub use trust_region::{
93    cauchy_point, dogleg_step, minimize_trust_exact, minimize_trust_krylov, minimize_trust_ncg,
94    solve_trust_subproblem, trust_region_minimize, TrustRegionConfig, TrustRegionResult,
95};
96
97/// Optimization methods for unconstrained minimization.
98#[derive(Debug, Clone, Copy)]
99pub enum Method {
100    /// Nelder-Mead simplex method
101    NelderMead,
102    /// Powell's method
103    Powell,
104    /// Conjugate gradient method
105    CG,
106    /// BFGS quasi-Newton method
107    BFGS,
108    /// SR1 quasi-Newton method
109    SR1,
110    /// DFP quasi-Newton method
111    DFP,
112    /// Limited-memory BFGS method
113    LBFGS,
114    /// Limited-memory BFGS method with bounds support
115    LBFGSB,
116    /// Newton's method with conjugate gradient solver
117    NewtonCG,
118    /// Trust-region Newton method with conjugate gradient solver
119    TrustNCG,
120    /// Trust-region method with Krylov subproblem solver
121    TrustKrylov,
122    /// Trust-region method with exact subproblem solver
123    TrustExact,
124    /// Truncated Newton method
125    TruncatedNewton,
126    /// Trust-region Newton method with truncated CG
127    TrustRegionNewton,
128    /// Trust-region dogleg method (Cauchy point + dogleg step)
129    TrustRegionDogleg,
130}
131
132/// Bounds for optimization variables
133#[derive(Debug, Clone)]
134pub struct Bounds {
135    /// Lower bounds
136    pub lower: Vec<Option<f64>>,
137    /// Upper bounds
138    pub upper: Vec<Option<f64>>,
139}
140
141/// Options for optimization algorithms
142#[derive(Debug, Clone)]
143pub struct Options {
144    /// Maximum number of iterations
145    pub max_iter: usize,
146    /// Maximum number of function evaluations
147    pub max_fev: Option<usize>,
148    /// Function tolerance for convergence
149    pub ftol: f64,
150    /// Change tolerance for convergence
151    pub xtol: f64,
152    /// Gradient tolerance for convergence
153    pub gtol: f64,
154    /// Initial step size
155    pub initial_step: Option<f64>,
156    /// Maximum step size for line search
157    pub maxstep: Option<f64>,
158    /// Whether to use finite differences for gradient
159    pub finite_diff: bool,
160    /// Finite difference step size
161    pub eps: f64,
162    /// Initial trust-region radius for trust-region methods
163    pub trust_radius: Option<f64>,
164    /// Maximum trust-region radius for trust-region methods
165    pub max_trust_radius: Option<f64>,
166    /// Minimum trust-region radius for trust-region methods
167    pub min_trust_radius: Option<f64>,
168    /// Tolerance for the trust-region subproblem
169    pub trust_tol: Option<f64>,
170    /// Maximum iterations for trust-region subproblem
171    pub trust_max_iter: Option<usize>,
172    /// Threshold for accepting a step in the trust-region method
173    pub trust_eta: Option<f64>,
174    /// Bounds constraints for variables
175    pub bounds: Option<Bounds>,
176}
177
178// Implement Display for Method
179impl fmt::Display for Method {
180    fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
181        match self {
182            Method::NelderMead => write!(f, "Nelder-Mead"),
183            Method::Powell => write!(f, "Powell"),
184            Method::CG => write!(f, "Conjugate Gradient"),
185            Method::BFGS => write!(f, "BFGS"),
186            Method::SR1 => write!(f, "SR1"),
187            Method::DFP => write!(f, "DFP"),
188            Method::LBFGS => write!(f, "L-BFGS"),
189            Method::LBFGSB => write!(f, "L-BFGS-B"),
190            Method::NewtonCG => write!(f, "Newton-CG"),
191            Method::TrustNCG => write!(f, "Trust-NCG"),
192            Method::TrustKrylov => write!(f, "Trust-Krylov"),
193            Method::TrustExact => write!(f, "Trust-Exact"),
194            Method::TruncatedNewton => write!(f, "Truncated Newton"),
195            Method::TrustRegionNewton => write!(f, "Trust-Region Newton"),
196            Method::TrustRegionDogleg => write!(f, "Trust-Region Dogleg"),
197        }
198    }
199}
200
201// Implement Default for Options
202impl Default for Options {
203    fn default() -> Self {
204        Options {
205            max_iter: 1000,
206            max_fev: None,
207            ftol: 1e-8,
208            xtol: 1e-8,
209            gtol: 1e-5,
210            initial_step: None,
211            maxstep: None,
212            finite_diff: false,
213            eps: 1.4901161193847656e-8,
214            trust_radius: Some(1.0),
215            max_trust_radius: Some(100.0),
216            min_trust_radius: Some(1e-10),
217            trust_tol: Some(1e-8),
218            trust_max_iter: Some(100),
219            trust_eta: Some(0.1),
220            bounds: None,
221        }
222    }
223}
224
225// Implement Bounds methods
226impl Bounds {
227    /// Create new bounds from arrays
228    pub fn new(bounds: &[(Option<f64>, Option<f64>)]) -> Self {
229        let (lower, upper): (Vec<_>, Vec<_>) = bounds.iter().cloned().unzip();
230        Self { lower, upper }
231    }
232
233    /// Create bounds from vectors
234    pub fn from_vecs(lb: Vec<Option<f64>>, ub: Vec<Option<f64>>) -> Result<Self, OptimizeError> {
235        if lb.len() != ub.len() {
236            return Err(OptimizeError::ValueError(
237                "Lower and upper bounds must have the same length".to_string(),
238            ));
239        }
240
241        for (l, u) in lb.iter().zip(ub.iter()) {
242            if let (Some(l_val), Some(u_val)) = (l, u) {
243                if l_val > u_val {
244                    return Err(OptimizeError::ValueError(
245                        "Lower bound must be less than or equal to upper bound".to_string(),
246                    ));
247                }
248            }
249        }
250
251        Ok(Self {
252            lower: lb,
253            upper: ub,
254        })
255    }
256
257    /// Check if point is feasible
258    pub fn is_feasible(&self, x: &[f64]) -> bool {
259        if x.len() != self.lower.len() {
260            return false;
261        }
262
263        for (&xi, (&lb, &ub)) in x.iter().zip(self.lower.iter().zip(self.upper.iter())) {
264            if let Some(l) = lb {
265                if xi < l {
266                    return false;
267                }
268            }
269            if let Some(u) = ub {
270                if xi > u {
271                    return false;
272                }
273            }
274        }
275        true
276    }
277
278    /// Project point onto feasible region
279    pub fn project(&self, x: &mut [f64]) {
280        for (xi, (&lb, &ub)) in x.iter_mut().zip(self.lower.iter().zip(self.upper.iter())) {
281            if let Some(l) = lb {
282                if *xi < l {
283                    *xi = l;
284                }
285            }
286            if let Some(u) = ub {
287                if *xi > u {
288                    *xi = u;
289                }
290            }
291        }
292    }
293
294    /// Check if bounds are active
295    pub fn has_bounds(&self) -> bool {
296        self.lower.iter().any(|b| b.is_some()) || self.upper.iter().any(|b| b.is_some())
297    }
298}
299
300/// Main minimize function for unconstrained optimization
301#[allow(dead_code)]
302pub fn minimize<F, S>(
303    fun: F,
304    x0: &[f64],
305    method: Method,
306    options: Option<Options>,
307) -> Result<OptimizeResult<S>, OptimizeError>
308where
309    F: FnMut(&ArrayView1<f64>) -> S + Clone,
310    S: Into<f64> + Clone + From<f64>,
311{
312    let options = &options.unwrap_or_default();
313    let x0 = Array1::from_vec(x0.to_vec());
314
315    // Check initial point feasibility if bounds are provided
316    if let Some(ref bounds) = options.bounds {
317        if !bounds.is_feasible(x0.as_slice().expect("Operation failed")) {
318            return Err(OptimizeError::ValueError(
319                "Initial point is not feasible".to_string(),
320            ));
321        }
322    }
323
324    match method {
325        Method::NelderMead => nelder_mead::minimize_nelder_mead(fun, x0, options),
326        Method::Powell => powell::minimize_powell(fun, x0, options),
327        Method::CG => conjugate_gradient::minimize_conjugate_gradient(
328            fun,
329            None::<fn(&ArrayView1<f64>) -> Array1<f64>>,
330            x0,
331            options,
332        ),
333        Method::BFGS => bfgs::minimize_bfgs(
334            fun,
335            None::<fn(&ArrayView1<f64>) -> Array1<f64>>,
336            x0,
337            options,
338        ),
339        Method::SR1 => quasi_newton::minimize_sr1(fun, x0, options),
340        Method::DFP => quasi_newton::minimize_dfp(fun, x0, options),
341        Method::LBFGS => lbfgs::minimize_lbfgs(
342            fun,
343            None::<fn(&ArrayView1<f64>) -> Array1<f64>>,
344            x0,
345            options,
346        ),
347        Method::LBFGSB => lbfgs::minimize_lbfgsb(
348            fun,
349            None::<fn(&ArrayView1<f64>) -> Array1<f64>>,
350            x0,
351            options,
352        ),
353        Method::NewtonCG => newton::minimize_newton_cg(fun, x0, options),
354        Method::TrustNCG => trust_region::minimize_trust_ncg(fun, x0, options),
355        Method::TrustKrylov => trust_region::minimize_trust_krylov(fun, x0, options),
356        Method::TrustExact => trust_region::minimize_trust_exact(fun, x0, options),
357        Method::TruncatedNewton => truncated_newton_wrapper(fun, x0, options),
358        Method::TrustRegionNewton => trust_region_newton_wrapper(fun, x0, options),
359        Method::TrustRegionDogleg => trust_region_dogleg_wrapper(fun, x0, options),
360    }
361}
362
363/// Wrapper function for truncated Newton method
364#[allow(dead_code)]
365fn truncated_newton_wrapper<F, S>(
366    mut fun: F,
367    x0: Array1<f64>,
368    options: &Options,
369) -> Result<OptimizeResult<S>, OptimizeError>
370where
371    F: FnMut(&ArrayView1<f64>) -> S + Clone,
372    S: Into<f64> + Clone + From<f64>,
373{
374    let fun_f64 = move |x: &ArrayView1<f64>| fun(x).into();
375
376    let truncated_options = TruncatedNewtonOptions {
377        max_iter: options.max_iter,
378        tol: options.gtol,
379        max_cg_iter: options.trust_max_iter.unwrap_or(100),
380        cg_tol: options.trust_tol.unwrap_or(0.1),
381        ..Default::default()
382    };
383
384    // Convert result back to generic type
385    let result = truncated_newton::minimize_truncated_newton(
386        fun_f64,
387        None::<fn(&ArrayView1<f64>) -> Array1<f64>>,
388        x0,
389        Some(truncated_options),
390    )?;
391
392    Ok(OptimizeResult {
393        x: result.x,
394        fun: result.fun.into(),
395        nit: result.nit,
396        func_evals: result.func_evals,
397        nfev: result.nfev,
398        jacobian: result.jacobian,
399        hessian: result.hessian,
400        success: result.success,
401        message: result.message,
402    })
403}
404
405/// Wrapper function for trust-region Newton method
406#[allow(dead_code)]
407fn trust_region_newton_wrapper<F, S>(
408    mut fun: F,
409    x0: Array1<f64>,
410    options: &Options,
411) -> Result<OptimizeResult<S>, OptimizeError>
412where
413    F: FnMut(&ArrayView1<f64>) -> S + Clone,
414    S: Into<f64> + Clone + From<f64>,
415{
416    let fun_f64 = move |x: &ArrayView1<f64>| fun(x).into();
417
418    let truncated_options = TruncatedNewtonOptions {
419        max_iter: options.max_iter,
420        tol: options.gtol,
421        max_cg_iter: options.trust_max_iter.unwrap_or(100),
422        cg_tol: options.trust_tol.unwrap_or(0.1),
423        trust_radius: options.trust_radius,
424        ..Default::default()
425    };
426
427    // Convert result back to generic type
428    let result = truncated_newton::minimize_trust_region_newton(
429        fun_f64,
430        None::<fn(&ArrayView1<f64>) -> Array1<f64>>,
431        x0,
432        Some(truncated_options),
433    )?;
434
435    Ok(OptimizeResult {
436        x: result.x,
437        fun: result.fun.into(),
438        nit: result.nit,
439        func_evals: result.func_evals,
440        nfev: result.nfev,
441        jacobian: result.jacobian,
442        hessian: result.hessian,
443        success: result.success,
444        message: result.message,
445    })
446}
447
448/// Wrapper function for trust-region dogleg method
449#[allow(dead_code)]
450fn trust_region_dogleg_wrapper<F, S>(
451    mut fun: F,
452    x0: Array1<f64>,
453    options: &Options,
454) -> Result<OptimizeResult<S>, OptimizeError>
455where
456    F: FnMut(&ArrayView1<f64>) -> S + Clone,
457    S: Into<f64> + Clone + From<f64>,
458{
459    let mut fun_f64 = move |x: &ArrayView1<f64>| fun(x).into();
460
461    let config = TrustRegionConfig {
462        initial_radius: options.trust_radius.unwrap_or(1.0),
463        max_radius: options.max_trust_radius.unwrap_or(100.0),
464        max_iter: options.max_iter,
465        tolerance: options.gtol,
466        ftol: options.ftol,
467        eps: options.eps,
468        min_radius: options.min_trust_radius.unwrap_or(1e-14),
469        ..Default::default()
470    };
471
472    let result = trust_region::trust_region_minimize(
473        &mut fun_f64,
474        None::<fn(&ArrayView1<f64>) -> Array1<f64>>,
475        None::<fn(&ArrayView1<f64>) -> Array2<f64>>,
476        x0,
477        Some(config),
478    )?;
479
480    Ok(OptimizeResult {
481        x: result.x,
482        fun: S::from(result.f_val),
483        nit: result.n_iter,
484        func_evals: result.n_fev,
485        nfev: result.n_fev,
486        jacobian: None,
487        hessian: None,
488        success: result.converged,
489        message: result.message,
490    })
491}
492
493#[cfg(test)]
494mod tests {
495    use super::*;
496    use approx::assert_abs_diff_eq;
497
498    #[test]
499    fn test_simple_quadratic() {
500        let quadratic = |x: &ArrayView1<f64>| -> f64 { x[0] * x[0] + x[1] * x[1] };
501
502        let x0 = vec![1.0, 1.0];
503        let result = minimize(quadratic, &x0, Method::BFGS, None);
504        assert!(result.is_ok());
505
506        let result = result.expect("Operation failed");
507        assert!(result.success);
508        assert_abs_diff_eq!(result.x[0], 0.0, epsilon = 1e-6);
509        assert_abs_diff_eq!(result.x[1], 0.0, epsilon = 1e-6);
510    }
511}