Skip to main content

apex_solver/optimizer/
dog_leg.rs

1//! Dog Leg trust region optimization algorithm implementation.
2//!
3//! The Dog Leg method is a robust trust region algorithm for solving nonlinear least squares problems:
4//!
5//! ```text
6//! min f(x) = ½||r(x)||² = ½Σᵢ rᵢ(x)²
7//! ```
8//!
9//! where `r: ℝⁿ → ℝᵐ` is the residual vector function.
10//!
11//! # Algorithm Overview
12//!
13//! Powell's Dog Leg method constructs a piecewise linear path within a spherical trust region
14//! of radius Δ, connecting three key points:
15//!
16//! 1. **Origin** (current position)
17//! 2. **Cauchy Point** p_c = -α·g (optimal steepest descent step)
18//! 3. **Gauss-Newton Point** h_gn (full Newton step solving J^T·J·h = -J^T·r)
19//!
20//! The "dog leg" path travels from the origin to the Cauchy point, then from the Cauchy point
21//! toward the Gauss-Newton point, stopping at the trust region boundary.
22//!
23//! ## Step Selection Strategy
24//!
25//! Given trust region radius Δ, the algorithm selects a step h based on three cases:
26//!
27//! **Case 1: GN step inside trust region** (`||h_gn|| ≤ Δ`)
28//! ```text
29//! h = h_gn  (full Gauss-Newton step)
30//! ```
31//!
32//! **Case 2: Even Cauchy point outside trust region** (`||p_c|| ≥ Δ`)
33//! ```text
34//! h = (Δ / ||g||) · (-g)  (scaled steepest descent to boundary)
35//! ```
36//!
37//! **Case 3: Dog leg interpolation** (`||p_c|| < Δ < ||h_gn||`)
38//! ```text
39//! h(β) = p_c + β·(h_gn - p_c),  where β ∈ [0,1] satisfies ||h(β)|| = Δ
40//! ```
41//!
42//! ## Cauchy Point Computation
43//!
44//! The Cauchy point is the optimal step along the steepest descent direction -g:
45//!
46//! ```text
47//! α = (g^T·g) / (g^T·H·g)  where H = J^T·J
48//! p_c = -α·g
49//! ```
50//!
51//! This minimizes the quadratic model along the gradient direction.
52//!
53//! ## Trust Region Management
54//!
55//! The trust region radius Δ adapts based on the gain ratio:
56//!
57//! ```text
58//! ρ = (actual reduction) / (predicted reduction)
59//! ```
60//!
61//! **Good step** (`ρ > 0.75`): Increase radius `Δ ← max(Δ, 3·||h||)`
62//! **Poor step** (`ρ < 0.25`): Decrease radius `Δ ← Δ/2`
63//! **Moderate step**: Keep radius unchanged
64//!
65//! # Advanced Features (Ceres Solver Enhancements)
66//!
67//! This implementation includes several improvements from Google's Ceres Solver:
68//!
69//! ## 1. Adaptive μ Regularization
70//!
71//! Instead of solving `J^T·J·h = -J^T·r` directly, we solve:
72//!
73//! ```text
74//! (J^T·J + μI)·h = -J^T·r
75//! ```
76//!
77//! where μ adapts to handle ill-conditioned Hessians:
78//! - **Increases** (×10) when linear solve fails
79//! - **Decreases** (÷5) when steps are accepted
80//! - **Bounded** between `min_mu` (1e-8) and `max_mu` (1.0)
81//!
82//! Default `initial_mu = 1e-4` provides good numerical stability.
83//!
84//! ## 2. Numerically Robust Beta Computation
85//!
86//! When computing β for dog leg interpolation, solves: `||p_c + β·v||² = Δ²`
87//!
88//! Uses two formulas to avoid catastrophic cancellation:
89//! ```text
90//! If b ≤ 0:  β = (-b + √(b²-ac)) / a    (standard formula)
91//! If b > 0:  β = -c / (b + √(b²-ac))    (alternative, avoids cancellation)
92//! ```
93//!
94//! ## 3. Step Reuse Mechanism
95//!
96//! When a step is rejected, caches the Gauss-Newton step, Cauchy point, and gradient.
97//! On the next iteration (with smaller Δ), reuses these cached values instead of
98//! recomputing them. This avoids expensive linear solves when trust region shrinks.
99//!
100//! **Safety limits:**
101//! - Maximum 5 consecutive reuses before forcing fresh computation
102//! - Cache invalidated when steps are accepted (parameters have moved)
103//!
104//! ## 4. Jacobi Scaling (Diagonal Preconditioning)
105//!
106//! Optionally applies column scaling to J before forming J^T·J:
107//! ```text
108//! S_ii = 1 / (1 + ||J_i||)  where J_i is column i
109//! ```
110//!
111//! This creates an **elliptical trust region** instead of spherical, improving
112//! convergence for problems with mixed parameter scales.
113//!
114//! # Mathematical Background
115//!
116//! ## Why Dog Leg Works
117//!
118//! The dog leg path is a cheap approximation to the optimal trust region step
119//! (which would require solving a constrained optimization problem). It:
120//!
121//! 1. Exploits the fact that optimal steps often lie near the 2D subspace
122//!    spanned by the gradient and Gauss-Newton directions
123//! 2. Provides global convergence guarantees (always finds descent direction)
124//! 3. Achieves local quadratic convergence (like Gauss-Newton near solution)
125//! 4. Requires only one linear solve per iteration (same as Gauss-Newton)
126//!
127//! ## Convergence Properties
128//!
129//! - **Global convergence**: Guaranteed descent at each iteration
130//! - **Local quadratic convergence**: Reduces to Gauss-Newton near solution
131//! - **Robustness**: Handles ill-conditioning via trust region + μ regularization
132//! - **Efficiency**: Comparable cost to Gauss-Newton with better reliability
133//!
134//! # When to Use
135//!
136//! Dog Leg is an excellent choice when:
137//! - You want explicit control over step size (via trust region radius)
138//! - The problem may be ill-conditioned
139//! - You need guaranteed descent at each iteration
140//! - Initial guess may be poor but you want reliable convergence
141//!
142//! Compared to alternatives:
143//! - **vs Gauss-Newton**: More robust but similar computational cost
144//! - **vs Levenberg-Marquardt**: Explicit trust region vs implicit damping
145//! - Both Dog Leg and LM are excellent general-purpose choices
146//!
147//! # Examples
148//!
149//! ## Basic usage
150//!
151//! ```no_run
152//! use apex_solver::optimizer::DogLeg;
153//! use apex_solver::core::problem::Problem;
154//! use apex_solver::JacobianMode;
155//! use std::collections::HashMap;
156//!
157//! # type TestResult = Result<(), Box<dyn std::error::Error>>;
158//! # fn main() -> TestResult {
159//! let mut problem = Problem::new(JacobianMode::Sparse);
160//! // ... add residual blocks (factors) to problem ...
161//!
162//! let initial_values = HashMap::new();
163//! // ... initialize parameters ...
164//!
165//! let mut solver = DogLeg::new();
166//! let result = solver.optimize(&problem, &initial_values)?;
167//! # Ok(())
168//! # }
169//! ```
170//!
171//! ## Advanced configuration with Ceres enhancements
172//!
173//! ```no_run
174//! use apex_solver::optimizer::dog_leg::{DogLegConfig, DogLeg};
175//! use apex_solver::linalg::LinearSolverType;
176//!
177//! # fn main() {
178//! let config = DogLegConfig::new()
179//!     .with_max_iterations(100)
180//!     .with_trust_region_radius(1e4)  // Large initial radius
181//!     .with_trust_region_bounds(1e-3, 1e6)  // Min/max radius
182//!     .with_mu_params(1e-4, 1e-8, 1.0, 10.0)  // Conservative regularization
183//!     .with_jacobi_scaling(true)  // Enable elliptical trust regions
184//!     .with_step_reuse(true);  // Enable Ceres-style caching
185//!
186//! let mut solver = DogLeg::with_config(config);
187//! # }
188//! ```
189//!
190//! # References
191//!
192//! - Powell, M. J. D. (1970). "A Hybrid Method for Nonlinear Equations". *Numerical Methods for Nonlinear Algebraic Equations*. Gordon and Breach.
193//! - Nocedal, J. & Wright, S. (2006). *Numerical Optimization* (2nd ed.). Springer. Chapter 4 (Trust Region Methods).
194//! - Madsen, K., Nielsen, H. B., & Tingleff, O. (2004). *Methods for Non-Linear Least Squares Problems* (2nd ed.). Chapter 6.
195//! - Conn, A. R., Gould, N. I. M., & Toint, P. L. (2000). *Trust-Region Methods*. SIAM.
196//! - Ceres Solver: http://ceres-solver.org/ - Google's C++ nonlinear least squares library
197
198use crate::{core::problem, error, linalg, optimizer};
199use apex_manifolds as manifold;
200
201use std::{collections, fmt, time};
202use tracing::debug;
203
204use crate::linalg::{
205    DenseCholeskySolver, DenseMode, DenseQRSolver, JacobianMode, LinearSolver, LinearSolverType,
206    SparseCholeskySolver, SparseMode, SparseQRSolver,
207};
208use crate::optimizer::{AssemblyBackend, IterationStats};
209
210/// Configuration parameters for the Dog Leg trust region optimizer.
211///
212/// Controls trust region management, convergence criteria, adaptive regularization,
213/// and Ceres Solver enhancements for the Dog Leg algorithm.
214///
215/// # Builder Pattern
216///
217/// All configuration options can be set using the builder pattern:
218///
219/// ```
220/// use apex_solver::optimizer::dog_leg::DogLegConfig;
221///
222/// let config = DogLegConfig::new()
223///     .with_max_iterations(100)
224///     .with_trust_region_radius(1e4)
225///     .with_mu_params(1e-4, 1e-8, 1.0, 10.0)
226///     .with_jacobi_scaling(true)
227///     .with_step_reuse(true);
228/// ```
229///
230/// # Trust Region Behavior
231///
232/// The trust region radius Δ controls the maximum allowed step size:
233///
234/// - **Initial radius** (`trust_region_radius`): Starting value (default: 1e4)
235/// - **Bounds** (`trust_region_min`, `trust_region_max`): Valid range (default: 1e-3 to 1e6)
236/// - **Adaptation**: Increases for good steps, decreases for poor steps
237///
238/// # Adaptive μ Regularization (Ceres Enhancement)
239///
240/// Controls the regularization parameter in `(J^T·J + μI)·h = -J^T·r`:
241///
242/// - `initial_mu`: Starting value (default: 1e-4 for numerical stability)
243/// - `min_mu`, `max_mu`: Bounds (default: 1e-8 to 1.0)
244/// - `mu_increase_factor`: Multiplier when solve fails (default: 10.0)
245///
246/// # Convergence Criteria
247///
248/// The optimizer terminates when ANY of the following conditions is met:
249///
250/// - **Cost tolerance**: `|cost_k - cost_{k-1}| < cost_tolerance`
251/// - **Parameter tolerance**: `||step|| < parameter_tolerance`
252/// - **Gradient tolerance**: `||J^T·r|| < gradient_tolerance`
253/// - **Maximum iterations**: `iteration >= max_iterations`
254/// - **Timeout**: `elapsed_time >= timeout`
255///
256/// # See Also
257///
258/// - [`DogLeg`] - The solver that uses this configuration
259/// - [`LevenbergMarquardtConfig`](crate::optimizer::LevenbergMarquardtConfig) - Alternative damping approach
260/// - [`GaussNewtonConfig`](crate::optimizer::GaussNewtonConfig) - Undamped variant
261#[derive(Clone)]
262pub struct DogLegConfig {
263    /// Type of linear solver for the linear systems
264    pub linear_solver_type: linalg::LinearSolverType,
265    /// Maximum number of iterations
266    pub max_iterations: usize,
267    /// Convergence tolerance for cost function
268    pub cost_tolerance: f64,
269    /// Convergence tolerance for parameter updates
270    pub parameter_tolerance: f64,
271    /// Convergence tolerance for gradient norm
272    pub gradient_tolerance: f64,
273    /// Timeout duration
274    pub timeout: Option<time::Duration>,
275    /// Initial trust region radius
276    pub trust_region_radius: f64,
277    /// Minimum trust region radius
278    pub trust_region_min: f64,
279    /// Maximum trust region radius
280    pub trust_region_max: f64,
281    /// Trust region increase factor (for good steps, rho > 0.75)
282    pub trust_region_increase_factor: f64,
283    /// Trust region decrease factor (for poor steps, rho < 0.25)
284    pub trust_region_decrease_factor: f64,
285    /// Minimum step quality for acceptance (typically 0.0)
286    pub min_step_quality: f64,
287    /// Good step quality threshold (typically 0.75)
288    pub good_step_quality: f64,
289    /// Poor step quality threshold (typically 0.25)
290    pub poor_step_quality: f64,
291    /// Use Jacobi column scaling (preconditioning)
292    pub use_jacobi_scaling: bool,
293
294    // Ceres-style adaptive mu regularization parameters
295    /// Initial mu regularization parameter for Gauss-Newton step
296    pub initial_mu: f64,
297    /// Minimum mu regularization parameter
298    pub min_mu: f64,
299    /// Maximum mu regularization parameter
300    pub max_mu: f64,
301    /// Factor to increase mu when linear solver fails
302    pub mu_increase_factor: f64,
303
304    // Ceres-style step reuse optimization
305    /// Enable step reuse after rejection (Ceres-style efficiency optimization)
306    pub enable_step_reuse: bool,
307
308    /// Minimum objective function cutoff (optional early termination)
309    ///
310    /// If set, optimization terminates when cost falls below this threshold.
311    /// Useful for early stopping when a "good enough" solution is acceptable.
312    ///
313    /// Default: None (disabled)
314    pub min_cost_threshold: Option<f64>,
315
316    /// Maximum condition number for Jacobian matrix (optional check)
317    ///
318    /// If set, the optimizer checks if condition_number(J^T*J) exceeds this
319    /// threshold and terminates with IllConditionedJacobian status.
320    /// Note: Computing condition number is expensive, so this is disabled by default.
321    ///
322    /// Default: None (disabled)
323    pub max_condition_number: Option<f64>,
324
325    /// Minimum relative cost decrease for step acceptance
326    ///
327    /// Used in computing step quality (rho = actual_reduction / predicted_reduction).
328    /// Steps with rho < min_relative_decrease are rejected. Matches Ceres Solver's
329    /// min_relative_decrease parameter.
330    ///
331    /// Default: 1e-3 (Ceres-compatible)
332    pub min_relative_decrease: f64,
333
334    /// Compute per-variable covariance matrices (uncertainty estimation)
335    ///
336    /// When enabled, computes covariance by inverting the Hessian matrix after
337    /// convergence. The full covariance matrix is extracted into per-variable
338    /// blocks stored in both Variable structs and optimizer::SolverResult.
339    ///
340    /// Default: false (to avoid performance overhead)
341    pub compute_covariances: bool,
342
343    /// Enable real-time visualization (graphical debugging).
344    ///
345    /// When enabled, optimization progress is logged to a Rerun viewer.
346    /// **Note:** Requires the `visualization` feature to be enabled in `Cargo.toml`.
347    ///
348    /// Default: false
349    #[cfg(feature = "visualization")]
350    pub enable_visualization: bool,
351}
352
353impl Default for DogLegConfig {
354    fn default() -> Self {
355        Self {
356            linear_solver_type: linalg::LinearSolverType::default(),
357            // Ceres Solver default: 50 (changed from 100 for compatibility)
358            max_iterations: 50,
359            // Ceres Solver default: 1e-6 (changed from 1e-8 for compatibility)
360            cost_tolerance: 1e-6,
361            // Ceres Solver default: 1e-8 (unchanged)
362            parameter_tolerance: 1e-8,
363            // Ceres Solver default: 1e-10 (changed from 1e-8 for compatibility)
364            gradient_tolerance: 1e-10,
365            timeout: None,
366            // Ceres-style: larger initial radius for better global convergence
367            trust_region_radius: 1e4,
368            trust_region_min: 1e-12,
369            trust_region_max: 1e12,
370            // Ceres uses adaptive increase (max(radius, 3*step_norm)),
371            // but we keep factor for simpler config
372            trust_region_increase_factor: 3.0,
373            trust_region_decrease_factor: 0.5,
374            min_step_quality: 0.0,
375            good_step_quality: 0.75,
376            poor_step_quality: 0.25,
377            // Ceres-style: Enable diagonal scaling by default for elliptical trust region
378            use_jacobi_scaling: true,
379
380            // Ceres-style adaptive mu regularization defaults
381            // Start with more conservative regularization to avoid singular Hessian
382            initial_mu: 1e-4,
383            min_mu: 1e-8,
384            max_mu: 1.0,
385            mu_increase_factor: 10.0,
386
387            // Ceres-style step reuse optimization
388            enable_step_reuse: true,
389
390            // New Ceres-compatible termination parameters
391            min_cost_threshold: None,
392            max_condition_number: None,
393            min_relative_decrease: 1e-3,
394
395            compute_covariances: false,
396            #[cfg(feature = "visualization")]
397            enable_visualization: false,
398        }
399    }
400}
401
402impl DogLegConfig {
403    /// Create a new Dog Leg configuration with default values.
404    pub fn new() -> Self {
405        Self::default()
406    }
407
408    /// Set the linear solver type
409    pub fn with_linear_solver_type(mut self, linear_solver_type: linalg::LinearSolverType) -> Self {
410        self.linear_solver_type = linear_solver_type;
411        self
412    }
413
414    /// Set the maximum number of iterations
415    pub fn with_max_iterations(mut self, max_iterations: usize) -> Self {
416        self.max_iterations = max_iterations;
417        self
418    }
419
420    /// Set the cost tolerance
421    pub fn with_cost_tolerance(mut self, cost_tolerance: f64) -> Self {
422        self.cost_tolerance = cost_tolerance;
423        self
424    }
425
426    /// Set the parameter tolerance
427    pub fn with_parameter_tolerance(mut self, parameter_tolerance: f64) -> Self {
428        self.parameter_tolerance = parameter_tolerance;
429        self
430    }
431
432    /// Set the gradient tolerance
433    pub fn with_gradient_tolerance(mut self, gradient_tolerance: f64) -> Self {
434        self.gradient_tolerance = gradient_tolerance;
435        self
436    }
437
438    /// Set the timeout duration
439    pub fn with_timeout(mut self, timeout: time::Duration) -> Self {
440        self.timeout = Some(timeout);
441        self
442    }
443
444    /// Set the initial trust region radius
445    pub fn with_trust_region_radius(mut self, radius: f64) -> Self {
446        self.trust_region_radius = radius;
447        self
448    }
449
450    /// Set the trust region radius bounds
451    pub fn with_trust_region_bounds(mut self, min: f64, max: f64) -> Self {
452        self.trust_region_min = min;
453        self.trust_region_max = max;
454        self
455    }
456
457    /// Set the trust region adjustment factors
458    pub fn with_trust_region_factors(mut self, increase: f64, decrease: f64) -> Self {
459        self.trust_region_increase_factor = increase;
460        self.trust_region_decrease_factor = decrease;
461        self
462    }
463
464    /// Set the trust region quality thresholds
465    pub fn with_step_quality_thresholds(
466        mut self,
467        min_quality: f64,
468        poor_quality: f64,
469        good_quality: f64,
470    ) -> Self {
471        self.min_step_quality = min_quality;
472        self.poor_step_quality = poor_quality;
473        self.good_step_quality = good_quality;
474        self
475    }
476
477    /// Enable or disable Jacobi column scaling (preconditioning)
478    pub fn with_jacobi_scaling(mut self, use_jacobi_scaling: bool) -> Self {
479        self.use_jacobi_scaling = use_jacobi_scaling;
480        self
481    }
482
483    /// Set adaptive mu regularization parameters (Ceres-style)
484    pub fn with_mu_params(
485        mut self,
486        initial_mu: f64,
487        min_mu: f64,
488        max_mu: f64,
489        increase_factor: f64,
490    ) -> Self {
491        self.initial_mu = initial_mu;
492        self.min_mu = min_mu;
493        self.max_mu = max_mu;
494        self.mu_increase_factor = increase_factor;
495        self
496    }
497
498    /// Enable or disable step reuse optimization (Ceres-style)
499    pub fn with_step_reuse(mut self, enable_step_reuse: bool) -> Self {
500        self.enable_step_reuse = enable_step_reuse;
501        self
502    }
503
504    /// Set minimum objective function cutoff for early termination.
505    ///
506    /// When set, optimization terminates with MinCostThresholdReached status
507    /// if the cost falls below this threshold. Useful for early stopping when
508    /// a "good enough" solution is acceptable.
509    pub fn with_min_cost_threshold(mut self, min_cost: f64) -> Self {
510        self.min_cost_threshold = Some(min_cost);
511        self
512    }
513
514    /// Set maximum condition number for Jacobian matrix.
515    ///
516    /// If set, the optimizer checks if condition_number(J^T*J) exceeds this
517    /// threshold and terminates with IllConditionedJacobian status.
518    /// Note: Computing condition number is expensive, disabled by default.
519    pub fn with_max_condition_number(mut self, max_cond: f64) -> Self {
520        self.max_condition_number = Some(max_cond);
521        self
522    }
523
524    /// Set minimum relative cost decrease for step acceptance.
525    ///
526    /// Steps with rho = (actual_reduction / predicted_reduction) below this
527    /// threshold are rejected. Default: 1e-3 (Ceres-compatible)
528    pub fn with_min_relative_decrease(mut self, min_decrease: f64) -> Self {
529        self.min_relative_decrease = min_decrease;
530        self
531    }
532
533    /// Enable or disable covariance computation (uncertainty estimation).
534    ///
535    /// When enabled, computes the full covariance matrix by inverting the Hessian
536    /// after convergence, then extracts per-variable covariance blocks.
537    pub fn with_compute_covariances(mut self, compute_covariances: bool) -> Self {
538        self.compute_covariances = compute_covariances;
539        self
540    }
541
542    /// Enable real-time visualization.
543    ///
544    /// **Note:** Requires the `visualization` feature to be enabled in `Cargo.toml`.
545    ///
546    /// # Arguments
547    ///
548    /// * `enable` - Whether to enable visualization
549    #[cfg(feature = "visualization")]
550    pub fn with_visualization(mut self, enable: bool) -> Self {
551        self.enable_visualization = enable;
552        self
553    }
554
555    /// Print configuration parameters (info level logging)
556    pub fn print_configuration(&self) {
557        debug!(
558            "Configuration:\n  Solver:        Dog-Leg\n  Linear solver: {:?}\n  Loss function: N/A\n\nConvergence Criteria:\n  Max iterations:      {}\n  Cost tolerance:      {:.2e}\n  Parameter tolerance: {:.2e}\n  Gradient tolerance:  {:.2e}\n  Timeout:             {:?}\n\nTrust Region:\n  Initial radius:      {:.2e}\n  Radius range:        [{:.2e}, {:.2e}]\n  Min step quality:    {:.2}\n  Good step quality:   {:.2}\n  Poor step quality:   {:.2}\n\nRegularization:\n  Initial mu:          {:.2e}\n  Mu range:            [{:.2e}, {:.2e}]\n  Mu increase factor:  {:.2}\n\nNumerical Settings:\n  Jacobi scaling:      {}\n  Step reuse:          {}\n  Compute covariances: {}",
559            self.linear_solver_type,
560            self.max_iterations,
561            self.cost_tolerance,
562            self.parameter_tolerance,
563            self.gradient_tolerance,
564            self.timeout,
565            self.trust_region_radius,
566            self.trust_region_min,
567            self.trust_region_max,
568            self.min_step_quality,
569            self.good_step_quality,
570            self.poor_step_quality,
571            self.initial_mu,
572            self.min_mu,
573            self.max_mu,
574            self.mu_increase_factor,
575            if self.use_jacobi_scaling {
576                "enabled"
577            } else {
578                "disabled"
579            },
580            if self.enable_step_reuse {
581                "enabled"
582            } else {
583                "disabled"
584            },
585            if self.compute_covariances {
586                "enabled"
587            } else {
588                "disabled"
589            }
590        );
591    }
592}
593
594/// Result from step computation
595struct StepResult {
596    step: faer::Mat<f64>,
597    gradient_norm: f64,
598    predicted_reduction: f64,
599}
600
601/// Type of step taken
602#[derive(Debug, Clone, Copy)]
603enum StepType {
604    /// Full Gauss-Newton step
605    GaussNewton,
606    /// Scaled steepest descent (Cauchy point)
607    SteepestDescent,
608    /// Dog leg interpolation
609    DogLeg,
610}
611
612impl fmt::Display for StepType {
613    fn fmt(&self, f: &mut fmt::Formatter<'_>) -> fmt::Result {
614        match self {
615            StepType::GaussNewton => write!(f, "GN"),
616            StepType::SteepestDescent => write!(f, "SD"),
617            StepType::DogLeg => write!(f, "DL"),
618        }
619    }
620}
621
622/// Result from step evaluation
623struct StepEvaluation {
624    accepted: bool,
625    cost_reduction: f64,
626    rho: f64,
627}
628
629/// Dog Leg trust region solver for nonlinear least squares optimization.
630///
631/// Implements Powell's Dog Leg algorithm with Ceres Solver enhancements including
632/// adaptive μ regularization, numerically robust beta computation, and step reuse caching.
633///
634/// # Algorithm
635///
636/// At each iteration k:
637/// 1. Compute residual `r(xₖ)` and Jacobian `J(xₖ)`
638/// 2. Solve for Gauss-Newton step: `(J^T·J + μI)·h_gn = -J^T·r`
639/// 3. Compute steepest descent direction: `-g` where `g = J^T·r`
640/// 4. Compute Cauchy point: `p_c = -α·g` (optimal along steepest descent)
641/// 5. Construct dog leg step based on trust region radius Δ:
642///    - If `||h_gn|| ≤ Δ`: Take full GN step
643///    - Else if `||p_c|| ≥ Δ`: Take scaled SD to boundary
644///    - Else: Interpolate `h = p_c + β·(h_gn - p_c)` where `||h|| = Δ`
645/// 6. Evaluate gain ratio: `ρ = (actual reduction) / (predicted reduction)`
646/// 7. Update trust region radius based on ρ
647/// 8. Accept/reject step and update parameters
648///
649/// # Ceres Solver Enhancements
650///
651/// This implementation includes four major improvements from Google's Ceres Solver:
652///
653/// **1. Adaptive μ Regularization:** Dynamically adjusts regularization parameter
654/// to handle ill-conditioned Hessians (increases on failure, decreases on success).
655///
656/// **2. Numerically Robust Beta:** Uses two formulas for computing dog leg
657/// interpolation parameter β to avoid catastrophic cancellation.
658///
659/// **3. Step Reuse Mechanism:** Caches GN step, Cauchy point, and gradient when
660/// steps are rejected. Limited to 5 consecutive reuses to prevent staleness.
661///
662/// **4. Jacobi Scaling:** Optional diagonal preconditioning creates elliptical
663/// trust regions for better handling of mixed-scale problems.
664///
665/// # Examples
666///
667/// ```no_run
668/// use apex_solver::optimizer::DogLeg;
669/// use apex_solver::core::problem::Problem;
670/// use apex_solver::JacobianMode;
671/// use std::collections::HashMap;
672///
673/// # type TestResult = Result<(), Box<dyn std::error::Error>>;
674/// # fn main() -> TestResult {
675/// let mut problem = Problem::new(JacobianMode::Sparse);
676/// // ... add factors to problem ...
677///
678/// let initial_values = HashMap::new();
679/// // ... initialize parameters ...
680///
681/// let mut solver = DogLeg::new();
682/// let result = solver.optimize(&problem, &initial_values)?;
683/// # Ok(())
684/// # }
685/// ```
686///
687/// # See Also
688///
689/// - [`DogLegConfig`] - Configuration options
690/// - [`LevenbergMarquardt`](crate::optimizer::LevenbergMarquardt) - Alternative adaptive damping
691/// - [`GaussNewton`](crate::optimizer::GaussNewton) - Undamped variant
692pub struct DogLeg {
693    config: DogLegConfig,
694    jacobi_scaling: Option<Vec<f64>>,
695    observers: optimizer::OptObserverVec,
696
697    // Adaptive mu regularization (Ceres-style)
698    mu: f64,
699    min_mu: f64,
700    max_mu: f64,
701    mu_increase_factor: f64,
702
703    // Step reuse mechanism (Ceres-style efficiency optimization)
704    reuse_step_on_rejection: bool,
705    cached_gn_step: Option<faer::Mat<f64>>,
706    cached_cauchy_point: Option<faer::Mat<f64>>,
707    cached_gradient: Option<faer::Mat<f64>>,
708    cached_alpha: Option<f64>,
709    cache_reuse_count: usize, // Track consecutive reuses to prevent staleness
710}
711
712impl Default for DogLeg {
713    fn default() -> Self {
714        Self::new()
715    }
716}
717
718impl DogLeg {
719    /// Create a new Dog Leg solver with default configuration.
720    pub fn new() -> Self {
721        Self::with_config(DogLegConfig::default())
722    }
723
724    /// Create a new Dog Leg solver with the given configuration.
725    pub fn with_config(config: DogLegConfig) -> Self {
726        Self {
727            // Initialize adaptive mu from config
728            mu: config.initial_mu,
729            min_mu: config.min_mu,
730            max_mu: config.max_mu,
731            mu_increase_factor: config.mu_increase_factor,
732
733            // Initialize step reuse mechanism (disabled initially, enabled after first rejection)
734            reuse_step_on_rejection: false,
735            cached_gn_step: None,
736            cached_cauchy_point: None,
737            cached_gradient: None,
738            cached_alpha: None,
739            cache_reuse_count: 0,
740
741            config,
742            jacobi_scaling: None,
743            observers: optimizer::OptObserverVec::new(),
744        }
745    }
746
747    /// Add an observer to the solver.
748    ///
749    /// Observers are notified at each iteration with the current variable values.
750    /// This enables real-time visualization, logging, metrics collection, etc.
751    ///
752    /// # Examples
753    ///
754    /// ```no_run
755    /// use apex_solver::optimizer::DogLeg;
756    /// # use apex_solver::optimizer::OptObserver;
757    /// # use std::collections::HashMap;
758    /// # use apex_solver::core::problem::VariableEnum;
759    ///
760    /// # struct MyObserver;
761    /// # impl OptObserver for MyObserver {
762    /// #     fn on_step(&self, _: &HashMap<String, VariableEnum>, _: usize) {}
763    /// # }
764    /// let mut solver = DogLeg::new();
765    /// solver.add_observer(MyObserver);
766    /// ```
767    pub fn add_observer(&mut self, observer: impl optimizer::OptObserver + 'static) {
768        self.observers.add(observer);
769    }
770
771    /// Compute Cauchy point and optimal step length for steepest descent (generic).
772    ///
773    /// Returns (alpha, cauchy_point) where:
774    /// - alpha: optimal step length α = ||g||² / (g^T H g)
775    /// - cauchy_point: p_c = -α * g (the Cauchy point)
776    fn compute_cauchy_point_and_alpha_generic<M: AssemblyBackend>(
777        &self,
778        gradient: &faer::Mat<f64>,
779        hessian: &M::Hessian,
780    ) -> (f64, faer::Mat<f64>) {
781        // Optimal step size along steepest descent: α = (g^T*g) / (g^T*H*g)
782        let g_norm_sq_mat = gradient.transpose() * gradient;
783        let g_norm_sq = g_norm_sq_mat[(0, 0)];
784
785        let h_g = M::hessian_vec_product(hessian, gradient);
786        let g_h_g_mat = gradient.transpose() * &h_g;
787        let g_h_g = g_h_g_mat[(0, 0)];
788
789        // Avoid division by zero
790        let alpha = if g_h_g.abs() > 1e-15 {
791            g_norm_sq / g_h_g
792        } else {
793            1.0
794        };
795
796        // Compute Cauchy point: p_c = -α * gradient
797        let mut cauchy_point = faer::Mat::zeros(gradient.nrows(), 1);
798        for i in 0..gradient.nrows() {
799            cauchy_point[(i, 0)] = -alpha * gradient[(i, 0)];
800        }
801
802        (alpha, cauchy_point)
803    }
804
805    /// Compute the dog leg step using Powell's Dog Leg method
806    ///
807    /// The dog leg path consists of two segments:
808    /// 1. From origin to Cauchy point (optimal along steepest descent)
809    /// 2. From Cauchy point to Gauss-Newton step
810    ///
811    /// Arguments:
812    /// - steepest_descent_dir: -gradient (steepest descent direction, not scaled)
813    /// - cauchy_point: p_c = α * (-gradient), the optimal steepest descent step
814    /// - h_gn: Gauss-Newton step
815    /// - delta: Trust region radius
816    ///
817    /// Returns: (step, step_type)
818    fn compute_dog_leg_step(
819        &self,
820        steepest_descent_dir: &faer::Mat<f64>,
821        cauchy_point: &faer::Mat<f64>,
822        h_gn: &faer::Mat<f64>,
823        delta: f64,
824    ) -> (faer::Mat<f64>, StepType) {
825        let gn_norm = h_gn.norm_l2();
826        let cauchy_norm = cauchy_point.norm_l2();
827        let sd_norm = steepest_descent_dir.norm_l2();
828
829        // Case 1: Full Gauss-Newton step fits in trust region
830        if gn_norm <= delta {
831            return (h_gn.clone(), StepType::GaussNewton);
832        }
833
834        // Case 2: Even Cauchy point is outside trust region
835        // Scale steepest descent direction to boundary: (delta / ||δ_sd||) * δ_sd
836        if cauchy_norm >= delta {
837            let mut scaled_sd = faer::Mat::zeros(steepest_descent_dir.nrows(), 1);
838            let scale = delta / sd_norm;
839            for i in 0..steepest_descent_dir.nrows() {
840                scaled_sd[(i, 0)] = steepest_descent_dir[(i, 0)] * scale;
841            }
842            return (scaled_sd, StepType::SteepestDescent);
843        }
844
845        // Case 3: Dog leg interpolation between Cauchy point and GN step
846        // Use Ceres-style numerically robust formula
847        //
848        // Following Ceres solver implementation for numerical stability:
849        // Compute intersection of trust region boundary with line from Cauchy point to GN step
850        //
851        // Let v = δ_gn - p_c
852        // Solve: ||p_c + β*v||² = delta²
853        // This gives: a*β² + 2*b*β + c = 0
854        // where:
855        //   a = v^T·v = ||v||²
856        //   b = p_c^T·v
857        //   c = p_c^T·p_c - delta² = ||p_c||² - delta²
858
859        let mut v = faer::Mat::zeros(cauchy_point.nrows(), 1);
860        for i in 0..cauchy_point.nrows() {
861            v[(i, 0)] = h_gn[(i, 0)] - cauchy_point[(i, 0)];
862        }
863
864        // Compute coefficients
865        let v_squared_norm = v.transpose() * &v;
866        let a = v_squared_norm[(0, 0)];
867
868        let pc_dot_v = cauchy_point.transpose() * &v;
869        let b = pc_dot_v[(0, 0)];
870
871        let c = cauchy_norm * cauchy_norm - delta * delta;
872
873        // Ceres-style numerically robust beta computation
874        // Uses two different formulas based on sign of b to avoid catastrophic cancellation
875        let d_squared = b * b - a * c;
876
877        let beta = if d_squared < 0.0 {
878            // Should not happen geometrically, but handle gracefully
879            1.0
880        } else if a.abs() < 1e-15 {
881            // Degenerate case: v is nearly zero
882            1.0
883        } else {
884            let d = d_squared.sqrt();
885
886            // Ceres formula: choose formula based on sign of b to avoid cancellation
887            // If b <= 0: beta = (-b + d) / a  (standard formula, no cancellation)
888            // If b > 0:  beta = -c / (b + d)  (alternative formula, avoids cancellation)
889            if b <= 0.0 { (-b + d) / a } else { -c / (b + d) }
890        };
891
892        // Clamp beta to [0, 1] for safety
893        let beta = beta.clamp(0.0, 1.0);
894
895        // Compute dog leg step: p_dl = p_c + β*(δ_gn - p_c)
896        let mut dog_leg = faer::Mat::zeros(cauchy_point.nrows(), 1);
897        for i in 0..cauchy_point.nrows() {
898            dog_leg[(i, 0)] = cauchy_point[(i, 0)] + beta * v[(i, 0)];
899        }
900
901        (dog_leg, StepType::DogLeg)
902    }
903
904    /// Update trust region radius based on step quality (Ceres-style)
905    fn update_trust_region(&mut self, rho: f64, step_norm: f64) -> bool {
906        if rho > self.config.good_step_quality {
907            // Good step, increase trust region (Ceres-style: max(radius, 3*step_norm))
908            let new_radius = self.config.trust_region_radius.max(3.0 * step_norm);
909            self.config.trust_region_radius = new_radius.min(self.config.trust_region_max);
910
911            // Decrease mu on successful step (Ceres-style adaptive regularization)
912            self.mu = (self.mu / (0.5 * self.mu_increase_factor)).max(self.min_mu);
913
914            // Clear reuse flag and invalidate cache on acceptance (parameters have moved)
915            self.reuse_step_on_rejection = false;
916            self.cached_gn_step = None;
917            self.cached_cauchy_point = None;
918            self.cached_gradient = None;
919            self.cached_alpha = None;
920            self.cache_reuse_count = 0;
921
922            true
923        } else if rho < self.config.poor_step_quality {
924            // Poor step, decrease trust region
925            self.config.trust_region_radius = (self.config.trust_region_radius
926                * self.config.trust_region_decrease_factor)
927                .max(self.config.trust_region_min);
928
929            // Enable step reuse flag for next iteration (Ceres-style)
930            self.reuse_step_on_rejection = self.config.enable_step_reuse;
931
932            false
933        } else {
934            // Moderate step, keep trust region unchanged
935            // Clear reuse flag and invalidate cache on acceptance (parameters have moved)
936            self.reuse_step_on_rejection = false;
937            self.cached_gn_step = None;
938            self.cached_cauchy_point = None;
939            self.cached_gradient = None;
940            self.cached_alpha = None;
941            self.cache_reuse_count = 0;
942
943            true
944        }
945    }
946
947    /// Compute predicted cost reduction from linear model (generic over assembly mode).
948    fn compute_predicted_reduction_generic<M: AssemblyBackend>(
949        &self,
950        step: &faer::Mat<f64>,
951        gradient: &faer::Mat<f64>,
952        hessian: &M::Hessian,
953    ) -> f64 {
954        // Dog Leg predicted reduction: -step^T * gradient - 0.5 * step^T * H * step
955        let linear_term = step.transpose() * gradient;
956        let hessian_step = M::hessian_vec_product(hessian, step);
957        let quadratic_term = step.transpose() * &hessian_step;
958
959        -linear_term[(0, 0)] - 0.5 * quadratic_term[(0, 0)]
960    }
961
962    /// Compute dog leg optimization step (generic over assembly mode).
963    fn compute_optimization_step_generic<M: AssemblyBackend>(
964        &mut self,
965        residuals: &faer::Mat<f64>,
966        scaled_jacobian: &M::Jacobian,
967        linear_solver: &mut dyn LinearSolver<M>,
968    ) -> Option<StepResult> {
969        // Check if we can reuse cached step (Ceres-style optimization)
970        // Safety limit: prevent excessive reuse that could lead to stale gradient/Hessian
971        const MAX_CACHE_REUSE: usize = 5;
972
973        if self.reuse_step_on_rejection
974            && self.config.enable_step_reuse
975            && self.cache_reuse_count < MAX_CACHE_REUSE
976            && let (Some(cached_gn), Some(cached_cauchy), Some(cached_grad), Some(_cached_a)) = (
977                &self.cached_gn_step,
978                &self.cached_cauchy_point,
979                &self.cached_gradient,
980                &self.cached_alpha,
981            )
982        {
983            // Increment reuse counter
984            self.cache_reuse_count += 1;
985
986            let gradient_norm = cached_grad.norm_l2();
987            let mut steepest_descent = faer::Mat::zeros(cached_grad.nrows(), 1);
988            for i in 0..cached_grad.nrows() {
989                steepest_descent[(i, 0)] = -cached_grad[(i, 0)];
990            }
991
992            let (scaled_step, _step_type) = self.compute_dog_leg_step(
993                &steepest_descent,
994                cached_cauchy,
995                cached_gn,
996                self.config.trust_region_radius,
997            );
998
999            let step = if self.config.use_jacobi_scaling {
1000                let scaling = self.jacobi_scaling.as_ref()?;
1001                M::apply_inverse_scaling(&scaled_step, scaling)
1002            } else {
1003                scaled_step.clone()
1004            };
1005
1006            // For cached reuse, we need the hessian for predicted reduction.
1007            // Use hessian from linear solver if available.
1008            let hessian = linear_solver.get_hessian()?;
1009            let predicted_reduction =
1010                self.compute_predicted_reduction_generic::<M>(&scaled_step, cached_grad, hessian);
1011
1012            return Some(StepResult {
1013                step,
1014                gradient_norm,
1015                predicted_reduction,
1016            });
1017        }
1018
1019        // Not reusing, compute fresh step
1020        // 1. Solve for Gauss-Newton step with adaptive mu regularization (Ceres-style)
1021        let residuals_owned = residuals.as_ref().to_owned();
1022        let mut scaled_gn_step = None;
1023        let mut mu_attempts = 0;
1024
1025        // Try to solve with current mu, increasing if necessary
1026        while mu_attempts < 10 && self.mu <= self.max_mu {
1027            let damping = self.mu;
1028
1029            if let Ok(step) =
1030                linear_solver.solve_augmented_equation(&residuals_owned, scaled_jacobian, damping)
1031            {
1032                scaled_gn_step = Some(step);
1033                break;
1034            }
1035
1036            // Increase mu (Ceres-style)
1037            self.mu = (self.mu * self.mu_increase_factor).min(self.max_mu);
1038            mu_attempts += 1;
1039        }
1040
1041        let scaled_gn_step = scaled_gn_step?;
1042
1043        // 2. Get gradient and Hessian (cached by solve_augmented_equation)
1044        let gradient = linear_solver.get_gradient()?;
1045        let hessian = linear_solver.get_hessian()?;
1046        let gradient_norm = gradient.norm_l2();
1047
1048        // 3. Compute steepest descent direction: δ_sd = -gradient
1049        let mut steepest_descent = faer::Mat::zeros(gradient.nrows(), 1);
1050        for i in 0..gradient.nrows() {
1051            steepest_descent[(i, 0)] = -gradient[(i, 0)];
1052        }
1053
1054        // 4. Compute Cauchy point and optimal step length
1055        let (alpha, cauchy_point) =
1056            self.compute_cauchy_point_and_alpha_generic::<M>(gradient, hessian);
1057
1058        // 5. Compute dog leg step based on trust region radius
1059        let (scaled_step, _step_type) = self.compute_dog_leg_step(
1060            &steepest_descent,
1061            &cauchy_point,
1062            &scaled_gn_step,
1063            self.config.trust_region_radius,
1064        );
1065
1066        // 6. Apply inverse Jacobi scaling if enabled
1067        let step = if self.config.use_jacobi_scaling {
1068            let scaling = self.jacobi_scaling.as_ref()?;
1069            M::apply_inverse_scaling(&scaled_step, scaling)
1070        } else {
1071            scaled_step.clone()
1072        };
1073
1074        // 7. Compute predicted reduction
1075        let predicted_reduction =
1076            self.compute_predicted_reduction_generic::<M>(&scaled_step, gradient, hessian);
1077
1078        // 8. Cache step components for potential reuse (Ceres-style)
1079        self.cached_gn_step = Some(scaled_gn_step.clone());
1080        self.cached_cauchy_point = Some(cauchy_point.clone());
1081        self.cached_gradient = Some(gradient.clone());
1082        self.cached_alpha = Some(alpha);
1083
1084        Some(StepResult {
1085            step,
1086            gradient_norm,
1087            predicted_reduction,
1088        })
1089    }
1090
1091    /// Evaluate and apply step
1092    fn evaluate_and_apply_step(
1093        &mut self,
1094        step_result: &StepResult,
1095        state: &mut optimizer::InitializedState,
1096        problem: &problem::Problem,
1097    ) -> error::ApexSolverResult<StepEvaluation> {
1098        // Apply parameter updates
1099        let step_norm = optimizer::apply_parameter_step(
1100            &mut state.variables,
1101            step_result.step.as_ref(),
1102            &state.sorted_vars,
1103        );
1104
1105        // Compute new cost (residual only, no Jacobian needed for step evaluation)
1106        let new_residual = problem.compute_residual_sparse(&state.variables)?;
1107        let new_cost = optimizer::compute_cost(&new_residual);
1108
1109        // Compute step quality
1110        let rho = optimizer::compute_step_quality(
1111            state.current_cost,
1112            new_cost,
1113            step_result.predicted_reduction,
1114        );
1115
1116        // Update trust region and decide acceptance
1117        // Filter out numerical noise with small threshold
1118        let accepted = rho > 1e-4;
1119        let _trust_region_updated = self.update_trust_region(rho, step_norm);
1120
1121        let cost_reduction = if accepted {
1122            let reduction = state.current_cost - new_cost;
1123            state.current_cost = new_cost;
1124            reduction
1125        } else {
1126            // Reject step - revert changes
1127            optimizer::apply_negative_parameter_step(
1128                &mut state.variables,
1129                step_result.step.as_ref(),
1130                &state.sorted_vars,
1131            );
1132            0.0
1133        };
1134
1135        Ok(StepEvaluation {
1136            accepted,
1137            cost_reduction,
1138            rho,
1139        })
1140    }
1141
1142    /// Run optimization using the specified assembly mode and linear solver.
1143    fn optimize_with_mode<M: AssemblyBackend>(
1144        &mut self,
1145        problem: &problem::Problem,
1146        initial_params: &collections::HashMap<
1147            String,
1148            (manifold::ManifoldType, nalgebra::DVector<f64>),
1149        >,
1150        linear_solver: &mut dyn LinearSolver<M>,
1151    ) -> Result<
1152        optimizer::SolverResult<collections::HashMap<String, problem::VariableEnum>>,
1153        error::ApexSolverError,
1154    > {
1155        let start_time = time::Instant::now();
1156        let mut iteration = 0;
1157        let mut cost_evaluations = 1;
1158        let mut jacobian_evaluations = 0;
1159        let mut successful_steps = 0;
1160        let mut unsuccessful_steps = 0;
1161
1162        let mut state = optimizer::initialize_optimization_state(problem, initial_params)?;
1163
1164        let mut max_gradient_norm: f64 = 0.0;
1165        let mut max_parameter_update_norm: f64 = 0.0;
1166        let mut total_cost_reduction = 0.0;
1167        let mut final_gradient_norm;
1168        let mut final_parameter_update_norm;
1169
1170        // Initialize iteration statistics tracking
1171        let mut iteration_stats = Vec::with_capacity(self.config.max_iterations);
1172        let mut previous_cost = state.current_cost;
1173
1174        // Print configuration and header if debug level is enabled
1175        if tracing::enabled!(tracing::Level::DEBUG) {
1176            self.config.print_configuration();
1177            IterationStats::print_header();
1178        }
1179
1180        loop {
1181            let iter_start = time::Instant::now();
1182
1183            // Evaluate residuals and Jacobian using the assembly mode
1184            let (residuals, jacobian) = M::assemble(
1185                problem,
1186                &state.variables,
1187                &state.variable_index_map,
1188                state.symbolic_structure.as_ref(),
1189                state.total_dof,
1190            )?;
1191            jacobian_evaluations += 1;
1192
1193            // Process Jacobian (apply scaling if enabled)
1194            let scaled_jacobian = if self.config.use_jacobi_scaling {
1195                optimizer::process_jacobian_generic::<M>(
1196                    &jacobian,
1197                    &mut self.jacobi_scaling,
1198                    iteration,
1199                )?
1200            } else {
1201                jacobian
1202            };
1203
1204            // Compute dog leg step
1205            let step_result = match self.compute_optimization_step_generic::<M>(
1206                &residuals,
1207                &scaled_jacobian,
1208                linear_solver,
1209            ) {
1210                Some(result) => result,
1211                None => {
1212                    return Err(optimizer::OptimizerError::LinearSolveFailed(
1213                        "Linear solver failed to solve system".to_string(),
1214                    )
1215                    .into());
1216                }
1217            };
1218
1219            // Update tracking
1220            max_gradient_norm = max_gradient_norm.max(step_result.gradient_norm);
1221            final_gradient_norm = step_result.gradient_norm;
1222            let step_norm = step_result.step.norm_l2();
1223            max_parameter_update_norm = max_parameter_update_norm.max(step_norm);
1224            final_parameter_update_norm = step_norm;
1225
1226            // Evaluate and apply step
1227            let step_eval = self.evaluate_and_apply_step(&step_result, &mut state, problem)?;
1228            cost_evaluations += 1;
1229
1230            if step_eval.accepted {
1231                successful_steps += 1;
1232                total_cost_reduction += step_eval.cost_reduction;
1233            } else {
1234                unsuccessful_steps += 1;
1235            }
1236
1237            // OPTIMIZATION: Only collect iteration statistics if debug level is enabled
1238            if tracing::enabled!(tracing::Level::DEBUG) {
1239                let iter_elapsed_ms = iter_start.elapsed().as_secs_f64() * 1000.0;
1240                let total_elapsed_ms = start_time.elapsed().as_secs_f64() * 1000.0;
1241
1242                let stats = IterationStats {
1243                    iteration,
1244                    cost: state.current_cost,
1245                    cost_change: previous_cost - state.current_cost,
1246                    gradient_norm: step_result.gradient_norm,
1247                    step_norm,
1248                    tr_ratio: step_eval.rho,
1249                    tr_radius: self.config.trust_region_radius,
1250                    ls_iter: 0,
1251                    iter_time_ms: iter_elapsed_ms,
1252                    total_time_ms: total_elapsed_ms,
1253                    accepted: step_eval.accepted,
1254                };
1255
1256                iteration_stats.push(stats.clone());
1257                stats.print_line();
1258            }
1259
1260            previous_cost = state.current_cost;
1261
1262            // Notify all observers with current state
1263            optimizer::notify_observers_generic::<M>(
1264                &mut self.observers,
1265                &state.variables,
1266                iteration,
1267                state.current_cost,
1268                step_result.gradient_norm,
1269                Some(self.config.trust_region_radius),
1270                step_norm,
1271                Some(step_eval.rho),
1272                linear_solver,
1273            );
1274
1275            // Check convergence
1276            let elapsed = start_time.elapsed();
1277            let parameter_norm = optimizer::compute_parameter_norm(&state.variables);
1278            let new_cost = state.current_cost;
1279
1280            let cost_before_step = if step_eval.accepted {
1281                state.current_cost + step_eval.cost_reduction
1282            } else {
1283                state.current_cost
1284            };
1285
1286            if let Some(status) = optimizer::check_convergence(&optimizer::ConvergenceParams {
1287                iteration,
1288                current_cost: cost_before_step,
1289                new_cost,
1290                parameter_norm,
1291                parameter_update_norm: step_norm,
1292                gradient_norm: step_result.gradient_norm,
1293                elapsed,
1294                step_accepted: step_eval.accepted,
1295                max_iterations: self.config.max_iterations,
1296                gradient_tolerance: self.config.gradient_tolerance,
1297                parameter_tolerance: self.config.parameter_tolerance,
1298                cost_tolerance: self.config.cost_tolerance,
1299                min_cost_threshold: self.config.min_cost_threshold,
1300                timeout: self.config.timeout,
1301                trust_region_radius: Some(self.config.trust_region_radius),
1302                min_trust_region_radius: Some(self.config.trust_region_min),
1303            }) {
1304                // Print summary only if debug level is enabled
1305                if tracing::enabled!(tracing::Level::DEBUG) {
1306                    let summary = optimizer::create_optimizer_summary(
1307                        "Dog-Leg",
1308                        state.initial_cost,
1309                        state.current_cost,
1310                        iteration + 1,
1311                        Some(successful_steps),
1312                        Some(unsuccessful_steps),
1313                        max_gradient_norm,
1314                        final_gradient_norm,
1315                        max_parameter_update_norm,
1316                        final_parameter_update_norm,
1317                        total_cost_reduction,
1318                        elapsed,
1319                        iteration_stats.clone(),
1320                        status.clone(),
1321                        None,
1322                        Some(self.config.trust_region_radius),
1323                        None,
1324                    );
1325                    debug!("{}", summary);
1326                }
1327
1328                // Compute covariances if enabled
1329                let covariances = if self.config.compute_covariances {
1330                    problem.compute_and_set_covariances_generic::<M>(
1331                        linear_solver,
1332                        &mut state.variables,
1333                        &state.variable_index_map,
1334                    )
1335                } else {
1336                    None
1337                };
1338
1339                return Ok(optimizer::build_solver_result(
1340                    status,
1341                    iteration + 1,
1342                    state,
1343                    elapsed,
1344                    final_gradient_norm,
1345                    final_parameter_update_norm,
1346                    cost_evaluations,
1347                    jacobian_evaluations,
1348                    covariances,
1349                ));
1350            }
1351
1352            iteration += 1;
1353        }
1354    }
1355
1356    /// Run optimization, automatically selecting sparse or dense path based on config.
1357    pub fn optimize(
1358        &mut self,
1359        problem: &problem::Problem,
1360        initial_params: &collections::HashMap<
1361            String,
1362            (manifold::ManifoldType, nalgebra::DVector<f64>),
1363        >,
1364    ) -> Result<
1365        optimizer::SolverResult<collections::HashMap<String, problem::VariableEnum>>,
1366        error::ApexSolverError,
1367    > {
1368        match problem.jacobian_mode {
1369            JacobianMode::Dense => match self.config.linear_solver_type {
1370                LinearSolverType::DenseQR => {
1371                    let mut solver = DenseQRSolver::new();
1372                    self.optimize_with_mode::<DenseMode>(problem, initial_params, &mut solver)
1373                }
1374                _ => {
1375                    let mut solver = DenseCholeskySolver::new();
1376                    self.optimize_with_mode::<DenseMode>(problem, initial_params, &mut solver)
1377                }
1378            },
1379            JacobianMode::Sparse => match self.config.linear_solver_type {
1380                linalg::LinearSolverType::SparseQR => {
1381                    let mut solver = SparseQRSolver::new();
1382                    self.optimize_with_mode::<SparseMode>(problem, initial_params, &mut solver)
1383                }
1384                _ => {
1385                    // SparseCholesky (default), SparseSchurComplement or DenseCholesky with
1386                    // sparse mode → SparseCholeskySolver
1387                    let mut solver = SparseCholeskySolver::new();
1388                    self.optimize_with_mode::<SparseMode>(problem, initial_params, &mut solver)
1389                }
1390            },
1391        }
1392    }
1393}
1394
1395impl optimizer::Optimizer for DogLeg {
1396    fn optimize(
1397        &mut self,
1398        problem: &problem::Problem,
1399        initial_params: &collections::HashMap<
1400            String,
1401            (manifold::ManifoldType, nalgebra::DVector<f64>),
1402        >,
1403    ) -> Result<
1404        optimizer::SolverResult<collections::HashMap<String, problem::VariableEnum>>,
1405        crate::error::ApexSolverError,
1406    > {
1407        self.optimize(problem, initial_params)
1408    }
1409}
1410
1411#[cfg(test)]
1412mod tests {
1413    use super::*;
1414    use crate::factors;
1415    use apex_manifolds as manifold;
1416    use nalgebra;
1417
1418    type TestResult = Result<(), Box<dyn std::error::Error>>;
1419
1420    /// Custom Rosenbrock Factor 1: r1 = 10(x2 - x1²)
1421    /// Demonstrates extensibility - custom factors can be defined outside of factors.rs
1422    #[derive(Debug, Clone)]
1423    struct RosenbrockFactor1;
1424
1425    impl factors::Factor for RosenbrockFactor1 {
1426        fn linearize(
1427            &self,
1428            params: &[nalgebra::DVector<f64>],
1429            compute_jacobian: bool,
1430        ) -> (nalgebra::DVector<f64>, Option<nalgebra::DMatrix<f64>>) {
1431            let x1 = params[0][0];
1432            let x2 = params[1][0];
1433
1434            // Residual: r1 = 10(x2 - x1²)
1435            let residual = nalgebra::dvector![10.0 * (x2 - x1 * x1)];
1436
1437            // Jacobian: ∂r1/∂x1 = -20*x1, ∂r1/∂x2 = 10
1438            let jacobian = if compute_jacobian {
1439                let mut jac = nalgebra::DMatrix::zeros(1, 2);
1440                jac[(0, 0)] = -20.0 * x1;
1441                jac[(0, 1)] = 10.0;
1442                Some(jac)
1443            } else {
1444                None
1445            };
1446
1447            (residual, jacobian)
1448        }
1449
1450        fn get_dimension(&self) -> usize {
1451            1
1452        }
1453    }
1454
1455    /// Custom Rosenbrock Factor 2: r2 = 1 - x1
1456    /// Demonstrates extensibility - custom factors can be defined outside of factors.rs
1457    #[derive(Debug, Clone)]
1458    struct RosenbrockFactor2;
1459
1460    impl factors::Factor for RosenbrockFactor2 {
1461        fn linearize(
1462            &self,
1463            params: &[nalgebra::DVector<f64>],
1464            compute_jacobian: bool,
1465        ) -> (nalgebra::DVector<f64>, Option<nalgebra::DMatrix<f64>>) {
1466            let x1 = params[0][0];
1467
1468            // Residual: r2 = 1 - x1
1469            let residual = nalgebra::dvector![1.0 - x1];
1470
1471            // Jacobian: ∂r2/∂x1 = -1
1472            let jacobian = if compute_jacobian {
1473                Some(nalgebra::DMatrix::from_element(1, 1, -1.0))
1474            } else {
1475                None
1476            };
1477
1478            (residual, jacobian)
1479        }
1480
1481        fn get_dimension(&self) -> usize {
1482            1
1483        }
1484    }
1485
1486    #[test]
1487    fn test_rosenbrock_optimization() -> TestResult {
1488        // Rosenbrock function test:
1489        // Minimize: r1² + r2² where
1490        //   r1 = 10(x2 - x1²)
1491        //   r2 = 1 - x1
1492        // Starting point: [-1.2, 1.0]
1493        // Expected minimum: [1.0, 1.0]
1494
1495        let mut problem = problem::Problem::new(JacobianMode::Sparse);
1496        let mut initial_values = collections::HashMap::new();
1497
1498        // Add variables using Rn manifold (Euclidean space)
1499        initial_values.insert(
1500            "x1".to_string(),
1501            (manifold::ManifoldType::RN, nalgebra::dvector![-1.2]),
1502        );
1503        initial_values.insert(
1504            "x2".to_string(),
1505            (manifold::ManifoldType::RN, nalgebra::dvector![1.0]),
1506        );
1507
1508        // Add custom factors (demonstrates extensibility!)
1509        problem.add_residual_block(&["x1", "x2"], Box::new(RosenbrockFactor1), None);
1510        problem.add_residual_block(&["x1"], Box::new(RosenbrockFactor2), None);
1511
1512        // Configure Dog Leg optimizer with appropriate trust region
1513        let config = DogLegConfig::new()
1514            .with_max_iterations(100)
1515            .with_cost_tolerance(1e-8)
1516            .with_parameter_tolerance(1e-8)
1517            .with_gradient_tolerance(1e-10)
1518            .with_trust_region_radius(10.0); // Start with larger trust region
1519
1520        let mut solver = DogLeg::with_config(config);
1521        let result = solver.optimize(&problem, &initial_values)?;
1522
1523        // Extract final values
1524        let x1_final = result
1525            .parameters
1526            .get("x1")
1527            .ok_or("x1 not found")?
1528            .to_vector()[0];
1529        let x2_final = result
1530            .parameters
1531            .get("x2")
1532            .ok_or("x2 not found")?
1533            .to_vector()[0];
1534
1535        // Verify convergence to [1.0, 1.0]
1536        assert!(
1537            matches!(
1538                result.status,
1539                optimizer::OptimizationStatus::Converged
1540                    | optimizer::OptimizationStatus::CostToleranceReached
1541                    | optimizer::OptimizationStatus::ParameterToleranceReached
1542                    | optimizer::OptimizationStatus::GradientToleranceReached
1543            ),
1544            "Optimization should converge"
1545        );
1546        assert!(
1547            (x1_final - 1.0).abs() < 1e-4,
1548            "x1 should converge to 1.0, got {}",
1549            x1_final
1550        );
1551        assert!(
1552            (x2_final - 1.0).abs() < 1e-4,
1553            "x2 should converge to 1.0, got {}",
1554            x2_final
1555        );
1556        assert!(
1557            result.final_cost < 1e-6,
1558            "Final cost should be near zero, got {}",
1559            result.final_cost
1560        );
1561        Ok(())
1562    }
1563
1564    /// Trivial factor: r = x - target, J = [[1.0]]
1565    struct LinearFactor {
1566        target: f64,
1567    }
1568
1569    impl factors::Factor for LinearFactor {
1570        fn linearize(
1571            &self,
1572            params: &[nalgebra::DVector<f64>],
1573            compute_jacobian: bool,
1574        ) -> (nalgebra::DVector<f64>, Option<nalgebra::DMatrix<f64>>) {
1575            let residual = nalgebra::dvector![params[0][0] - self.target];
1576            let jacobian = if compute_jacobian {
1577                Some(nalgebra::DMatrix::from_element(1, 1, 1.0))
1578            } else {
1579                None
1580            };
1581            (residual, jacobian)
1582        }
1583
1584        fn get_dimension(&self) -> usize {
1585            1
1586        }
1587    }
1588
1589    fn rosenbrock_problem() -> (
1590        problem::Problem,
1591        collections::HashMap<String, (manifold::ManifoldType, nalgebra::DVector<f64>)>,
1592    ) {
1593        let mut prob = problem::Problem::new(JacobianMode::Sparse);
1594        let mut init = collections::HashMap::new();
1595        init.insert(
1596            "x1".to_string(),
1597            (manifold::ManifoldType::RN, nalgebra::dvector![-1.2]),
1598        );
1599        init.insert(
1600            "x2".to_string(),
1601            (manifold::ManifoldType::RN, nalgebra::dvector![1.0]),
1602        );
1603        prob.add_residual_block(&["x1", "x2"], Box::new(RosenbrockFactor1), None);
1604        prob.add_residual_block(&["x1"], Box::new(RosenbrockFactor2), None);
1605        (prob, init)
1606    }
1607
1608    fn linear_problem(
1609        start: f64,
1610    ) -> (
1611        problem::Problem,
1612        collections::HashMap<String, (manifold::ManifoldType, nalgebra::DVector<f64>)>,
1613    ) {
1614        let mut prob = problem::Problem::new(JacobianMode::Sparse);
1615        let mut init = collections::HashMap::new();
1616        init.insert(
1617            "x".to_string(),
1618            (manifold::ManifoldType::RN, nalgebra::dvector![start]),
1619        );
1620        prob.add_residual_block(&["x"], Box::new(LinearFactor { target: 0.0 }), None);
1621        (prob, init)
1622    }
1623
1624    // -------------------------------------------------------------------------
1625    // DogLegConfig builder tests
1626    // -------------------------------------------------------------------------
1627
1628    #[test]
1629    fn test_dl_config_default() {
1630        let cfg = DogLegConfig::default();
1631        assert_eq!(cfg.max_iterations, 50);
1632        assert!((cfg.cost_tolerance - 1e-6).abs() < 1e-15);
1633        assert!(cfg.use_jacobi_scaling); // DogLeg enables Jacobi scaling by default
1634        assert!(!cfg.compute_covariances);
1635        assert!(cfg.enable_step_reuse);
1636    }
1637
1638    #[test]
1639    fn test_dl_config_builders() {
1640        use crate::linalg::LinearSolverType;
1641        let cfg = DogLegConfig::new()
1642            .with_max_iterations(20)
1643            .with_cost_tolerance(1e-5)
1644            .with_parameter_tolerance(1e-6)
1645            .with_gradient_tolerance(1e-7)
1646            .with_trust_region_radius(200.0)
1647            .with_trust_region_bounds(1e-15, 1e15)
1648            .with_trust_region_factors(4.0, 0.3)
1649            .with_step_quality_thresholds(0.01, 0.3, 0.7)
1650            .with_jacobi_scaling(false)
1651            .with_mu_params(1e-3, 1e-10, 2.0, 5.0)
1652            .with_step_reuse(false)
1653            .with_min_cost_threshold(1e-9)
1654            .with_compute_covariances(true)
1655            .with_linear_solver_type(LinearSolverType::SparseQR);
1656        assert_eq!(cfg.max_iterations, 20);
1657        assert!((cfg.trust_region_radius - 200.0).abs() < 1e-10);
1658        assert!(!cfg.use_jacobi_scaling);
1659        assert!(!cfg.enable_step_reuse);
1660        assert!(cfg.min_cost_threshold.is_some());
1661        assert!(cfg.compute_covariances);
1662        assert!((cfg.min_step_quality - 0.01).abs() < 1e-15);
1663        assert!((cfg.poor_step_quality - 0.3).abs() < 1e-15);
1664        assert!((cfg.good_step_quality - 0.7).abs() < 1e-15);
1665        assert!((cfg.initial_mu - 1e-3).abs() < 1e-15);
1666    }
1667
1668    #[test]
1669    fn test_dl_print_configuration_no_panic() {
1670        DogLegConfig::default().print_configuration();
1671    }
1672
1673    #[test]
1674    fn test_dl_default_equals_new() {
1675        let _a = DogLeg::new();
1676        let _b = DogLeg::default();
1677    }
1678
1679    #[test]
1680    fn test_dl_with_config_method() {
1681        let cfg = DogLegConfig::new().with_max_iterations(5);
1682        let _solver = DogLeg::with_config(cfg);
1683    }
1684
1685    // -------------------------------------------------------------------------
1686    // Convergence termination paths
1687    // -------------------------------------------------------------------------
1688
1689    #[test]
1690    fn test_dl_max_iterations_termination() -> TestResult {
1691        let (problem, initial_values) = rosenbrock_problem();
1692        let cfg = DogLegConfig::new().with_max_iterations(2);
1693        let mut solver = DogLeg::with_config(cfg);
1694        let result = solver.optimize(&problem, &initial_values)?;
1695        assert_eq!(
1696            result.status,
1697            optimizer::OptimizationStatus::MaxIterationsReached
1698        );
1699        assert!(result.iterations <= 3, "iterations={}", result.iterations);
1700        Ok(())
1701    }
1702
1703    #[test]
1704    fn test_dl_gradient_tolerance_convergence() -> TestResult {
1705        let (problem, initial_values) = linear_problem(1.0);
1706        let cfg = DogLegConfig::new()
1707            .with_gradient_tolerance(1e3)
1708            .with_cost_tolerance(1e-20)
1709            .with_parameter_tolerance(1e-20);
1710        let mut solver = DogLeg::with_config(cfg);
1711        let result = solver.optimize(&problem, &initial_values)?;
1712        assert_eq!(
1713            result.status,
1714            optimizer::OptimizationStatus::GradientToleranceReached
1715        );
1716        Ok(())
1717    }
1718
1719    #[test]
1720    fn test_dl_cost_tolerance_convergence() -> TestResult {
1721        let (problem, initial_values) = rosenbrock_problem();
1722        let cfg = DogLegConfig::new()
1723            .with_cost_tolerance(1e2) // very loose
1724            .with_gradient_tolerance(1e-20)
1725            .with_parameter_tolerance(1e-20);
1726        let mut solver = DogLeg::with_config(cfg);
1727        let result = solver.optimize(&problem, &initial_values)?;
1728        assert!(matches!(
1729            result.status,
1730            optimizer::OptimizationStatus::CostToleranceReached
1731                | optimizer::OptimizationStatus::GradientToleranceReached
1732                | optimizer::OptimizationStatus::ParameterToleranceReached
1733                | optimizer::OptimizationStatus::Converged
1734        ));
1735        Ok(())
1736    }
1737
1738    #[test]
1739    fn test_dl_qr_solver() -> TestResult {
1740        use crate::linalg::LinearSolverType;
1741        let (problem, initial_values) = rosenbrock_problem();
1742        let cfg = DogLegConfig::new()
1743            .with_linear_solver_type(LinearSolverType::SparseQR)
1744            .with_max_iterations(100);
1745        let mut solver = DogLeg::with_config(cfg);
1746        let result = solver.optimize(&problem, &initial_values)?;
1747        assert!(result.final_cost < 1e-6);
1748        Ok(())
1749    }
1750
1751    #[test]
1752    fn test_dl_jacobi_scaling_disabled() -> TestResult {
1753        // DogLeg has Jacobi scaling ON by default; test with it explicitly disabled
1754        let (problem, initial_values) = rosenbrock_problem();
1755        let cfg = DogLegConfig::new()
1756            .with_jacobi_scaling(false)
1757            .with_max_iterations(100);
1758        let mut solver = DogLeg::with_config(cfg);
1759        let result = solver.optimize(&problem, &initial_values)?;
1760        assert!(result.final_cost < 1e-6);
1761        Ok(())
1762    }
1763
1764    #[test]
1765    fn test_dl_min_cost_threshold() -> TestResult {
1766        let (problem, initial_values) = rosenbrock_problem();
1767        let cfg = DogLegConfig::new()
1768            .with_min_cost_threshold(1e10)
1769            .with_cost_tolerance(1e-20)
1770            .with_gradient_tolerance(1e-20)
1771            .with_parameter_tolerance(1e-20);
1772        let mut solver = DogLeg::with_config(cfg);
1773        let result = solver.optimize(&problem, &initial_values)?;
1774        assert_eq!(
1775            result.status,
1776            optimizer::OptimizationStatus::MinCostThresholdReached
1777        );
1778        Ok(())
1779    }
1780
1781    #[test]
1782    fn test_dl_trust_region_radius_in_config() -> TestResult {
1783        let (problem, initial_values) = rosenbrock_problem();
1784        let cfg = DogLegConfig::new()
1785            .with_trust_region_radius(0.1) // small initial radius
1786            .with_max_iterations(200);
1787        let mut solver = DogLeg::with_config(cfg);
1788        let result = solver.optimize(&problem, &initial_values)?;
1789        assert!(result.iterations > 0);
1790        Ok(())
1791    }
1792
1793    #[test]
1794    fn test_dl_step_reuse_config() -> TestResult {
1795        let (problem, initial_values) = rosenbrock_problem();
1796        let cfg = DogLegConfig::new()
1797            .with_step_reuse(true)
1798            .with_max_iterations(100);
1799        let mut solver = DogLeg::with_config(cfg);
1800        let result = solver.optimize(&problem, &initial_values)?;
1801        assert!(result.final_cost < 1e-6);
1802        Ok(())
1803    }
1804
1805    #[test]
1806    fn test_dl_mu_params_config() {
1807        let cfg = DogLegConfig::new().with_mu_params(1e-2, 1e-9, 10.0, 20.0);
1808        assert!((cfg.initial_mu - 1e-2).abs() < 1e-15);
1809        assert!((cfg.min_mu - 1e-9).abs() < 1e-20);
1810        assert!((cfg.max_mu - 10.0).abs() < 1e-12);
1811        assert!((cfg.mu_increase_factor - 20.0).abs() < 1e-12);
1812    }
1813
1814    #[test]
1815    fn test_dl_step_quality_thresholds() {
1816        let cfg = DogLegConfig::new().with_step_quality_thresholds(0.05, 0.25, 0.8);
1817        assert!((cfg.min_step_quality - 0.05).abs() < 1e-15);
1818        assert!((cfg.poor_step_quality - 0.25).abs() < 1e-15);
1819        assert!((cfg.good_step_quality - 0.8).abs() < 1e-15);
1820    }
1821
1822    #[test]
1823    fn test_dl_result_fields() -> TestResult {
1824        let (problem, initial_values) = rosenbrock_problem();
1825        let mut solver = DogLeg::new();
1826        let result = solver.optimize(&problem, &initial_values)?;
1827        assert!(result.initial_cost > result.final_cost);
1828        assert!(result.iterations > 0);
1829        assert!(result.convergence_info.is_some());
1830        Ok(())
1831    }
1832
1833    #[test]
1834    fn test_dl_timeout_config() {
1835        let cfg = DogLegConfig::new().with_timeout(time::Duration::from_secs(60));
1836        assert!(cfg.timeout.is_some());
1837    }
1838}