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 std::collections::HashMap;
155//!
156//! # fn main() -> Result<(), Box<dyn std::error::Error>> {
157//! let mut problem = Problem::new();
158//! // ... add residual blocks (factors) to problem ...
159//!
160//! let initial_values = HashMap::new();
161//! // ... initialize parameters ...
162//!
163//! let mut solver = DogLeg::new();
164//! let result = solver.optimize(&problem, &initial_values)?;
165//! # Ok(())
166//! # }
167//! ```
168//!
169//! ## Advanced configuration with Ceres enhancements
170//!
171//! ```no_run
172//! use apex_solver::optimizer::dog_leg::{DogLegConfig, DogLeg};
173//! use apex_solver::linalg::LinearSolverType;
174//!
175//! # fn main() {
176//! let config = DogLegConfig::new()
177//!     .with_max_iterations(100)
178//!     .with_trust_region_radius(1e4)  // Large initial radius
179//!     .with_trust_region_bounds(1e-3, 1e6)  // Min/max radius
180//!     .with_mu_params(1e-4, 1e-8, 1.0, 10.0)  // Conservative regularization
181//!     .with_jacobi_scaling(true)  // Enable elliptical trust regions
182//!     .with_step_reuse(true);  // Enable Ceres-style caching
183//!
184//! let mut solver = DogLeg::with_config(config);
185//! # }
186//! ```
187//!
188//! # References
189//!
190//! - Powell, M. J. D. (1970). "A Hybrid Method for Nonlinear Equations". *Numerical Methods for Nonlinear Algebraic Equations*. Gordon and Breach.
191//! - Nocedal, J. & Wright, S. (2006). *Numerical Optimization* (2nd ed.). Springer. Chapter 4 (Trust Region Methods).
192//! - Madsen, K., Nielsen, H. B., & Tingleff, O. (2004). *Methods for Non-Linear Least Squares Problems* (2nd ed.). Chapter 6.
193//! - Conn, A. R., Gould, N. I. M., & Toint, P. L. (2000). *Trust-Region Methods*. SIAM.
194//! - Ceres Solver: http://ceres-solver.org/ - Google's C++ nonlinear least squares library
195
196use crate::{core::problem, error, linalg, optimizer};
197use apex_manifolds as manifold;
198
199use faer::sparse;
200use std::{collections, fmt, time};
201use tracing::debug;
202
203use crate::optimizer::IterationStats;
204
205/// Configuration parameters for the Dog Leg trust region optimizer.
206///
207/// Controls trust region management, convergence criteria, adaptive regularization,
208/// and Ceres Solver enhancements for the Dog Leg algorithm.
209///
210/// # Builder Pattern
211///
212/// All configuration options can be set using the builder pattern:
213///
214/// ```
215/// use apex_solver::optimizer::dog_leg::DogLegConfig;
216///
217/// let config = DogLegConfig::new()
218///     .with_max_iterations(100)
219///     .with_trust_region_radius(1e4)
220///     .with_mu_params(1e-4, 1e-8, 1.0, 10.0)
221///     .with_jacobi_scaling(true)
222///     .with_step_reuse(true);
223/// ```
224///
225/// # Trust Region Behavior
226///
227/// The trust region radius Δ controls the maximum allowed step size:
228///
229/// - **Initial radius** (`trust_region_radius`): Starting value (default: 1e4)
230/// - **Bounds** (`trust_region_min`, `trust_region_max`): Valid range (default: 1e-3 to 1e6)
231/// - **Adaptation**: Increases for good steps, decreases for poor steps
232///
233/// # Adaptive μ Regularization (Ceres Enhancement)
234///
235/// Controls the regularization parameter in `(J^T·J + μI)·h = -J^T·r`:
236///
237/// - `initial_mu`: Starting value (default: 1e-4 for numerical stability)
238/// - `min_mu`, `max_mu`: Bounds (default: 1e-8 to 1.0)
239/// - `mu_increase_factor`: Multiplier when solve fails (default: 10.0)
240///
241/// # Convergence Criteria
242///
243/// The optimizer terminates when ANY of the following conditions is met:
244///
245/// - **Cost tolerance**: `|cost_k - cost_{k-1}| < cost_tolerance`
246/// - **Parameter tolerance**: `||step|| < parameter_tolerance`
247/// - **Gradient tolerance**: `||J^T·r|| < gradient_tolerance`
248/// - **Maximum iterations**: `iteration >= max_iterations`
249/// - **Timeout**: `elapsed_time >= timeout`
250///
251/// # See Also
252///
253/// - [`DogLeg`] - The solver that uses this configuration
254/// - [`LevenbergMarquardtConfig`](crate::optimizer::LevenbergMarquardtConfig) - Alternative damping approach
255/// - [`GaussNewtonConfig`](crate::optimizer::GaussNewtonConfig) - Undamped variant
256#[derive(Clone)]
257pub struct DogLegConfig {
258    /// Type of linear solver for the linear systems
259    pub linear_solver_type: linalg::LinearSolverType,
260    /// Maximum number of iterations
261    pub max_iterations: usize,
262    /// Convergence tolerance for cost function
263    pub cost_tolerance: f64,
264    /// Convergence tolerance for parameter updates
265    pub parameter_tolerance: f64,
266    /// Convergence tolerance for gradient norm
267    pub gradient_tolerance: f64,
268    /// Timeout duration
269    pub timeout: Option<time::Duration>,
270    /// Initial trust region radius
271    pub trust_region_radius: f64,
272    /// Minimum trust region radius
273    pub trust_region_min: f64,
274    /// Maximum trust region radius
275    pub trust_region_max: f64,
276    /// Trust region increase factor (for good steps, rho > 0.75)
277    pub trust_region_increase_factor: f64,
278    /// Trust region decrease factor (for poor steps, rho < 0.25)
279    pub trust_region_decrease_factor: f64,
280    /// Minimum step quality for acceptance (typically 0.0)
281    pub min_step_quality: f64,
282    /// Good step quality threshold (typically 0.75)
283    pub good_step_quality: f64,
284    /// Poor step quality threshold (typically 0.25)
285    pub poor_step_quality: f64,
286    /// Use Jacobi column scaling (preconditioning)
287    pub use_jacobi_scaling: bool,
288
289    // Ceres-style adaptive mu regularization parameters
290    /// Initial mu regularization parameter for Gauss-Newton step
291    pub initial_mu: f64,
292    /// Minimum mu regularization parameter
293    pub min_mu: f64,
294    /// Maximum mu regularization parameter
295    pub max_mu: f64,
296    /// Factor to increase mu when linear solver fails
297    pub mu_increase_factor: f64,
298
299    // Ceres-style step reuse optimization
300    /// Enable step reuse after rejection (Ceres-style efficiency optimization)
301    pub enable_step_reuse: bool,
302
303    /// Minimum objective function cutoff (optional early termination)
304    ///
305    /// If set, optimization terminates when cost falls below this threshold.
306    /// Useful for early stopping when a "good enough" solution is acceptable.
307    ///
308    /// Default: None (disabled)
309    pub min_cost_threshold: Option<f64>,
310
311    /// Maximum condition number for Jacobian matrix (optional check)
312    ///
313    /// If set, the optimizer checks if condition_number(J^T*J) exceeds this
314    /// threshold and terminates with IllConditionedJacobian status.
315    /// Note: Computing condition number is expensive, so this is disabled by default.
316    ///
317    /// Default: None (disabled)
318    pub max_condition_number: Option<f64>,
319
320    /// Minimum relative cost decrease for step acceptance
321    ///
322    /// Used in computing step quality (rho = actual_reduction / predicted_reduction).
323    /// Steps with rho < min_relative_decrease are rejected. Matches Ceres Solver's
324    /// min_relative_decrease parameter.
325    ///
326    /// Default: 1e-3 (Ceres-compatible)
327    pub min_relative_decrease: f64,
328
329    /// Compute per-variable covariance matrices (uncertainty estimation)
330    ///
331    /// When enabled, computes covariance by inverting the Hessian matrix after
332    /// convergence. The full covariance matrix is extracted into per-variable
333    /// blocks stored in both Variable structs and optimizer::SolverResult.
334    ///
335    /// Default: false (to avoid performance overhead)
336    pub compute_covariances: bool,
337
338    /// Enable real-time visualization (graphical debugging).
339    ///
340    /// When enabled, optimization progress is logged to a Rerun viewer.
341    /// **Note:** Requires the `visualization` feature to be enabled in `Cargo.toml`.
342    ///
343    /// Default: false
344    #[cfg(feature = "visualization")]
345    pub enable_visualization: bool,
346}
347
348impl Default for DogLegConfig {
349    fn default() -> Self {
350        Self {
351            linear_solver_type: linalg::LinearSolverType::default(),
352            // Ceres Solver default: 50 (changed from 100 for compatibility)
353            max_iterations: 50,
354            // Ceres Solver default: 1e-6 (changed from 1e-8 for compatibility)
355            cost_tolerance: 1e-6,
356            // Ceres Solver default: 1e-8 (unchanged)
357            parameter_tolerance: 1e-8,
358            // Ceres Solver default: 1e-10 (changed from 1e-8 for compatibility)
359            gradient_tolerance: 1e-10,
360            timeout: None,
361            // Ceres-style: larger initial radius for better global convergence
362            trust_region_radius: 1e4,
363            trust_region_min: 1e-12,
364            trust_region_max: 1e12,
365            // Ceres uses adaptive increase (max(radius, 3*step_norm)),
366            // but we keep factor for simpler config
367            trust_region_increase_factor: 3.0,
368            trust_region_decrease_factor: 0.5,
369            min_step_quality: 0.0,
370            good_step_quality: 0.75,
371            poor_step_quality: 0.25,
372            // Ceres-style: Enable diagonal scaling by default for elliptical trust region
373            use_jacobi_scaling: true,
374
375            // Ceres-style adaptive mu regularization defaults
376            // Start with more conservative regularization to avoid singular Hessian
377            initial_mu: 1e-4,
378            min_mu: 1e-8,
379            max_mu: 1.0,
380            mu_increase_factor: 10.0,
381
382            // Ceres-style step reuse optimization
383            enable_step_reuse: true,
384
385            // New Ceres-compatible termination parameters
386            min_cost_threshold: None,
387            max_condition_number: None,
388            min_relative_decrease: 1e-3,
389
390            compute_covariances: false,
391            #[cfg(feature = "visualization")]
392            enable_visualization: false,
393        }
394    }
395}
396
397impl DogLegConfig {
398    /// Create a new Dog Leg configuration with default values.
399    pub fn new() -> Self {
400        Self::default()
401    }
402
403    /// Set the linear solver type
404    pub fn with_linear_solver_type(mut self, linear_solver_type: linalg::LinearSolverType) -> Self {
405        self.linear_solver_type = linear_solver_type;
406        self
407    }
408
409    /// Set the maximum number of iterations
410    pub fn with_max_iterations(mut self, max_iterations: usize) -> Self {
411        self.max_iterations = max_iterations;
412        self
413    }
414
415    /// Set the cost tolerance
416    pub fn with_cost_tolerance(mut self, cost_tolerance: f64) -> Self {
417        self.cost_tolerance = cost_tolerance;
418        self
419    }
420
421    /// Set the parameter tolerance
422    pub fn with_parameter_tolerance(mut self, parameter_tolerance: f64) -> Self {
423        self.parameter_tolerance = parameter_tolerance;
424        self
425    }
426
427    /// Set the gradient tolerance
428    pub fn with_gradient_tolerance(mut self, gradient_tolerance: f64) -> Self {
429        self.gradient_tolerance = gradient_tolerance;
430        self
431    }
432
433    /// Set the timeout duration
434    pub fn with_timeout(mut self, timeout: time::Duration) -> Self {
435        self.timeout = Some(timeout);
436        self
437    }
438
439    /// Set the initial trust region radius
440    pub fn with_trust_region_radius(mut self, radius: f64) -> Self {
441        self.trust_region_radius = radius;
442        self
443    }
444
445    /// Set the trust region radius bounds
446    pub fn with_trust_region_bounds(mut self, min: f64, max: f64) -> Self {
447        self.trust_region_min = min;
448        self.trust_region_max = max;
449        self
450    }
451
452    /// Set the trust region adjustment factors
453    pub fn with_trust_region_factors(mut self, increase: f64, decrease: f64) -> Self {
454        self.trust_region_increase_factor = increase;
455        self.trust_region_decrease_factor = decrease;
456        self
457    }
458
459    /// Set the trust region quality thresholds
460    pub fn with_step_quality_thresholds(
461        mut self,
462        min_quality: f64,
463        poor_quality: f64,
464        good_quality: f64,
465    ) -> Self {
466        self.min_step_quality = min_quality;
467        self.poor_step_quality = poor_quality;
468        self.good_step_quality = good_quality;
469        self
470    }
471
472    /// Enable or disable Jacobi column scaling (preconditioning)
473    pub fn with_jacobi_scaling(mut self, use_jacobi_scaling: bool) -> Self {
474        self.use_jacobi_scaling = use_jacobi_scaling;
475        self
476    }
477
478    /// Set adaptive mu regularization parameters (Ceres-style)
479    pub fn with_mu_params(
480        mut self,
481        initial_mu: f64,
482        min_mu: f64,
483        max_mu: f64,
484        increase_factor: f64,
485    ) -> Self {
486        self.initial_mu = initial_mu;
487        self.min_mu = min_mu;
488        self.max_mu = max_mu;
489        self.mu_increase_factor = increase_factor;
490        self
491    }
492
493    /// Enable or disable step reuse optimization (Ceres-style)
494    pub fn with_step_reuse(mut self, enable_step_reuse: bool) -> Self {
495        self.enable_step_reuse = enable_step_reuse;
496        self
497    }
498
499    /// Set minimum objective function cutoff for early termination.
500    ///
501    /// When set, optimization terminates with MinCostThresholdReached status
502    /// if the cost falls below this threshold. Useful for early stopping when
503    /// a "good enough" solution is acceptable.
504    pub fn with_min_cost_threshold(mut self, min_cost: f64) -> Self {
505        self.min_cost_threshold = Some(min_cost);
506        self
507    }
508
509    /// Set maximum condition number for Jacobian matrix.
510    ///
511    /// If set, the optimizer checks if condition_number(J^T*J) exceeds this
512    /// threshold and terminates with IllConditionedJacobian status.
513    /// Note: Computing condition number is expensive, disabled by default.
514    pub fn with_max_condition_number(mut self, max_cond: f64) -> Self {
515        self.max_condition_number = Some(max_cond);
516        self
517    }
518
519    /// Set minimum relative cost decrease for step acceptance.
520    ///
521    /// Steps with rho = (actual_reduction / predicted_reduction) below this
522    /// threshold are rejected. Default: 1e-3 (Ceres-compatible)
523    pub fn with_min_relative_decrease(mut self, min_decrease: f64) -> Self {
524        self.min_relative_decrease = min_decrease;
525        self
526    }
527
528    /// Enable or disable covariance computation (uncertainty estimation).
529    ///
530    /// When enabled, computes the full covariance matrix by inverting the Hessian
531    /// after convergence, then extracts per-variable covariance blocks.
532    pub fn with_compute_covariances(mut self, compute_covariances: bool) -> Self {
533        self.compute_covariances = compute_covariances;
534        self
535    }
536
537    /// Enable real-time visualization.
538    ///
539    /// **Note:** Requires the `visualization` feature to be enabled in `Cargo.toml`.
540    ///
541    /// # Arguments
542    ///
543    /// * `enable` - Whether to enable visualization
544    #[cfg(feature = "visualization")]
545    pub fn with_visualization(mut self, enable: bool) -> Self {
546        self.enable_visualization = enable;
547        self
548    }
549
550    /// Print configuration parameters (info level logging)
551    pub fn print_configuration(&self) {
552        debug!(
553            "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: {}",
554            self.linear_solver_type,
555            self.max_iterations,
556            self.cost_tolerance,
557            self.parameter_tolerance,
558            self.gradient_tolerance,
559            self.timeout,
560            self.trust_region_radius,
561            self.trust_region_min,
562            self.trust_region_max,
563            self.min_step_quality,
564            self.good_step_quality,
565            self.poor_step_quality,
566            self.initial_mu,
567            self.min_mu,
568            self.max_mu,
569            self.mu_increase_factor,
570            if self.use_jacobi_scaling {
571                "enabled"
572            } else {
573                "disabled"
574            },
575            if self.enable_step_reuse {
576                "enabled"
577            } else {
578                "disabled"
579            },
580            if self.compute_covariances {
581                "enabled"
582            } else {
583                "disabled"
584            }
585        );
586    }
587}
588
589/// Result from step computation
590struct StepResult {
591    step: faer::Mat<f64>,
592    gradient_norm: f64,
593    predicted_reduction: f64,
594}
595
596/// Type of step taken
597#[derive(Debug, Clone, Copy)]
598enum StepType {
599    /// Full Gauss-Newton step
600    GaussNewton,
601    /// Scaled steepest descent (Cauchy point)
602    SteepestDescent,
603    /// Dog leg interpolation
604    DogLeg,
605}
606
607impl fmt::Display for StepType {
608    fn fmt(&self, f: &mut fmt::Formatter<'_>) -> fmt::Result {
609        match self {
610            StepType::GaussNewton => write!(f, "GN"),
611            StepType::SteepestDescent => write!(f, "SD"),
612            StepType::DogLeg => write!(f, "DL"),
613        }
614    }
615}
616
617/// Result from step evaluation
618struct StepEvaluation {
619    accepted: bool,
620    cost_reduction: f64,
621    rho: f64,
622}
623
624/// Dog Leg trust region solver for nonlinear least squares optimization.
625///
626/// Implements Powell's Dog Leg algorithm with Ceres Solver enhancements including
627/// adaptive μ regularization, numerically robust beta computation, and step reuse caching.
628///
629/// # Algorithm
630///
631/// At each iteration k:
632/// 1. Compute residual `r(xₖ)` and Jacobian `J(xₖ)`
633/// 2. Solve for Gauss-Newton step: `(J^T·J + μI)·h_gn = -J^T·r`
634/// 3. Compute steepest descent direction: `-g` where `g = J^T·r`
635/// 4. Compute Cauchy point: `p_c = -α·g` (optimal along steepest descent)
636/// 5. Construct dog leg step based on trust region radius Δ:
637///    - If `||h_gn|| ≤ Δ`: Take full GN step
638///    - Else if `||p_c|| ≥ Δ`: Take scaled SD to boundary
639///    - Else: Interpolate `h = p_c + β·(h_gn - p_c)` where `||h|| = Δ`
640/// 6. Evaluate gain ratio: `ρ = (actual reduction) / (predicted reduction)`
641/// 7. Update trust region radius based on ρ
642/// 8. Accept/reject step and update parameters
643///
644/// # Ceres Solver Enhancements
645///
646/// This implementation includes four major improvements from Google's Ceres Solver:
647///
648/// **1. Adaptive μ Regularization:** Dynamically adjusts regularization parameter
649/// to handle ill-conditioned Hessians (increases on failure, decreases on success).
650///
651/// **2. Numerically Robust Beta:** Uses two formulas for computing dog leg
652/// interpolation parameter β to avoid catastrophic cancellation.
653///
654/// **3. Step Reuse Mechanism:** Caches GN step, Cauchy point, and gradient when
655/// steps are rejected. Limited to 5 consecutive reuses to prevent staleness.
656///
657/// **4. Jacobi Scaling:** Optional diagonal preconditioning creates elliptical
658/// trust regions for better handling of mixed-scale problems.
659///
660/// # Examples
661///
662/// ```no_run
663/// use apex_solver::optimizer::DogLeg;
664/// use apex_solver::core::problem::Problem;
665/// use std::collections::HashMap;
666///
667/// # fn main() -> Result<(), Box<dyn std::error::Error>> {
668/// let mut problem = Problem::new();
669/// // ... add factors to problem ...
670///
671/// let initial_values = HashMap::new();
672/// // ... initialize parameters ...
673///
674/// let mut solver = DogLeg::new();
675/// let result = solver.optimize(&problem, &initial_values)?;
676/// # Ok(())
677/// # }
678/// ```
679///
680/// # See Also
681///
682/// - [`DogLegConfig`] - Configuration options
683/// - [`LevenbergMarquardt`](crate::optimizer::LevenbergMarquardt) - Alternative adaptive damping
684/// - [`GaussNewton`](crate::optimizer::GaussNewton) - Undamped variant
685pub struct DogLeg {
686    config: DogLegConfig,
687    jacobi_scaling: Option<sparse::SparseColMat<usize, f64>>,
688    observers: optimizer::OptObserverVec,
689
690    // Adaptive mu regularization (Ceres-style)
691    mu: f64,
692    min_mu: f64,
693    max_mu: f64,
694    mu_increase_factor: f64,
695
696    // Step reuse mechanism (Ceres-style efficiency optimization)
697    reuse_step_on_rejection: bool,
698    cached_gn_step: Option<faer::Mat<f64>>,
699    cached_cauchy_point: Option<faer::Mat<f64>>,
700    cached_gradient: Option<faer::Mat<f64>>,
701    cached_alpha: Option<f64>,
702    cache_reuse_count: usize, // Track consecutive reuses to prevent staleness
703}
704
705impl Default for DogLeg {
706    fn default() -> Self {
707        Self::new()
708    }
709}
710
711impl DogLeg {
712    /// Create a new Dog Leg solver with default configuration.
713    pub fn new() -> Self {
714        Self::with_config(DogLegConfig::default())
715    }
716
717    /// Create a new Dog Leg solver with the given configuration.
718    pub fn with_config(config: DogLegConfig) -> Self {
719        Self {
720            // Initialize adaptive mu from config
721            mu: config.initial_mu,
722            min_mu: config.min_mu,
723            max_mu: config.max_mu,
724            mu_increase_factor: config.mu_increase_factor,
725
726            // Initialize step reuse mechanism (disabled initially, enabled after first rejection)
727            reuse_step_on_rejection: false,
728            cached_gn_step: None,
729            cached_cauchy_point: None,
730            cached_gradient: None,
731            cached_alpha: None,
732            cache_reuse_count: 0,
733
734            config,
735            jacobi_scaling: None,
736            observers: optimizer::OptObserverVec::new(),
737        }
738    }
739
740    /// Add an observer to the solver.
741    ///
742    /// Observers are notified at each iteration with the current variable values.
743    /// This enables real-time visualization, logging, metrics collection, etc.
744    ///
745    /// # Examples
746    ///
747    /// ```no_run
748    /// use apex_solver::optimizer::DogLeg;
749    /// # use apex_solver::optimizer::OptObserver;
750    /// # use std::collections::HashMap;
751    /// # use apex_solver::core::problem::VariableEnum;
752    ///
753    /// # struct MyObserver;
754    /// # impl OptObserver for MyObserver {
755    /// #     fn on_step(&self, _: &HashMap<String, VariableEnum>, _: usize) {}
756    /// # }
757    /// let mut solver = DogLeg::new();
758    /// solver.add_observer(MyObserver);
759    /// ```
760    pub fn add_observer(&mut self, observer: impl optimizer::OptObserver + 'static) {
761        self.observers.add(observer);
762    }
763
764    /// Compute the Cauchy point (steepest descent step)
765    /// Returns the optimal step along the negative gradient direction
766    /// Compute Cauchy point and optimal step length for steepest descent
767    ///
768    /// Returns (alpha, cauchy_point) where:
769    /// - alpha: optimal step length α = ||g||² / (g^T H g)
770    /// - cauchy_point: p_c = -α * g (the Cauchy point)
771    ///
772    /// This is the optimal point along the steepest descent direction within
773    /// the quadratic approximation of the objective function.
774    fn compute_cauchy_point_and_alpha(
775        &self,
776        gradient: &faer::Mat<f64>,
777        hessian: &sparse::SparseColMat<usize, f64>,
778    ) -> (f64, faer::Mat<f64>) {
779        // Optimal step size along steepest descent: α = (g^T*g) / (g^T*H*g)
780        let g_norm_sq_mat = gradient.transpose() * gradient;
781        let g_norm_sq = g_norm_sq_mat[(0, 0)];
782
783        let h_g = hessian * gradient;
784        let g_h_g_mat = gradient.transpose() * &h_g;
785        let g_h_g = g_h_g_mat[(0, 0)];
786
787        // Avoid division by zero
788        let alpha = if g_h_g.abs() > 1e-15 {
789            g_norm_sq / g_h_g
790        } else {
791            1.0
792        };
793
794        // Compute Cauchy point: p_c = -α * gradient
795        let mut cauchy_point = faer::Mat::zeros(gradient.nrows(), 1);
796        for i in 0..gradient.nrows() {
797            cauchy_point[(i, 0)] = -alpha * gradient[(i, 0)];
798        }
799
800        (alpha, cauchy_point)
801    }
802
803    /// Compute the dog leg step using Powell's Dog Leg method
804    ///
805    /// The dog leg path consists of two segments:
806    /// 1. From origin to Cauchy point (optimal along steepest descent)
807    /// 2. From Cauchy point to Gauss-Newton step
808    ///
809    /// Arguments:
810    /// - steepest_descent_dir: -gradient (steepest descent direction, not scaled)
811    /// - cauchy_point: p_c = α * (-gradient), the optimal steepest descent step
812    /// - h_gn: Gauss-Newton step
813    /// - delta: Trust region radius
814    ///
815    /// Returns: (step, step_type)
816    fn compute_dog_leg_step(
817        &self,
818        steepest_descent_dir: &faer::Mat<f64>,
819        cauchy_point: &faer::Mat<f64>,
820        h_gn: &faer::Mat<f64>,
821        delta: f64,
822    ) -> (faer::Mat<f64>, StepType) {
823        let gn_norm = h_gn.norm_l2();
824        let cauchy_norm = cauchy_point.norm_l2();
825        let sd_norm = steepest_descent_dir.norm_l2();
826
827        // Case 1: Full Gauss-Newton step fits in trust region
828        if gn_norm <= delta {
829            return (h_gn.clone(), StepType::GaussNewton);
830        }
831
832        // Case 2: Even Cauchy point is outside trust region
833        // Scale steepest descent direction to boundary: (delta / ||δ_sd||) * δ_sd
834        if cauchy_norm >= delta {
835            let mut scaled_sd = faer::Mat::zeros(steepest_descent_dir.nrows(), 1);
836            let scale = delta / sd_norm;
837            for i in 0..steepest_descent_dir.nrows() {
838                scaled_sd[(i, 0)] = steepest_descent_dir[(i, 0)] * scale;
839            }
840            return (scaled_sd, StepType::SteepestDescent);
841        }
842
843        // Case 3: Dog leg interpolation between Cauchy point and GN step
844        // Use Ceres-style numerically robust formula
845        //
846        // Following Ceres solver implementation for numerical stability:
847        // Compute intersection of trust region boundary with line from Cauchy point to GN step
848        //
849        // Let v = δ_gn - p_c
850        // Solve: ||p_c + β*v||² = delta²
851        // This gives: a*β² + 2*b*β + c = 0
852        // where:
853        //   a = v^T·v = ||v||²
854        //   b = p_c^T·v
855        //   c = p_c^T·p_c - delta² = ||p_c||² - delta²
856
857        let mut v = faer::Mat::zeros(cauchy_point.nrows(), 1);
858        for i in 0..cauchy_point.nrows() {
859            v[(i, 0)] = h_gn[(i, 0)] - cauchy_point[(i, 0)];
860        }
861
862        // Compute coefficients
863        let v_squared_norm = v.transpose() * &v;
864        let a = v_squared_norm[(0, 0)];
865
866        let pc_dot_v = cauchy_point.transpose() * &v;
867        let b = pc_dot_v[(0, 0)];
868
869        let c = cauchy_norm * cauchy_norm - delta * delta;
870
871        // Ceres-style numerically robust beta computation
872        // Uses two different formulas based on sign of b to avoid catastrophic cancellation
873        let d_squared = b * b - a * c;
874
875        let beta = if d_squared < 0.0 {
876            // Should not happen geometrically, but handle gracefully
877            1.0
878        } else if a.abs() < 1e-15 {
879            // Degenerate case: v is nearly zero
880            1.0
881        } else {
882            let d = d_squared.sqrt();
883
884            // Ceres formula: choose formula based on sign of b to avoid cancellation
885            // If b <= 0: beta = (-b + d) / a  (standard formula, no cancellation)
886            // If b > 0:  beta = -c / (b + d)  (alternative formula, avoids cancellation)
887            if b <= 0.0 { (-b + d) / a } else { -c / (b + d) }
888        };
889
890        // Clamp beta to [0, 1] for safety
891        let beta = beta.clamp(0.0, 1.0);
892
893        // Compute dog leg step: p_dl = p_c + β*(δ_gn - p_c)
894        let mut dog_leg = faer::Mat::zeros(cauchy_point.nrows(), 1);
895        for i in 0..cauchy_point.nrows() {
896            dog_leg[(i, 0)] = cauchy_point[(i, 0)] + beta * v[(i, 0)];
897        }
898
899        (dog_leg, StepType::DogLeg)
900    }
901
902    /// Update trust region radius based on step quality (Ceres-style)
903    fn update_trust_region(&mut self, rho: f64, step_norm: f64) -> bool {
904        if rho > self.config.good_step_quality {
905            // Good step, increase trust region (Ceres-style: max(radius, 3*step_norm))
906            let new_radius = self.config.trust_region_radius.max(3.0 * step_norm);
907            self.config.trust_region_radius = new_radius.min(self.config.trust_region_max);
908
909            // Decrease mu on successful step (Ceres-style adaptive regularization)
910            self.mu = (self.mu / (0.5 * self.mu_increase_factor)).max(self.min_mu);
911
912            // Clear reuse flag and invalidate cache on acceptance (parameters have moved)
913            self.reuse_step_on_rejection = false;
914            self.cached_gn_step = None;
915            self.cached_cauchy_point = None;
916            self.cached_gradient = None;
917            self.cached_alpha = None;
918            self.cache_reuse_count = 0;
919
920            true
921        } else if rho < self.config.poor_step_quality {
922            // Poor step, decrease trust region
923            self.config.trust_region_radius = (self.config.trust_region_radius
924                * self.config.trust_region_decrease_factor)
925                .max(self.config.trust_region_min);
926
927            // Enable step reuse flag for next iteration (Ceres-style)
928            self.reuse_step_on_rejection = self.config.enable_step_reuse;
929
930            false
931        } else {
932            // Moderate step, keep trust region unchanged
933            // Clear reuse flag and invalidate cache on acceptance (parameters have moved)
934            self.reuse_step_on_rejection = false;
935            self.cached_gn_step = None;
936            self.cached_cauchy_point = None;
937            self.cached_gradient = None;
938            self.cached_alpha = None;
939            self.cache_reuse_count = 0;
940
941            true
942        }
943    }
944
945    /// Compute predicted cost reduction from linear model
946    fn compute_predicted_reduction(
947        &self,
948        step: &faer::Mat<f64>,
949        gradient: &faer::Mat<f64>,
950        hessian: &sparse::SparseColMat<usize, f64>,
951    ) -> f64 {
952        // Dog Leg predicted reduction: -step^T * gradient - 0.5 * step^T * H * step
953        let linear_term = step.transpose() * gradient;
954        let hessian_step = hessian * step;
955        let quadratic_term = step.transpose() * &hessian_step;
956
957        -linear_term[(0, 0)] - 0.5 * quadratic_term[(0, 0)]
958    }
959
960    /// Compute dog leg optimization step
961    fn compute_optimization_step(
962        &mut self,
963        residuals: &faer::Mat<f64>,
964        scaled_jacobian: &sparse::SparseColMat<usize, f64>,
965        linear_solver: &mut Box<dyn linalg::SparseLinearSolver>,
966    ) -> Option<StepResult> {
967        // Check if we can reuse cached step (Ceres-style optimization)
968        // Safety limit: prevent excessive reuse that could lead to stale gradient/Hessian
969        const MAX_CACHE_REUSE: usize = 5;
970
971        if self.reuse_step_on_rejection
972            && self.config.enable_step_reuse
973            && self.cache_reuse_count < MAX_CACHE_REUSE
974            && let (Some(cached_gn), Some(cached_cauchy), Some(cached_grad), Some(_cached_a)) = (
975                &self.cached_gn_step,
976                &self.cached_cauchy_point,
977                &self.cached_gradient,
978                &self.cached_alpha,
979            )
980        {
981            // Increment reuse counter
982            self.cache_reuse_count += 1;
983
984            let gradient_norm = cached_grad.norm_l2();
985            let mut steepest_descent = faer::Mat::zeros(cached_grad.nrows(), 1);
986            for i in 0..cached_grad.nrows() {
987                steepest_descent[(i, 0)] = -cached_grad[(i, 0)];
988            }
989
990            let (scaled_step, _step_type) = self.compute_dog_leg_step(
991                &steepest_descent,
992                cached_cauchy,
993                cached_gn,
994                self.config.trust_region_radius,
995            );
996
997            let step = if self.config.use_jacobi_scaling {
998                let scaling = self.jacobi_scaling.as_ref()?;
999                scaling * &scaled_step
1000            } else {
1001                scaled_step.clone()
1002            };
1003
1004            let hessian = linear_solver.get_hessian()?;
1005            let predicted_reduction =
1006                self.compute_predicted_reduction(&scaled_step, cached_grad, hessian);
1007
1008            return Some(StepResult {
1009                step,
1010                gradient_norm,
1011                predicted_reduction,
1012            });
1013        }
1014
1015        // Not reusing, compute fresh step
1016        // 1. Solve for Gauss-Newton step with adaptive mu regularization (Ceres-style)
1017        let residuals_owned = residuals.as_ref().to_owned();
1018        let mut scaled_gn_step = None;
1019        let mut mu_attempts = 0;
1020
1021        // Try to solve with current mu, increasing if necessary
1022        while mu_attempts < 10 && self.mu <= self.max_mu {
1023            let damping = self.mu;
1024
1025            if let Ok(step) =
1026                linear_solver.solve_augmented_equation(&residuals_owned, scaled_jacobian, damping)
1027            {
1028                scaled_gn_step = Some(step);
1029                break;
1030            }
1031
1032            // Increase mu (Ceres-style)
1033            self.mu = (self.mu * self.mu_increase_factor).min(self.max_mu);
1034            mu_attempts += 1;
1035        }
1036
1037        let scaled_gn_step = scaled_gn_step?;
1038
1039        // 2. Get gradient and Hessian (cached by solve_augmented_equation)
1040        let gradient = linear_solver.get_gradient()?;
1041        let hessian = linear_solver.get_hessian()?;
1042        let gradient_norm = gradient.norm_l2();
1043
1044        // 3. Compute steepest descent direction: δ_sd = -gradient
1045        let mut steepest_descent = faer::Mat::zeros(gradient.nrows(), 1);
1046        for i in 0..gradient.nrows() {
1047            steepest_descent[(i, 0)] = -gradient[(i, 0)];
1048        }
1049
1050        // 4. Compute Cauchy point and optimal step length
1051        let (alpha, cauchy_point) = self.compute_cauchy_point_and_alpha(gradient, hessian);
1052
1053        // 5. Compute dog leg step based on trust region radius
1054        let (scaled_step, _step_type) = self.compute_dog_leg_step(
1055            &steepest_descent,
1056            &cauchy_point,
1057            &scaled_gn_step,
1058            self.config.trust_region_radius,
1059        );
1060
1061        // 6. Apply inverse Jacobi scaling if enabled
1062        let step = if self.config.use_jacobi_scaling {
1063            let scaling = self.jacobi_scaling.as_ref()?;
1064            scaling * &scaled_step
1065        } else {
1066            scaled_step.clone()
1067        };
1068
1069        // 7. Compute predicted reduction
1070        let predicted_reduction = self.compute_predicted_reduction(&scaled_step, gradient, hessian);
1071
1072        // 8. Cache step components for potential reuse (Ceres-style)
1073        self.cached_gn_step = Some(scaled_gn_step.clone());
1074        self.cached_cauchy_point = Some(cauchy_point.clone());
1075        self.cached_gradient = Some(gradient.clone());
1076        self.cached_alpha = Some(alpha);
1077
1078        Some(StepResult {
1079            step,
1080            gradient_norm,
1081            predicted_reduction,
1082        })
1083    }
1084
1085    /// Evaluate and apply step
1086    fn evaluate_and_apply_step(
1087        &mut self,
1088        step_result: &StepResult,
1089        state: &mut optimizer::InitializedState,
1090        problem: &problem::Problem,
1091    ) -> error::ApexSolverResult<StepEvaluation> {
1092        // Apply parameter updates
1093        let step_norm = optimizer::apply_parameter_step(
1094            &mut state.variables,
1095            step_result.step.as_ref(),
1096            &state.sorted_vars,
1097        );
1098
1099        // Compute new cost (residual only, no Jacobian needed for step evaluation)
1100        let new_residual = problem.compute_residual_sparse(&state.variables)?;
1101        let new_residual_norm = new_residual.norm_l2();
1102        let new_cost = new_residual_norm * new_residual_norm;
1103
1104        // Compute step quality
1105        let rho = optimizer::compute_step_quality(
1106            state.current_cost,
1107            new_cost,
1108            step_result.predicted_reduction,
1109        );
1110
1111        // Update trust region and decide acceptance
1112        // Filter out numerical noise with small threshold
1113        let accepted = rho > 1e-4;
1114        let _trust_region_updated = self.update_trust_region(rho, step_norm);
1115
1116        let cost_reduction = if accepted {
1117            let reduction = state.current_cost - new_cost;
1118            state.current_cost = new_cost;
1119            reduction
1120        } else {
1121            // Reject step - revert changes
1122            optimizer::apply_negative_parameter_step(
1123                &mut state.variables,
1124                step_result.step.as_ref(),
1125                &state.sorted_vars,
1126            );
1127            0.0
1128        };
1129
1130        Ok(StepEvaluation {
1131            accepted,
1132            cost_reduction,
1133            rho,
1134        })
1135    }
1136
1137    /// Minimize the optimization problem using Dog Leg algorithm
1138    pub fn optimize(
1139        &mut self,
1140        problem: &problem::Problem,
1141        initial_params: &collections::HashMap<
1142            String,
1143            (manifold::ManifoldType, nalgebra::DVector<f64>),
1144        >,
1145    ) -> Result<
1146        optimizer::SolverResult<collections::HashMap<String, problem::VariableEnum>>,
1147        error::ApexSolverError,
1148    > {
1149        let start_time = time::Instant::now();
1150        let mut iteration = 0;
1151        let mut cost_evaluations = 1;
1152        let mut jacobian_evaluations = 0;
1153        let mut successful_steps = 0;
1154        let mut unsuccessful_steps = 0;
1155
1156        let mut state = optimizer::initialize_optimization_state(problem, initial_params)?;
1157        let mut linear_solver = optimizer::create_linear_solver(&self.config.linear_solver_type);
1158
1159        let mut max_gradient_norm: f64 = 0.0;
1160        let mut max_parameter_update_norm: f64 = 0.0;
1161        let mut total_cost_reduction = 0.0;
1162        let mut final_gradient_norm;
1163        let mut final_parameter_update_norm;
1164
1165        // Initialize iteration statistics tracking
1166        let mut iteration_stats = Vec::with_capacity(self.config.max_iterations);
1167        let mut previous_cost = state.current_cost;
1168
1169        // Print configuration and header if info/debug level is enabled
1170        if tracing::enabled!(tracing::Level::DEBUG) {
1171            self.config.print_configuration();
1172            IterationStats::print_header();
1173        }
1174
1175        loop {
1176            let iter_start = time::Instant::now();
1177            // Evaluate residuals and Jacobian
1178            let (residuals, jacobian) = problem.compute_residual_and_jacobian_sparse(
1179                &state.variables,
1180                &state.variable_index_map,
1181                &state.symbolic_structure,
1182            )?;
1183            jacobian_evaluations += 1;
1184
1185            // Process Jacobian (apply scaling if enabled)
1186            let scaled_jacobian = if self.config.use_jacobi_scaling {
1187                optimizer::process_jacobian(&jacobian, &mut self.jacobi_scaling, iteration)?
1188            } else {
1189                jacobian
1190            };
1191
1192            // Compute dog leg step
1193            let step_result = match self.compute_optimization_step(
1194                &residuals,
1195                &scaled_jacobian,
1196                &mut linear_solver,
1197            ) {
1198                Some(result) => result,
1199                None => {
1200                    return Err(optimizer::OptimizerError::LinearSolveFailed(
1201                        "Linear solver failed to solve system".to_string(),
1202                    )
1203                    .into());
1204                }
1205            };
1206
1207            // Update tracking
1208            max_gradient_norm = max_gradient_norm.max(step_result.gradient_norm);
1209            final_gradient_norm = step_result.gradient_norm;
1210            let step_norm = step_result.step.norm_l2();
1211            max_parameter_update_norm = max_parameter_update_norm.max(step_norm);
1212            final_parameter_update_norm = step_norm;
1213
1214            // Evaluate and apply step
1215            let step_eval = self.evaluate_and_apply_step(&step_result, &mut state, problem)?;
1216            cost_evaluations += 1;
1217
1218            if step_eval.accepted {
1219                successful_steps += 1;
1220                total_cost_reduction += step_eval.cost_reduction;
1221            } else {
1222                unsuccessful_steps += 1;
1223            }
1224
1225            // OPTIMIZATION: Only collect iteration statistics if debug level is enabled
1226            // This eliminates ~2-5ms overhead per iteration for non-debug optimization
1227            if tracing::enabled!(tracing::Level::DEBUG) {
1228                let iter_elapsed_ms = iter_start.elapsed().as_secs_f64() * 1000.0;
1229                let total_elapsed_ms = start_time.elapsed().as_secs_f64() * 1000.0;
1230
1231                let stats = IterationStats {
1232                    iteration,
1233                    cost: state.current_cost,
1234                    cost_change: previous_cost - state.current_cost,
1235                    gradient_norm: step_result.gradient_norm,
1236                    step_norm,
1237                    tr_ratio: step_eval.rho,
1238                    tr_radius: self.config.trust_region_radius,
1239                    ls_iter: 0, // Direct solver (Cholesky) has no iterations
1240                    iter_time_ms: iter_elapsed_ms,
1241                    total_time_ms: total_elapsed_ms,
1242                    accepted: step_eval.accepted,
1243                };
1244
1245                iteration_stats.push(stats.clone());
1246                stats.print_line();
1247            }
1248
1249            previous_cost = state.current_cost;
1250
1251            // Notify all observers with current state
1252            optimizer::notify_observers(
1253                &mut self.observers,
1254                &state.variables,
1255                iteration,
1256                state.current_cost,
1257                step_result.gradient_norm,
1258                Some(self.config.trust_region_radius),
1259                step_norm,
1260                Some(step_eval.rho),
1261                linear_solver.as_ref(),
1262            );
1263
1264            // Check convergence
1265            let elapsed = start_time.elapsed();
1266
1267            // Compute parameter norm for relative parameter tolerance check
1268            let parameter_norm = optimizer::compute_parameter_norm(&state.variables);
1269
1270            // Compute costs for convergence check
1271            let new_cost = state.current_cost;
1272
1273            let cost_before_step = if step_eval.accepted {
1274                state.current_cost + step_eval.cost_reduction
1275            } else {
1276                state.current_cost
1277            };
1278
1279            if let Some(status) = optimizer::check_convergence(&optimizer::ConvergenceParams {
1280                iteration,
1281                current_cost: cost_before_step,
1282                new_cost,
1283                parameter_norm,
1284                parameter_update_norm: step_norm,
1285                gradient_norm: step_result.gradient_norm,
1286                elapsed,
1287                step_accepted: step_eval.accepted,
1288                max_iterations: self.config.max_iterations,
1289                gradient_tolerance: self.config.gradient_tolerance,
1290                parameter_tolerance: self.config.parameter_tolerance,
1291                cost_tolerance: self.config.cost_tolerance,
1292                min_cost_threshold: self.config.min_cost_threshold,
1293                timeout: self.config.timeout,
1294                trust_region_radius: Some(self.config.trust_region_radius),
1295                min_trust_region_radius: Some(self.config.trust_region_min),
1296            }) {
1297                // Print summary only if debug level is enabled
1298                if tracing::enabled!(tracing::Level::DEBUG) {
1299                    let summary = optimizer::create_optimizer_summary(
1300                        "Dog-Leg",
1301                        state.initial_cost,
1302                        state.current_cost,
1303                        iteration + 1,
1304                        Some(successful_steps),
1305                        Some(unsuccessful_steps),
1306                        max_gradient_norm,
1307                        final_gradient_norm,
1308                        max_parameter_update_norm,
1309                        final_parameter_update_norm,
1310                        total_cost_reduction,
1311                        elapsed,
1312                        iteration_stats.clone(),
1313                        status.clone(),
1314                        None,
1315                        Some(self.config.trust_region_radius),
1316                        None,
1317                    );
1318                    debug!("{}", summary);
1319                }
1320
1321                // Compute covariances if enabled
1322                let covariances = if self.config.compute_covariances {
1323                    problem.compute_and_set_covariances(
1324                        &mut linear_solver,
1325                        &mut state.variables,
1326                        &state.variable_index_map,
1327                    )
1328                } else {
1329                    None
1330                };
1331
1332                return Ok(optimizer::build_solver_result(
1333                    status,
1334                    iteration + 1,
1335                    state,
1336                    elapsed,
1337                    final_gradient_norm,
1338                    final_parameter_update_norm,
1339                    cost_evaluations,
1340                    jacobian_evaluations,
1341                    covariances,
1342                ));
1343            }
1344
1345            // Note: Max iterations and timeout checks are now handled inside check_convergence()
1346
1347            iteration += 1;
1348        }
1349    }
1350}
1351
1352impl optimizer::Solver for DogLeg {
1353    type Config = DogLegConfig;
1354    type Error = error::ApexSolverError;
1355
1356    fn new() -> Self {
1357        Self::default()
1358    }
1359
1360    fn optimize(
1361        &mut self,
1362        problem: &problem::Problem,
1363        initial_params: &collections::HashMap<
1364            String,
1365            (manifold::ManifoldType, nalgebra::DVector<f64>),
1366        >,
1367    ) -> Result<
1368        optimizer::SolverResult<collections::HashMap<String, problem::VariableEnum>>,
1369        Self::Error,
1370    > {
1371        self.optimize(problem, initial_params)
1372    }
1373}
1374
1375#[cfg(test)]
1376mod tests {
1377    use super::*;
1378    use crate::factors;
1379    use apex_manifolds as manifold;
1380    use nalgebra;
1381
1382    /// Custom Rosenbrock Factor 1: r1 = 10(x2 - x1²)
1383    /// Demonstrates extensibility - custom factors can be defined outside of factors.rs
1384    #[derive(Debug, Clone)]
1385    struct RosenbrockFactor1;
1386
1387    impl factors::Factor for RosenbrockFactor1 {
1388        fn linearize(
1389            &self,
1390            params: &[nalgebra::DVector<f64>],
1391            compute_jacobian: bool,
1392        ) -> (nalgebra::DVector<f64>, Option<nalgebra::DMatrix<f64>>) {
1393            let x1 = params[0][0];
1394            let x2 = params[1][0];
1395
1396            // Residual: r1 = 10(x2 - x1²)
1397            let residual = nalgebra::dvector![10.0 * (x2 - x1 * x1)];
1398
1399            // Jacobian: ∂r1/∂x1 = -20*x1, ∂r1/∂x2 = 10
1400            let jacobian = if compute_jacobian {
1401                let mut jac = nalgebra::DMatrix::zeros(1, 2);
1402                jac[(0, 0)] = -20.0 * x1;
1403                jac[(0, 1)] = 10.0;
1404                Some(jac)
1405            } else {
1406                None
1407            };
1408
1409            (residual, jacobian)
1410        }
1411
1412        fn get_dimension(&self) -> usize {
1413            1
1414        }
1415    }
1416
1417    /// Custom Rosenbrock Factor 2: r2 = 1 - x1
1418    /// Demonstrates extensibility - custom factors can be defined outside of factors.rs
1419    #[derive(Debug, Clone)]
1420    struct RosenbrockFactor2;
1421
1422    impl factors::Factor for RosenbrockFactor2 {
1423        fn linearize(
1424            &self,
1425            params: &[nalgebra::DVector<f64>],
1426            compute_jacobian: bool,
1427        ) -> (nalgebra::DVector<f64>, Option<nalgebra::DMatrix<f64>>) {
1428            let x1 = params[0][0];
1429
1430            // Residual: r2 = 1 - x1
1431            let residual = nalgebra::dvector![1.0 - x1];
1432
1433            // Jacobian: ∂r2/∂x1 = -1
1434            let jacobian = if compute_jacobian {
1435                Some(nalgebra::DMatrix::from_element(1, 1, -1.0))
1436            } else {
1437                None
1438            };
1439
1440            (residual, jacobian)
1441        }
1442
1443        fn get_dimension(&self) -> usize {
1444            1
1445        }
1446    }
1447
1448    #[test]
1449    fn test_rosenbrock_optimization() -> Result<(), Box<dyn std::error::Error>> {
1450        // Rosenbrock function test:
1451        // Minimize: r1² + r2² where
1452        //   r1 = 10(x2 - x1²)
1453        //   r2 = 1 - x1
1454        // Starting point: [-1.2, 1.0]
1455        // Expected minimum: [1.0, 1.0]
1456
1457        let mut problem = problem::Problem::new();
1458        let mut initial_values = collections::HashMap::new();
1459
1460        // Add variables using Rn manifold (Euclidean space)
1461        initial_values.insert(
1462            "x1".to_string(),
1463            (manifold::ManifoldType::RN, nalgebra::dvector![-1.2]),
1464        );
1465        initial_values.insert(
1466            "x2".to_string(),
1467            (manifold::ManifoldType::RN, nalgebra::dvector![1.0]),
1468        );
1469
1470        // Add custom factors (demonstrates extensibility!)
1471        problem.add_residual_block(&["x1", "x2"], Box::new(RosenbrockFactor1), None);
1472        problem.add_residual_block(&["x1"], Box::new(RosenbrockFactor2), None);
1473
1474        // Configure Dog Leg optimizer with appropriate trust region
1475        let config = DogLegConfig::new()
1476            .with_max_iterations(100)
1477            .with_cost_tolerance(1e-8)
1478            .with_parameter_tolerance(1e-8)
1479            .with_gradient_tolerance(1e-10)
1480            .with_trust_region_radius(10.0); // Start with larger trust region
1481
1482        let mut solver = DogLeg::with_config(config);
1483        let result = solver.optimize(&problem, &initial_values)?;
1484
1485        // Extract final values
1486        let x1_final = result
1487            .parameters
1488            .get("x1")
1489            .ok_or("x1 not found")?
1490            .to_vector()[0];
1491        let x2_final = result
1492            .parameters
1493            .get("x2")
1494            .ok_or("x2 not found")?
1495            .to_vector()[0];
1496
1497        // Verify convergence to [1.0, 1.0]
1498        assert!(
1499            matches!(
1500                result.status,
1501                optimizer::OptimizationStatus::Converged
1502                    | optimizer::OptimizationStatus::CostToleranceReached
1503                    | optimizer::OptimizationStatus::ParameterToleranceReached
1504                    | optimizer::OptimizationStatus::GradientToleranceReached
1505            ),
1506            "Optimization should converge"
1507        );
1508        assert!(
1509            (x1_final - 1.0).abs() < 1e-4,
1510            "x1 should converge to 1.0, got {}",
1511            x1_final
1512        );
1513        assert!(
1514            (x2_final - 1.0).abs() < 1e-4,
1515            "x2 should converge to 1.0, got {}",
1516            x2_final
1517        );
1518        assert!(
1519            result.final_cost < 1e-6,
1520            "Final cost should be near zero, got {}",
1521            result.final_cost
1522        );
1523        Ok(())
1524    }
1525}