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::GaussNewton) if well-conditioned with good initialization
76//! - [`DogLeg`](crate::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::LevenbergMarquardt;
104//! use apex_solver::core::problem::Problem;
105//! use apex_solver::JacobianMode;
106//! use std::collections::HashMap;
107//!
108//! # type TestResult = Result<(), Box<dyn std::error::Error>>;
109//! # fn main() -> TestResult {
110//! let mut problem = Problem::new(JacobianMode::Sparse);
111//! // ... add residual blocks (factors) to problem ...
112//!
113//! let initial_values = HashMap::new();
114//! // ... initialize parameters ...
115//!
116//! let mut solver = LevenbergMarquardt::new();
117//! let result = solver.optimize(&problem, &initial_values)?;
118//!
119//! # Ok(())
120//! # }
121//! ```
122//!
123//! ## Advanced configuration
124//!
125//! ```no_run
126//! use apex_solver::optimizer::levenberg_marquardt::{LevenbergMarquardtConfig, LevenbergMarquardt};
127//! use apex_solver::linalg::LinearSolverType;
128//!
129//! # fn main() {
130//! let config = LevenbergMarquardtConfig::new()
131//! .with_max_iterations(100)
132//! .with_cost_tolerance(1e-6)
133//! .with_damping(1e-3) // Initial damping
134//! .with_damping_bounds(1e-12, 1e12) // Min/max damping
135//! .with_jacobi_scaling(true); // Improve conditioning
136//!
137//! let mut solver = LevenbergMarquardt::with_config(config);
138//! # }
139//! ```
140//!
141//! # References
142//!
143//! - Levenberg, K. (1944). "A Method for the Solution of Certain Non-Linear Problems in Least Squares". *Quarterly of Applied Mathematics*.
144//! - Marquardt, D. W. (1963). "An Algorithm for Least-Squares Estimation of Nonlinear Parameters". *Journal of the Society for Industrial and Applied Mathematics*.
145//! - Madsen, K., Nielsen, H. B., & Tingleff, O. (2004). *Methods for Non-Linear Least Squares Problems* (2nd ed.). Chapter 3.
146//! - Nocedal, J. & Wright, S. (2006). *Numerical Optimization* (2nd ed.). Springer. Chapter 10.
147//! - Nielsen, H. B. (1999). "Damping Parameter in Marquardt's Method". Technical Report IMM-REP-1999-05.
148
149use crate::core::problem::{Problem, VariableEnum};
150use crate::error;
151use crate::error::ErrorLogging;
152use crate::linalg::{
153 DenseCholeskySolver, DenseMode, DenseQRSolver, JacobianMode, LinearSolver, LinearSolverType,
154 SchurPreconditioner, SchurVariant, SparseCholeskySolver, SparseMode, SparseQRSolver,
155 SparseSchurComplementSolver, StructureAware,
156};
157use crate::optimizer::{
158 AssemblyBackend, ConvergenceParams, InitializedState, IterationStats, OptObserverVec,
159 OptimizerError, SolverResult, apply_negative_parameter_step, apply_parameter_step,
160 compute_cost,
161};
162use apex_manifolds::ManifoldType;
163
164use faer::Mat;
165use nalgebra::DVector;
166use std::collections::HashMap;
167use std::time::{Duration, Instant};
168use tracing::debug;
169
170/// Configuration parameters for the Levenberg-Marquardt optimizer.
171///
172/// Controls the adaptive damping strategy, convergence criteria, and numerical stability
173/// enhancements for the Levenberg-Marquardt algorithm.
174///
175/// # Builder Pattern
176///
177/// All configuration options can be set using the builder pattern:
178///
179/// ```
180/// use apex_solver::optimizer::levenberg_marquardt::LevenbergMarquardtConfig;
181///
182/// let config = LevenbergMarquardtConfig::new()
183/// .with_max_iterations(100)
184/// .with_damping(1e-3)
185/// .with_damping_bounds(1e-12, 1e12)
186/// .with_jacobi_scaling(true);
187/// ```
188///
189/// # Damping Parameter Behavior
190///
191/// The damping parameter λ controls the trade-off between Gauss-Newton and gradient descent:
192///
193/// - **Initial damping** (`damping`): Starting value (default: 1e-4)
194/// - **Damping bounds** (`damping_min`, `damping_max`): Valid range (default: 1e-12 to 1e12)
195/// - **Adaptation**: Automatically adjusted based on step quality using Nielsen's formula
196///
197/// # Convergence Criteria
198///
199/// The optimizer terminates when ANY of the following conditions is met:
200///
201/// - **Cost tolerance**: `|cost_k - cost_{k-1}| < cost_tolerance`
202/// - **Parameter tolerance**: `||step|| < parameter_tolerance`
203/// - **Gradient tolerance**: `||J^T·r|| < gradient_tolerance`
204/// - **Maximum iterations**: `iteration >= max_iterations`
205/// - **Timeout**: `elapsed_time >= timeout`
206///
207/// # See Also
208///
209/// - [`LevenbergMarquardt`] - The solver that uses this configuration
210/// - [`GaussNewtonConfig`](crate::GaussNewtonConfig) - Undamped variant
211/// - [`DogLegConfig`](crate::DogLegConfig) - Trust region alternative
212#[derive(Clone)]
213pub struct LevenbergMarquardtConfig {
214 /// Type of linear solver for the linear systems
215 pub linear_solver_type: LinearSolverType,
216 /// Maximum number of iterations
217 pub max_iterations: usize,
218 /// Convergence tolerance for cost function
219 pub cost_tolerance: f64,
220 /// Convergence tolerance for parameter updates
221 pub parameter_tolerance: f64,
222 /// Convergence tolerance for gradient norm
223 pub gradient_tolerance: f64,
224 /// Timeout duration
225 pub timeout: Option<Duration>,
226 /// Initial damping parameter
227 pub damping: f64,
228 /// Minimum damping parameter
229 pub damping_min: f64,
230 /// Maximum damping parameter
231 pub damping_max: f64,
232 /// Damping increase factor (when step rejected)
233 pub damping_increase_factor: f64,
234 /// Damping decrease factor (when step accepted)
235 pub damping_decrease_factor: f64,
236 /// Damping nu parameter
237 pub damping_nu: f64,
238 /// Trust region radius
239 pub trust_region_radius: f64,
240 /// Minimum step quality for acceptance
241 pub min_step_quality: f64,
242 /// Good step quality threshold
243 pub good_step_quality: f64,
244 /// Minimum diagonal value for regularization
245 pub min_diagonal: f64,
246 /// Maximum diagonal value for regularization
247 pub max_diagonal: f64,
248 /// Minimum objective function cutoff (optional early termination)
249 ///
250 /// If set, optimization terminates when cost falls below this threshold.
251 /// Useful for early stopping when a "good enough" solution is acceptable.
252 ///
253 /// Default: None (disabled)
254 pub min_cost_threshold: Option<f64>,
255 /// Minimum trust region radius before termination
256 ///
257 /// When the trust region radius falls below this value, the optimizer
258 /// terminates as it indicates the search has converged or the problem
259 /// is ill-conditioned. Matches Ceres Solver's min_trust_region_radius.
260 ///
261 /// Default: 1e-32 (Ceres-compatible)
262 pub min_trust_region_radius: f64,
263 /// Maximum condition number for Jacobian matrix (optional check)
264 ///
265 /// If set, the optimizer checks if condition_number(J^T*J) exceeds this
266 /// threshold and terminates with IllConditionedJacobian status.
267 /// Note: Computing condition number is expensive, so this is disabled by default.
268 ///
269 /// Default: None (disabled)
270 pub max_condition_number: Option<f64>,
271 /// Minimum relative cost decrease for step acceptance
272 ///
273 /// Used in computing step quality (rho = actual_reduction / predicted_reduction).
274 /// Steps with rho < min_relative_decrease are rejected. Matches Ceres Solver's
275 /// min_relative_decrease parameter.
276 ///
277 /// Default: 1e-3 (Ceres-compatible)
278 pub min_relative_decrease: f64,
279 /// Use Jacobi column scaling (preconditioning)
280 ///
281 /// When enabled, normalizes Jacobian columns by their L2 norm before solving.
282 /// This can improve convergence for problems with mixed parameter scales
283 /// (e.g., positions in meters + angles in radians) but adds ~5-10% overhead.
284 ///
285 /// Default: false (to avoid performance overhead and faster convergence)
286 pub use_jacobi_scaling: bool,
287 /// Compute per-variable covariance matrices (uncertainty estimation)
288 ///
289 /// When enabled, computes covariance by inverting the Hessian matrix after
290 /// convergence. The full covariance matrix is extracted into per-variable
291 /// blocks stored in both Variable structs and SolverResult.
292 ///
293 /// Default: false (to avoid performance overhead)
294 pub compute_covariances: bool,
295 /// Schur complement solver variant (for bundle adjustment problems)
296 ///
297 /// When using LinearSolverType::SparseSchurComplement, this determines which
298 /// variant of the Schur complement method to use:
299 /// - Sparse: Direct sparse Cholesky factorization (most accurate, moderate speed)
300 /// - Iterative: Preconditioned Conjugate Gradients (memory efficient, good for large problems)
301 /// - PowerSeries: Power series approximation (fastest, less accurate)
302 ///
303 /// Default: Sparse
304 pub schur_variant: SchurVariant,
305 /// Schur complement preconditioner type
306 ///
307 /// Determines the preconditioning strategy for iterative Schur methods:
308 /// - Diagonal: Simple diagonal preconditioner (fast, less effective)
309 /// - BlockDiagonal: Block-diagonal preconditioner (balanced)
310 /// - IncompleteCholesky: Incomplete Cholesky factorization (slower, more effective)
311 ///
312 /// Default: Diagonal
313 pub schur_preconditioner: SchurPreconditioner,
314 // Note: Visualization is now handled via the observer pattern.
315 // Use `solver.add_observer(RerunObserver::new(true)?)` to enable visualization.
316 // This provides cleaner separation of concerns and allows multiple observers.
317}
318
319impl Default for LevenbergMarquardtConfig {
320 fn default() -> Self {
321 Self {
322 linear_solver_type: LinearSolverType::default(),
323 // Ceres Solver default: 50 (changed from 100 for compatibility)
324 max_iterations: 50,
325 // Ceres Solver default: 1e-6 (changed from 1e-8 for compatibility)
326 cost_tolerance: 1e-6,
327 // Ceres Solver default: 1e-8 (unchanged)
328 parameter_tolerance: 1e-8,
329 // Ceres Solver default: 1e-10 (changed from 1e-8 for compatibility)
330 // Note: Typically should be 1e-4 * cost_tolerance per Ceres docs
331 gradient_tolerance: 1e-10,
332 timeout: None,
333 damping: 1e-3, // Increased from 1e-4 for better initial convergence on BA
334 damping_min: 1e-12,
335 damping_max: 1e12,
336 damping_increase_factor: 10.0,
337 damping_decrease_factor: 0.3,
338 damping_nu: 2.0,
339 trust_region_radius: 1e4,
340 min_step_quality: 0.0,
341 good_step_quality: 0.75,
342 min_diagonal: 1e-6,
343 max_diagonal: 1e32,
344 // New Ceres-compatible parameters
345 min_cost_threshold: None,
346 min_trust_region_radius: 1e-32,
347 max_condition_number: None,
348 min_relative_decrease: 1e-3,
349 // Existing parameters
350 // Jacobi scaling disabled by default for Schur solvers (incompatible with block structure)
351 // Enable manually for Cholesky/QR solvers on mixed-scale problems
352 use_jacobi_scaling: false,
353 compute_covariances: false,
354 // Schur complement parameters
355 schur_variant: SchurVariant::default(),
356 schur_preconditioner: SchurPreconditioner::default(),
357 }
358 }
359}
360
361impl LevenbergMarquardtConfig {
362 /// Create a new Levenberg-Marquardt configuration with default values.
363 pub fn new() -> Self {
364 Self::default()
365 }
366
367 /// Set the linear solver type
368 pub fn with_linear_solver_type(mut self, linear_solver_type: LinearSolverType) -> Self {
369 self.linear_solver_type = linear_solver_type;
370 self
371 }
372
373 /// Set the maximum number of iterations
374 pub fn with_max_iterations(mut self, max_iterations: usize) -> Self {
375 self.max_iterations = max_iterations;
376 self
377 }
378
379 /// Set the cost tolerance
380 pub fn with_cost_tolerance(mut self, cost_tolerance: f64) -> Self {
381 self.cost_tolerance = cost_tolerance;
382 self
383 }
384
385 /// Set the parameter tolerance
386 pub fn with_parameter_tolerance(mut self, parameter_tolerance: f64) -> Self {
387 self.parameter_tolerance = parameter_tolerance;
388 self
389 }
390
391 /// Set the gradient tolerance
392 pub fn with_gradient_tolerance(mut self, gradient_tolerance: f64) -> Self {
393 self.gradient_tolerance = gradient_tolerance;
394 self
395 }
396
397 /// Set the timeout duration
398 pub fn with_timeout(mut self, timeout: Duration) -> Self {
399 self.timeout = Some(timeout);
400 self
401 }
402
403 /// Set the initial damping parameter.
404 pub fn with_damping(mut self, damping: f64) -> Self {
405 self.damping = damping;
406 self
407 }
408
409 /// Set the damping parameter bounds.
410 pub fn with_damping_bounds(mut self, min: f64, max: f64) -> Self {
411 self.damping_min = min;
412 self.damping_max = max;
413 self
414 }
415
416 /// Set the damping adjustment factors.
417 pub fn with_damping_factors(mut self, increase: f64, decrease: f64) -> Self {
418 self.damping_increase_factor = increase;
419 self.damping_decrease_factor = decrease;
420 self
421 }
422
423 /// Set the trust region parameters.
424 pub fn with_trust_region(mut self, radius: f64, min_quality: f64, good_quality: f64) -> Self {
425 self.trust_region_radius = radius;
426 self.min_step_quality = min_quality;
427 self.good_step_quality = good_quality;
428 self
429 }
430
431 /// Set minimum objective function cutoff for early termination.
432 ///
433 /// When set, optimization terminates with MinCostThresholdReached status
434 /// if the cost falls below this threshold. Useful for early stopping when
435 /// a "good enough" solution is acceptable.
436 pub fn with_min_cost_threshold(mut self, min_cost: f64) -> Self {
437 self.min_cost_threshold = Some(min_cost);
438 self
439 }
440
441 /// Set minimum trust region radius before termination.
442 ///
443 /// When the trust region radius falls below this value, optimization
444 /// terminates with TrustRegionRadiusTooSmall status.
445 /// Default: 1e-32 (Ceres-compatible)
446 pub fn with_min_trust_region_radius(mut self, min_radius: f64) -> Self {
447 self.min_trust_region_radius = min_radius;
448 self
449 }
450
451 /// Set maximum condition number for Jacobian matrix.
452 ///
453 /// If set, the optimizer checks if condition_number(J^T*J) exceeds this
454 /// threshold and terminates with IllConditionedJacobian status.
455 /// Note: Computing condition number is expensive, disabled by default.
456 pub fn with_max_condition_number(mut self, max_cond: f64) -> Self {
457 self.max_condition_number = Some(max_cond);
458 self
459 }
460
461 /// Set minimum relative cost decrease for step acceptance.
462 ///
463 /// Steps with rho = (actual_reduction / predicted_reduction) below this
464 /// threshold are rejected. Default: 1e-3 (Ceres-compatible)
465 pub fn with_min_relative_decrease(mut self, min_decrease: f64) -> Self {
466 self.min_relative_decrease = min_decrease;
467 self
468 }
469
470 /// Enable or disable Jacobi column scaling (preconditioning).
471 ///
472 /// When enabled, normalizes Jacobian columns by their L2 norm before solving.
473 /// Can improve convergence for mixed-scale problems but adds ~5-10% overhead.
474 pub fn with_jacobi_scaling(mut self, use_jacobi_scaling: bool) -> Self {
475 self.use_jacobi_scaling = use_jacobi_scaling;
476 self
477 }
478
479 /// Enable or disable covariance computation (uncertainty estimation).
480 ///
481 /// When enabled, computes the full covariance matrix by inverting the Hessian
482 /// after convergence, then extracts per-variable covariance blocks.
483 pub fn with_compute_covariances(mut self, compute_covariances: bool) -> Self {
484 self.compute_covariances = compute_covariances;
485 self
486 }
487
488 /// Set Schur complement solver variant
489 pub fn with_schur_variant(mut self, variant: SchurVariant) -> Self {
490 self.schur_variant = variant;
491 self
492 }
493
494 /// Set Schur complement preconditioner
495 pub fn with_schur_preconditioner(mut self, preconditioner: SchurPreconditioner) -> Self {
496 self.schur_preconditioner = preconditioner;
497 self
498 }
499
500 /// Configuration optimized for bundle adjustment problems.
501 ///
502 /// This preset uses settings tuned for large-scale bundle adjustment:
503 /// - **Schur complement solver** with iterative PCG (memory efficient)
504 /// - **Schur-Jacobi preconditioner** (Ceres-style, best PCG convergence)
505 /// - **Moderate initial damping** (1e-3) - not too aggressive
506 /// - **200 max iterations** (BA often needs more iterations for full convergence)
507 /// - **Very tight tolerances** matching Ceres Solver for accurate reconstruction
508 ///
509 /// This configuration matches Ceres Solver's recommended BA settings and
510 /// should achieve similar convergence quality.
511 ///
512 /// # Example
513 ///
514 /// ```
515 /// use apex_solver::optimizer::levenberg_marquardt::LevenbergMarquardtConfig;
516 ///
517 /// let config = LevenbergMarquardtConfig::for_bundle_adjustment();
518 /// ```
519 pub fn for_bundle_adjustment() -> Self {
520 Self::default()
521 .with_linear_solver_type(LinearSolverType::SparseSchurComplement)
522 .with_schur_variant(SchurVariant::Iterative)
523 .with_schur_preconditioner(SchurPreconditioner::SchurJacobi)
524 .with_damping(1e-3) // Moderate initial damping (Ceres default)
525 .with_max_iterations(20) // Reduced for early stop when RMSE < 1px
526 // Match Ceres tolerances for faster convergence
527 .with_cost_tolerance(1e-6) // Ceres function_tolerance (was 1e-12)
528 .with_parameter_tolerance(1e-8) // Ceres parameter_tolerance (was 1e-14)
529 .with_gradient_tolerance(1e-10) // Relaxed (was 1e-16)
530 }
531
532 /// Enable real-time visualization (graphical debugging).
533 ///
534 /// When enabled, optimization progress is logged to a Rerun viewer with:
535 /// - Time series plots of cost, gradient norm, damping, step quality
536 /// - Sparse Hessian matrix visualization as heat map
537 /// - Gradient vector visualization
538 /// - Real-time manifold state updates (for SE2/SE3 problems)
539 ///
540 /// **Note:** Requires the `visualization` feature to be enabled in `Cargo.toml`.
541 /// Use `verbose` for terminal logging.
542 ///
543 /// # Arguments
544 ///
545 /// * `enable` - Whether to enable visualization
546 // Note: with_visualization() method has been removed.
547 // Use the observer pattern instead:
548 // let mut solver = LevenbergMarquardt::with_config(config);
549 // solver.add_observer(RerunObserver::new(true)?);
550 // This provides cleaner separation and allows multiple observers.
551 /// Print configuration parameters (verbose mode only)
552 pub fn print_configuration(&self) {
553 debug!(
554 "Configuration:\n Solver: Levenberg-Marquardt\n Linear solver: {:?}\n Convergence Criteria:\n Max iterations: {}\n Cost tolerance: {:.2e}\n Parameter tolerance: {:.2e}\n Gradient tolerance: {:.2e}\n Timeout: {:?}\n Damping Parameters:\n Initial damping: {:.2e}\n Damping range: [{:.2e}, {:.2e}]\n Increase factor: {:.2}\n Decrease factor: {:.2}\n Trust Region:\n Initial radius: {:.2e}\n Min step quality: {:.2}\n Good step quality: {:.2}\n Numerical Settings:\n Jacobi scaling: {}\n Compute covariances: {}",
555 self.linear_solver_type,
556 self.max_iterations,
557 self.cost_tolerance,
558 self.parameter_tolerance,
559 self.gradient_tolerance,
560 self.timeout,
561 self.damping,
562 self.damping_min,
563 self.damping_max,
564 self.damping_increase_factor,
565 self.damping_decrease_factor,
566 self.trust_region_radius,
567 self.min_step_quality,
568 self.good_step_quality,
569 if self.use_jacobi_scaling {
570 "enabled"
571 } else {
572 "disabled"
573 },
574 if self.compute_covariances {
575 "enabled"
576 } else {
577 "disabled"
578 }
579 );
580 }
581}
582
583/// Result from step computation
584struct StepResult {
585 step: Mat<f64>,
586 gradient_norm: f64,
587 predicted_reduction: f64,
588}
589
590/// Result from step evaluation
591struct StepEvaluation {
592 accepted: bool,
593 cost_reduction: f64,
594 rho: f64,
595}
596
597/// Levenberg-Marquardt solver for nonlinear least squares optimization.
598///
599/// Implements the damped Gauss-Newton method with adaptive damping parameter λ that
600/// interpolates between Gauss-Newton and gradient descent based on step quality.
601///
602/// # Algorithm
603///
604/// At each iteration k:
605/// 1. Compute residual `r(xₖ)` and Jacobian `J(xₖ)`
606/// 2. Solve augmented system: `(J^T·J + λI)·h = -J^T·r`
607/// 3. Evaluate step quality: `ρ = (actual reduction) / (predicted reduction)`
608/// 4. If `ρ > 0`: Accept step and update `xₖ₊₁ = xₖ ⊕ h`, decrease λ
609/// 5. If `ρ ≤ 0`: Reject step (keep `xₖ₊₁ = xₖ`), increase λ
610/// 6. Check convergence criteria
611///
612/// The damping parameter λ is updated using Nielsen's smooth formula:
613/// `λₖ₊₁ = λₖ · max(1/3, 1 - (2ρ - 1)³)` for accepted steps,
614/// or `λₖ₊₁ = λₖ · ν` (with increasing ν) for rejected steps.
615///
616/// # Examples
617///
618/// ```no_run
619/// use apex_solver::optimizer::levenberg_marquardt::{LevenbergMarquardtConfig, LevenbergMarquardt};
620/// use apex_solver::linalg::LinearSolverType;
621///
622/// # fn main() {
623/// let config = LevenbergMarquardtConfig::new()
624/// .with_max_iterations(100)
625/// .with_damping(1e-3)
626/// .with_damping_bounds(1e-12, 1e12)
627/// .with_jacobi_scaling(true);
628///
629/// let mut solver = LevenbergMarquardt::with_config(config);
630/// # }
631/// ```
632///
633/// # See Also
634///
635/// - [`LevenbergMarquardtConfig`] - Configuration options
636/// - [`GaussNewton`](crate::GaussNewton) - Undamped variant (faster but less robust)
637/// - [`DogLeg`](crate::DogLeg) - Alternative trust region method
638pub struct LevenbergMarquardt {
639 config: LevenbergMarquardtConfig,
640 jacobi_scaling: Option<Vec<f64>>,
641 observers: OptObserverVec,
642}
643
644impl Default for LevenbergMarquardt {
645 fn default() -> Self {
646 Self::new()
647 }
648}
649
650impl LevenbergMarquardt {
651 /// Create a new Levenberg-Marquardt solver with default configuration.
652 pub fn new() -> Self {
653 Self::with_config(LevenbergMarquardtConfig::default())
654 }
655
656 /// Create a new Levenberg-Marquardt solver with the given configuration.
657 pub fn with_config(config: LevenbergMarquardtConfig) -> Self {
658 Self {
659 config,
660 jacobi_scaling: None,
661 observers: OptObserverVec::new(),
662 }
663 }
664
665 /// Add an observer to monitor optimization progress.
666 ///
667 /// Observers are notified at each iteration with the current variable values.
668 /// This enables real-time visualization, logging, metrics collection, etc.
669 ///
670 /// # Arguments
671 ///
672 /// * `observer` - Any type implementing `OptObserver`
673 ///
674 /// # Examples
675 ///
676 /// ```no_run
677 /// use apex_solver::{LevenbergMarquardt, LevenbergMarquardtConfig};
678 /// # use apex_solver::core::problem::Problem;
679 /// # use std::collections::HashMap;
680 ///
681 /// # fn main() -> Result<(), Box<dyn std::error::Error>> {
682 /// let mut solver = LevenbergMarquardt::new();
683 ///
684 /// #[cfg(feature = "visualization")]
685 /// {
686 /// use apex_solver::observers::RerunObserver;
687 /// let rerun_observer = RerunObserver::new(true)?;
688 /// solver.add_observer(rerun_observer);
689 /// }
690 ///
691 /// // ... optimize ...
692 /// # Ok(())
693 /// # }
694 /// ```
695 pub fn add_observer(&mut self, observer: impl crate::optimizer::OptObserver + 'static) {
696 self.observers.add(observer);
697 }
698
699 /// Update damping parameter based on step quality using trust region approach
700 /// Reference: Introduction to Optimization and Data Fitting
701 /// Algorithm 6.18
702 fn update_damping(&mut self, rho: f64) -> bool {
703 if rho > 0.0 {
704 // Step accepted - decrease damping
705 let coff = 2.0 * rho - 1.0;
706 self.config.damping *= (1.0_f64 / 3.0).max(1.0 - coff * coff * coff);
707 self.config.damping = self.config.damping.max(self.config.damping_min);
708 self.config.damping_nu = 2.0;
709 true
710 } else {
711 // Step rejected - increase damping
712 self.config.damping *= self.config.damping_nu;
713 self.config.damping_nu *= 2.0;
714 self.config.damping = self.config.damping.min(self.config.damping_max);
715 false
716 }
717 }
718
719 /// Compute predicted cost reduction from linear model
720 /// Standard LM formula: 0.5 * step^T * (damping * step - gradient)
721 fn compute_predicted_reduction(&self, step: &Mat<f64>, gradient: &Mat<f64>) -> f64 {
722 // Standard Levenberg-Marquardt predicted reduction formula
723 // predicted_reduction = -step^T * gradient - 0.5 * step^T * H * step
724 // = 0.5 * step^T * (damping * step - gradient)
725 let diff = self.config.damping * step - gradient;
726 (0.5 * step.transpose() * &diff)[(0, 0)]
727 }
728
729 /// Compute optimization step by solving the augmented system (generic over assembly mode).
730 fn compute_step_generic<M: AssemblyBackend>(
731 &self,
732 residuals: &Mat<f64>,
733 scaled_jacobian: &M::Jacobian,
734 linear_solver: &mut dyn LinearSolver<M>,
735 ) -> Result<StepResult, OptimizerError> {
736 // Solve augmented equation: (J_scaled^T * J_scaled + λI) * dx_scaled = -J_scaled^T * r
737 let residuals_owned = residuals.as_ref().to_owned();
738 let scaled_step = linear_solver
739 .solve_augmented_equation(&residuals_owned, scaled_jacobian, self.config.damping)
740 .map_err(|e| OptimizerError::LinearSolveFailed(e.to_string()).log_with_source(e))?;
741
742 // Get cached gradient from the solver
743 let gradient = linear_solver.get_gradient().ok_or_else(|| {
744 OptimizerError::NumericalInstability("Gradient not available".into()).log()
745 })?;
746 let gradient_norm = gradient.norm_l2();
747
748 // Apply inverse Jacobi scaling to get final step (if enabled)
749 let step = if self.config.use_jacobi_scaling {
750 let scaling = self
751 .jacobi_scaling
752 .as_ref()
753 .ok_or_else(|| OptimizerError::JacobiScalingNotInitialized.log())?;
754 M::apply_inverse_scaling(&scaled_step, scaling)
755 } else {
756 scaled_step
757 };
758
759 // Compute predicted reduction using scaled values
760 let predicted_reduction = self.compute_predicted_reduction(&step, gradient);
761
762 Ok(StepResult {
763 step,
764 gradient_norm,
765 predicted_reduction,
766 })
767 }
768
769 /// Evaluate and apply step, handling acceptance/rejection based on step quality
770 fn evaluate_and_apply_step(
771 &mut self,
772 step_result: &StepResult,
773 state: &mut InitializedState,
774 problem: &Problem,
775 ) -> error::ApexSolverResult<StepEvaluation> {
776 // Apply parameter updates using manifold operations
777 let _step_norm = apply_parameter_step(
778 &mut state.variables,
779 step_result.step.as_ref(),
780 &state.sorted_vars,
781 );
782
783 // Compute new cost (residual only, no Jacobian needed for step evaluation)
784 let new_residual = problem.compute_residual_sparse(&state.variables)?;
785 let new_cost = compute_cost(&new_residual);
786
787 // Compute step quality
788 let rho = crate::optimizer::compute_step_quality(
789 state.current_cost,
790 new_cost,
791 step_result.predicted_reduction,
792 );
793
794 // Update damping and decide whether to accept step
795 let accepted = self.update_damping(rho);
796
797 let cost_reduction = if accepted {
798 // Accept the step - parameters already updated
799 let reduction = state.current_cost - new_cost;
800 state.current_cost = new_cost;
801 reduction
802 } else {
803 // Reject the step - revert parameter changes
804 apply_negative_parameter_step(
805 &mut state.variables,
806 step_result.step.as_ref(),
807 &state.sorted_vars,
808 );
809 0.0
810 };
811
812 Ok(StepEvaluation {
813 accepted,
814 cost_reduction,
815 rho,
816 })
817 }
818
819 /// Run optimization using the specified assembly mode and linear solver.
820 ///
821 /// This is the core generic optimization loop. The public `optimize()` method
822 /// dispatches to this based on `LinearSolverType`.
823 fn optimize_with_mode<M: AssemblyBackend>(
824 &mut self,
825 problem: &Problem,
826 initial_params: &HashMap<String, (ManifoldType, DVector<f64>)>,
827 linear_solver: &mut dyn LinearSolver<M>,
828 ) -> Result<SolverResult<HashMap<String, VariableEnum>>, error::ApexSolverError> {
829 let start_time = Instant::now();
830 let mut iteration = 0;
831 let mut cost_evaluations = 1;
832 let mut jacobian_evaluations = 0;
833 let mut successful_steps = 0;
834 let mut unsuccessful_steps = 0;
835
836 // Initialize optimization state
837 let mut state = crate::optimizer::initialize_optimization_state(problem, initial_params)?;
838
839 // Initialize summary tracking variables
840 let mut max_gradient_norm: f64 = 0.0;
841 let mut max_parameter_update_norm: f64 = 0.0;
842 let mut total_cost_reduction = 0.0;
843 let mut final_gradient_norm;
844 let mut final_parameter_update_norm;
845
846 // Initialize iteration statistics tracking
847 let mut iteration_stats = Vec::with_capacity(self.config.max_iterations);
848 let mut previous_cost = state.current_cost;
849
850 // Print configuration and header if debug level is enabled
851 if tracing::enabled!(tracing::Level::DEBUG) {
852 self.config.print_configuration();
853 IterationStats::print_header();
854 }
855
856 // Main optimization loop
857 loop {
858 let iter_start = Instant::now();
859
860 // Evaluate residuals and Jacobian using the assembly mode
861 let (residuals, jacobian) = M::assemble(
862 problem,
863 &state.variables,
864 &state.variable_index_map,
865 state.symbolic_structure.as_ref(),
866 state.total_dof,
867 )?;
868 jacobian_evaluations += 1;
869
870 // Process Jacobian (apply scaling if enabled)
871 let scaled_jacobian = if self.config.use_jacobi_scaling {
872 crate::optimizer::process_jacobian_generic::<M>(
873 &jacobian,
874 &mut self.jacobi_scaling,
875 iteration,
876 )?
877 } else {
878 jacobian
879 };
880
881 // Compute optimization step
882 let step_result =
883 self.compute_step_generic::<M>(&residuals, &scaled_jacobian, linear_solver)?;
884
885 // Update tracking variables
886 max_gradient_norm = max_gradient_norm.max(step_result.gradient_norm);
887 final_gradient_norm = step_result.gradient_norm;
888 let step_norm = step_result.step.norm_l2();
889 max_parameter_update_norm = max_parameter_update_norm.max(step_norm);
890 final_parameter_update_norm = step_norm;
891
892 // Evaluate and apply step (handles accept/reject)
893 let step_eval = self.evaluate_and_apply_step(&step_result, &mut state, problem)?;
894 cost_evaluations += 1;
895
896 // Update counters based on acceptance
897 if step_eval.accepted {
898 successful_steps += 1;
899 total_cost_reduction += step_eval.cost_reduction;
900 } else {
901 unsuccessful_steps += 1;
902 }
903
904 // OPTIMIZATION: Only collect iteration statistics if debug level is enabled
905 if tracing::enabled!(tracing::Level::DEBUG) {
906 let iter_elapsed_ms = iter_start.elapsed().as_secs_f64() * 1000.0;
907 let total_elapsed_ms = start_time.elapsed().as_secs_f64() * 1000.0;
908
909 let stats = IterationStats {
910 iteration,
911 cost: state.current_cost,
912 cost_change: previous_cost - state.current_cost,
913 gradient_norm: step_result.gradient_norm,
914 step_norm,
915 tr_ratio: step_eval.rho,
916 tr_radius: self.config.damping,
917 ls_iter: 0,
918 iter_time_ms: iter_elapsed_ms,
919 total_time_ms: total_elapsed_ms,
920 accepted: step_eval.accepted,
921 };
922
923 iteration_stats.push(stats.clone());
924 stats.print_line();
925 }
926
927 previous_cost = state.current_cost;
928
929 // Notify all observers with current state
930 crate::optimizer::notify_observers_generic::<M>(
931 &mut self.observers,
932 &state.variables,
933 iteration,
934 state.current_cost,
935 step_result.gradient_norm,
936 Some(self.config.damping),
937 step_norm,
938 Some(step_eval.rho),
939 linear_solver,
940 );
941
942 // Check convergence
943 let elapsed = start_time.elapsed();
944 let parameter_norm = crate::optimizer::compute_parameter_norm(&state.variables);
945 let new_cost = state.current_cost;
946 let cost_before_step = if step_eval.accepted {
947 state.current_cost + step_eval.cost_reduction
948 } else {
949 state.current_cost
950 };
951
952 if let Some(status) = crate::optimizer::check_convergence(&ConvergenceParams {
953 iteration,
954 current_cost: cost_before_step,
955 new_cost,
956 parameter_norm,
957 parameter_update_norm: step_norm,
958 gradient_norm: step_result.gradient_norm,
959 elapsed,
960 step_accepted: step_eval.accepted,
961 max_iterations: self.config.max_iterations,
962 gradient_tolerance: self.config.gradient_tolerance,
963 parameter_tolerance: self.config.parameter_tolerance,
964 cost_tolerance: self.config.cost_tolerance,
965 min_cost_threshold: self.config.min_cost_threshold,
966 timeout: self.config.timeout,
967 trust_region_radius: Some(self.config.trust_region_radius),
968 min_trust_region_radius: Some(self.config.min_trust_region_radius),
969 }) {
970 if tracing::enabled!(tracing::Level::DEBUG) {
971 let summary = crate::optimizer::create_optimizer_summary(
972 "Levenberg-Marquardt",
973 state.initial_cost,
974 state.current_cost,
975 iteration + 1,
976 Some(successful_steps),
977 Some(unsuccessful_steps),
978 max_gradient_norm,
979 final_gradient_norm,
980 max_parameter_update_norm,
981 final_parameter_update_norm,
982 total_cost_reduction,
983 elapsed,
984 iteration_stats.clone(),
985 status.clone(),
986 Some(self.config.damping),
987 None,
988 Some(step_eval.rho),
989 );
990 debug!("{}", summary);
991 }
992
993 // Compute covariances if enabled
994 let covariances = if self.config.compute_covariances {
995 problem.compute_and_set_covariances_generic::<M>(
996 linear_solver,
997 &mut state.variables,
998 &state.variable_index_map,
999 )
1000 } else {
1001 None
1002 };
1003
1004 // Notify observers that optimization is complete
1005 let final_parameters: HashMap<String, VariableEnum> = state
1006 .variables
1007 .iter()
1008 .map(|(k, v)| (k.clone(), v.clone()))
1009 .collect();
1010 self.observers
1011 .notify_complete(&final_parameters, iteration + 1);
1012
1013 return Ok(crate::optimizer::build_solver_result(
1014 status,
1015 iteration + 1,
1016 state,
1017 elapsed,
1018 final_gradient_norm,
1019 final_parameter_update_norm,
1020 cost_evaluations,
1021 jacobian_evaluations,
1022 covariances,
1023 ));
1024 }
1025
1026 iteration += 1;
1027 }
1028 }
1029
1030 /// Run optimization, dispatching based on `problem.jacobian_mode`.
1031 ///
1032 /// - `JacobianMode::Dense` → always uses `DenseCholeskySolver`
1033 /// - `JacobianMode::Sparse` → uses the solver selected by `config.linear_solver_type`
1034 pub fn optimize(
1035 &mut self,
1036 problem: &Problem,
1037 initial_params: &HashMap<String, (ManifoldType, DVector<f64>)>,
1038 ) -> Result<SolverResult<HashMap<String, VariableEnum>>, error::ApexSolverError> {
1039 match problem.jacobian_mode {
1040 JacobianMode::Dense => match self.config.linear_solver_type {
1041 LinearSolverType::DenseQR => {
1042 let mut solver = DenseQRSolver::new();
1043 self.optimize_with_mode::<DenseMode>(problem, initial_params, &mut solver)
1044 }
1045 _ => {
1046 let mut solver = DenseCholeskySolver::new();
1047 self.optimize_with_mode::<DenseMode>(problem, initial_params, &mut solver)
1048 }
1049 },
1050 JacobianMode::Sparse => match self.config.linear_solver_type {
1051 LinearSolverType::SparseQR => {
1052 let mut solver = SparseQRSolver::new();
1053 self.optimize_with_mode::<SparseMode>(problem, initial_params, &mut solver)
1054 }
1055 LinearSolverType::SparseSchurComplement => {
1056 // Schur complement needs variable structure before the first solve.
1057 // Initialize state once to get variables, then create and prepare the solver.
1058 // optimize_with_mode will re-initialize state internally, which is acceptable
1059 // since Schur structure initialization is cheap.
1060 let state =
1061 crate::optimizer::initialize_optimization_state(problem, initial_params)?;
1062 let mut solver = SparseSchurComplementSolver::new()
1063 .with_variant(self.config.schur_variant)
1064 .with_preconditioner(self.config.schur_preconditioner);
1065 solver
1066 .initialize_structure(&state.variables, &state.variable_index_map)
1067 .map_err(|e| {
1068 OptimizerError::LinearSolveFailed(format!(
1069 "Failed to initialize Schur solver: {}",
1070 e
1071 ))
1072 .log()
1073 })?;
1074 self.optimize_with_mode::<SparseMode>(problem, initial_params, &mut solver)
1075 }
1076 _ => {
1077 // SparseCholesky (default) or DenseCholesky with sparse mode → SparseCholeskySolver
1078 let mut solver = SparseCholeskySolver::new();
1079 self.optimize_with_mode::<SparseMode>(problem, initial_params, &mut solver)
1080 }
1081 },
1082 }
1083 }
1084}
1085impl crate::optimizer::Optimizer for LevenbergMarquardt {
1086 fn optimize(
1087 &mut self,
1088 problem: &Problem,
1089 initial_params: &HashMap<String, (ManifoldType, DVector<f64>)>,
1090 ) -> Result<SolverResult<HashMap<String, VariableEnum>>, crate::error::ApexSolverError> {
1091 self.optimize(problem, initial_params)
1092 }
1093}
1094
1095#[cfg(test)]
1096mod tests {
1097 use super::*;
1098 use crate::factors::Factor;
1099 use crate::optimizer::OptimizationStatus;
1100 use nalgebra::{DMatrix, dvector};
1101
1102 type TestResult = Result<(), Box<dyn std::error::Error>>;
1103 /// Custom Rosenbrock Factor 1: r1 = 10(x2 - x1²)
1104 /// Demonstrates extensibility - custom factors can be defined outside of factors.rs
1105 #[derive(Debug, Clone)]
1106 struct RosenbrockFactor1;
1107
1108 impl Factor for RosenbrockFactor1 {
1109 fn linearize(
1110 &self,
1111 params: &[DVector<f64>],
1112 compute_jacobian: bool,
1113 ) -> (DVector<f64>, Option<DMatrix<f64>>) {
1114 let x1 = params[0][0];
1115 let x2 = params[1][0];
1116
1117 // Residual: r1 = 10(x2 - x1²)
1118 let residual = dvector![10.0 * (x2 - x1 * x1)];
1119
1120 let jacobian = if compute_jacobian {
1121 // Jacobian: ∂r1/∂x1 = -20*x1, ∂r1/∂x2 = 10
1122 let mut jacobian = DMatrix::zeros(1, 2);
1123 jacobian[(0, 0)] = -20.0 * x1;
1124 jacobian[(0, 1)] = 10.0;
1125
1126 Some(jacobian)
1127 } else {
1128 None
1129 };
1130
1131 (residual, jacobian)
1132 }
1133
1134 fn get_dimension(&self) -> usize {
1135 1
1136 }
1137 }
1138
1139 /// Custom Rosenbrock Factor 2: r2 = 1 - x1
1140 /// Demonstrates extensibility - custom factors can be defined outside of factors.rs
1141 #[derive(Debug, Clone)]
1142 struct RosenbrockFactor2;
1143
1144 impl Factor for RosenbrockFactor2 {
1145 fn linearize(
1146 &self,
1147 params: &[DVector<f64>],
1148 compute_jacobian: bool,
1149 ) -> (DVector<f64>, Option<DMatrix<f64>>) {
1150 let x1 = params[0][0];
1151
1152 // Residual: r2 = 1 - x1
1153 let residual = dvector![1.0 - x1];
1154
1155 let jacobian = if compute_jacobian {
1156 // Jacobian: ∂r2/∂x1 = -1
1157 Some(DMatrix::from_element(1, 1, -1.0))
1158 } else {
1159 None
1160 };
1161
1162 (residual, jacobian)
1163 }
1164
1165 fn get_dimension(&self) -> usize {
1166 1
1167 }
1168 }
1169
1170 #[test]
1171 fn test_rosenbrock_optimization() -> TestResult {
1172 // Rosenbrock function test:
1173 // Minimize: r1² + r2² where
1174 // r1 = 10(x2 - x1²)
1175 // r2 = 1 - x1
1176 // Starting point: [-1.2, 1.0]
1177 // Expected minimum: [1.0, 1.0]
1178
1179 let mut problem = Problem::new(JacobianMode::Sparse);
1180 let mut initial_values = HashMap::new();
1181
1182 // Add variables using Rn manifold (Euclidean space)
1183 initial_values.insert("x1".to_string(), (ManifoldType::RN, dvector![-1.2]));
1184 initial_values.insert("x2".to_string(), (ManifoldType::RN, dvector![1.0]));
1185
1186 // Add custom factors (demonstrates extensibility!)
1187 problem.add_residual_block(&["x1", "x2"], Box::new(RosenbrockFactor1), None);
1188 problem.add_residual_block(&["x1"], Box::new(RosenbrockFactor2), None);
1189
1190 // Configure Levenberg-Marquardt optimizer
1191 let config = LevenbergMarquardtConfig::new()
1192 .with_max_iterations(100)
1193 .with_cost_tolerance(1e-8)
1194 .with_parameter_tolerance(1e-8)
1195 .with_gradient_tolerance(1e-10);
1196
1197 let mut solver = LevenbergMarquardt::with_config(config);
1198 let result = solver.optimize(&problem, &initial_values)?;
1199
1200 // Extract final values
1201 let x1_final = result
1202 .parameters
1203 .get("x1")
1204 .ok_or("x1 not found")?
1205 .to_vector()[0];
1206 let x2_final = result
1207 .parameters
1208 .get("x2")
1209 .ok_or("x2 not found")?
1210 .to_vector()[0];
1211
1212 // Verify convergence to [1.0, 1.0]
1213 assert!(
1214 matches!(
1215 result.status,
1216 OptimizationStatus::Converged
1217 | OptimizationStatus::CostToleranceReached
1218 | OptimizationStatus::ParameterToleranceReached
1219 | OptimizationStatus::GradientToleranceReached
1220 ),
1221 "Optimization should converge"
1222 );
1223 assert!(
1224 (x1_final - 1.0).abs() < 1e-4,
1225 "x1 should converge to 1.0, got {}",
1226 x1_final
1227 );
1228 assert!(
1229 (x2_final - 1.0).abs() < 1e-4,
1230 "x2 should converge to 1.0, got {}",
1231 x2_final
1232 );
1233 assert!(
1234 result.final_cost < 1e-6,
1235 "Final cost should be near zero, got {}",
1236 result.final_cost
1237 );
1238 Ok(())
1239 }
1240
1241 /// Trivial factor: r = x - target, J = [[1.0]]
1242 struct LinearFactor {
1243 target: f64,
1244 }
1245
1246 impl Factor for LinearFactor {
1247 fn linearize(
1248 &self,
1249 params: &[nalgebra::DVector<f64>],
1250 compute_jacobian: bool,
1251 ) -> (nalgebra::DVector<f64>, Option<nalgebra::DMatrix<f64>>) {
1252 let residual = dvector![params[0][0] - self.target];
1253 let jacobian = if compute_jacobian {
1254 Some(nalgebra::DMatrix::from_element(1, 1, 1.0))
1255 } else {
1256 None
1257 };
1258 (residual, jacobian)
1259 }
1260
1261 fn get_dimension(&self) -> usize {
1262 1
1263 }
1264 }
1265
1266 fn rosenbrock_problem() -> (
1267 Problem,
1268 HashMap<String, (apex_manifolds::ManifoldType, nalgebra::DVector<f64>)>,
1269 ) {
1270 let mut problem = Problem::new(JacobianMode::Sparse);
1271 let mut initial_values = HashMap::new();
1272 initial_values.insert("x1".to_string(), (ManifoldType::RN, dvector![-1.2]));
1273 initial_values.insert("x2".to_string(), (ManifoldType::RN, dvector![1.0]));
1274 problem.add_residual_block(&["x1", "x2"], Box::new(RosenbrockFactor1), None);
1275 problem.add_residual_block(&["x1"], Box::new(RosenbrockFactor2), None);
1276 (problem, initial_values)
1277 }
1278
1279 fn linear_problem(
1280 start: f64,
1281 ) -> (
1282 Problem,
1283 HashMap<String, (apex_manifolds::ManifoldType, nalgebra::DVector<f64>)>,
1284 ) {
1285 let mut problem = Problem::new(JacobianMode::Sparse);
1286 let mut initial_values = HashMap::new();
1287 initial_values.insert("x".to_string(), (ManifoldType::RN, dvector![start]));
1288 problem.add_residual_block(&["x"], Box::new(LinearFactor { target: 0.0 }), None);
1289 (problem, initial_values)
1290 }
1291
1292 // -------------------------------------------------------------------------
1293 // Config builder tests
1294 // -------------------------------------------------------------------------
1295
1296 #[test]
1297 fn test_lm_config_default() {
1298 let cfg = LevenbergMarquardtConfig::default();
1299 assert_eq!(cfg.max_iterations, 50);
1300 assert!((cfg.cost_tolerance - 1e-6).abs() < 1e-15);
1301 assert!((cfg.damping - 1e-3).abs() < 1e-15);
1302 assert!(!cfg.use_jacobi_scaling);
1303 assert!(!cfg.compute_covariances);
1304 }
1305
1306 #[test]
1307 fn test_lm_config_builders() {
1308 let cfg = LevenbergMarquardtConfig::new()
1309 .with_max_iterations(42)
1310 .with_cost_tolerance(1e-4)
1311 .with_parameter_tolerance(1e-5)
1312 .with_gradient_tolerance(1e-6)
1313 .with_damping(1e-2)
1314 .with_damping_bounds(1e-15, 1e15)
1315 .with_damping_factors(8.0, 0.2)
1316 .with_trust_region(500.0, 0.1, 0.8)
1317 .with_min_cost_threshold(1e-12)
1318 .with_min_trust_region_radius(1e-35)
1319 .with_jacobi_scaling(true)
1320 .with_compute_covariances(true)
1321 .with_linear_solver_type(LinearSolverType::SparseQR);
1322 assert_eq!(cfg.max_iterations, 42);
1323 assert!((cfg.cost_tolerance - 1e-4).abs() < 1e-20);
1324 assert!((cfg.parameter_tolerance - 1e-5).abs() < 1e-20);
1325 assert!((cfg.gradient_tolerance - 1e-6).abs() < 1e-20);
1326 assert!((cfg.damping - 1e-2).abs() < 1e-15);
1327 assert!((cfg.damping_min - 1e-15).abs() < 1e-25);
1328 assert!((cfg.damping_max - 1e15).abs() < 1.0);
1329 assert!((cfg.damping_increase_factor - 8.0).abs() < 1e-12);
1330 assert!((cfg.damping_decrease_factor - 0.2).abs() < 1e-12);
1331 assert!((cfg.trust_region_radius - 500.0).abs() < 1e-10);
1332 assert!(cfg.min_cost_threshold.is_some());
1333 assert!(cfg.use_jacobi_scaling);
1334 assert!(cfg.compute_covariances);
1335 assert!(matches!(cfg.linear_solver_type, LinearSolverType::SparseQR));
1336 }
1337
1338 #[test]
1339 fn test_lm_for_bundle_adjustment() {
1340 let cfg = LevenbergMarquardtConfig::for_bundle_adjustment();
1341 assert!(matches!(
1342 cfg.linear_solver_type,
1343 LinearSolverType::SparseSchurComplement
1344 ));
1345 assert_eq!(cfg.max_iterations, 20);
1346 }
1347
1348 #[test]
1349 fn test_lm_print_configuration_no_panic() {
1350 LevenbergMarquardtConfig::default().print_configuration();
1351 }
1352
1353 #[test]
1354 fn test_lm_default_equals_new() {
1355 let a = LevenbergMarquardt::new();
1356 let b = LevenbergMarquardt::default();
1357 // Both should solve the same problem identically (smoke check)
1358 drop(a);
1359 drop(b);
1360 }
1361
1362 #[test]
1363 fn test_lm_with_config_method() {
1364 let cfg = LevenbergMarquardtConfig::new().with_max_iterations(7);
1365 let solver = LevenbergMarquardt::with_config(cfg);
1366 drop(solver);
1367 }
1368
1369 // -------------------------------------------------------------------------
1370 // Convergence termination paths
1371 // -------------------------------------------------------------------------
1372
1373 #[test]
1374 fn test_lm_max_iterations_termination() -> TestResult {
1375 let (problem, initial_values) = rosenbrock_problem();
1376 let cfg = LevenbergMarquardtConfig::new().with_max_iterations(2);
1377 let mut solver = LevenbergMarquardt::with_config(cfg);
1378 let result = solver.optimize(&problem, &initial_values)?;
1379 assert_eq!(result.status, OptimizationStatus::MaxIterationsReached);
1380 assert!(result.iterations <= 3, "iterations={}", result.iterations);
1381 Ok(())
1382 }
1383
1384 #[test]
1385 fn test_lm_gradient_tolerance_convergence() -> TestResult {
1386 let (problem, initial_values) = linear_problem(1.0);
1387 // Very loose gradient tolerance → triggers after first accepted step
1388 let cfg = LevenbergMarquardtConfig::new()
1389 .with_gradient_tolerance(1e3)
1390 .with_cost_tolerance(1e-20)
1391 .with_parameter_tolerance(1e-20);
1392 let mut solver = LevenbergMarquardt::with_config(cfg);
1393 let result = solver.optimize(&problem, &initial_values)?;
1394 assert_eq!(result.status, OptimizationStatus::GradientToleranceReached);
1395 Ok(())
1396 }
1397
1398 #[test]
1399 fn test_lm_min_cost_threshold() -> TestResult {
1400 let (problem, initial_values) = rosenbrock_problem();
1401 // Set threshold very high so even initial cost triggers it
1402 let cfg = LevenbergMarquardtConfig::new()
1403 .with_min_cost_threshold(1e10)
1404 .with_cost_tolerance(1e-20)
1405 .with_gradient_tolerance(1e-20)
1406 .with_parameter_tolerance(1e-20);
1407 let mut solver = LevenbergMarquardt::with_config(cfg);
1408 let result = solver.optimize(&problem, &initial_values)?;
1409 assert_eq!(result.status, OptimizationStatus::MinCostThresholdReached);
1410 Ok(())
1411 }
1412
1413 #[test]
1414 fn test_lm_qr_solver() -> TestResult {
1415 let (problem, initial_values) = rosenbrock_problem();
1416 let cfg = LevenbergMarquardtConfig::new()
1417 .with_linear_solver_type(LinearSolverType::SparseQR)
1418 .with_max_iterations(100);
1419 let mut solver = LevenbergMarquardt::with_config(cfg);
1420 let result = solver.optimize(&problem, &initial_values)?;
1421 assert!(result.final_cost < 1e-6);
1422 Ok(())
1423 }
1424
1425 #[test]
1426 fn test_lm_jacobi_scaling_enabled() -> TestResult {
1427 let (problem, initial_values) = rosenbrock_problem();
1428 let cfg = LevenbergMarquardtConfig::new()
1429 .with_jacobi_scaling(true)
1430 .with_max_iterations(100);
1431 let mut solver = LevenbergMarquardt::with_config(cfg);
1432 let result = solver.optimize(&problem, &initial_values)?;
1433 assert!(result.final_cost < 1e-6);
1434 Ok(())
1435 }
1436
1437 #[test]
1438 fn test_lm_result_initial_cost_greater_than_final() -> TestResult {
1439 let (problem, initial_values) = rosenbrock_problem();
1440 let mut solver = LevenbergMarquardt::new();
1441 let result = solver.optimize(&problem, &initial_values)?;
1442 assert!(
1443 result.initial_cost > result.final_cost,
1444 "initial={} final={}",
1445 result.initial_cost,
1446 result.final_cost
1447 );
1448 Ok(())
1449 }
1450
1451 #[test]
1452 fn test_lm_convergence_info_populated() -> TestResult {
1453 let (problem, initial_values) = rosenbrock_problem();
1454 let mut solver = LevenbergMarquardt::new();
1455 let result = solver.optimize(&problem, &initial_values)?;
1456 assert!(result.convergence_info.is_some());
1457 Ok(())
1458 }
1459
1460 #[test]
1461 fn test_lm_iterations_positive() -> TestResult {
1462 let (problem, initial_values) = rosenbrock_problem();
1463 let mut solver = LevenbergMarquardt::new();
1464 let result = solver.optimize(&problem, &initial_values)?;
1465 assert!(result.iterations > 0);
1466 Ok(())
1467 }
1468
1469 #[test]
1470 fn test_lm_timeout_config() {
1471 let cfg = LevenbergMarquardtConfig::new().with_timeout(Duration::from_secs(30));
1472 assert!(cfg.timeout.is_some());
1473 }
1474
1475 #[test]
1476 fn test_lm_config_schur_variant_and_preconditioner() {
1477 use crate::linalg::{SchurPreconditioner, SchurVariant};
1478 let cfg = LevenbergMarquardtConfig::new()
1479 .with_schur_variant(SchurVariant::Iterative)
1480 .with_schur_preconditioner(SchurPreconditioner::BlockDiagonal);
1481 assert!(matches!(cfg.schur_variant, SchurVariant::Iterative));
1482 assert!(matches!(
1483 cfg.schur_preconditioner,
1484 SchurPreconditioner::BlockDiagonal
1485 ));
1486 }
1487
1488 // -------------------------------------------------------------------------
1489 // Dense Jacobian mode dispatch
1490 // -------------------------------------------------------------------------
1491
1492 /// Exercises the `JacobianMode::Dense + _ => DenseCholeskySolver` arm of `optimize()`.
1493 /// All existing tests use `JacobianMode::Sparse`, so this branch was previously uncovered.
1494 #[test]
1495 fn test_lm_dense_cholesky_solver() -> TestResult {
1496 let mut problem = Problem::new(JacobianMode::Dense);
1497 let mut initial_values = HashMap::new();
1498 initial_values.insert("x1".to_string(), (ManifoldType::RN, dvector![-1.2]));
1499 initial_values.insert("x2".to_string(), (ManifoldType::RN, dvector![1.0]));
1500 problem.add_residual_block(&["x1", "x2"], Box::new(RosenbrockFactor1), None);
1501 problem.add_residual_block(&["x1"], Box::new(RosenbrockFactor2), None);
1502
1503 // Default linear solver type (SparseCholesky) with Dense mode → DenseCholeskySolver
1504 let cfg = LevenbergMarquardtConfig::new().with_max_iterations(100);
1505 let mut solver = LevenbergMarquardt::with_config(cfg);
1506 let result = solver.optimize(&problem, &initial_values)?;
1507 assert!(
1508 result.final_cost < 1e-6,
1509 "Dense Cholesky mode should converge Rosenbrock, got cost={}",
1510 result.final_cost
1511 );
1512 Ok(())
1513 }
1514
1515 /// Exercises the `JacobianMode::Dense + DenseQR` arm of `optimize()`.
1516 #[test]
1517 fn test_lm_dense_qr_solver() -> TestResult {
1518 let mut problem = Problem::new(JacobianMode::Dense);
1519 let mut initial_values = HashMap::new();
1520 initial_values.insert("x1".to_string(), (ManifoldType::RN, dvector![-1.2]));
1521 initial_values.insert("x2".to_string(), (ManifoldType::RN, dvector![1.0]));
1522 problem.add_residual_block(&["x1", "x2"], Box::new(RosenbrockFactor1), None);
1523 problem.add_residual_block(&["x1"], Box::new(RosenbrockFactor2), None);
1524
1525 let cfg = LevenbergMarquardtConfig::new()
1526 .with_linear_solver_type(LinearSolverType::DenseQR)
1527 .with_max_iterations(100);
1528 let mut solver = LevenbergMarquardt::with_config(cfg);
1529 let result = solver.optimize(&problem, &initial_values)?;
1530 assert!(
1531 result.final_cost < 1e-6,
1532 "Dense QR mode should converge Rosenbrock, got cost={}",
1533 result.final_cost
1534 );
1535 Ok(())
1536 }
1537
1538 // -------------------------------------------------------------------------
1539 // Covariance computation
1540 // -------------------------------------------------------------------------
1541
1542 /// Exercises the `if self.config.compute_covariances { ... }` block at convergence.
1543 /// This block was completely unreachable in prior tests.
1544 #[test]
1545 fn test_lm_compute_covariances_enabled() -> TestResult {
1546 let (problem, initial_values) = rosenbrock_problem();
1547 let cfg = LevenbergMarquardtConfig::new()
1548 .with_max_iterations(100)
1549 .with_compute_covariances(true);
1550 let mut solver = LevenbergMarquardt::with_config(cfg);
1551 let result = solver.optimize(&problem, &initial_values)?;
1552 // result.covariances is Option<HashMap<String, Mat<f64>>>
1553 assert!(
1554 result.covariances.is_some(),
1555 "compute_covariances=true should populate result.covariances"
1556 );
1557 Ok(())
1558 }
1559
1560 // -------------------------------------------------------------------------
1561 // update_damping() direct unit tests
1562 // -------------------------------------------------------------------------
1563
1564 /// `update_damping(rho > 0)` should accept the step, decrease damping, and reset nu.
1565 #[test]
1566 fn test_update_damping_accepted_step() {
1567 let cfg = LevenbergMarquardtConfig::new()
1568 .with_damping(1e-2)
1569 .with_damping_bounds(1e-15, 1e15);
1570 let mut solver = LevenbergMarquardt::with_config(cfg);
1571 let initial_damping = solver.config.damping;
1572
1573 // rho = 0.8 > 0 → accepted branch
1574 let accepted = solver.update_damping(0.8);
1575
1576 assert!(accepted, "rho > 0 should return true (step accepted)");
1577 assert!(
1578 solver.config.damping < initial_damping,
1579 "accepted step should decrease damping: {} < {}",
1580 solver.config.damping,
1581 initial_damping
1582 );
1583 // damping_nu should be reset to 2.0 on acceptance
1584 assert!(
1585 (solver.config.damping_nu - 2.0).abs() < 1e-15,
1586 "damping_nu should be reset to 2.0 after accepted step, got {}",
1587 solver.config.damping_nu
1588 );
1589 }
1590
1591 /// `update_damping(rho <= 0)` should reject the step, increase damping, and double nu.
1592 #[test]
1593 fn test_update_damping_rejected_step() {
1594 let cfg = LevenbergMarquardtConfig::new()
1595 .with_damping(1e-2)
1596 .with_damping_bounds(1e-15, 1e15);
1597 let initial_nu = cfg.damping_nu; // default 2.0
1598 let mut solver = LevenbergMarquardt::with_config(cfg);
1599 let initial_damping = solver.config.damping;
1600
1601 // rho = -0.5 <= 0 → rejected branch
1602 let rejected = solver.update_damping(-0.5);
1603
1604 assert!(!rejected, "rho <= 0 should return false (step rejected)");
1605 assert!(
1606 solver.config.damping > initial_damping,
1607 "rejected step should increase damping: {} > {}",
1608 solver.config.damping,
1609 initial_damping
1610 );
1611 // damping_nu doubles on rejection
1612 assert!(
1613 (solver.config.damping_nu - initial_nu * 2.0).abs() < 1e-15,
1614 "damping_nu should double on rejected step: expected {}, got {}",
1615 initial_nu * 2.0,
1616 solver.config.damping_nu
1617 );
1618 }
1619
1620 // -------------------------------------------------------------------------
1621 // Untested config builder methods
1622 // -------------------------------------------------------------------------
1623
1624 /// Verifies `with_max_condition_number` and `with_min_relative_decrease` builder methods.
1625 #[test]
1626 fn test_lm_config_condition_number_and_relative_decrease() -> TestResult {
1627 let cfg = LevenbergMarquardtConfig::new()
1628 .with_max_condition_number(1e8)
1629 .with_min_relative_decrease(1e-4);
1630 let max_cond = cfg
1631 .max_condition_number
1632 .ok_or("max_condition_number should be Some")?;
1633 assert!((max_cond - 1e8).abs() < 1.0);
1634 assert!((cfg.min_relative_decrease - 1e-4).abs() < 1e-20);
1635 Ok(())
1636 }
1637
1638 // -------------------------------------------------------------------------
1639 // Observer integration
1640 // -------------------------------------------------------------------------
1641
1642 /// Verifies that `add_observer` registers an observer and `notify_complete` is called
1643 /// exactly once after optimization finishes.
1644 #[test]
1645 fn test_lm_add_observer_called_on_completion() -> TestResult {
1646 use crate::optimizer::OptObserver;
1647 use std::sync::{Arc, Mutex};
1648
1649 struct CountObserver {
1650 complete_calls: Arc<Mutex<usize>>,
1651 }
1652
1653 impl OptObserver for CountObserver {
1654 fn on_step(&self, _values: &HashMap<String, VariableEnum>, _iteration: usize) {}
1655
1656 fn on_optimization_complete(
1657 &self,
1658 _values: &HashMap<String, VariableEnum>,
1659 _iterations: usize,
1660 ) {
1661 if let Ok(mut guard) = self.complete_calls.lock() {
1662 *guard += 1;
1663 }
1664 }
1665 }
1666
1667 let call_count = Arc::new(Mutex::new(0usize));
1668 let observer = CountObserver {
1669 complete_calls: Arc::clone(&call_count),
1670 };
1671
1672 let (problem, initial_values) = rosenbrock_problem();
1673 let mut solver = LevenbergMarquardt::new();
1674 solver.add_observer(observer);
1675 let _ = solver.optimize(&problem, &initial_values)?;
1676
1677 let count = *call_count
1678 .lock()
1679 .map_err(|e| format!("mutex poisoned: {e}"))?;
1680 assert_eq!(
1681 count, 1,
1682 "on_optimization_complete should be called exactly once"
1683 );
1684 Ok(())
1685 }
1686}