Skip to main content

numra_optim/
auto.rs

1//! Automatic solver selection for optimization problems.
2//!
3//! Author: Moussa Leblouba
4//! Date: 5 March 2026
5//! Modified: 2 May 2026
6
7use std::sync::Arc;
8
9use numra_core::Scalar;
10
11use crate::augmented_lagrangian::augmented_lagrangian_minimize;
12use crate::bfgs::bfgs_minimize;
13use crate::error::OptimError;
14use crate::lbfgs::{lbfgs_minimize, LbfgsOptions};
15use crate::lbfgsb::{lbfgsb_minimize, LbfgsBOptions};
16use crate::levenberg_marquardt::{lm_minimize, LmOptions};
17use crate::problem::{
18    finite_diff_gradient, finite_diff_jacobian, ConstraintKind, ObjectiveKind, OptimProblem,
19    VarType,
20};
21use crate::types::OptimResult;
22
23/// Explicit solver choice for `OptimProblem::solve_with`.
24#[derive(Clone, Copy, Debug, PartialEq, Eq)]
25pub enum SolverChoice {
26    Simplex,
27    ActiveSetQP,
28    MixedIntegerLP,
29    DifferentialEvolution,
30    NsgaII,
31    Bfgs,
32    Lbfgs,
33    LbfgsB,
34    LevenbergMarquardt,
35    AugmentedLagrangian,
36    NelderMead,
37    Powell,
38    CmaEs,
39    Sqp,
40}
41
42/// Select the best solver for the given problem structure.
43pub fn select_solver<S: Scalar>(problem: &OptimProblem<S>) -> SolverChoice {
44    if problem.is_multi_objective() {
45        return SolverChoice::NsgaII;
46    }
47    if problem.global {
48        return SolverChoice::DifferentialEvolution;
49    }
50
51    // Structural variants: exact coefficient access available
52    if problem.is_linear() && problem.has_integer_vars() {
53        return SolverChoice::MixedIntegerLP;
54    }
55    if problem.is_linear() {
56        return SolverChoice::Simplex;
57    }
58    if problem.is_quadratic() {
59        return SolverChoice::ActiveSetQP;
60    }
61
62    // Hint-based dispatch: closure is known to be linear/quadratic but
63    // coefficients aren't available, so use NLP solvers that work well
64    // on these problem classes.
65    if problem.is_hint_linear() || problem.is_hint_quadratic() {
66        if problem.has_constraints() {
67            return SolverChoice::Sqp;
68        }
69        return SolverChoice::LbfgsB;
70    }
71
72    if problem.is_least_squares() && !problem.has_constraints() && !problem.has_bounds() {
73        SolverChoice::LevenbergMarquardt
74    } else if problem.has_constraints() {
75        SolverChoice::AugmentedLagrangian
76    } else if problem.has_bounds() {
77        SolverChoice::LbfgsB
78    } else {
79        SolverChoice::Lbfgs
80    }
81}
82
83/// Round integer and binary variables to their nearest feasible integer values.
84///
85/// Note: `result.f` is NOT re-evaluated (the objective closure may have been consumed).
86/// The caller (bridge) can re-evaluate if needed.
87fn round_integer_vars<S: Scalar>(x: &mut [S], var_types: &[VarType], bounds: &[Option<(S, S)>]) {
88    for (i, vt) in var_types.iter().enumerate() {
89        match vt {
90            VarType::Integer => {
91                let rounded = x[i].round();
92                x[i] = match bounds.get(i).copied().flatten() {
93                    Some((lo, hi)) => rounded.max(lo).min(hi),
94                    None => rounded,
95                };
96            }
97            VarType::Binary => {
98                x[i] = if x[i] > S::from_f64(0.5) {
99                    S::ONE
100                } else {
101                    S::ZERO
102                };
103            }
104            VarType::Continuous => {}
105        }
106    }
107}
108
109/// Solve a problem using automatic solver selection.
110pub fn auto_minimize<
111    S: Scalar + faer::SimpleEntity + faer::Conjugate<Canonical = S> + faer::ComplexField,
112>(
113    problem: OptimProblem<S>,
114) -> Result<OptimResult<S>, OptimError> {
115    let negate = problem.negate_result;
116    let has_int = problem.has_integer_vars();
117    let var_types = problem.var_types.clone();
118    let bounds = problem.bounds.clone();
119    let choice = select_solver(&problem);
120    let mut result = dispatch(problem, choice)?;
121
122    // Post-solve integer rounding when NLP solver was used for integer problems
123    if has_int && choice != SolverChoice::MixedIntegerLP {
124        round_integer_vars(&mut result.x, &var_types, &bounds);
125        result
126            .message
127            .push_str(" (integer vars rounded post-solve; f is pre-round value)");
128    }
129
130    if negate {
131        result.f = -result.f;
132        for gi in result.grad.iter_mut() {
133            *gi = -*gi;
134        }
135        // Also negate objective in history records
136        for rec in result.history.iter_mut() {
137            rec.objective = -rec.objective;
138        }
139    }
140    Ok(result)
141}
142
143/// Dispatch to a specific solver.
144#[allow(clippy::type_complexity)]
145pub fn dispatch<
146    S: Scalar + faer::SimpleEntity + faer::Conjugate<Canonical = S> + faer::ComplexField,
147>(
148    problem: OptimProblem<S>,
149    choice: SolverChoice,
150) -> Result<OptimResult<S>, OptimError> {
151    // Augmented Lagrangian consumes the whole problem
152    if choice == SolverChoice::AugmentedLagrangian {
153        let aug_opts = problem.aug_lag_options.clone().unwrap_or_default();
154        return augmented_lagrangian_minimize(problem, &aug_opts);
155    }
156
157    // SQP consumes the whole problem (it needs the constraints)
158    if choice == SolverChoice::Sqp {
159        let x0 = problem.x0.ok_or(OptimError::NoInitialPoint)?;
160        let objective = problem.objective.ok_or(OptimError::NoObjective)?;
161        let (func, grad) = match objective {
162            ObjectiveKind::Minimize { func, grad } => (func, grad),
163            _ => return Err(OptimError::Other("SQP requires scalar objective".into())),
164        };
165        let func_arc: Arc<dyn Fn(&[S]) -> S> = Arc::from(func);
166        let grad_fn: Box<dyn Fn(&[S], &mut [S])> = match grad {
167            Some(g) => g,
168            None => {
169                let fa = Arc::clone(&func_arc);
170                Box::new(move |x: &[S], g: &mut [S]| {
171                    finite_diff_gradient(&*fa, x, g);
172                })
173            }
174        };
175        let sqp_opts = crate::sqp::SqpOptions {
176            max_iter: problem.options.max_iter,
177            tol: problem.options.gtol,
178            ..Default::default()
179        };
180        return crate::sqp::sqp_minimize(
181            &*func_arc,
182            &*grad_fn,
183            &problem.constraints,
184            &x0,
185            &sqp_opts,
186        );
187    }
188
189    // Destructure problem
190    let OptimProblem {
191        x0,
192        bounds,
193        var_types,
194        objective,
195        options: opts,
196        linear_constraints,
197        global: _global,
198        de_options,
199        ..
200    } = problem;
201
202    // DE and NSGA-II don't require x0 -- handle them before the x0 check
203    if choice == SolverChoice::DifferentialEvolution {
204        let objective = objective.ok_or(OptimError::NoObjective)?;
205        let func = match objective {
206            ObjectiveKind::Minimize { func, .. } => func,
207            _ => return Err(OptimError::Other("DE requires scalar objective".into())),
208        };
209        let de_bounds: Vec<(S, S)> = bounds
210            .iter()
211            .map(|b: &Option<(S, S)>| {
212                b.ok_or(OptimError::Other(
213                    "DE requires all variables to be bounded".into(),
214                ))
215            })
216            .collect::<Result<_, _>>()?;
217        let de_opts = de_options.unwrap_or_default();
218        return crate::global::de_minimize(&*func, &de_bounds, &de_opts);
219    }
220
221    if choice == SolverChoice::NsgaII {
222        let objective = objective.ok_or(OptimError::NoObjective)?;
223        let funcs = match objective {
224            ObjectiveKind::MultiObjective { funcs } => funcs,
225            _ => return Err(OptimError::Other("NSGA-II requires multi-objective".into())),
226        };
227        let nsga_bounds: Vec<(S, S)> = bounds
228            .iter()
229            .map(|b: &Option<(S, S)>| {
230                b.ok_or(OptimError::Other(
231                    "NSGA-II requires all variables to be bounded".into(),
232                ))
233            })
234            .collect::<Result<_, _>>()?;
235        let obj_refs: Vec<&dyn Fn(&[S]) -> S> =
236            funcs.iter().map(|f| &**f as &dyn Fn(&[S]) -> S).collect();
237        let nsga_opts = crate::multiobjective::NsgaIIOptions::default();
238        return crate::multiobjective::nsga2_optimize(&obj_refs, &nsga_bounds, &nsga_opts);
239    }
240
241    // LP, QP, and MILP solvers don't require x0 -- handle them before the x0 check
242    if choice == SolverChoice::MixedIntegerLP {
243        let objective = objective.ok_or(OptimError::NoObjective)?;
244        let c = match objective {
245            ObjectiveKind::Linear { c } => c,
246            _ => return Err(OptimError::Other("MILP requires linear objective".into())),
247        };
248
249        let mut a_ineq: Vec<Vec<S>> = Vec::new();
250        let mut b_ineq_vec: Vec<S> = Vec::new();
251        let mut a_eq: Vec<Vec<S>> = Vec::new();
252        let mut b_eq_vec: Vec<S> = Vec::new();
253
254        for lc in &linear_constraints {
255            match lc.kind {
256                ConstraintKind::Inequality => {
257                    a_ineq.push(lc.a.clone());
258                    b_ineq_vec.push(lc.b);
259                }
260                ConstraintKind::Equality => {
261                    a_eq.push(lc.a.clone());
262                    b_eq_vec.push(lc.b);
263                }
264            }
265        }
266
267        let min_lp_tol = S::from_f64(1e-10);
268        let milp_opts = crate::milp::MILPOptions {
269            max_nodes: opts.max_iter,
270            lp_tol: if opts.gtol > min_lp_tol {
271                opts.gtol
272            } else {
273                min_lp_tol
274            },
275            verbose: opts.verbose,
276            ..Default::default()
277        };
278        return crate::milp::milp_solve(
279            &c,
280            &a_ineq,
281            &b_ineq_vec,
282            &a_eq,
283            &b_eq_vec,
284            &var_types,
285            &bounds,
286            &milp_opts,
287        );
288    }
289
290    if choice == SolverChoice::Simplex {
291        let objective = objective.ok_or(OptimError::NoObjective)?;
292        let c = match objective {
293            ObjectiveKind::Linear { c } => c,
294            _ => {
295                return Err(OptimError::Other(
296                    "Simplex requires linear objective".into(),
297                ))
298            }
299        };
300
301        let mut a_ineq: Vec<Vec<S>> = Vec::new();
302        let mut b_ineq_vec: Vec<S> = Vec::new();
303        let mut a_eq: Vec<Vec<S>> = Vec::new();
304        let mut b_eq_vec: Vec<S> = Vec::new();
305
306        for lc in &linear_constraints {
307            match lc.kind {
308                ConstraintKind::Inequality => {
309                    a_ineq.push(lc.a.clone());
310                    b_ineq_vec.push(lc.b);
311                }
312                ConstraintKind::Equality => {
313                    a_eq.push(lc.a.clone());
314                    b_eq_vec.push(lc.b);
315                }
316            }
317        }
318
319        let min_lp_tol = S::from_f64(1e-10);
320        let lp_opts = crate::lp::LPOptions {
321            max_iter: opts.max_iter,
322            tol: if opts.gtol > min_lp_tol {
323                opts.gtol
324            } else {
325                min_lp_tol
326            },
327            verbose: opts.verbose,
328        };
329        return crate::lp::simplex_solve(&c, &a_ineq, &b_ineq_vec, &a_eq, &b_eq_vec, &lp_opts);
330    }
331
332    if choice == SolverChoice::ActiveSetQP {
333        let objective = objective.ok_or(OptimError::NoObjective)?;
334        let (h_row_major, qp_c, qp_n) = match objective {
335            ObjectiveKind::Quadratic { h_row_major, c, n } => (h_row_major, c, n),
336            _ => {
337                return Err(OptimError::Other(
338                    "ActiveSetQP requires quadratic objective".into(),
339                ))
340            }
341        };
342
343        let mut a_ineq: Vec<Vec<S>> = Vec::new();
344        let mut b_ineq_vec: Vec<S> = Vec::new();
345        let mut a_eq: Vec<Vec<S>> = Vec::new();
346        let mut b_eq_vec: Vec<S> = Vec::new();
347
348        for lc in &linear_constraints {
349            match lc.kind {
350                ConstraintKind::Inequality => {
351                    a_ineq.push(lc.a.clone());
352                    b_ineq_vec.push(lc.b);
353                }
354                ConstraintKind::Equality => {
355                    a_eq.push(lc.a.clone());
356                    b_eq_vec.push(lc.b);
357                }
358            }
359        }
360
361        let min_qp_tol = S::from_f64(1e-10);
362        let qp_opts = crate::qp::QPOptions {
363            max_iter: opts.max_iter,
364            tol: if opts.gtol > min_qp_tol {
365                opts.gtol
366            } else {
367                min_qp_tol
368            },
369            verbose: opts.verbose,
370        };
371        return crate::qp::active_set_qp_solve(
372            &h_row_major,
373            &qp_c,
374            qp_n,
375            &a_ineq,
376            &b_ineq_vec,
377            &a_eq,
378            &b_eq_vec,
379            &bounds,
380            &qp_opts,
381        );
382    }
383
384    // Derivative-free and CMA-ES require x0 but handle objective extraction differently
385    if choice == SolverChoice::NelderMead {
386        let x0 = x0.ok_or(OptimError::NoInitialPoint)?;
387        let objective = objective.ok_or(OptimError::NoObjective)?;
388        let func = match objective {
389            ObjectiveKind::Minimize { func, .. } => func,
390            _ => {
391                return Err(OptimError::Other(
392                    "Nelder-Mead requires scalar objective".into(),
393                ))
394            }
395        };
396        let nm_opts = crate::derivative_free::NelderMeadOptions {
397            max_iter: opts.max_iter,
398            ..Default::default()
399        };
400        return crate::derivative_free::nelder_mead(&*func, &x0, &nm_opts);
401    }
402
403    if choice == SolverChoice::Powell {
404        let x0 = x0.ok_or(OptimError::NoInitialPoint)?;
405        let objective = objective.ok_or(OptimError::NoObjective)?;
406        let func = match objective {
407            ObjectiveKind::Minimize { func, .. } => func,
408            _ => return Err(OptimError::Other("Powell requires scalar objective".into())),
409        };
410        let p_opts = crate::derivative_free::PowellOptions {
411            max_iter: opts.max_iter,
412            ..Default::default()
413        };
414        return crate::derivative_free::powell(&*func, &x0, &p_opts);
415    }
416
417    if choice == SolverChoice::CmaEs {
418        let x0 = x0.ok_or(OptimError::NoInitialPoint)?;
419        let objective = objective.ok_or(OptimError::NoObjective)?;
420        let func = match objective {
421            ObjectiveKind::Minimize { func, .. } => func,
422            _ => return Err(OptimError::Other("CMA-ES requires scalar objective".into())),
423        };
424        let cma_opts = crate::cmaes::CmaEsOptions {
425            max_iter: opts.max_iter,
426            ..Default::default()
427        };
428        return crate::cmaes::cmaes_minimize(&*func, &x0, &cma_opts);
429    }
430
431    let x0 = x0.ok_or(OptimError::NoInitialPoint)?;
432    let objective = objective.ok_or(OptimError::NoObjective)?;
433
434    match choice {
435        SolverChoice::LevenbergMarquardt => {
436            let (residual, jacobian, m) = match objective {
437                ObjectiveKind::LeastSquares {
438                    residual,
439                    jacobian,
440                    n_residuals,
441                } => (residual, jacobian, n_residuals),
442                _ => {
443                    return Err(OptimError::Other(
444                        "LM requires least-squares objective".into(),
445                    ))
446                }
447            };
448
449            let res_arc: Arc<dyn Fn(&[S], &mut [S])> = Arc::from(residual);
450            let jac_fn: Box<dyn Fn(&[S], &mut [S])> = match jacobian {
451                Some(j) => j,
452                None => {
453                    let r = Arc::clone(&res_arc);
454                    Box::new(move |x: &[S], jac: &mut [S]| {
455                        finite_diff_jacobian(&*r, x, m, jac);
456                    })
457                }
458            };
459
460            let lm_opts = LmOptions {
461                max_iter: opts.max_iter,
462                gtol: opts.gtol,
463                xtol: opts.xtol,
464                ftol: opts.ftol,
465                ..LmOptions::default()
466            };
467            lm_minimize(&*res_arc, &*jac_fn, &x0, m, &lm_opts)
468        }
469
470        SolverChoice::Bfgs => {
471            let (func, grad) = match objective {
472                ObjectiveKind::Minimize { func, grad } => (func, grad),
473                _ => return Err(OptimError::Other("BFGS requires scalar objective".into())),
474            };
475
476            let func_arc: Arc<dyn Fn(&[S]) -> S> = Arc::from(func);
477            let grad_fn: Box<dyn Fn(&[S], &mut [S])> = match grad {
478                Some(g) => g,
479                None => {
480                    let f = Arc::clone(&func_arc);
481                    Box::new(move |x: &[S], g: &mut [S]| {
482                        finite_diff_gradient(&*f, x, g);
483                    })
484                }
485            };
486
487            bfgs_minimize(&*func_arc, &*grad_fn, &x0, &opts)
488        }
489
490        SolverChoice::Lbfgs => {
491            let (func, grad) = match objective {
492                ObjectiveKind::Minimize { func, grad } => (func, grad),
493                _ => return Err(OptimError::Other("L-BFGS requires scalar objective".into())),
494            };
495
496            let func_arc: Arc<dyn Fn(&[S]) -> S> = Arc::from(func);
497            let grad_fn: Box<dyn Fn(&[S], &mut [S])> = match grad {
498                Some(g) => g,
499                None => {
500                    let f = Arc::clone(&func_arc);
501                    Box::new(move |x: &[S], g: &mut [S]| {
502                        finite_diff_gradient(&*f, x, g);
503                    })
504                }
505            };
506
507            let lbfgs_opts = LbfgsOptions {
508                base: opts,
509                memory: 10,
510            };
511            lbfgs_minimize(&*func_arc, &*grad_fn, &x0, &lbfgs_opts)
512        }
513
514        SolverChoice::LbfgsB => {
515            let (func, grad) = match objective {
516                ObjectiveKind::Minimize { func, grad } => (func, grad),
517                _ => {
518                    return Err(OptimError::Other(
519                        "L-BFGS-B requires scalar objective".into(),
520                    ))
521                }
522            };
523
524            let func_arc: Arc<dyn Fn(&[S]) -> S> = Arc::from(func);
525            let grad_fn: Box<dyn Fn(&[S], &mut [S])> = match grad {
526                Some(g) => g,
527                None => {
528                    let f = Arc::clone(&func_arc);
529                    Box::new(move |x: &[S], g: &mut [S]| {
530                        finite_diff_gradient(&*f, x, g);
531                    })
532                }
533            };
534
535            let lbfgsb_opts = LbfgsBOptions {
536                base: opts,
537                memory: 10,
538            };
539            lbfgsb_minimize(&*func_arc, &*grad_fn, &x0, &bounds, &lbfgsb_opts)
540        }
541
542        SolverChoice::AugmentedLagrangian
543        | SolverChoice::Simplex
544        | SolverChoice::ActiveSetQP
545        | SolverChoice::MixedIntegerLP
546        | SolverChoice::DifferentialEvolution
547        | SolverChoice::NsgaII
548        | SolverChoice::NelderMead
549        | SolverChoice::Powell
550        | SolverChoice::CmaEs
551        | SolverChoice::Sqp => unreachable!(),
552    }
553}
554
555#[cfg(test)]
556mod tests {
557    use super::*;
558    use crate::problem::OptimProblem;
559
560    #[test]
561    fn test_select_solver_unconstrained() {
562        let p = OptimProblem::new(2)
563            .x0(&[0.0, 0.0])
564            .objective(|x: &[f64]| x[0] * x[0] + x[1] * x[1]);
565        assert_eq!(select_solver(&p), SolverChoice::Lbfgs);
566    }
567
568    #[test]
569    fn test_select_solver_bounded() {
570        let p = OptimProblem::new(2)
571            .x0(&[0.0, 0.0])
572            .objective(|x: &[f64]| x[0] * x[0] + x[1] * x[1])
573            .bounds(0, (0.0, 1.0));
574        assert_eq!(select_solver(&p), SolverChoice::LbfgsB);
575    }
576
577    #[test]
578    fn test_select_solver_constrained() {
579        let p = OptimProblem::new(2)
580            .x0(&[0.0, 0.0])
581            .objective(|x: &[f64]| x[0] + x[1])
582            .constraint_eq(|x: &[f64]| x[0] * x[0] + x[1] * x[1] - 1.0);
583        assert_eq!(select_solver(&p), SolverChoice::AugmentedLagrangian);
584    }
585
586    #[test]
587    fn test_select_solver_least_squares() {
588        let p =
589            OptimProblem::new(2)
590                .x0(&[0.0, 0.0])
591                .least_squares(3, |_x: &[f64], r: &mut [f64]| {
592                    r[0] = 0.0;
593                    r[1] = 0.0;
594                    r[2] = 0.0;
595                });
596        assert_eq!(select_solver(&p), SolverChoice::LevenbergMarquardt);
597    }
598
599    #[test]
600    fn test_select_solver_lp() {
601        let p = OptimProblem::new(2)
602            .x0(&[0.0, 0.0])
603            .linear_objective(&[1.0, 2.0])
604            .linear_constraint_ineq(&[1.0, 1.0], 10.0);
605        assert_eq!(select_solver(&p), SolverChoice::Simplex);
606    }
607
608    #[test]
609    fn test_select_solver_qp() {
610        let p = OptimProblem::new(2)
611            .x0(&[0.0, 0.0])
612            .quadratic_objective_dense(&[2.0, 0.0, 0.0, 2.0], &[0.0, 0.0]);
613        assert_eq!(select_solver(&p), SolverChoice::ActiveSetQP);
614    }
615
616    #[test]
617    fn test_auto_unconstrained_solve() {
618        let result = OptimProblem::new(2)
619            .x0(&[5.0, 3.0])
620            .objective(|x: &[f64]| x[0] * x[0] + 4.0 * x[1] * x[1])
621            .gradient(|x: &[f64], g: &mut [f64]| {
622                g[0] = 2.0 * x[0];
623                g[1] = 8.0 * x[1];
624            })
625            .solve()
626            .unwrap();
627
628        assert!(result.converged, "did not converge: {}", result.message);
629        assert!(result.x[0].abs() < 1e-6);
630        assert!(result.x[1].abs() < 1e-6);
631    }
632
633    #[test]
634    fn test_auto_bounded_solve() {
635        let result = OptimProblem::new(2)
636            .x0(&[0.5, 0.5])
637            .objective(|x: &[f64]| (x[0] - 2.0).powi(2) + (x[1] - 2.0).powi(2))
638            .gradient(|x: &[f64], g: &mut [f64]| {
639                g[0] = 2.0 * (x[0] - 2.0);
640                g[1] = 2.0 * (x[1] - 2.0);
641            })
642            .all_bounds(&[(0.0, 1.0), (0.0, 3.0)])
643            .solve()
644            .unwrap();
645
646        assert!(result.converged, "did not converge: {}", result.message);
647        assert!((result.x[0] - 1.0).abs() < 1e-4, "x0={}", result.x[0]);
648        assert!((result.x[1] - 2.0).abs() < 1e-4, "x1={}", result.x[1]);
649    }
650
651    #[test]
652    fn test_auto_lp_solve() {
653        let result = OptimProblem::new(2)
654            .x0(&[0.0, 0.0])
655            .linear_objective(&[-1.0, -1.0])
656            .linear_constraint_ineq(&[1.0, 1.0], 4.0)
657            .linear_constraint_ineq(&[1.0, 0.0], 3.0)
658            .linear_constraint_ineq(&[0.0, 1.0], 3.0)
659            .solve()
660            .unwrap();
661        assert!(result.converged, "LP did not converge: {}", result.message);
662        assert!((result.f - (-4.0)).abs() < 1e-6, "f={}", result.f);
663    }
664
665    #[test]
666    fn test_solve_with_explicit_simplex() {
667        let result = OptimProblem::new(2)
668            .x0(&[0.0, 0.0])
669            .linear_objective(&[1.0, 2.0])
670            .linear_constraint_eq(&[1.0, 1.0], 3.0)
671            .solve_with(SolverChoice::Simplex)
672            .unwrap();
673        assert!(result.converged);
674        assert!((result.f - 3.0).abs() < 1e-6, "f={}", result.f);
675    }
676
677    #[test]
678    fn test_auto_qp_solve() {
679        let result = OptimProblem::new(2)
680            .x0(&[0.0, 0.0])
681            .quadratic_objective_dense(&[1.0, 0.0, 0.0, 1.0], &[-3.0, -3.0])
682            .all_bounds(&[(0.0, 1.0), (0.0, 1.0)])
683            .solve()
684            .unwrap();
685        assert!(result.converged, "QP did not converge: {}", result.message);
686        assert!((result.x[0] - 1.0).abs() < 1e-4, "x0={}", result.x[0]);
687        assert!((result.x[1] - 1.0).abs() < 1e-4, "x1={}", result.x[1]);
688    }
689
690    #[test]
691    fn test_auto_qp_equality() {
692        let result = OptimProblem::new(2)
693            .x0(&[0.0, 0.0])
694            .quadratic_objective_dense(&[1.0, 0.0, 0.0, 1.0], &[0.0, 0.0])
695            .linear_constraint_eq(&[1.0, 1.0], 1.0)
696            .solve()
697            .unwrap();
698        assert!(result.converged, "QP did not converge: {}", result.message);
699        assert!((result.x[0] - 0.5).abs() < 1e-4, "x0={}", result.x[0]);
700        assert!((result.x[1] - 0.5).abs() < 1e-4, "x1={}", result.x[1]);
701    }
702
703    #[test]
704    fn test_solve_with_explicit_active_set() {
705        let result = OptimProblem::new(2)
706            .x0(&[0.0, 0.0])
707            .quadratic_objective_dense(&[2.0, 0.0, 0.0, 2.0], &[-2.0, -4.0])
708            .solve_with(SolverChoice::ActiveSetQP)
709            .unwrap();
710        assert!(result.converged);
711        assert!((result.x[0] - 1.0).abs() < 1e-6, "x0={}", result.x[0]);
712        assert!((result.x[1] - 2.0).abs() < 1e-6, "x1={}", result.x[1]);
713    }
714
715    #[test]
716    fn test_select_solver_milp() {
717        let p = OptimProblem::new(2)
718            .linear_objective(&[1.0, 2.0])
719            .integer_var(0)
720            .linear_constraint_ineq(&[1.0, 1.0], 10.0);
721        assert_eq!(select_solver(&p), SolverChoice::MixedIntegerLP);
722    }
723
724    #[test]
725    fn test_auto_milp_solve() {
726        // minimize -3x1 - 5x2 s.t. x1+2x2<=6, 2x1+x2<=8, integer
727        let result = OptimProblem::new(2)
728            .linear_objective(&[-3.0, -5.0])
729            .linear_constraint_ineq(&[1.0, 2.0], 6.0)
730            .linear_constraint_ineq(&[2.0, 1.0], 8.0)
731            .integer_var(0)
732            .integer_var(1)
733            .solve()
734            .unwrap();
735        assert!(result.converged);
736        assert!((result.x[0] - result.x[0].round()).abs() < 1e-6);
737        assert!((result.x[1] - result.x[1].round()).abs() < 1e-6);
738    }
739
740    #[test]
741    fn test_milp_via_solve_with() {
742        let result = OptimProblem::new(2)
743            .linear_objective(&[-1.0, -1.0])
744            .linear_constraint_ineq(&[1.0, 1.0], 3.5)
745            .integer_var(0)
746            .solve_with(SolverChoice::MixedIntegerLP)
747            .unwrap();
748        assert!(result.converged);
749        assert!((result.x[0] - result.x[0].round()).abs() < 1e-6);
750    }
751
752    #[test]
753    fn test_select_solver_global() {
754        let p = OptimProblem::new(2)
755            .x0(&[0.0, 0.0])
756            .objective(|x: &[f64]| x[0] * x[0] + x[1] * x[1])
757            .all_bounds(&[(-5.0, 5.0), (-5.0, 5.0)])
758            .global(true);
759        assert_eq!(select_solver(&p), SolverChoice::DifferentialEvolution);
760    }
761
762    #[test]
763    fn test_auto_de_solve() {
764        let result = OptimProblem::new(2)
765            .objective(|x: &[f64]| x[0] * x[0] + x[1] * x[1])
766            .all_bounds(&[(-5.0, 5.0), (-5.0, 5.0)])
767            .global(true)
768            .solve()
769            .unwrap();
770        assert!(result.f < 0.1, "f={}", result.f);
771    }
772
773    #[test]
774    fn test_select_solver_multi_objective() {
775        let p = OptimProblem::new(2)
776            .multi_objective(vec![
777                Box::new(|x: &[f64]| x[0] * x[0]),
778                Box::new(|x: &[f64]| (x[0] - 2.0).powi(2)),
779            ])
780            .all_bounds(&[(0.0, 4.0), (0.0, 4.0)]);
781        assert_eq!(select_solver(&p), SolverChoice::NsgaII);
782    }
783
784    #[test]
785    fn test_auto_nsga2_solve() {
786        let result = OptimProblem::new(1)
787            .multi_objective(vec![
788                Box::new(|x: &[f64]| x[0] * x[0]),
789                Box::new(|x: &[f64]| (x[0] - 2.0).powi(2)),
790            ])
791            .all_bounds(&[(0.0, 4.0)])
792            .solve()
793            .unwrap();
794        assert!(result.pareto.is_some());
795        let pareto = result.pareto.unwrap();
796        assert!(!pareto.points.is_empty());
797    }
798
799    #[test]
800    fn test_maximize_quadratic() {
801        // maximize -(x0^2 + x1^2) => maximum at (0, 0) with value 0
802        let result = OptimProblem::new(2)
803            .x0(&[5.0, 3.0])
804            .maximize(|x: &[f64]| -(x[0] * x[0] + x[1] * x[1]))
805            .solve()
806            .unwrap();
807        assert!(result.converged, "did not converge: {}", result.message);
808        assert!(result.x[0].abs() < 1e-6);
809        assert!(result.x[1].abs() < 1e-6);
810        assert!(result.f.abs() < 1e-10, "f={}", result.f); // Maximum value is 0
811    }
812
813    #[test]
814    fn test_maximize_with_bounds() {
815        // maximize x0 + x1 subject to 0 <= xi <= 1
816        // Maximum at (1, 1) with value 2
817        let result = OptimProblem::new(2)
818            .x0(&[0.5, 0.5])
819            .maximize(|x: &[f64]| x[0] + x[1])
820            .maximize_gradient(|_x: &[f64], g: &mut [f64]| {
821                g[0] = 1.0;
822                g[1] = 1.0;
823            })
824            .all_bounds(&[(0.0, 1.0), (0.0, 1.0)])
825            .solve()
826            .unwrap();
827        assert!(result.converged, "did not converge: {}", result.message);
828        assert!((result.x[0] - 1.0).abs() < 1e-4, "x0={}", result.x[0]);
829        assert!((result.x[1] - 1.0).abs() < 1e-4, "x1={}", result.x[1]);
830        assert!((result.f - 2.0).abs() < 1e-4, "f={}", result.f);
831    }
832}