apex_solver/optimizer/levenberg_marquardt.rs
1//! Levenberg-Marquardt algorithm implementation.
2//!
3//! The Levenberg-Marquardt (LM) method is a robust and widely-used algorithm for solving
4//! nonlinear least squares problems of the form:
5//!
6//! ```text
7//! min f(x) = ½||r(x)||² = ½Σᵢ rᵢ(x)²
8//! ```
9//!
10//! where `r: ℝⁿ → ℝᵐ` is the residual vector function.
11//!
12//! # Algorithm Overview
13//!
14//! The Levenberg-Marquardt method solves the damped normal equations at each iteration:
15//!
16//! ```text
17//! (J^T·J + λI)·h = -J^T·r
18//! ```
19//!
20//! where:
21//! - `J` is the Jacobian matrix (m × n)
22//! - `r` is the residual vector (m × 1)
23//! - `h` is the step vector (n × 1)
24//! - `λ` is the adaptive damping parameter (scalar)
25//! - `I` is the identity matrix (or diagonal scaling matrix)
26//!
27//! ## Damping Parameter Strategy
28//!
29//! The damping parameter λ adapts based on step quality:
30//!
31//! - **λ → 0** (small damping): Behaves like Gauss-Newton with fast quadratic convergence
32//! - **λ → ∞** (large damping): Behaves like gradient descent with guaranteed descent direction
33//!
34//! This interpolation between Newton and gradient descent provides excellent robustness
35//! while maintaining fast convergence near the solution.
36//!
37//! ## Step Acceptance and Damping Update
38//!
39//! The algorithm evaluates each proposed step using the gain ratio:
40//!
41//! ```text
42//! ρ = (actual reduction) / (predicted reduction)
43//! = [f(xₖ) - f(xₖ + h)] / [f(xₖ) - L(h)]
44//! ```
45//!
46//! where `L(h) = f(xₖ) + h^T·g + ½h^T·H·h` is the local quadratic model.
47//!
48//! **Step acceptance:**
49//! - If `ρ > 0`: Accept step (cost decreased), decrease λ to trust the model more
50//! - If `ρ ≤ 0`: Reject step (cost increased), increase λ to be more conservative
51//!
52//! **Damping update** (Nielsen's formula):
53//! ```text
54//! λₖ₊₁ = λₖ · max(1/3, 1 - (2ρ - 1)³)
55//! ```
56//!
57//! This provides smooth, data-driven adaptation of the damping parameter.
58//!
59//! ## Convergence Properties
60//!
61//! - **Global convergence**: Guaranteed to find a stationary point from any starting guess
62//! - **Local quadratic convergence**: Near the solution, behaves like Gauss-Newton
63//! - **Robust to poor initialization**: Adaptive damping prevents divergence
64//! - **Handles ill-conditioning**: Large λ stabilizes nearly singular Hessian
65//!
66//! ## When to Use
67//!
68//! Levenberg-Marquardt is the best general-purpose choice when:
69//! - Initial parameter guess may be far from the optimum
70//! - Problem conditioning is unknown
71//! - Robustness is prioritized over raw speed
72//! - You want reliable convergence across diverse problem types
73//!
74//! For problems with specific structure, consider:
75//! - [`GaussNewton`](crate::optimizer::GaussNewton) if well-conditioned with good initialization
76//! - [`DogLeg`](crate::optimizer::DogLeg) for explicit trust region control
77//!
78//! # Implementation Features
79//!
80//! - **Sparse matrix support**: Efficient handling of large-scale problems via `faer` sparse library
81//! - **Adaptive damping**: Nielsen's formula for smooth parameter adaptation
82//! - **Robust linear solvers**: Cholesky (fast) or QR (stable) factorization
83//! - **Jacobi scaling**: Optional diagonal preconditioning for mixed-scale problems
84//! - **Covariance computation**: Optional uncertainty quantification after convergence
85//! - **Manifold operations**: Native support for optimization on Lie groups (SE2, SE3, SO2, SO3)
86//! - **Comprehensive diagnostics**: Detailed summaries of convergence and performance
87//!
88//! # Mathematical Background
89//!
90//! The augmented Hessian `J^T·J + λI` combines two beneficial properties:
91//!
92//! 1. **Positive definiteness**: Always solvable even when `J^T·J` is singular
93//! 2. **Regularization**: Prevents taking steps in poorly-determined directions
94//!
95//! The trust region interpretation: λ controls an implicit spherical trust region where
96//! larger λ restricts step size, ensuring the linear model remains valid.
97//!
98//! # Examples
99//!
100//! ## Basic usage
101//!
102//! ```no_run
103//! use apex_solver::optimizer::LevenbergMarquardt;
104//! use apex_solver::core::problem::Problem;
105//! use std::collections::HashMap;
106//!
107//! # fn main() -> Result<(), Box<dyn std::error::Error>> {
108//! let mut problem = Problem::new();
109//! // ... add residual blocks (factors) to problem ...
110//!
111//! let initial_values = HashMap::new();
112//! // ... initialize parameters ...
113//!
114//! let mut solver = LevenbergMarquardt::new();
115//! let result = solver.optimize(&problem, &initial_values)?;
116//!
117//! println!("Status: {:?}", result.status);
118//! println!("Final cost: {:.6e}", result.final_cost);
119//! println!("Iterations: {}", result.iterations);
120//! # Ok(())
121//! # }
122//! ```
123//!
124//! ## Advanced configuration
125//!
126//! ```no_run
127//! use apex_solver::optimizer::levenberg_marquardt::{LevenbergMarquardtConfig, LevenbergMarquardt};
128//! use apex_solver::linalg::LinearSolverType;
129//!
130//! # fn main() {
131//! let config = LevenbergMarquardtConfig::new()
132//! .with_max_iterations(100)
133//! .with_cost_tolerance(1e-6)
134//! .with_damping(1e-3) // Initial damping
135//! .with_damping_bounds(1e-12, 1e12) // Min/max damping
136//! .with_jacobi_scaling(true) // Improve conditioning
137//! .with_verbose(true);
138//!
139//! let mut solver = LevenbergMarquardt::with_config(config);
140//! # }
141//! ```
142//!
143//! # References
144//!
145//! - Levenberg, K. (1944). "A Method for the Solution of Certain Non-Linear Problems in Least Squares". *Quarterly of Applied Mathematics*.
146//! - Marquardt, D. W. (1963). "An Algorithm for Least-Squares Estimation of Nonlinear Parameters". *Journal of the Society for Industrial and Applied Mathematics*.
147//! - Madsen, K., Nielsen, H. B., & Tingleff, O. (2004). *Methods for Non-Linear Least Squares Problems* (2nd ed.). Chapter 3.
148//! - Nocedal, J. & Wright, S. (2006). *Numerical Optimization* (2nd ed.). Springer. Chapter 10.
149//! - Nielsen, H. B. (1999). "Damping Parameter in Marquardt's Method". Technical Report IMM-REP-1999-05.
150
151use crate::core::problem;
152use crate::error;
153use crate::linalg;
154use crate::manifold;
155use crate::optimizer;
156use faer::sparse;
157use nalgebra;
158use std::collections;
159use std::fmt;
160use std::time;
161
162/// Summary statistics for the Levenberg-Marquardt optimization process.
163#[derive(Debug, Clone)]
164pub struct LevenbergMarquardtSummary {
165 /// Initial cost value
166 pub initial_cost: f64,
167 /// Final cost value
168 pub final_cost: f64,
169 /// Total number of iterations performed
170 pub iterations: usize,
171 /// Number of successful steps (cost decreased)
172 pub successful_steps: usize,
173 /// Number of unsuccessful steps (cost increased, damping increased)
174 pub unsuccessful_steps: usize,
175 /// Final damping parameter value
176 pub final_damping: f64,
177 /// Average cost reduction per iteration
178 pub average_cost_reduction: f64,
179 /// Maximum gradient norm encountered
180 pub max_gradient_norm: f64,
181 /// Final gradient norm
182 pub final_gradient_norm: f64,
183 /// Maximum parameter update norm
184 pub max_parameter_update_norm: f64,
185 /// Final parameter update norm
186 pub final_parameter_update_norm: f64,
187 /// Total time elapsed
188 pub total_time: std::time::Duration,
189 /// Average time per iteration
190 pub average_time_per_iteration: std::time::Duration,
191}
192
193impl fmt::Display for LevenbergMarquardtSummary {
194 fn fmt(&self, f: &mut fmt::Formatter<'_>) -> fmt::Result {
195 writeln!(f, "Levenberg-Marquardt Optimization Summary")?;
196 writeln!(f, "Initial cost: {:.6e}", self.initial_cost)?;
197 writeln!(f, "Final cost: {:.6e}", self.final_cost)?;
198 writeln!(
199 f,
200 "Cost reduction: {:.6e} ({:.2}%)",
201 self.initial_cost - self.final_cost,
202 100.0 * (self.initial_cost - self.final_cost) / self.initial_cost.max(1e-12)
203 )?;
204 writeln!(f, "Total iterations: {}", self.iterations)?;
205 writeln!(
206 f,
207 "Successful steps: {} ({:.1}%)",
208 self.successful_steps,
209 100.0 * self.successful_steps as f64 / self.iterations.max(1) as f64
210 )?;
211 writeln!(
212 f,
213 "Unsuccessful steps: {} ({:.1}%)",
214 self.unsuccessful_steps,
215 100.0 * self.unsuccessful_steps as f64 / self.iterations.max(1) as f64
216 )?;
217 writeln!(f, "Final damping parameter: {:.6e}", self.final_damping)?;
218 writeln!(
219 f,
220 "Average cost reduction: {:.6e}",
221 self.average_cost_reduction
222 )?;
223 writeln!(
224 f,
225 "Max gradient norm: {:.6e}",
226 self.max_gradient_norm
227 )?;
228 writeln!(
229 f,
230 "Final gradient norm: {:.6e}",
231 self.final_gradient_norm
232 )?;
233 writeln!(
234 f,
235 "Max parameter update norm: {:.6e}",
236 self.max_parameter_update_norm
237 )?;
238 writeln!(
239 f,
240 "Final param update norm: {:.6e}",
241 self.final_parameter_update_norm
242 )?;
243 writeln!(f, "Total time: {:?}", self.total_time)?;
244 writeln!(
245 f,
246 "Average time per iteration: {:?}",
247 self.average_time_per_iteration
248 )?;
249 Ok(())
250 }
251}
252
253/// Configuration parameters for the Levenberg-Marquardt optimizer.
254///
255/// Controls the adaptive damping strategy, convergence criteria, and numerical stability
256/// enhancements for the Levenberg-Marquardt algorithm.
257///
258/// # Builder Pattern
259///
260/// All configuration options can be set using the builder pattern:
261///
262/// ```
263/// use apex_solver::optimizer::levenberg_marquardt::LevenbergMarquardtConfig;
264///
265/// let config = LevenbergMarquardtConfig::new()
266/// .with_max_iterations(100)
267/// .with_damping(1e-3)
268/// .with_damping_bounds(1e-12, 1e12)
269/// .with_jacobi_scaling(true)
270/// .with_verbose(true);
271/// ```
272///
273/// # Damping Parameter Behavior
274///
275/// The damping parameter λ controls the trade-off between Gauss-Newton and gradient descent:
276///
277/// - **Initial damping** (`damping`): Starting value (default: 1e-4)
278/// - **Damping bounds** (`damping_min`, `damping_max`): Valid range (default: 1e-12 to 1e12)
279/// - **Adaptation**: Automatically adjusted based on step quality using Nielsen's formula
280///
281/// # Convergence Criteria
282///
283/// The optimizer terminates when ANY of the following conditions is met:
284///
285/// - **Cost tolerance**: `|cost_k - cost_{k-1}| < cost_tolerance`
286/// - **Parameter tolerance**: `||step|| < parameter_tolerance`
287/// - **Gradient tolerance**: `||J^T·r|| < gradient_tolerance`
288/// - **Maximum iterations**: `iteration >= max_iterations`
289/// - **Timeout**: `elapsed_time >= timeout`
290///
291/// # See Also
292///
293/// - [`LevenbergMarquardt`] - The solver that uses this configuration
294/// - [`GaussNewtonConfig`](crate::optimizer::GaussNewtonConfig) - Undamped variant
295/// - [`DogLegConfig`](crate::optimizer::DogLegConfig) - Trust region alternative
296#[derive(Clone)]
297pub struct LevenbergMarquardtConfig {
298 /// Type of linear solver for the linear systems
299 pub linear_solver_type: linalg::LinearSolverType,
300 /// Maximum number of iterations
301 pub max_iterations: usize,
302 /// Convergence tolerance for cost function
303 pub cost_tolerance: f64,
304 /// Convergence tolerance for parameter updates
305 pub parameter_tolerance: f64,
306 /// Convergence tolerance for gradient norm
307 pub gradient_tolerance: f64,
308 /// Timeout duration
309 pub timeout: Option<time::Duration>,
310 /// Enable detailed logging
311 pub verbose: bool,
312 /// Initial damping parameter
313 pub damping: f64,
314 /// Minimum damping parameter
315 pub damping_min: f64,
316 /// Maximum damping parameter
317 pub damping_max: f64,
318 /// Damping increase factor (when step rejected)
319 pub damping_increase_factor: f64,
320 /// Damping decrease factor (when step accepted)
321 pub damping_decrease_factor: f64,
322 /// Damping nu parameter
323 pub damping_nu: f64,
324 /// Trust region radius
325 pub trust_region_radius: f64,
326 /// Minimum step quality for acceptance
327 pub min_step_quality: f64,
328 /// Good step quality threshold
329 pub good_step_quality: f64,
330 /// Minimum diagonal value for regularization
331 pub min_diagonal: f64,
332 /// Maximum diagonal value for regularization
333 pub max_diagonal: f64,
334 /// Minimum objective function cutoff (optional early termination)
335 ///
336 /// If set, optimization terminates when cost falls below this threshold.
337 /// Useful for early stopping when a "good enough" solution is acceptable.
338 ///
339 /// Default: None (disabled)
340 pub min_cost_threshold: Option<f64>,
341 /// Minimum trust region radius before termination
342 ///
343 /// When the trust region radius falls below this value, the optimizer
344 /// terminates as it indicates the search has converged or the problem
345 /// is ill-conditioned. Matches Ceres Solver's min_trust_region_radius.
346 ///
347 /// Default: 1e-32 (Ceres-compatible)
348 pub min_trust_region_radius: f64,
349 /// Maximum condition number for Jacobian matrix (optional check)
350 ///
351 /// If set, the optimizer checks if condition_number(J^T*J) exceeds this
352 /// threshold and terminates with IllConditionedJacobian status.
353 /// Note: Computing condition number is expensive, so this is disabled by default.
354 ///
355 /// Default: None (disabled)
356 pub max_condition_number: Option<f64>,
357 /// Minimum relative cost decrease for step acceptance
358 ///
359 /// Used in computing step quality (rho = actual_reduction / predicted_reduction).
360 /// Steps with rho < min_relative_decrease are rejected. Matches Ceres Solver's
361 /// min_relative_decrease parameter.
362 ///
363 /// Default: 1e-3 (Ceres-compatible)
364 pub min_relative_decrease: f64,
365 /// Use Jacobi column scaling (preconditioning)
366 ///
367 /// When enabled, normalizes Jacobian columns by their L2 norm before solving.
368 /// This can improve convergence for problems with mixed parameter scales
369 /// (e.g., positions in meters + angles in radians) but adds ~5-10% overhead.
370 ///
371 /// Default: false (to avoid performance overhead and faster convergence)
372 pub use_jacobi_scaling: bool,
373 /// Compute per-variable covariance matrices (uncertainty estimation)
374 ///
375 /// When enabled, computes covariance by inverting the Hessian matrix after
376 /// convergence. The full covariance matrix is extracted into per-variable
377 /// blocks stored in both Variable structs and SolverResult.
378 ///
379 /// Default: false (to avoid performance overhead)
380 pub compute_covariances: bool,
381 /// Enable real-time Rerun visualization (graphical debugging)
382 ///
383 /// When enabled, logs optimization progress to Rerun viewer including:
384 /// - Time series plots (cost, gradient norm, damping, step quality)
385 /// - Sparse Hessian heat map visualization
386 /// - Gradient vector visualization
387 /// - Manifold state updates (for SE2/SE3 problems)
388 ///
389 /// This is automatically disabled in release builds (zero overhead).
390 /// Use `verbose` for terminal output; this is for graphical visualization.
391 ///
392 /// Default: false
393 pub enable_visualization: bool,
394}
395
396impl Default for LevenbergMarquardtConfig {
397 fn default() -> Self {
398 Self {
399 linear_solver_type: linalg::LinearSolverType::default(),
400 // Ceres Solver default: 50 (changed from 100 for compatibility)
401 max_iterations: 50,
402 // Ceres Solver default: 1e-6 (changed from 1e-8 for compatibility)
403 cost_tolerance: 1e-6,
404 // Ceres Solver default: 1e-8 (unchanged)
405 parameter_tolerance: 1e-8,
406 // Ceres Solver default: 1e-10 (changed from 1e-8 for compatibility)
407 // Note: Typically should be 1e-4 * cost_tolerance per Ceres docs
408 gradient_tolerance: 1e-10,
409 timeout: None,
410 verbose: false,
411 damping: 1e-4,
412 damping_min: 1e-12,
413 damping_max: 1e12,
414 damping_increase_factor: 10.0,
415 damping_decrease_factor: 0.3,
416 damping_nu: 2.0,
417 trust_region_radius: 1e4,
418 min_step_quality: 0.0,
419 good_step_quality: 0.75,
420 min_diagonal: 1e-6,
421 max_diagonal: 1e32,
422 // New Ceres-compatible parameters
423 min_cost_threshold: None,
424 min_trust_region_radius: 1e-32,
425 max_condition_number: None,
426 min_relative_decrease: 1e-3,
427 // Existing parameters
428 use_jacobi_scaling: false,
429 compute_covariances: false,
430 enable_visualization: false,
431 }
432 }
433}
434
435impl LevenbergMarquardtConfig {
436 /// Create a new Levenberg-Marquardt configuration with default values.
437 pub fn new() -> Self {
438 Self::default()
439 }
440
441 /// Set the linear solver type
442 pub fn with_linear_solver_type(mut self, linear_solver_type: linalg::LinearSolverType) -> Self {
443 self.linear_solver_type = linear_solver_type;
444 self
445 }
446
447 /// Set the maximum number of iterations
448 pub fn with_max_iterations(mut self, max_iterations: usize) -> Self {
449 self.max_iterations = max_iterations;
450 self
451 }
452
453 /// Set the cost tolerance
454 pub fn with_cost_tolerance(mut self, cost_tolerance: f64) -> Self {
455 self.cost_tolerance = cost_tolerance;
456 self
457 }
458
459 /// Set the parameter tolerance
460 pub fn with_parameter_tolerance(mut self, parameter_tolerance: f64) -> Self {
461 self.parameter_tolerance = parameter_tolerance;
462 self
463 }
464
465 /// Set the gradient tolerance
466 pub fn with_gradient_tolerance(mut self, gradient_tolerance: f64) -> Self {
467 self.gradient_tolerance = gradient_tolerance;
468 self
469 }
470
471 /// Set the timeout duration
472 pub fn with_timeout(mut self, timeout: time::Duration) -> Self {
473 self.timeout = Some(timeout);
474 self
475 }
476
477 /// Enable or disable verbose logging
478 pub fn with_verbose(mut self, verbose: bool) -> Self {
479 self.verbose = verbose;
480 self
481 }
482
483 /// Set the initial damping parameter.
484 pub fn with_damping(mut self, damping: f64) -> Self {
485 self.damping = damping;
486 self
487 }
488
489 /// Set the damping parameter bounds.
490 pub fn with_damping_bounds(mut self, min: f64, max: f64) -> Self {
491 self.damping_min = min;
492 self.damping_max = max;
493 self
494 }
495
496 /// Set the damping adjustment factors.
497 pub fn with_damping_factors(mut self, increase: f64, decrease: f64) -> Self {
498 self.damping_increase_factor = increase;
499 self.damping_decrease_factor = decrease;
500 self
501 }
502
503 /// Set the trust region parameters.
504 pub fn with_trust_region(mut self, radius: f64, min_quality: f64, good_quality: f64) -> Self {
505 self.trust_region_radius = radius;
506 self.min_step_quality = min_quality;
507 self.good_step_quality = good_quality;
508 self
509 }
510
511 /// Set minimum objective function cutoff for early termination.
512 ///
513 /// When set, optimization terminates with MinCostThresholdReached status
514 /// if the cost falls below this threshold. Useful for early stopping when
515 /// a "good enough" solution is acceptable.
516 pub fn with_min_cost_threshold(mut self, min_cost: f64) -> Self {
517 self.min_cost_threshold = Some(min_cost);
518 self
519 }
520
521 /// Set minimum trust region radius before termination.
522 ///
523 /// When the trust region radius falls below this value, optimization
524 /// terminates with TrustRegionRadiusTooSmall status.
525 /// Default: 1e-32 (Ceres-compatible)
526 pub fn with_min_trust_region_radius(mut self, min_radius: f64) -> Self {
527 self.min_trust_region_radius = min_radius;
528 self
529 }
530
531 /// Set maximum condition number for Jacobian matrix.
532 ///
533 /// If set, the optimizer checks if condition_number(J^T*J) exceeds this
534 /// threshold and terminates with IllConditionedJacobian status.
535 /// Note: Computing condition number is expensive, disabled by default.
536 pub fn with_max_condition_number(mut self, max_cond: f64) -> Self {
537 self.max_condition_number = Some(max_cond);
538 self
539 }
540
541 /// Set minimum relative cost decrease for step acceptance.
542 ///
543 /// Steps with rho = (actual_reduction / predicted_reduction) below this
544 /// threshold are rejected. Default: 1e-3 (Ceres-compatible)
545 pub fn with_min_relative_decrease(mut self, min_decrease: f64) -> Self {
546 self.min_relative_decrease = min_decrease;
547 self
548 }
549
550 /// Enable or disable Jacobi column scaling (preconditioning).
551 ///
552 /// When enabled, normalizes Jacobian columns by their L2 norm before solving.
553 /// Can improve convergence for mixed-scale problems but adds ~5-10% overhead.
554 pub fn with_jacobi_scaling(mut self, use_jacobi_scaling: bool) -> Self {
555 self.use_jacobi_scaling = use_jacobi_scaling;
556 self
557 }
558
559 /// Enable or disable covariance computation (uncertainty estimation).
560 ///
561 /// When enabled, computes the full covariance matrix by inverting the Hessian
562 /// after convergence, then extracts per-variable covariance blocks.
563 pub fn with_compute_covariances(mut self, compute_covariances: bool) -> Self {
564 self.compute_covariances = compute_covariances;
565 self
566 }
567
568 /// Enable real-time visualization (graphical debugging).
569 ///
570 /// When enabled, optimization progress is logged to a Rerun viewer with:
571 /// - Time series plots of cost, gradient norm, damping, step quality
572 /// - Sparse Hessian matrix visualization as heat map
573 /// - Gradient vector visualization
574 /// - Real-time manifold state updates (for SE2/SE3 problems)
575 ///
576 /// Note: Has zero overhead when disabled. Use `verbose` for terminal logging.
577 ///
578 /// # Arguments
579 ///
580 /// * `enable` - Whether to enable visualization
581 pub fn with_visualization(mut self, enable: bool) -> Self {
582 self.enable_visualization = enable;
583 self
584 }
585}
586
587/// State for optimization iteration
588struct LinearizerResult {
589 variables: collections::HashMap<String, problem::VariableEnum>,
590 variable_index_map: collections::HashMap<String, usize>,
591 sorted_vars: Vec<String>,
592 symbolic_structure: problem::SymbolicStructure,
593 current_cost: f64,
594 initial_cost: f64,
595}
596
597/// Result from step computation
598struct StepResult {
599 step: faer::Mat<f64>,
600 gradient_norm: f64,
601 predicted_reduction: f64,
602}
603
604/// Result from step evaluation
605struct StepEvaluation {
606 accepted: bool,
607 cost_reduction: f64,
608 rho: f64,
609}
610
611/// Levenberg-Marquardt solver for nonlinear least squares optimization.
612///
613/// Implements the damped Gauss-Newton method with adaptive damping parameter λ that
614/// interpolates between Gauss-Newton and gradient descent based on step quality.
615///
616/// # Algorithm
617///
618/// At each iteration k:
619/// 1. Compute residual `r(xₖ)` and Jacobian `J(xₖ)`
620/// 2. Solve augmented system: `(J^T·J + λI)·h = -J^T·r`
621/// 3. Evaluate step quality: `ρ = (actual reduction) / (predicted reduction)`
622/// 4. If `ρ > 0`: Accept step and update `xₖ₊₁ = xₖ ⊕ h`, decrease λ
623/// 5. If `ρ ≤ 0`: Reject step (keep `xₖ₊₁ = xₖ`), increase λ
624/// 6. Check convergence criteria
625///
626/// The damping parameter λ is updated using Nielsen's smooth formula:
627/// `λₖ₊₁ = λₖ · max(1/3, 1 - (2ρ - 1)³)` for accepted steps,
628/// or `λₖ₊₁ = λₖ · ν` (with increasing ν) for rejected steps.
629///
630/// # Examples
631///
632/// ```no_run
633/// use apex_solver::optimizer::LevenbergMarquardt;
634/// use apex_solver::core::problem::Problem;
635/// use std::collections::HashMap;
636///
637/// # fn main() -> Result<(), Box<dyn std::error::Error>> {
638/// let mut problem = Problem::new();
639/// // ... add factors to problem ...
640///
641/// let initial_values = HashMap::new();
642/// // ... initialize parameters ...
643///
644/// let mut solver = LevenbergMarquardt::new();
645/// let result = solver.optimize(&problem, &initial_values)?;
646///
647/// println!("Status: {:?}", result.status);
648/// println!("Final cost: {}", result.final_cost);
649/// println!("Iterations: {}", result.iterations);
650/// # Ok(())
651/// # }
652/// ```
653///
654/// # See Also
655///
656/// - [`LevenbergMarquardtConfig`] - Configuration options
657/// - [`GaussNewton`](crate::optimizer::GaussNewton) - Undamped variant (faster but less robust)
658/// - [`DogLeg`](crate::optimizer::DogLeg) - Alternative trust region method
659pub struct LevenbergMarquardt {
660 config: LevenbergMarquardtConfig,
661 jacobi_scaling: Option<sparse::SparseColMat<usize, f64>>,
662 visualizer: Option<optimizer::OptimizationVisualizer>,
663}
664
665impl Default for LevenbergMarquardt {
666 fn default() -> Self {
667 Self::new()
668 }
669}
670
671impl LevenbergMarquardt {
672 /// Create a new Levenberg-Marquardt solver with default configuration.
673 pub fn new() -> Self {
674 Self::with_config(LevenbergMarquardtConfig::default())
675 }
676
677 /// Create a new Levenberg-Marquardt solver with the given configuration.
678 pub fn with_config(config: LevenbergMarquardtConfig) -> Self {
679 // Create visualizer if enabled (zero overhead when disabled)
680 let visualizer = if config.enable_visualization {
681 match optimizer::OptimizationVisualizer::new(true) {
682 Ok(vis) => Some(vis),
683 Err(e) => {
684 eprintln!("[WARNING] Failed to create visualizer: {}", e);
685 None
686 }
687 }
688 } else {
689 None
690 };
691
692 Self {
693 config,
694 jacobi_scaling: None,
695 visualizer,
696 }
697 }
698
699 /// Create the appropriate linear solver based on configuration
700 fn create_linear_solver(&self) -> Box<dyn linalg::SparseLinearSolver> {
701 match self.config.linear_solver_type {
702 linalg::LinearSolverType::SparseCholesky => {
703 Box::new(linalg::SparseCholeskySolver::new())
704 }
705 linalg::LinearSolverType::SparseQR => Box::new(linalg::SparseQRSolver::new()),
706 }
707 }
708
709 /// Update damping parameter based on step quality using trust region approach
710 /// Reference: Introduction to Optimization and Data Fitting
711 /// Algorithm 6.18
712 fn update_damping(&mut self, rho: f64) -> bool {
713 if rho > 0.0 {
714 // Step accepted - decrease damping
715 let coff = 2.0 * rho - 1.0;
716 self.config.damping *= (1.0_f64 / 3.0).max(1.0 - coff * coff * coff);
717 self.config.damping = self.config.damping.max(self.config.damping_min);
718 self.config.damping_nu = 2.0;
719 true
720 } else {
721 // Step rejected - increase damping
722 self.config.damping *= self.config.damping_nu;
723 self.config.damping_nu *= 2.0;
724 self.config.damping = self.config.damping.min(self.config.damping_max);
725 false
726 }
727 }
728
729 /// Compute step quality ratio (actual vs predicted reduction)
730 /// Reference: Introduction to Optimization and Data Fitting
731 /// Reference: Damping parameter in marquardt's method
732 /// Formula 2.2
733 fn compute_step_quality(
734 &self,
735 current_cost: f64,
736 new_cost: f64,
737 predicted_reduction: f64,
738 ) -> f64 {
739 let actual_reduction = current_cost - new_cost;
740 if predicted_reduction.abs() < 1e-15 {
741 if actual_reduction > 0.0 { 1.0 } else { 0.0 }
742 } else {
743 actual_reduction / predicted_reduction
744 }
745 }
746
747 /// Compute predicted cost reduction from linear model
748 /// Standard LM formula: 0.5 * step^T * (damping * step - gradient)
749 fn compute_predicted_reduction(&self, step: &faer::Mat<f64>, gradient: &faer::Mat<f64>) -> f64 {
750 // Standard Levenberg-Marquardt predicted reduction formula
751 // predicted_reduction = -step^T * gradient - 0.5 * step^T * H * step
752 // = 0.5 * step^T * (damping * step - gradient)
753 let diff = self.config.damping * step - gradient;
754 (0.5 * step.transpose() * &diff)[(0, 0)]
755 }
756
757 /// Check convergence criteria
758 /// Check convergence using comprehensive termination criteria.
759 ///
760 /// Implements 9 termination criteria following Ceres Solver standards:
761 ///
762 /// 1. **Gradient Norm (First-Order Optimality)**: ||g||∞ ≤ gradient_tolerance
763 /// 2. **Parameter Change Tolerance**: ||h|| ≤ parameter_tolerance * (||x|| + parameter_tolerance)
764 /// 3. **Function Value Change Tolerance**: |ΔF| < cost_tolerance * F
765 /// 4. **Objective Function Cutoff**: F_new < min_cost_threshold (optional)
766 /// 5. **Trust Region Radius**: radius < min_trust_region_radius
767 /// 6. **Singular/Ill-Conditioned Jacobian**: Detected during linear solve
768 /// 7. **Invalid Numerical Values**: NaN or Inf in cost or parameters
769 /// 8. **Maximum Iterations**: iteration >= max_iterations
770 /// 9. **Timeout**: elapsed >= timeout
771 ///
772 /// # Arguments
773 ///
774 /// * `iteration` - Current iteration number
775 /// * `current_cost` - Cost before applying the step
776 /// * `new_cost` - Cost after applying the step
777 /// * `parameter_norm` - L2 norm of current parameter vector ||x||
778 /// * `parameter_update_norm` - L2 norm of parameter update step ||h||
779 /// * `gradient_norm` - Infinity norm of gradient ||g||∞
780 /// * `trust_region_radius` - Current trust region radius
781 /// * `elapsed` - Elapsed time since optimization start
782 /// * `step_accepted` - Whether the current step was accepted
783 ///
784 /// # Returns
785 ///
786 /// `Some(OptimizationStatus)` if any termination criterion is satisfied, `None` otherwise.
787 #[allow(clippy::too_many_arguments)]
788 fn check_convergence(
789 &self,
790 iteration: usize,
791 current_cost: f64,
792 new_cost: f64,
793 parameter_norm: f64,
794 parameter_update_norm: f64,
795 gradient_norm: f64,
796 trust_region_radius: f64,
797 elapsed: std::time::Duration,
798 step_accepted: bool,
799 ) -> Option<optimizer::OptimizationStatus> {
800 // ========================================================================
801 // CRITICAL SAFETY CHECKS (perform first, before convergence checks)
802 // ========================================================================
803
804 // CRITERION 7: Invalid Numerical Values (NaN/Inf)
805 // Always check for numerical instability first
806 if !new_cost.is_finite() || !parameter_update_norm.is_finite() || !gradient_norm.is_finite()
807 {
808 return Some(optimizer::OptimizationStatus::InvalidNumericalValues);
809 }
810
811 // CRITERION 9: Timeout
812 // Check wall-clock time limit
813 if let Some(timeout) = self.config.timeout
814 && elapsed >= timeout
815 {
816 return Some(optimizer::OptimizationStatus::Timeout);
817 }
818
819 // CRITERION 8: Maximum Iterations
820 // Check iteration count limit
821 if iteration >= self.config.max_iterations {
822 return Some(optimizer::OptimizationStatus::MaxIterationsReached);
823 }
824
825 // ========================================================================
826 // CONVERGENCE CRITERIA (only check after successful steps)
827 // ========================================================================
828
829 // Only check convergence criteria after accepted steps
830 // (rejected steps don't indicate convergence)
831 if !step_accepted {
832 return None;
833 }
834
835 // CRITERION 1: Gradient Norm (First-Order Optimality)
836 // Check if gradient infinity norm is below threshold
837 // This indicates we're at a critical point (local minimum, saddle, or maximum)
838 if gradient_norm < self.config.gradient_tolerance {
839 return Some(optimizer::OptimizationStatus::GradientToleranceReached);
840 }
841
842 // Only check parameter and cost criteria after first iteration
843 if iteration > 0 {
844 // CRITERION 2: Parameter Change Tolerance (xtol)
845 // Ceres formula: ||h|| ≤ ε_param * (||x|| + ε_param)
846 // This is a relative measure that scales with parameter magnitude
847 let relative_step_tolerance = self.config.parameter_tolerance
848 * (parameter_norm + self.config.parameter_tolerance);
849
850 if parameter_update_norm <= relative_step_tolerance {
851 return Some(optimizer::OptimizationStatus::ParameterToleranceReached);
852 }
853
854 // CRITERION 3: Function Value Change Tolerance (ftol)
855 // Ceres formula: |ΔF| < ε_cost * F
856 // Check relative cost change (not absolute)
857 let cost_change = (current_cost - new_cost).abs();
858 let relative_cost_change = cost_change / current_cost.max(1e-10); // Avoid division by zero
859
860 if relative_cost_change < self.config.cost_tolerance {
861 return Some(optimizer::OptimizationStatus::CostToleranceReached);
862 }
863 }
864
865 // CRITERION 4: Objective Function Cutoff (optional early stopping)
866 // Useful for "good enough" solutions
867 if let Some(min_cost) = self.config.min_cost_threshold
868 && new_cost < min_cost
869 {
870 return Some(optimizer::OptimizationStatus::MinCostThresholdReached);
871 }
872
873 // CRITERION 5: Trust Region Radius
874 // If trust region has collapsed, optimization has converged or problem is ill-conditioned
875 if trust_region_radius < self.config.min_trust_region_radius {
876 return Some(optimizer::OptimizationStatus::TrustRegionRadiusTooSmall);
877 }
878
879 // CRITERION 6: Singular/Ill-Conditioned Jacobian
880 // Note: This is typically detected during the linear solve and handled there
881 // The max_condition_number check would be expensive to compute here
882 // If linear solve fails, it returns an error that's converted to NumericalFailure
883
884 // No termination criterion satisfied
885 None
886 }
887
888 /// Compute total parameter vector norm ||x||.
889 ///
890 /// Computes the L2 norm of all parameter vectors concatenated together.
891 /// This is used in the relative parameter tolerance check.
892 ///
893 /// # Arguments
894 ///
895 /// * `variables` - Map of variable names to their current values
896 ///
897 /// # Returns
898 ///
899 /// The L2 norm of the concatenated parameter vector
900 fn compute_parameter_norm(
901 variables: &collections::HashMap<String, problem::VariableEnum>,
902 ) -> f64 {
903 variables
904 .values()
905 .map(|v| {
906 let vec = v.to_vector();
907 vec.norm_squared()
908 })
909 .sum::<f64>()
910 .sqrt()
911 }
912
913 /// Create Jacobi scaling matrix from Jacobian
914 fn create_jacobi_scaling(
915 &self,
916 jacobian: &sparse::SparseColMat<usize, f64>,
917 ) -> sparse::SparseColMat<usize, f64> {
918 use faer::sparse::Triplet;
919
920 let cols = jacobian.ncols();
921 let jacobi_scaling_vec: Vec<Triplet<usize, usize, f64>> = (0..cols)
922 .map(|c| {
923 // Compute column norm: sqrt(sum(J_col^2))
924 let col_norm_squared: f64 = jacobian
925 .triplet_iter()
926 .filter(|t| t.col == c)
927 .map(|t| t.val * t.val)
928 .sum();
929 let col_norm = col_norm_squared.sqrt();
930 // Scaling factor: 1.0 / (1.0 + col_norm)
931 let scaling = 1.0 / (1.0 + col_norm);
932 Triplet::new(c, c, scaling)
933 })
934 .collect();
935
936 sparse::SparseColMat::try_new_from_triplets(cols, cols, &jacobi_scaling_vec)
937 .expect("Failed to create Jacobi scaling matrix")
938 }
939
940 /// Initialize optimization state from problem and initial parameters
941 fn initialize_optimization_state(
942 &self,
943 problem: &problem::Problem,
944 initial_params: &collections::HashMap<
945 String,
946 (manifold::ManifoldType, nalgebra::DVector<f64>),
947 >,
948 ) -> Result<LinearizerResult, error::ApexError> {
949 // Initialize variables from initial values
950 let variables = problem.initialize_variables(initial_params);
951
952 // Create column mapping for variables
953 let mut variable_index_map = collections::HashMap::new();
954 let mut col_offset = 0;
955 let mut sorted_vars: Vec<String> = variables.keys().cloned().collect();
956 sorted_vars.sort();
957
958 for var_name in &sorted_vars {
959 variable_index_map.insert(var_name.clone(), col_offset);
960 col_offset += variables[var_name].get_size();
961 }
962
963 // Build symbolic structure for sparse operations
964 let symbolic_structure =
965 problem.build_symbolic_structure(&variables, &variable_index_map, col_offset)?;
966
967 // Initial cost evaluation (residual only, no Jacobian needed)
968 let residual = problem.compute_residual_sparse(&variables)?;
969 let current_cost = optimizer::compute_cost(&residual);
970 let initial_cost = current_cost;
971
972 if self.config.verbose {
973 println!(
974 "Starting Levenberg-Marquardt optimization with {} max iterations",
975 self.config.max_iterations
976 );
977 println!(
978 "Initial cost: {:.6e}, initial damping: {:.6e}",
979 current_cost, self.config.damping
980 );
981 }
982
983 Ok(LinearizerResult {
984 variables,
985 variable_index_map,
986 sorted_vars,
987 symbolic_structure,
988 current_cost,
989 initial_cost,
990 })
991 }
992
993 /// Process Jacobian by creating and applying Jacobi scaling if enabled
994 fn process_jacobian(
995 &mut self,
996 jacobian: &sparse::SparseColMat<usize, f64>,
997 iteration: usize,
998 ) -> sparse::SparseColMat<usize, f64> {
999 // Create Jacobi scaling on first iteration if enabled
1000 if iteration == 0 {
1001 let scaling = self.create_jacobi_scaling(jacobian);
1002
1003 if self.config.verbose {
1004 println!("Jacobi scaling computed for {} columns", scaling.ncols());
1005 }
1006
1007 self.jacobi_scaling = Some(scaling);
1008 }
1009 jacobian * self.jacobi_scaling.as_ref().unwrap()
1010 }
1011
1012 /// Compute optimization step by solving the augmented system
1013 fn compute_levenberg_marquardt_step(
1014 &self,
1015 residuals: &faer::Mat<f64>,
1016 scaled_jacobian: &sparse::SparseColMat<usize, f64>,
1017 linear_solver: &mut Box<dyn linalg::SparseLinearSolver>,
1018 ) -> Option<StepResult> {
1019 // Solve augmented equation: (J_scaled^T * J_scaled + λI) * dx_scaled = -J_scaled^T * r
1020 let residuals_owned = residuals.as_ref().to_owned();
1021 let scaled_step = linear_solver
1022 .solve_augmented_equation(&residuals_owned, scaled_jacobian, self.config.damping)
1023 .ok()?;
1024
1025 // Get cached gradient and Hessian from the solver
1026 let gradient = linear_solver.get_gradient()?;
1027 let hessian = linear_solver.get_hessian()?;
1028 let gradient_norm = gradient.norm_l2();
1029
1030 if self.config.verbose {
1031 println!("Gradient (J^T*r) norm: {:.12e}", gradient_norm);
1032 println!("Hessian shape: ({}, {})", hessian.nrows(), hessian.ncols());
1033 println!("Damping parameter: {:.12e}", self.config.damping);
1034 }
1035
1036 // Apply inverse Jacobi scaling to get final step (if enabled)
1037 let step = if self.config.use_jacobi_scaling {
1038 &scaled_step * self.jacobi_scaling.as_ref().unwrap()
1039 } else {
1040 scaled_step
1041 };
1042
1043 if self.config.verbose {
1044 println!("Final step norm: {:.12e}", step.norm_l2());
1045 }
1046
1047 // Compute predicted reduction using scaled values
1048 let predicted_reduction = self.compute_predicted_reduction(&step, gradient);
1049
1050 if self.config.verbose {
1051 println!("Predicted reduction: {:.12e}", predicted_reduction);
1052 }
1053
1054 Some(StepResult {
1055 step,
1056 gradient_norm,
1057 predicted_reduction,
1058 })
1059 }
1060
1061 /// Evaluate and apply step, handling acceptance/rejection based on step quality
1062 fn evaluate_and_apply_step(
1063 &mut self,
1064 step_result: &StepResult,
1065 state: &mut LinearizerResult,
1066 problem: &problem::Problem,
1067 ) -> error::ApexResult<StepEvaluation> {
1068 // Apply parameter updates using manifold operations
1069 let _step_norm = optimizer::apply_parameter_step(
1070 &mut state.variables,
1071 step_result.step.as_ref(),
1072 &state.sorted_vars,
1073 );
1074
1075 // Compute new cost (residual only, no Jacobian needed for step evaluation)
1076 let new_residual = problem.compute_residual_sparse(&state.variables)?;
1077 let new_cost = optimizer::compute_cost(&new_residual);
1078
1079 // Compute step quality
1080 let rho = self.compute_step_quality(
1081 state.current_cost,
1082 new_cost,
1083 step_result.predicted_reduction,
1084 );
1085
1086 if self.config.verbose {
1087 println!("RHO (Gain Factor) CALCULATION DETAILS");
1088 println!("Old cost: {:.12e}", state.current_cost);
1089 println!("New cost: {:.12e}", new_cost);
1090 let actual_reduction = state.current_cost - new_cost;
1091 println!("Actual cost reduction: {:.12e}", actual_reduction);
1092 println!(
1093 "Predicted cost reduction: {:.12e}",
1094 step_result.predicted_reduction
1095 );
1096 println!("Rho (actual/predicted): {:.12e}", rho);
1097 }
1098
1099 // Update damping and decide whether to accept step
1100 let accepted = self.update_damping(rho);
1101
1102 let cost_reduction = if accepted {
1103 // Accept the step - parameters already updated
1104 let reduction = state.current_cost - new_cost;
1105 state.current_cost = new_cost;
1106 reduction
1107 } else {
1108 // Reject the step - revert parameter changes
1109 optimizer::apply_negative_parameter_step(
1110 &mut state.variables,
1111 step_result.step.as_ref(),
1112 &state.sorted_vars,
1113 );
1114 0.0
1115 };
1116
1117 Ok(StepEvaluation {
1118 accepted,
1119 cost_reduction,
1120 rho,
1121 })
1122 }
1123
1124 /// Create optimization summary
1125 #[allow(clippy::too_many_arguments)]
1126 fn create_summary(
1127 &self,
1128 initial_cost: f64,
1129 final_cost: f64,
1130 iterations: usize,
1131 successful_steps: usize,
1132 unsuccessful_steps: usize,
1133 max_gradient_norm: f64,
1134 final_gradient_norm: f64,
1135 max_parameter_update_norm: f64,
1136 final_parameter_update_norm: f64,
1137 total_cost_reduction: f64,
1138 total_time: std::time::Duration,
1139 ) -> LevenbergMarquardtSummary {
1140 LevenbergMarquardtSummary {
1141 initial_cost,
1142 final_cost,
1143 iterations,
1144 successful_steps,
1145 unsuccessful_steps,
1146 final_damping: self.config.damping,
1147 average_cost_reduction: if iterations > 0 {
1148 total_cost_reduction / iterations as f64
1149 } else {
1150 0.0
1151 },
1152 max_gradient_norm,
1153 final_gradient_norm,
1154 max_parameter_update_norm,
1155 final_parameter_update_norm,
1156 total_time,
1157 average_time_per_iteration: if iterations > 0 {
1158 total_time / iterations as u32
1159 } else {
1160 std::time::Duration::from_secs(0)
1161 },
1162 }
1163 }
1164
1165 pub fn optimize(
1166 &mut self,
1167 problem: &problem::Problem,
1168 initial_params: &collections::HashMap<
1169 String,
1170 (manifold::ManifoldType, nalgebra::DVector<f64>),
1171 >,
1172 ) -> Result<
1173 optimizer::SolverResult<collections::HashMap<String, problem::VariableEnum>>,
1174 error::ApexError,
1175 > {
1176 let start_time = time::Instant::now();
1177 let mut iteration = 0;
1178 let mut cost_evaluations = 1; // Initial cost evaluation
1179 let mut jacobian_evaluations = 0;
1180 let mut successful_steps = 0;
1181 let mut unsuccessful_steps = 0;
1182
1183 // Initialize optimization state
1184 let mut state = self.initialize_optimization_state(problem, initial_params)?;
1185
1186 // Create linear solver
1187 let mut linear_solver = self.create_linear_solver();
1188
1189 // Initialize summary tracking variables
1190 let mut max_gradient_norm: f64 = 0.0;
1191 let mut max_parameter_update_norm: f64 = 0.0;
1192 let mut total_cost_reduction = 0.0;
1193 let mut final_gradient_norm;
1194 let mut final_parameter_update_norm;
1195
1196 // Main optimization loop
1197 loop {
1198 // Evaluate residuals and Jacobian
1199 let (residuals, jacobian) = problem.compute_residual_and_jacobian_sparse(
1200 &state.variables,
1201 &state.variable_index_map,
1202 &state.symbolic_structure,
1203 )?;
1204 jacobian_evaluations += 1;
1205
1206 if self.config.verbose {
1207 println!("APEX-SOLVER DEBUG ITERATION {}", iteration);
1208 println!("Current cost: {:.12e}", state.current_cost);
1209 println!(
1210 "Residuals shape: ({}, {})",
1211 residuals.nrows(),
1212 residuals.ncols()
1213 );
1214 println!("Residuals norm: {:.12e}", residuals.norm_l2());
1215 println!(
1216 "Jacobian shape: ({}, {})",
1217 jacobian.nrows(),
1218 jacobian.ncols()
1219 );
1220 }
1221
1222 // Process Jacobian (apply scaling if enabled)
1223 let scaled_jacobian = if self.config.use_jacobi_scaling {
1224 self.process_jacobian(&jacobian, iteration)
1225 } else {
1226 jacobian
1227 };
1228
1229 if self.config.verbose {
1230 println!(
1231 "Scaled Jacobian shape: ({}, {})",
1232 scaled_jacobian.nrows(),
1233 scaled_jacobian.ncols()
1234 );
1235 }
1236
1237 // Compute optimization step
1238 let step_result = match self.compute_levenberg_marquardt_step(
1239 &residuals,
1240 &scaled_jacobian,
1241 &mut linear_solver,
1242 ) {
1243 Some(result) => result,
1244 None => {
1245 return Err(error::ApexError::Solver(
1246 "Linear solver failed to solve augmented system".to_string(),
1247 ));
1248 }
1249 };
1250
1251 // Update tracking variables
1252 max_gradient_norm = max_gradient_norm.max(step_result.gradient_norm);
1253 final_gradient_norm = step_result.gradient_norm;
1254 let step_norm = step_result.step.norm_l2();
1255 max_parameter_update_norm = max_parameter_update_norm.max(step_norm);
1256 final_parameter_update_norm = step_norm;
1257
1258 // Evaluate and apply step (handles accept/reject)
1259 let step_eval = self.evaluate_and_apply_step(&step_result, &mut state, problem)?;
1260 cost_evaluations += 1;
1261
1262 // Update counters based on acceptance
1263 if step_eval.accepted {
1264 successful_steps += 1;
1265 total_cost_reduction += step_eval.cost_reduction;
1266 } else {
1267 unsuccessful_steps += 1;
1268 }
1269
1270 // Rerun visualization
1271 if let Some(ref vis) = self.visualizer {
1272 if let Err(e) = vis.log_scalars(
1273 iteration,
1274 state.current_cost,
1275 step_result.gradient_norm,
1276 self.config.damping,
1277 step_norm,
1278 Some(step_eval.rho),
1279 ) {
1280 eprintln!("[WARNING] Failed to log scalars: {}", e);
1281 }
1282
1283 // Log expensive visualizations (Hessian, gradient, manifolds)
1284 if let Err(e) = vis.log_hessian(linear_solver.get_hessian(), iteration) {
1285 eprintln!("[WARNING] Failed to log Hessian: {}", e);
1286 }
1287 if let Err(e) = vis.log_gradient(linear_solver.get_gradient(), iteration) {
1288 eprintln!("[WARNING] Failed to log gradient: {}", e);
1289 }
1290 if let Err(e) = vis.log_manifolds(&state.variables, iteration) {
1291 eprintln!("[WARNING] Failed to log manifolds: {}", e);
1292 }
1293 }
1294
1295 // Check convergence
1296 let elapsed = start_time.elapsed();
1297
1298 // Compute parameter norm for relative parameter tolerance check
1299 let parameter_norm = Self::compute_parameter_norm(&state.variables);
1300
1301 // Compute new cost for convergence check (state may already have new cost if step accepted)
1302 let new_cost = if step_eval.accepted {
1303 state.current_cost
1304 } else {
1305 // Use cost before step application
1306 state.current_cost
1307 };
1308
1309 // Cost before this step (need to add back reduction if step was accepted)
1310 let cost_before_step = if step_eval.accepted {
1311 state.current_cost + step_eval.cost_reduction
1312 } else {
1313 state.current_cost
1314 };
1315
1316 if let Some(status) = self.check_convergence(
1317 iteration,
1318 cost_before_step,
1319 new_cost,
1320 parameter_norm,
1321 step_norm,
1322 step_result.gradient_norm,
1323 self.config.trust_region_radius,
1324 elapsed,
1325 step_eval.accepted,
1326 ) {
1327 let summary = self.create_summary(
1328 state.initial_cost,
1329 state.current_cost,
1330 iteration + 1,
1331 successful_steps,
1332 unsuccessful_steps,
1333 max_gradient_norm,
1334 final_gradient_norm,
1335 max_parameter_update_norm,
1336 final_parameter_update_norm,
1337 total_cost_reduction,
1338 elapsed,
1339 );
1340
1341 if self.config.verbose {
1342 println!("{}", summary);
1343 }
1344
1345 // Log convergence to Rerun
1346 if let Some(ref vis) = self.visualizer {
1347 let _ = vis.log_convergence(&format!("Converged: {}", status));
1348 }
1349
1350 // Compute covariances if enabled
1351 let covariances = if self.config.compute_covariances {
1352 problem.compute_and_set_covariances(
1353 &mut linear_solver,
1354 &mut state.variables,
1355 &state.variable_index_map,
1356 )
1357 } else {
1358 None
1359 };
1360
1361 return Ok(optimizer::SolverResult {
1362 status,
1363 iterations: iteration + 1,
1364 initial_cost: state.initial_cost,
1365 final_cost: state.current_cost,
1366 parameters: state.variables.into_iter().collect(),
1367 elapsed_time: elapsed,
1368 convergence_info: Some(optimizer::ConvergenceInfo {
1369 final_gradient_norm,
1370 final_parameter_update_norm,
1371 cost_evaluations,
1372 jacobian_evaluations,
1373 }),
1374 covariances,
1375 });
1376 }
1377
1378 // Note: Max iterations and timeout checks are now handled inside check_convergence()
1379
1380 iteration += 1;
1381 }
1382 }
1383}
1384// Implement Solver trait
1385impl optimizer::Solver for LevenbergMarquardt {
1386 type Config = LevenbergMarquardtConfig;
1387 type Error = error::ApexError;
1388
1389 fn new() -> Self {
1390 Self::default()
1391 }
1392
1393 fn optimize(
1394 &mut self,
1395 problem: &problem::Problem,
1396 initial_params: &collections::HashMap<
1397 String,
1398 (manifold::ManifoldType, nalgebra::DVector<f64>),
1399 >,
1400 ) -> Result<
1401 optimizer::SolverResult<collections::HashMap<String, problem::VariableEnum>>,
1402 Self::Error,
1403 > {
1404 self.optimize(problem, initial_params)
1405 }
1406}
1407
1408#[cfg(test)]
1409mod tests {
1410 use super::*;
1411 use crate::{core::factors, manifold};
1412 use nalgebra::dvector;
1413 /// Custom Rosenbrock Factor 1: r1 = 10(x2 - x1²)
1414 /// Demonstrates extensibility - custom factors can be defined outside of factors.rs
1415 #[derive(Debug, Clone)]
1416 struct RosenbrockFactor1;
1417
1418 impl factors::Factor for RosenbrockFactor1 {
1419 fn linearize(
1420 &self,
1421 params: &[nalgebra::DVector<f64>],
1422 compute_jacobian: bool,
1423 ) -> (nalgebra::DVector<f64>, Option<nalgebra::DMatrix<f64>>) {
1424 let x1 = params[0][0];
1425 let x2 = params[1][0];
1426
1427 // Residual: r1 = 10(x2 - x1²)
1428 let residual = dvector![10.0 * (x2 - x1 * x1)];
1429
1430 let jacobian = if compute_jacobian {
1431 // Jacobian: ∂r1/∂x1 = -20*x1, ∂r1/∂x2 = 10
1432 let mut jacobian = nalgebra::DMatrix::zeros(1, 2);
1433 jacobian[(0, 0)] = -20.0 * x1;
1434 jacobian[(0, 1)] = 10.0;
1435
1436 Some(jacobian)
1437 } else {
1438 None
1439 };
1440
1441 (residual, jacobian)
1442 }
1443
1444 fn get_dimension(&self) -> usize {
1445 1
1446 }
1447 }
1448
1449 /// Custom Rosenbrock Factor 2: r2 = 1 - x1
1450 /// Demonstrates extensibility - custom factors can be defined outside of factors.rs
1451 #[derive(Debug, Clone)]
1452 struct RosenbrockFactor2;
1453
1454 impl factors::Factor for RosenbrockFactor2 {
1455 fn linearize(
1456 &self,
1457 params: &[nalgebra::DVector<f64>],
1458 compute_jacobian: bool,
1459 ) -> (nalgebra::DVector<f64>, Option<nalgebra::DMatrix<f64>>) {
1460 let x1 = params[0][0];
1461
1462 // Residual: r2 = 1 - x1
1463 let residual = dvector![1.0 - x1];
1464
1465 let jacobian = if compute_jacobian {
1466 // Jacobian: ∂r2/∂x1 = -1
1467 Some(nalgebra::DMatrix::from_element(1, 1, -1.0))
1468 } else {
1469 None
1470 };
1471
1472 (residual, jacobian)
1473 }
1474
1475 fn get_dimension(&self) -> usize {
1476 1
1477 }
1478 }
1479
1480 #[test]
1481 fn test_rosenbrock_optimization() {
1482 // Rosenbrock function test:
1483 // Minimize: r1² + r2² where
1484 // r1 = 10(x2 - x1²)
1485 // r2 = 1 - x1
1486 // Starting point: [-1.2, 1.0]
1487 // Expected minimum: [1.0, 1.0]
1488
1489 let mut problem = problem::Problem::new();
1490 let mut initial_values = collections::HashMap::new();
1491
1492 // Add variables using Rn manifold (Euclidean space)
1493 initial_values.insert(
1494 "x1".to_string(),
1495 (manifold::ManifoldType::RN, dvector![-1.2]),
1496 );
1497 initial_values.insert(
1498 "x2".to_string(),
1499 (manifold::ManifoldType::RN, dvector![1.0]),
1500 );
1501
1502 // Add custom factors (demonstrates extensibility!)
1503 problem.add_residual_block(&["x1", "x2"], Box::new(RosenbrockFactor1), None);
1504 problem.add_residual_block(&["x1"], Box::new(RosenbrockFactor2), None);
1505
1506 // Configure Levenberg-Marquardt optimizer
1507 let config = LevenbergMarquardtConfig::new()
1508 .with_max_iterations(100)
1509 .with_cost_tolerance(1e-8)
1510 .with_parameter_tolerance(1e-8)
1511 .with_gradient_tolerance(1e-10);
1512
1513 let mut solver = LevenbergMarquardt::with_config(config);
1514 let result = solver.optimize(&problem, &initial_values).unwrap();
1515
1516 // Extract final values
1517 let x1_final = result.parameters.get("x1").unwrap().to_vector()[0];
1518 let x2_final = result.parameters.get("x2").unwrap().to_vector()[0];
1519
1520 println!("Rosenbrock optimization result:");
1521 println!(" Status: {:?}", result.status);
1522 println!(" Initial cost: {:.6e}", result.initial_cost);
1523 println!(" Final cost: {:.6e}", result.final_cost);
1524 println!(" Iterations: {}", result.iterations);
1525 println!(" x1: {:.6} (expected 1.0)", x1_final);
1526 println!(" x2: {:.6} (expected 1.0)", x2_final);
1527
1528 // Verify convergence to [1.0, 1.0]
1529 assert!(
1530 matches!(
1531 result.status,
1532 optimizer::OptimizationStatus::Converged
1533 | optimizer::OptimizationStatus::CostToleranceReached
1534 | optimizer::OptimizationStatus::ParameterToleranceReached
1535 | optimizer::OptimizationStatus::GradientToleranceReached
1536 ),
1537 "Optimization should converge"
1538 );
1539 assert!(
1540 (x1_final - 1.0).abs() < 1e-4,
1541 "x1 should converge to 1.0, got {}",
1542 x1_final
1543 );
1544 assert!(
1545 (x2_final - 1.0).abs() < 1e-4,
1546 "x2 should converge to 1.0, got {}",
1547 x2_final
1548 );
1549 assert!(
1550 result.final_cost < 1e-6,
1551 "Final cost should be near zero, got {}",
1552 result.final_cost
1553 );
1554 }
1555}