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, manifold, optimizer};
197
198// Note: Visualization support via observer pattern will be added in a future update
199// For now, use LevenbergMarquardt optimizer for visualization support
200
201use faer::sparse;
202use std::{collections, fmt, time};
203use tracing::debug;
204
205/// Summary statistics for the Dog Leg optimization process.
206#[derive(Debug, Clone)]
207pub struct DogLegSummary {
208 /// Initial cost value
209 pub initial_cost: f64,
210 /// Final cost value
211 pub final_cost: f64,
212 /// Total number of iterations performed
213 pub iterations: usize,
214 /// Number of successful steps (cost decreased)
215 pub successful_steps: usize,
216 /// Number of unsuccessful steps (cost increased, step rejected)
217 pub unsuccessful_steps: usize,
218 /// Final trust region radius
219 pub final_trust_region_radius: f64,
220 /// Average cost reduction per iteration
221 pub average_cost_reduction: f64,
222 /// Maximum gradient norm encountered
223 pub max_gradient_norm: f64,
224 /// Final gradient norm
225 pub final_gradient_norm: f64,
226 /// Maximum parameter update norm
227 pub max_parameter_update_norm: f64,
228 /// Final parameter update norm
229 pub final_parameter_update_norm: f64,
230 /// Total time elapsed
231 pub total_time: time::Duration,
232 /// Average time per iteration
233 pub average_time_per_iteration: time::Duration,
234 /// Detailed per-iteration statistics history
235 pub iteration_history: Vec<IterationStats>,
236 /// Convergence status
237 pub convergence_status: optimizer::OptimizationStatus,
238}
239
240/// Per-iteration statistics for detailed logging (Ceres-style output).
241///
242/// Captures all relevant metrics for each optimization iteration, enabling
243/// detailed analysis and debugging of the optimization process.
244#[derive(Debug, Clone)]
245pub struct IterationStats {
246 /// Iteration number (0-indexed)
247 pub iteration: usize,
248 /// Cost function value at this iteration
249 pub cost: f64,
250 /// Change in cost from previous iteration
251 pub cost_change: f64,
252 /// L2 norm of the gradient (||J^T·r||)
253 pub gradient_norm: f64,
254 /// L2 norm of the parameter update step (||Δx||)
255 pub step_norm: f64,
256 /// Trust region ratio (ρ = actual_reduction / predicted_reduction)
257 pub tr_ratio: f64,
258 /// Trust region radius (Δ)
259 pub tr_radius: f64,
260 /// Linear solver iterations (0 for direct solvers like Cholesky)
261 pub ls_iter: usize,
262 /// Time taken for this iteration in milliseconds
263 pub iter_time_ms: f64,
264 /// Total elapsed time since optimization started in milliseconds
265 pub total_time_ms: f64,
266 /// Whether the step was accepted (true) or rejected (false)
267 pub accepted: bool,
268}
269
270impl IterationStats {
271 /// Print table header in Ceres-style format
272 pub fn print_header() {
273 debug!(
274 "{:>4} {:>13} {:>13} {:>13} {:>13} {:>11} {:>11} {:>7} {:>11} {:>13} {:>6}",
275 "iter",
276 "cost",
277 "cost_change",
278 "|gradient|",
279 "|step|",
280 "tr_ratio",
281 "tr_radius",
282 "ls_iter",
283 "iter_time",
284 "total_time",
285 "status"
286 );
287 }
288
289 /// Print single iteration line in Ceres-style format with scientific notation
290 pub fn print_line(&self) {
291 let status = if self.iteration == 0 {
292 "-"
293 } else if self.accepted {
294 "✓"
295 } else {
296 "✗"
297 };
298
299 debug!(
300 "{:>4} {:>13.6e} {:>13.2e} {:>13.2e} {:>13.2e} {:>11.2e} {:>11.2e} {:>7} {:>9.2}ms {:>11.2}ms {:>6}",
301 self.iteration,
302 self.cost,
303 self.cost_change,
304 self.gradient_norm,
305 self.step_norm,
306 self.tr_ratio,
307 self.tr_radius,
308 self.ls_iter,
309 self.iter_time_ms,
310 self.total_time_ms,
311 status
312 );
313 }
314}
315
316impl fmt::Display for DogLegSummary {
317 fn fmt(&self, f: &mut fmt::Formatter<'_>) -> fmt::Result {
318 // Determine if converged
319 let converged = matches!(
320 self.convergence_status,
321 optimizer::OptimizationStatus::Converged
322 | optimizer::OptimizationStatus::CostToleranceReached
323 | optimizer::OptimizationStatus::GradientToleranceReached
324 | optimizer::OptimizationStatus::ParameterToleranceReached
325 );
326
327 writeln!(f, "Dog-Leg Final Result")?;
328
329 // Title with convergence status
330 if converged {
331 writeln!(f, "CONVERGED ({:?})", self.convergence_status)?;
332 } else {
333 writeln!(f, "DIVERGED ({:?})", self.convergence_status)?;
334 }
335
336 writeln!(f)?;
337 writeln!(f, "Cost:")?;
338 writeln!(f, " Initial: {:.6e}", self.initial_cost)?;
339 writeln!(f, " Final: {:.6e}", self.final_cost)?;
340 writeln!(
341 f,
342 " Reduction: {:.6e} ({:.2}%)",
343 self.initial_cost - self.final_cost,
344 100.0 * (self.initial_cost - self.final_cost) / self.initial_cost.max(1e-12)
345 )?;
346 writeln!(f)?;
347 writeln!(f, "Iterations:")?;
348 writeln!(f, " Total: {}", self.iterations)?;
349 writeln!(
350 f,
351 " Successful steps: {} ({:.1}%)",
352 self.successful_steps,
353 100.0 * self.successful_steps as f64 / self.iterations.max(1) as f64
354 )?;
355 writeln!(
356 f,
357 " Unsuccessful steps: {} ({:.1}%)",
358 self.unsuccessful_steps,
359 100.0 * self.unsuccessful_steps as f64 / self.iterations.max(1) as f64
360 )?;
361 writeln!(f)?;
362 writeln!(f, "Trust Region:")?;
363 writeln!(f, " Final radius: {:.6e}", self.final_trust_region_radius)?;
364 writeln!(f)?;
365 writeln!(f, "Gradient:")?;
366 writeln!(f, " Max norm: {:.2e}", self.max_gradient_norm)?;
367 writeln!(f, " Final norm: {:.2e}", self.final_gradient_norm)?;
368 writeln!(f)?;
369 writeln!(f, "Parameter Update:")?;
370 writeln!(f, " Max norm: {:.2e}", self.max_parameter_update_norm)?;
371 writeln!(f, " Final norm: {:.2e}", self.final_parameter_update_norm)?;
372 writeln!(f)?;
373 writeln!(f, "Performance:")?;
374 writeln!(
375 f,
376 " Total time: {:.2}ms",
377 self.total_time.as_secs_f64() * 1000.0
378 )?;
379 writeln!(
380 f,
381 " Average per iteration: {:.2}ms",
382 self.average_time_per_iteration.as_secs_f64() * 1000.0
383 )?;
384
385 Ok(())
386 }
387}
388
389/// Configuration parameters for the Dog Leg trust region optimizer.
390///
391/// Controls trust region management, convergence criteria, adaptive regularization,
392/// and Ceres Solver enhancements for the Dog Leg algorithm.
393///
394/// # Builder Pattern
395///
396/// All configuration options can be set using the builder pattern:
397///
398/// ```
399/// use apex_solver::optimizer::dog_leg::DogLegConfig;
400///
401/// let config = DogLegConfig::new()
402/// .with_max_iterations(100)
403/// .with_trust_region_radius(1e4)
404/// .with_mu_params(1e-4, 1e-8, 1.0, 10.0)
405/// .with_jacobi_scaling(true)
406/// .with_step_reuse(true);
407/// ```
408///
409/// # Trust Region Behavior
410///
411/// The trust region radius Δ controls the maximum allowed step size:
412///
413/// - **Initial radius** (`trust_region_radius`): Starting value (default: 1e4)
414/// - **Bounds** (`trust_region_min`, `trust_region_max`): Valid range (default: 1e-3 to 1e6)
415/// - **Adaptation**: Increases for good steps, decreases for poor steps
416///
417/// # Adaptive μ Regularization (Ceres Enhancement)
418///
419/// Controls the regularization parameter in `(J^T·J + μI)·h = -J^T·r`:
420///
421/// - `initial_mu`: Starting value (default: 1e-4 for numerical stability)
422/// - `min_mu`, `max_mu`: Bounds (default: 1e-8 to 1.0)
423/// - `mu_increase_factor`: Multiplier when solve fails (default: 10.0)
424///
425/// # Convergence Criteria
426///
427/// The optimizer terminates when ANY of the following conditions is met:
428///
429/// - **Cost tolerance**: `|cost_k - cost_{k-1}| < cost_tolerance`
430/// - **Parameter tolerance**: `||step|| < parameter_tolerance`
431/// - **Gradient tolerance**: `||J^T·r|| < gradient_tolerance`
432/// - **Maximum iterations**: `iteration >= max_iterations`
433/// - **Timeout**: `elapsed_time >= timeout`
434///
435/// # See Also
436///
437/// - [`DogLeg`] - The solver that uses this configuration
438/// - [`LevenbergMarquardtConfig`](crate::optimizer::LevenbergMarquardtConfig) - Alternative damping approach
439/// - [`GaussNewtonConfig`](crate::optimizer::GaussNewtonConfig) - Undamped variant
440#[derive(Clone)]
441pub struct DogLegConfig {
442 /// Type of linear solver for the linear systems
443 pub linear_solver_type: linalg::LinearSolverType,
444 /// Maximum number of iterations
445 pub max_iterations: usize,
446 /// Convergence tolerance for cost function
447 pub cost_tolerance: f64,
448 /// Convergence tolerance for parameter updates
449 pub parameter_tolerance: f64,
450 /// Convergence tolerance for gradient norm
451 pub gradient_tolerance: f64,
452 /// Timeout duration
453 pub timeout: Option<time::Duration>,
454 /// Initial trust region radius
455 pub trust_region_radius: f64,
456 /// Minimum trust region radius
457 pub trust_region_min: f64,
458 /// Maximum trust region radius
459 pub trust_region_max: f64,
460 /// Trust region increase factor (for good steps, rho > 0.75)
461 pub trust_region_increase_factor: f64,
462 /// Trust region decrease factor (for poor steps, rho < 0.25)
463 pub trust_region_decrease_factor: f64,
464 /// Minimum step quality for acceptance (typically 0.0)
465 pub min_step_quality: f64,
466 /// Good step quality threshold (typically 0.75)
467 pub good_step_quality: f64,
468 /// Poor step quality threshold (typically 0.25)
469 pub poor_step_quality: f64,
470 /// Use Jacobi column scaling (preconditioning)
471 pub use_jacobi_scaling: bool,
472
473 // Ceres-style adaptive mu regularization parameters
474 /// Initial mu regularization parameter for Gauss-Newton step
475 pub initial_mu: f64,
476 /// Minimum mu regularization parameter
477 pub min_mu: f64,
478 /// Maximum mu regularization parameter
479 pub max_mu: f64,
480 /// Factor to increase mu when linear solver fails
481 pub mu_increase_factor: f64,
482
483 // Ceres-style step reuse optimization
484 /// Enable step reuse after rejection (Ceres-style efficiency optimization)
485 pub enable_step_reuse: bool,
486
487 /// Minimum objective function cutoff (optional early termination)
488 ///
489 /// If set, optimization terminates when cost falls below this threshold.
490 /// Useful for early stopping when a "good enough" solution is acceptable.
491 ///
492 /// Default: None (disabled)
493 pub min_cost_threshold: Option<f64>,
494
495 /// Maximum condition number for Jacobian matrix (optional check)
496 ///
497 /// If set, the optimizer checks if condition_number(J^T*J) exceeds this
498 /// threshold and terminates with IllConditionedJacobian status.
499 /// Note: Computing condition number is expensive, so this is disabled by default.
500 ///
501 /// Default: None (disabled)
502 pub max_condition_number: Option<f64>,
503
504 /// Minimum relative cost decrease for step acceptance
505 ///
506 /// Used in computing step quality (rho = actual_reduction / predicted_reduction).
507 /// Steps with rho < min_relative_decrease are rejected. Matches Ceres Solver's
508 /// min_relative_decrease parameter.
509 ///
510 /// Default: 1e-3 (Ceres-compatible)
511 pub min_relative_decrease: f64,
512
513 /// Compute per-variable covariance matrices (uncertainty estimation)
514 ///
515 /// When enabled, computes covariance by inverting the Hessian matrix after
516 /// convergence. The full covariance matrix is extracted into per-variable
517 /// blocks stored in both Variable structs and optimizer::SolverResult.
518 ///
519 /// Default: false (to avoid performance overhead)
520 pub compute_covariances: bool,
521
522 /// Enable real-time visualization (graphical debugging).
523 ///
524 /// When enabled, optimization progress is logged to a Rerun viewer.
525 /// **Note:** Requires the `visualization` feature to be enabled in `Cargo.toml`.
526 ///
527 /// Default: false
528 #[cfg(feature = "visualization")]
529 pub enable_visualization: bool,
530}
531
532impl Default for DogLegConfig {
533 fn default() -> Self {
534 Self {
535 linear_solver_type: linalg::LinearSolverType::default(),
536 // Ceres Solver default: 50 (changed from 100 for compatibility)
537 max_iterations: 50,
538 // Ceres Solver default: 1e-6 (changed from 1e-8 for compatibility)
539 cost_tolerance: 1e-6,
540 // Ceres Solver default: 1e-8 (unchanged)
541 parameter_tolerance: 1e-8,
542 // Ceres Solver default: 1e-10 (changed from 1e-8 for compatibility)
543 gradient_tolerance: 1e-10,
544 timeout: None,
545 // Ceres-style: larger initial radius for better global convergence
546 trust_region_radius: 1e4,
547 trust_region_min: 1e-12,
548 trust_region_max: 1e12,
549 // Ceres uses adaptive increase (max(radius, 3*step_norm)),
550 // but we keep factor for simpler config
551 trust_region_increase_factor: 3.0,
552 trust_region_decrease_factor: 0.5,
553 min_step_quality: 0.0,
554 good_step_quality: 0.75,
555 poor_step_quality: 0.25,
556 // Ceres-style: Enable diagonal scaling by default for elliptical trust region
557 use_jacobi_scaling: true,
558
559 // Ceres-style adaptive mu regularization defaults
560 // Start with more conservative regularization to avoid singular Hessian
561 initial_mu: 1e-4,
562 min_mu: 1e-8,
563 max_mu: 1.0,
564 mu_increase_factor: 10.0,
565
566 // Ceres-style step reuse optimization
567 enable_step_reuse: true,
568
569 // New Ceres-compatible termination parameters
570 min_cost_threshold: None,
571 max_condition_number: None,
572 min_relative_decrease: 1e-3,
573
574 compute_covariances: false,
575 #[cfg(feature = "visualization")]
576 enable_visualization: false,
577 }
578 }
579}
580
581impl DogLegConfig {
582 /// Create a new Dog Leg configuration with default values.
583 pub fn new() -> Self {
584 Self::default()
585 }
586
587 /// Set the linear solver type
588 pub fn with_linear_solver_type(mut self, linear_solver_type: linalg::LinearSolverType) -> Self {
589 self.linear_solver_type = linear_solver_type;
590 self
591 }
592
593 /// Set the maximum number of iterations
594 pub fn with_max_iterations(mut self, max_iterations: usize) -> Self {
595 self.max_iterations = max_iterations;
596 self
597 }
598
599 /// Set the cost tolerance
600 pub fn with_cost_tolerance(mut self, cost_tolerance: f64) -> Self {
601 self.cost_tolerance = cost_tolerance;
602 self
603 }
604
605 /// Set the parameter tolerance
606 pub fn with_parameter_tolerance(mut self, parameter_tolerance: f64) -> Self {
607 self.parameter_tolerance = parameter_tolerance;
608 self
609 }
610
611 /// Set the gradient tolerance
612 pub fn with_gradient_tolerance(mut self, gradient_tolerance: f64) -> Self {
613 self.gradient_tolerance = gradient_tolerance;
614 self
615 }
616
617 /// Set the timeout duration
618 pub fn with_timeout(mut self, timeout: time::Duration) -> Self {
619 self.timeout = Some(timeout);
620 self
621 }
622
623 /// Set the initial trust region radius
624 pub fn with_trust_region_radius(mut self, radius: f64) -> Self {
625 self.trust_region_radius = radius;
626 self
627 }
628
629 /// Set the trust region radius bounds
630 pub fn with_trust_region_bounds(mut self, min: f64, max: f64) -> Self {
631 self.trust_region_min = min;
632 self.trust_region_max = max;
633 self
634 }
635
636 /// Set the trust region adjustment factors
637 pub fn with_trust_region_factors(mut self, increase: f64, decrease: f64) -> Self {
638 self.trust_region_increase_factor = increase;
639 self.trust_region_decrease_factor = decrease;
640 self
641 }
642
643 /// Set the trust region quality thresholds
644 pub fn with_step_quality_thresholds(
645 mut self,
646 min_quality: f64,
647 poor_quality: f64,
648 good_quality: f64,
649 ) -> Self {
650 self.min_step_quality = min_quality;
651 self.poor_step_quality = poor_quality;
652 self.good_step_quality = good_quality;
653 self
654 }
655
656 /// Enable or disable Jacobi column scaling (preconditioning)
657 pub fn with_jacobi_scaling(mut self, use_jacobi_scaling: bool) -> Self {
658 self.use_jacobi_scaling = use_jacobi_scaling;
659 self
660 }
661
662 /// Set adaptive mu regularization parameters (Ceres-style)
663 pub fn with_mu_params(
664 mut self,
665 initial_mu: f64,
666 min_mu: f64,
667 max_mu: f64,
668 increase_factor: f64,
669 ) -> Self {
670 self.initial_mu = initial_mu;
671 self.min_mu = min_mu;
672 self.max_mu = max_mu;
673 self.mu_increase_factor = increase_factor;
674 self
675 }
676
677 /// Enable or disable step reuse optimization (Ceres-style)
678 pub fn with_step_reuse(mut self, enable_step_reuse: bool) -> Self {
679 self.enable_step_reuse = enable_step_reuse;
680 self
681 }
682
683 /// Set minimum objective function cutoff for early termination.
684 ///
685 /// When set, optimization terminates with MinCostThresholdReached status
686 /// if the cost falls below this threshold. Useful for early stopping when
687 /// a "good enough" solution is acceptable.
688 pub fn with_min_cost_threshold(mut self, min_cost: f64) -> Self {
689 self.min_cost_threshold = Some(min_cost);
690 self
691 }
692
693 /// Set maximum condition number for Jacobian matrix.
694 ///
695 /// If set, the optimizer checks if condition_number(J^T*J) exceeds this
696 /// threshold and terminates with IllConditionedJacobian status.
697 /// Note: Computing condition number is expensive, disabled by default.
698 pub fn with_max_condition_number(mut self, max_cond: f64) -> Self {
699 self.max_condition_number = Some(max_cond);
700 self
701 }
702
703 /// Set minimum relative cost decrease for step acceptance.
704 ///
705 /// Steps with rho = (actual_reduction / predicted_reduction) below this
706 /// threshold are rejected. Default: 1e-3 (Ceres-compatible)
707 pub fn with_min_relative_decrease(mut self, min_decrease: f64) -> Self {
708 self.min_relative_decrease = min_decrease;
709 self
710 }
711
712 /// Enable or disable covariance computation (uncertainty estimation).
713 ///
714 /// When enabled, computes the full covariance matrix by inverting the Hessian
715 /// after convergence, then extracts per-variable covariance blocks.
716 pub fn with_compute_covariances(mut self, compute_covariances: bool) -> Self {
717 self.compute_covariances = compute_covariances;
718 self
719 }
720
721 /// Enable real-time visualization.
722 ///
723 /// **Note:** Requires the `visualization` feature to be enabled in `Cargo.toml`.
724 ///
725 /// # Arguments
726 ///
727 /// * `enable` - Whether to enable visualization
728 #[cfg(feature = "visualization")]
729 pub fn with_visualization(mut self, enable: bool) -> Self {
730 self.enable_visualization = enable;
731 self
732 }
733
734 /// Print configuration parameters (info level logging)
735 pub fn print_configuration(&self) {
736 debug!(
737 "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: {}",
738 self.linear_solver_type,
739 self.max_iterations,
740 self.cost_tolerance,
741 self.parameter_tolerance,
742 self.gradient_tolerance,
743 self.timeout,
744 self.trust_region_radius,
745 self.trust_region_min,
746 self.trust_region_max,
747 self.min_step_quality,
748 self.good_step_quality,
749 self.poor_step_quality,
750 self.initial_mu,
751 self.min_mu,
752 self.max_mu,
753 self.mu_increase_factor,
754 if self.use_jacobi_scaling {
755 "enabled"
756 } else {
757 "disabled"
758 },
759 if self.enable_step_reuse {
760 "enabled"
761 } else {
762 "disabled"
763 },
764 if self.compute_covariances {
765 "enabled"
766 } else {
767 "disabled"
768 }
769 );
770 }
771}
772
773/// State for optimization iteration
774struct OptimizationState {
775 variables: collections::HashMap<String, problem::VariableEnum>,
776 variable_index_map: collections::HashMap<String, usize>,
777 sorted_vars: Vec<String>,
778 symbolic_structure: problem::SymbolicStructure,
779 current_cost: f64,
780 initial_cost: f64,
781}
782
783/// Result from step computation
784struct StepResult {
785 step: faer::Mat<f64>,
786 gradient_norm: f64,
787 predicted_reduction: f64,
788}
789
790/// Type of step taken
791#[derive(Debug, Clone, Copy)]
792enum StepType {
793 /// Full Gauss-Newton step
794 GaussNewton,
795 /// Scaled steepest descent (Cauchy point)
796 SteepestDescent,
797 /// Dog leg interpolation
798 DogLeg,
799}
800
801impl fmt::Display for StepType {
802 fn fmt(&self, f: &mut fmt::Formatter<'_>) -> fmt::Result {
803 match self {
804 StepType::GaussNewton => write!(f, "GN"),
805 StepType::SteepestDescent => write!(f, "SD"),
806 StepType::DogLeg => write!(f, "DL"),
807 }
808 }
809}
810
811/// Result from step evaluation
812struct StepEvaluation {
813 accepted: bool,
814 cost_reduction: f64,
815 rho: f64,
816}
817
818/// Dog Leg trust region solver for nonlinear least squares optimization.
819///
820/// Implements Powell's Dog Leg algorithm with Ceres Solver enhancements including
821/// adaptive μ regularization, numerically robust beta computation, and step reuse caching.
822///
823/// # Algorithm
824///
825/// At each iteration k:
826/// 1. Compute residual `r(xₖ)` and Jacobian `J(xₖ)`
827/// 2. Solve for Gauss-Newton step: `(J^T·J + μI)·h_gn = -J^T·r`
828/// 3. Compute steepest descent direction: `-g` where `g = J^T·r`
829/// 4. Compute Cauchy point: `p_c = -α·g` (optimal along steepest descent)
830/// 5. Construct dog leg step based on trust region radius Δ:
831/// - If `||h_gn|| ≤ Δ`: Take full GN step
832/// - Else if `||p_c|| ≥ Δ`: Take scaled SD to boundary
833/// - Else: Interpolate `h = p_c + β·(h_gn - p_c)` where `||h|| = Δ`
834/// 6. Evaluate gain ratio: `ρ = (actual reduction) / (predicted reduction)`
835/// 7. Update trust region radius based on ρ
836/// 8. Accept/reject step and update parameters
837///
838/// # Ceres Solver Enhancements
839///
840/// This implementation includes four major improvements from Google's Ceres Solver:
841///
842/// **1. Adaptive μ Regularization:** Dynamically adjusts regularization parameter
843/// to handle ill-conditioned Hessians (increases on failure, decreases on success).
844///
845/// **2. Numerically Robust Beta:** Uses two formulas for computing dog leg
846/// interpolation parameter β to avoid catastrophic cancellation.
847///
848/// **3. Step Reuse Mechanism:** Caches GN step, Cauchy point, and gradient when
849/// steps are rejected. Limited to 5 consecutive reuses to prevent staleness.
850///
851/// **4. Jacobi Scaling:** Optional diagonal preconditioning creates elliptical
852/// trust regions for better handling of mixed-scale problems.
853///
854/// # Examples
855///
856/// ```no_run
857/// use apex_solver::optimizer::DogLeg;
858/// use apex_solver::core::problem::Problem;
859/// use std::collections::HashMap;
860///
861/// # fn main() -> Result<(), Box<dyn std::error::Error>> {
862/// let mut problem = Problem::new();
863/// // ... add factors to problem ...
864///
865/// let initial_values = HashMap::new();
866/// // ... initialize parameters ...
867///
868/// let mut solver = DogLeg::new();
869/// let result = solver.optimize(&problem, &initial_values)?;
870/// # Ok(())
871/// # }
872/// ```
873///
874/// # See Also
875///
876/// - [`DogLegConfig`] - Configuration options
877/// - [`LevenbergMarquardt`](crate::optimizer::LevenbergMarquardt) - Alternative adaptive damping
878/// - [`GaussNewton`](crate::optimizer::GaussNewton) - Undamped variant
879pub struct DogLeg {
880 config: DogLegConfig,
881 jacobi_scaling: Option<sparse::SparseColMat<usize, f64>>,
882 observers: optimizer::OptObserverVec,
883
884 // Adaptive mu regularization (Ceres-style)
885 mu: f64,
886 min_mu: f64,
887 max_mu: f64,
888 mu_increase_factor: f64,
889
890 // Step reuse mechanism (Ceres-style efficiency optimization)
891 reuse_step_on_rejection: bool,
892 cached_gn_step: Option<faer::Mat<f64>>,
893 cached_cauchy_point: Option<faer::Mat<f64>>,
894 cached_gradient: Option<faer::Mat<f64>>,
895 cached_alpha: Option<f64>,
896 cache_reuse_count: usize, // Track consecutive reuses to prevent staleness
897}
898
899impl Default for DogLeg {
900 fn default() -> Self {
901 Self::new()
902 }
903}
904
905impl DogLeg {
906 /// Create a new Dog Leg solver with default configuration.
907 pub fn new() -> Self {
908 Self::with_config(DogLegConfig::default())
909 }
910
911 /// Create a new Dog Leg solver with the given configuration.
912 pub fn with_config(config: DogLegConfig) -> Self {
913 Self {
914 // Initialize adaptive mu from config
915 mu: config.initial_mu,
916 min_mu: config.min_mu,
917 max_mu: config.max_mu,
918 mu_increase_factor: config.mu_increase_factor,
919
920 // Initialize step reuse mechanism (disabled initially, enabled after first rejection)
921 reuse_step_on_rejection: false,
922 cached_gn_step: None,
923 cached_cauchy_point: None,
924 cached_gradient: None,
925 cached_alpha: None,
926 cache_reuse_count: 0,
927
928 config,
929 jacobi_scaling: None,
930 observers: optimizer::OptObserverVec::new(),
931 }
932 }
933
934 /// Add an observer to the solver.
935 ///
936 /// Observers are notified at each iteration with the current variable values.
937 /// This enables real-time visualization, logging, metrics collection, etc.
938 ///
939 /// # Examples
940 ///
941 /// ```no_run
942 /// use apex_solver::optimizer::DogLeg;
943 /// # use apex_solver::optimizer::OptObserver;
944 /// # use std::collections::HashMap;
945 /// # use apex_solver::core::problem::VariableEnum;
946 ///
947 /// # struct MyObserver;
948 /// # impl OptObserver for MyObserver {
949 /// # fn on_step(&self, _: &HashMap<String, VariableEnum>, _: usize) {}
950 /// # }
951 /// let mut solver = DogLeg::new();
952 /// solver.add_observer(MyObserver);
953 /// ```
954 pub fn add_observer(&mut self, observer: impl optimizer::OptObserver + 'static) {
955 self.observers.add(observer);
956 }
957
958 /// Create the appropriate linear solver based on configuration
959 fn create_linear_solver(&self) -> Box<dyn linalg::SparseLinearSolver> {
960 match self.config.linear_solver_type {
961 linalg::LinearSolverType::SparseCholesky => {
962 Box::new(linalg::SparseCholeskySolver::new())
963 }
964 linalg::LinearSolverType::SparseQR => Box::new(linalg::SparseQRSolver::new()),
965 linalg::LinearSolverType::SparseSchurComplement => {
966 // Schur complement solver requires special handling - fallback to Cholesky
967 Box::new(linalg::SparseCholeskySolver::new())
968 }
969 }
970 }
971
972 /// Compute the Cauchy point (steepest descent step)
973 /// Returns the optimal step along the negative gradient direction
974 /// Compute Cauchy point and optimal step length for steepest descent
975 ///
976 /// Returns (alpha, cauchy_point) where:
977 /// - alpha: optimal step length α = ||g||² / (g^T H g)
978 /// - cauchy_point: p_c = -α * g (the Cauchy point)
979 ///
980 /// This is the optimal point along the steepest descent direction within
981 /// the quadratic approximation of the objective function.
982 fn compute_cauchy_point_and_alpha(
983 &self,
984 gradient: &faer::Mat<f64>,
985 hessian: &sparse::SparseColMat<usize, f64>,
986 ) -> (f64, faer::Mat<f64>) {
987 // Optimal step size along steepest descent: α = (g^T*g) / (g^T*H*g)
988 let g_norm_sq_mat = gradient.transpose() * gradient;
989 let g_norm_sq = g_norm_sq_mat[(0, 0)];
990
991 let h_g = hessian * gradient;
992 let g_h_g_mat = gradient.transpose() * &h_g;
993 let g_h_g = g_h_g_mat[(0, 0)];
994
995 // Avoid division by zero
996 let alpha = if g_h_g.abs() > 1e-15 {
997 g_norm_sq / g_h_g
998 } else {
999 1.0
1000 };
1001
1002 // Compute Cauchy point: p_c = -α * gradient
1003 let mut cauchy_point = faer::Mat::zeros(gradient.nrows(), 1);
1004 for i in 0..gradient.nrows() {
1005 cauchy_point[(i, 0)] = -alpha * gradient[(i, 0)];
1006 }
1007
1008 (alpha, cauchy_point)
1009 }
1010
1011 /// Compute the dog leg step using Powell's Dog Leg method
1012 ///
1013 /// The dog leg path consists of two segments:
1014 /// 1. From origin to Cauchy point (optimal along steepest descent)
1015 /// 2. From Cauchy point to Gauss-Newton step
1016 ///
1017 /// Arguments:
1018 /// - steepest_descent_dir: -gradient (steepest descent direction, not scaled)
1019 /// - cauchy_point: p_c = α * (-gradient), the optimal steepest descent step
1020 /// - h_gn: Gauss-Newton step
1021 /// - delta: Trust region radius
1022 ///
1023 /// Returns: (step, step_type)
1024 fn compute_dog_leg_step(
1025 &self,
1026 steepest_descent_dir: &faer::Mat<f64>,
1027 cauchy_point: &faer::Mat<f64>,
1028 h_gn: &faer::Mat<f64>,
1029 delta: f64,
1030 ) -> (faer::Mat<f64>, StepType) {
1031 let gn_norm = h_gn.norm_l2();
1032 let cauchy_norm = cauchy_point.norm_l2();
1033 let sd_norm = steepest_descent_dir.norm_l2();
1034
1035 // Case 1: Full Gauss-Newton step fits in trust region
1036 if gn_norm <= delta {
1037 return (h_gn.clone(), StepType::GaussNewton);
1038 }
1039
1040 // Case 2: Even Cauchy point is outside trust region
1041 // Scale steepest descent direction to boundary: (delta / ||δ_sd||) * δ_sd
1042 if cauchy_norm >= delta {
1043 let mut scaled_sd = faer::Mat::zeros(steepest_descent_dir.nrows(), 1);
1044 let scale = delta / sd_norm;
1045 for i in 0..steepest_descent_dir.nrows() {
1046 scaled_sd[(i, 0)] = steepest_descent_dir[(i, 0)] * scale;
1047 }
1048 return (scaled_sd, StepType::SteepestDescent);
1049 }
1050
1051 // Case 3: Dog leg interpolation between Cauchy point and GN step
1052 // Use Ceres-style numerically robust formula
1053 //
1054 // Following Ceres solver implementation for numerical stability:
1055 // Compute intersection of trust region boundary with line from Cauchy point to GN step
1056 //
1057 // Let v = δ_gn - p_c
1058 // Solve: ||p_c + β*v||² = delta²
1059 // This gives: a*β² + 2*b*β + c = 0
1060 // where:
1061 // a = v^T·v = ||v||²
1062 // b = p_c^T·v
1063 // c = p_c^T·p_c - delta² = ||p_c||² - delta²
1064
1065 let mut v = faer::Mat::zeros(cauchy_point.nrows(), 1);
1066 for i in 0..cauchy_point.nrows() {
1067 v[(i, 0)] = h_gn[(i, 0)] - cauchy_point[(i, 0)];
1068 }
1069
1070 // Compute coefficients
1071 let v_squared_norm = v.transpose() * &v;
1072 let a = v_squared_norm[(0, 0)];
1073
1074 let pc_dot_v = cauchy_point.transpose() * &v;
1075 let b = pc_dot_v[(0, 0)];
1076
1077 let c = cauchy_norm * cauchy_norm - delta * delta;
1078
1079 // Ceres-style numerically robust beta computation
1080 // Uses two different formulas based on sign of b to avoid catastrophic cancellation
1081 let d_squared = b * b - a * c;
1082
1083 let beta = if d_squared < 0.0 {
1084 // Should not happen geometrically, but handle gracefully
1085 1.0
1086 } else if a.abs() < 1e-15 {
1087 // Degenerate case: v is nearly zero
1088 1.0
1089 } else {
1090 let d = d_squared.sqrt();
1091
1092 // Ceres formula: choose formula based on sign of b to avoid cancellation
1093 // If b <= 0: beta = (-b + d) / a (standard formula, no cancellation)
1094 // If b > 0: beta = -c / (b + d) (alternative formula, avoids cancellation)
1095 if b <= 0.0 { (-b + d) / a } else { -c / (b + d) }
1096 };
1097
1098 // Clamp beta to [0, 1] for safety
1099 let beta = beta.clamp(0.0, 1.0);
1100
1101 // Compute dog leg step: p_dl = p_c + β*(δ_gn - p_c)
1102 let mut dog_leg = faer::Mat::zeros(cauchy_point.nrows(), 1);
1103 for i in 0..cauchy_point.nrows() {
1104 dog_leg[(i, 0)] = cauchy_point[(i, 0)] + beta * v[(i, 0)];
1105 }
1106
1107 (dog_leg, StepType::DogLeg)
1108 }
1109
1110 /// Update trust region radius based on step quality (Ceres-style)
1111 fn update_trust_region(&mut self, rho: f64, step_norm: f64) -> bool {
1112 if rho > self.config.good_step_quality {
1113 // Good step, increase trust region (Ceres-style: max(radius, 3*step_norm))
1114 let new_radius = self.config.trust_region_radius.max(3.0 * step_norm);
1115 self.config.trust_region_radius = new_radius.min(self.config.trust_region_max);
1116
1117 // Decrease mu on successful step (Ceres-style adaptive regularization)
1118 self.mu = (self.mu / (0.5 * self.mu_increase_factor)).max(self.min_mu);
1119
1120 // Clear reuse flag and invalidate cache on acceptance (parameters have moved)
1121 self.reuse_step_on_rejection = false;
1122 self.cached_gn_step = None;
1123 self.cached_cauchy_point = None;
1124 self.cached_gradient = None;
1125 self.cached_alpha = None;
1126 self.cache_reuse_count = 0;
1127
1128 true
1129 } else if rho < self.config.poor_step_quality {
1130 // Poor step, decrease trust region
1131 self.config.trust_region_radius = (self.config.trust_region_radius
1132 * self.config.trust_region_decrease_factor)
1133 .max(self.config.trust_region_min);
1134
1135 // Enable step reuse flag for next iteration (Ceres-style)
1136 self.reuse_step_on_rejection = self.config.enable_step_reuse;
1137
1138 false
1139 } else {
1140 // Moderate step, keep trust region unchanged
1141 // Clear reuse flag and invalidate cache on acceptance (parameters have moved)
1142 self.reuse_step_on_rejection = false;
1143 self.cached_gn_step = None;
1144 self.cached_cauchy_point = None;
1145 self.cached_gradient = None;
1146 self.cached_alpha = None;
1147 self.cache_reuse_count = 0;
1148
1149 true
1150 }
1151 }
1152
1153 /// Compute step quality ratio (actual vs predicted reduction)
1154 fn compute_step_quality(
1155 &self,
1156 current_cost: f64,
1157 new_cost: f64,
1158 predicted_reduction: f64,
1159 ) -> f64 {
1160 let actual_reduction = current_cost - new_cost;
1161 if predicted_reduction.abs() < 1e-15 {
1162 if actual_reduction > 0.0 { 1.0 } else { 0.0 }
1163 } else {
1164 actual_reduction / predicted_reduction
1165 }
1166 }
1167
1168 /// Compute predicted cost reduction from linear model
1169 fn compute_predicted_reduction(
1170 &self,
1171 step: &faer::Mat<f64>,
1172 gradient: &faer::Mat<f64>,
1173 hessian: &sparse::SparseColMat<usize, f64>,
1174 ) -> f64 {
1175 // Dog Leg predicted reduction: -step^T * gradient - 0.5 * step^T * H * step
1176 let linear_term = step.transpose() * gradient;
1177 let hessian_step = hessian * step;
1178 let quadratic_term = step.transpose() * &hessian_step;
1179
1180 -linear_term[(0, 0)] - 0.5 * quadratic_term[(0, 0)]
1181 }
1182
1183 /// Check convergence criteria
1184 /// Check convergence using comprehensive termination criteria.
1185 ///
1186 /// Implements 9 termination criteria following Ceres Solver standards for Dog Leg:
1187 ///
1188 /// 1. **Gradient Norm (First-Order Optimality)**: ||g||∞ ≤ gradient_tolerance
1189 /// 2. **Parameter Change Tolerance**: ||h|| ≤ parameter_tolerance * (||x|| + parameter_tolerance)
1190 /// 3. **Function Value Change Tolerance**: |ΔF| < cost_tolerance * F
1191 /// 4. **Objective Function Cutoff**: F_new < min_cost_threshold (optional)
1192 /// 5. **Trust Region Radius**: radius < trust_region_min
1193 /// 6. **Singular/Ill-Conditioned Jacobian**: Detected during linear solve
1194 /// 7. **Invalid Numerical Values**: NaN or Inf in cost or parameters
1195 /// 8. **Maximum Iterations**: iteration >= max_iterations
1196 /// 9. **Timeout**: elapsed >= timeout
1197 ///
1198 /// # Arguments
1199 ///
1200 /// * `iteration` - Current iteration number
1201 /// * `current_cost` - Cost before applying the step
1202 /// * `new_cost` - Cost after applying the step
1203 /// * `parameter_norm` - L2 norm of current parameter vector ||x||
1204 /// * `parameter_update_norm` - L2 norm of parameter update step ||h||
1205 /// * `gradient_norm` - Infinity norm of gradient ||g||∞
1206 /// * `trust_region_radius` - Current trust region radius
1207 /// * `elapsed` - Elapsed time since optimization start
1208 /// * `step_accepted` - Whether the current step was accepted
1209 ///
1210 /// # Returns
1211 ///
1212 /// `Some(OptimizationStatus)` if any termination criterion is satisfied, `None` otherwise.
1213 ///
1214 /// # Design Note
1215 /// This internal convergence checker requires multiple scalar metrics for comprehensive
1216 /// termination criteria. Grouping into a struct would reduce clarity without performance benefit.
1217 #[allow(clippy::too_many_arguments)]
1218 fn check_convergence(
1219 &self,
1220 iteration: usize,
1221 current_cost: f64,
1222 new_cost: f64,
1223 parameter_norm: f64,
1224 parameter_update_norm: f64,
1225 gradient_norm: f64,
1226 trust_region_radius: f64,
1227 elapsed: time::Duration,
1228 step_accepted: bool,
1229 ) -> Option<optimizer::OptimizationStatus> {
1230 // ========================================================================
1231 // CRITICAL SAFETY CHECKS (perform first, before convergence checks)
1232 // ========================================================================
1233
1234 // CRITERION 7: Invalid Numerical Values (NaN/Inf)
1235 // Always check for numerical instability first
1236 if !new_cost.is_finite() || !parameter_update_norm.is_finite() || !gradient_norm.is_finite()
1237 {
1238 return Some(optimizer::OptimizationStatus::InvalidNumericalValues);
1239 }
1240
1241 // CRITERION 9: Timeout
1242 // Check wall-clock time limit
1243 if let Some(timeout) = self.config.timeout
1244 && elapsed >= timeout
1245 {
1246 return Some(optimizer::OptimizationStatus::Timeout);
1247 }
1248
1249 // CRITERION 8: Maximum Iterations
1250 // Check iteration count limit
1251 if iteration >= self.config.max_iterations {
1252 return Some(optimizer::OptimizationStatus::MaxIterationsReached);
1253 }
1254
1255 // ========================================================================
1256 // CONVERGENCE CRITERIA (only check after successful steps)
1257 // ========================================================================
1258
1259 // Only check convergence criteria after accepted steps
1260 // (rejected steps don't indicate convergence)
1261 if !step_accepted {
1262 return None;
1263 }
1264
1265 // CRITERION 1: Gradient Norm (First-Order Optimality)
1266 // Check if gradient infinity norm is below threshold
1267 // This indicates we're at a critical point (local minimum, saddle, or maximum)
1268 if gradient_norm < self.config.gradient_tolerance {
1269 return Some(optimizer::OptimizationStatus::GradientToleranceReached);
1270 }
1271
1272 // Only check parameter and cost criteria after first iteration
1273 if iteration > 0 {
1274 // CRITERION 2: Parameter Change Tolerance (xtol)
1275 // Ceres formula: ||h|| ≤ ε_param * (||x|| + ε_param)
1276 // This is a relative measure that scales with parameter magnitude
1277 let relative_step_tolerance = self.config.parameter_tolerance
1278 * (parameter_norm + self.config.parameter_tolerance);
1279
1280 if parameter_update_norm <= relative_step_tolerance {
1281 return Some(optimizer::OptimizationStatus::ParameterToleranceReached);
1282 }
1283
1284 // CRITERION 3: Function Value Change Tolerance (ftol)
1285 // Ceres formula: |ΔF| < ε_cost * F
1286 // Check relative cost change (not absolute)
1287 let cost_change = (current_cost - new_cost).abs();
1288 let relative_cost_change = cost_change / current_cost.max(1e-10); // Avoid division by zero
1289
1290 if relative_cost_change < self.config.cost_tolerance {
1291 return Some(optimizer::OptimizationStatus::CostToleranceReached);
1292 }
1293 }
1294
1295 // CRITERION 4: Objective Function Cutoff (optional early stopping)
1296 // Useful for "good enough" solutions
1297 if let Some(min_cost) = self.config.min_cost_threshold
1298 && new_cost < min_cost
1299 {
1300 return Some(optimizer::OptimizationStatus::MinCostThresholdReached);
1301 }
1302
1303 // CRITERION 5: Trust Region Radius
1304 // If trust region has collapsed, optimization has converged or problem is ill-conditioned
1305 if trust_region_radius < self.config.trust_region_min {
1306 return Some(optimizer::OptimizationStatus::TrustRegionRadiusTooSmall);
1307 }
1308
1309 // CRITERION 6: Singular/Ill-Conditioned Jacobian
1310 // Note: This is typically detected during the linear solve and handled there
1311 // The max_condition_number check would be expensive to compute here
1312 // If linear solve fails, it returns an error that's converted to NumericalFailure
1313
1314 // No termination criterion satisfied
1315 None
1316 }
1317
1318 /// Compute total parameter vector norm ||x||.
1319 ///
1320 /// Computes the L2 norm of all parameter vectors concatenated together.
1321 /// This is used in the relative parameter tolerance check.
1322 ///
1323 /// # Arguments
1324 ///
1325 /// * `variables` - Map of variable names to their current values
1326 ///
1327 /// # Returns
1328 ///
1329 /// The L2 norm of the concatenated parameter vector
1330 fn compute_parameter_norm(
1331 variables: &collections::HashMap<String, problem::VariableEnum>,
1332 ) -> f64 {
1333 variables
1334 .values()
1335 .map(|v| {
1336 let vec = v.to_vector();
1337 vec.norm_squared()
1338 })
1339 .sum::<f64>()
1340 .sqrt()
1341 }
1342
1343 /// Create Jacobi scaling matrix from Jacobian
1344 fn create_jacobi_scaling(
1345 &self,
1346 jacobian: &sparse::SparseColMat<usize, f64>,
1347 ) -> Result<sparse::SparseColMat<usize, f64>, optimizer::OptimizerError> {
1348 use faer::sparse::Triplet;
1349
1350 let cols = jacobian.ncols();
1351 let jacobi_scaling_vec: Vec<Triplet<usize, usize, f64>> = (0..cols)
1352 .map(|c| {
1353 let col_norm_squared: f64 = jacobian
1354 .triplet_iter()
1355 .filter(|t| t.col == c)
1356 .map(|t| t.val * t.val)
1357 .sum();
1358 let col_norm = col_norm_squared.sqrt();
1359 let scaling = 1.0 / (1.0 + col_norm);
1360 Triplet::new(c, c, scaling)
1361 })
1362 .collect();
1363
1364 sparse::SparseColMat::try_new_from_triplets(cols, cols, &jacobi_scaling_vec).map_err(|e| {
1365 optimizer::OptimizerError::JacobiScalingCreation(e.to_string()).log_with_source(e)
1366 })
1367 }
1368
1369 /// Initialize optimization state
1370 fn initialize_optimization_state(
1371 &self,
1372 problem: &problem::Problem,
1373 initial_params: &collections::HashMap<
1374 String,
1375 (manifold::ManifoldType, nalgebra::DVector<f64>),
1376 >,
1377 ) -> Result<OptimizationState, error::ApexSolverError> {
1378 let variables = problem.initialize_variables(initial_params);
1379
1380 let mut variable_index_map = collections::HashMap::new();
1381 let mut col_offset = 0;
1382 let mut sorted_vars: Vec<String> = variables.keys().cloned().collect();
1383 sorted_vars.sort();
1384
1385 for var_name in &sorted_vars {
1386 variable_index_map.insert(var_name.clone(), col_offset);
1387 col_offset += variables[var_name].get_size();
1388 }
1389
1390 let symbolic_structure =
1391 problem.build_symbolic_structure(&variables, &variable_index_map, col_offset)?;
1392
1393 // Initial cost evaluation (residual only, no Jacobian needed)
1394 let residual = problem.compute_residual_sparse(&variables)?;
1395
1396 let residual_norm = residual.norm_l2();
1397 let current_cost = residual_norm * residual_norm;
1398 let initial_cost = current_cost;
1399
1400 Ok(OptimizationState {
1401 variables,
1402 variable_index_map,
1403 sorted_vars,
1404 symbolic_structure,
1405 current_cost,
1406 initial_cost,
1407 })
1408 }
1409
1410 /// Process Jacobian by creating and applying Jacobi scaling if enabled
1411 fn process_jacobian(
1412 &mut self,
1413 jacobian: &sparse::SparseColMat<usize, f64>,
1414 iteration: usize,
1415 ) -> Result<sparse::SparseColMat<usize, f64>, optimizer::OptimizerError> {
1416 // Create Jacobi scaling on first iteration if enabled
1417 if iteration == 0 {
1418 let scaling = self.create_jacobi_scaling(jacobian)?;
1419 self.jacobi_scaling = Some(scaling);
1420 }
1421 let scaling = self
1422 .jacobi_scaling
1423 .as_ref()
1424 .ok_or_else(|| optimizer::OptimizerError::JacobiScalingNotInitialized.log())?;
1425 Ok(jacobian * scaling)
1426 }
1427
1428 /// Compute dog leg optimization step
1429 fn compute_optimization_step(
1430 &mut self,
1431 residuals: &faer::Mat<f64>,
1432 scaled_jacobian: &sparse::SparseColMat<usize, f64>,
1433 linear_solver: &mut Box<dyn linalg::SparseLinearSolver>,
1434 ) -> Option<StepResult> {
1435 // Check if we can reuse cached step (Ceres-style optimization)
1436 // Safety limit: prevent excessive reuse that could lead to stale gradient/Hessian
1437 const MAX_CACHE_REUSE: usize = 5;
1438
1439 if self.reuse_step_on_rejection
1440 && self.config.enable_step_reuse
1441 && self.cache_reuse_count < MAX_CACHE_REUSE
1442 && let (Some(cached_gn), Some(cached_cauchy), Some(cached_grad), Some(_cached_a)) = (
1443 &self.cached_gn_step,
1444 &self.cached_cauchy_point,
1445 &self.cached_gradient,
1446 &self.cached_alpha,
1447 )
1448 {
1449 // Increment reuse counter
1450 self.cache_reuse_count += 1;
1451
1452 let gradient_norm = cached_grad.norm_l2();
1453 let mut steepest_descent = faer::Mat::zeros(cached_grad.nrows(), 1);
1454 for i in 0..cached_grad.nrows() {
1455 steepest_descent[(i, 0)] = -cached_grad[(i, 0)];
1456 }
1457
1458 let (scaled_step, _step_type) = self.compute_dog_leg_step(
1459 &steepest_descent,
1460 cached_cauchy,
1461 cached_gn,
1462 self.config.trust_region_radius,
1463 );
1464
1465 let step = if self.config.use_jacobi_scaling {
1466 let scaling = self.jacobi_scaling.as_ref()?;
1467 scaling * &scaled_step
1468 } else {
1469 scaled_step.clone()
1470 };
1471
1472 let hessian = linear_solver.get_hessian()?;
1473 let predicted_reduction =
1474 self.compute_predicted_reduction(&scaled_step, cached_grad, hessian);
1475
1476 return Some(StepResult {
1477 step,
1478 gradient_norm,
1479 predicted_reduction,
1480 });
1481 }
1482
1483 // Not reusing, compute fresh step
1484 // 1. Solve for Gauss-Newton step with adaptive mu regularization (Ceres-style)
1485 let residuals_owned = residuals.as_ref().to_owned();
1486 let mut scaled_gn_step = None;
1487 let mut mu_attempts = 0;
1488
1489 // Try to solve with current mu, increasing if necessary
1490 while mu_attempts < 10 && self.mu <= self.max_mu {
1491 let damping = self.mu;
1492
1493 if let Ok(step) =
1494 linear_solver.solve_augmented_equation(&residuals_owned, scaled_jacobian, damping)
1495 {
1496 scaled_gn_step = Some(step);
1497 break;
1498 }
1499
1500 // Increase mu (Ceres-style)
1501 self.mu = (self.mu * self.mu_increase_factor).min(self.max_mu);
1502 mu_attempts += 1;
1503 }
1504
1505 let scaled_gn_step = scaled_gn_step?;
1506
1507 // 2. Get gradient and Hessian (cached by solve_augmented_equation)
1508 let gradient = linear_solver.get_gradient()?;
1509 let hessian = linear_solver.get_hessian()?;
1510 let gradient_norm = gradient.norm_l2();
1511
1512 // 3. Compute steepest descent direction: δ_sd = -gradient
1513 let mut steepest_descent = faer::Mat::zeros(gradient.nrows(), 1);
1514 for i in 0..gradient.nrows() {
1515 steepest_descent[(i, 0)] = -gradient[(i, 0)];
1516 }
1517
1518 // 4. Compute Cauchy point and optimal step length
1519 let (alpha, cauchy_point) = self.compute_cauchy_point_and_alpha(gradient, hessian);
1520
1521 // 5. Compute dog leg step based on trust region radius
1522 let (scaled_step, _step_type) = self.compute_dog_leg_step(
1523 &steepest_descent,
1524 &cauchy_point,
1525 &scaled_gn_step,
1526 self.config.trust_region_radius,
1527 );
1528
1529 // 6. Apply inverse Jacobi scaling if enabled
1530 let step = if self.config.use_jacobi_scaling {
1531 let scaling = self.jacobi_scaling.as_ref()?;
1532 scaling * &scaled_step
1533 } else {
1534 scaled_step.clone()
1535 };
1536
1537 // 7. Compute predicted reduction
1538 let predicted_reduction = self.compute_predicted_reduction(&scaled_step, gradient, hessian);
1539
1540 // 8. Cache step components for potential reuse (Ceres-style)
1541 self.cached_gn_step = Some(scaled_gn_step.clone());
1542 self.cached_cauchy_point = Some(cauchy_point.clone());
1543 self.cached_gradient = Some(gradient.clone());
1544 self.cached_alpha = Some(alpha);
1545
1546 Some(StepResult {
1547 step,
1548 gradient_norm,
1549 predicted_reduction,
1550 })
1551 }
1552
1553 /// Evaluate and apply step
1554 fn evaluate_and_apply_step(
1555 &mut self,
1556 step_result: &StepResult,
1557 state: &mut OptimizationState,
1558 problem: &problem::Problem,
1559 ) -> error::ApexSolverResult<StepEvaluation> {
1560 // Apply parameter updates
1561 let step_norm = optimizer::apply_parameter_step(
1562 &mut state.variables,
1563 step_result.step.as_ref(),
1564 &state.sorted_vars,
1565 );
1566
1567 // Compute new cost (residual only, no Jacobian needed for step evaluation)
1568 let new_residual = problem.compute_residual_sparse(&state.variables)?;
1569 let new_residual_norm = new_residual.norm_l2();
1570 let new_cost = new_residual_norm * new_residual_norm;
1571
1572 // Compute step quality
1573 let rho = self.compute_step_quality(
1574 state.current_cost,
1575 new_cost,
1576 step_result.predicted_reduction,
1577 );
1578
1579 // Update trust region and decide acceptance
1580 // Filter out numerical noise with small threshold
1581 let accepted = rho > 1e-4;
1582 let _trust_region_updated = self.update_trust_region(rho, step_norm);
1583
1584 let cost_reduction = if accepted {
1585 let reduction = state.current_cost - new_cost;
1586 state.current_cost = new_cost;
1587 reduction
1588 } else {
1589 // Reject step - revert changes
1590 optimizer::apply_negative_parameter_step(
1591 &mut state.variables,
1592 step_result.step.as_ref(),
1593 &state.sorted_vars,
1594 );
1595 0.0
1596 };
1597
1598 Ok(StepEvaluation {
1599 accepted,
1600 cost_reduction,
1601 rho,
1602 })
1603 }
1604
1605 /// Create optimization summary
1606 ///
1607 /// # Design Note
1608 /// This internal summary builder accepts individual scalar results from the optimization loop.
1609 /// A struct parameter would not improve readability for this private method.
1610 #[allow(clippy::too_many_arguments)]
1611 fn create_summary(
1612 &self,
1613 initial_cost: f64,
1614 final_cost: f64,
1615 iterations: usize,
1616 successful_steps: usize,
1617 unsuccessful_steps: usize,
1618 max_gradient_norm: f64,
1619 final_gradient_norm: f64,
1620 max_parameter_update_norm: f64,
1621 final_parameter_update_norm: f64,
1622 total_cost_reduction: f64,
1623 total_time: time::Duration,
1624 iteration_history: Vec<IterationStats>,
1625 convergence_status: optimizer::OptimizationStatus,
1626 ) -> DogLegSummary {
1627 DogLegSummary {
1628 initial_cost,
1629 final_cost,
1630 iterations,
1631 successful_steps,
1632 unsuccessful_steps,
1633 final_trust_region_radius: self.config.trust_region_radius,
1634 average_cost_reduction: if iterations > 0 {
1635 total_cost_reduction / iterations as f64
1636 } else {
1637 0.0
1638 },
1639 max_gradient_norm,
1640 final_gradient_norm,
1641 max_parameter_update_norm,
1642 final_parameter_update_norm,
1643 total_time,
1644 average_time_per_iteration: if iterations > 0 {
1645 total_time / iterations as u32
1646 } else {
1647 time::Duration::from_secs(0)
1648 },
1649 iteration_history,
1650 convergence_status,
1651 }
1652 }
1653
1654 /// Minimize the optimization problem using Dog Leg algorithm
1655 pub fn optimize(
1656 &mut self,
1657 problem: &problem::Problem,
1658 initial_params: &collections::HashMap<
1659 String,
1660 (manifold::ManifoldType, nalgebra::DVector<f64>),
1661 >,
1662 ) -> Result<
1663 optimizer::SolverResult<collections::HashMap<String, problem::VariableEnum>>,
1664 error::ApexSolverError,
1665 > {
1666 let start_time = time::Instant::now();
1667 let mut iteration = 0;
1668 let mut cost_evaluations = 1;
1669 let mut jacobian_evaluations = 0;
1670 let mut successful_steps = 0;
1671 let mut unsuccessful_steps = 0;
1672
1673 let mut state = self.initialize_optimization_state(problem, initial_params)?;
1674 let mut linear_solver = self.create_linear_solver();
1675
1676 let mut max_gradient_norm: f64 = 0.0;
1677 let mut max_parameter_update_norm: f64 = 0.0;
1678 let mut total_cost_reduction = 0.0;
1679 let mut final_gradient_norm;
1680 let mut final_parameter_update_norm;
1681
1682 // Initialize iteration statistics tracking
1683 let mut iteration_stats = Vec::with_capacity(self.config.max_iterations);
1684 let mut previous_cost = state.current_cost;
1685
1686 // Print configuration and header if info/debug level is enabled
1687 if tracing::enabled!(tracing::Level::DEBUG) {
1688 self.config.print_configuration();
1689 IterationStats::print_header();
1690 }
1691
1692 loop {
1693 let iter_start = time::Instant::now();
1694 // Evaluate residuals and Jacobian
1695 let (residuals, jacobian) = problem.compute_residual_and_jacobian_sparse(
1696 &state.variables,
1697 &state.variable_index_map,
1698 &state.symbolic_structure,
1699 )?;
1700 jacobian_evaluations += 1;
1701
1702 // Process Jacobian (apply scaling if enabled)
1703 let scaled_jacobian = if self.config.use_jacobi_scaling {
1704 self.process_jacobian(&jacobian, iteration)?
1705 } else {
1706 jacobian
1707 };
1708
1709 // Compute dog leg step
1710 let step_result = match self.compute_optimization_step(
1711 &residuals,
1712 &scaled_jacobian,
1713 &mut linear_solver,
1714 ) {
1715 Some(result) => result,
1716 None => {
1717 return Err(optimizer::OptimizerError::LinearSolveFailed(
1718 "Linear solver failed to solve system".to_string(),
1719 )
1720 .into());
1721 }
1722 };
1723
1724 // Update tracking
1725 max_gradient_norm = max_gradient_norm.max(step_result.gradient_norm);
1726 final_gradient_norm = step_result.gradient_norm;
1727 let step_norm = step_result.step.norm_l2();
1728 max_parameter_update_norm = max_parameter_update_norm.max(step_norm);
1729 final_parameter_update_norm = step_norm;
1730
1731 // Evaluate and apply step
1732 let step_eval = self.evaluate_and_apply_step(&step_result, &mut state, problem)?;
1733 cost_evaluations += 1;
1734
1735 if step_eval.accepted {
1736 successful_steps += 1;
1737 total_cost_reduction += step_eval.cost_reduction;
1738 } else {
1739 unsuccessful_steps += 1;
1740 }
1741
1742 // OPTIMIZATION: Only collect iteration statistics if debug level is enabled
1743 // This eliminates ~2-5ms overhead per iteration for non-debug optimization
1744 if tracing::enabled!(tracing::Level::DEBUG) {
1745 let iter_elapsed_ms = iter_start.elapsed().as_secs_f64() * 1000.0;
1746 let total_elapsed_ms = start_time.elapsed().as_secs_f64() * 1000.0;
1747
1748 let stats = IterationStats {
1749 iteration,
1750 cost: state.current_cost,
1751 cost_change: previous_cost - state.current_cost,
1752 gradient_norm: step_result.gradient_norm,
1753 step_norm,
1754 tr_ratio: step_eval.rho,
1755 tr_radius: self.config.trust_region_radius,
1756 ls_iter: 0, // Direct solver (Cholesky) has no iterations
1757 iter_time_ms: iter_elapsed_ms,
1758 total_time_ms: total_elapsed_ms,
1759 accepted: step_eval.accepted,
1760 };
1761
1762 iteration_stats.push(stats.clone());
1763 stats.print_line();
1764 }
1765
1766 previous_cost = state.current_cost;
1767
1768 // Notify all observers with current state
1769 // First set metrics data, then notify observers
1770 self.observers.set_iteration_metrics(
1771 state.current_cost,
1772 step_result.gradient_norm,
1773 Some(self.config.trust_region_radius), // Dog Leg uses trust region radius
1774 step_norm,
1775 Some(step_eval.rho),
1776 );
1777
1778 // Set matrix data if available and there are observers
1779 if !self.observers.is_empty()
1780 && let (Some(hessian), Some(gradient)) =
1781 (linear_solver.get_hessian(), linear_solver.get_gradient())
1782 {
1783 self.observers
1784 .set_matrix_data(Some(hessian.clone()), Some(gradient.clone()));
1785 }
1786
1787 // Notify observers with current variable values and iteration number
1788 self.observers.notify(&state.variables, iteration);
1789
1790 // Check convergence
1791 let elapsed = start_time.elapsed();
1792
1793 // Compute parameter norm for relative parameter tolerance check
1794 let parameter_norm = Self::compute_parameter_norm(&state.variables);
1795
1796 // Compute costs for convergence check
1797 let new_cost = state.current_cost;
1798
1799 let cost_before_step = if step_eval.accepted {
1800 state.current_cost + step_eval.cost_reduction
1801 } else {
1802 state.current_cost
1803 };
1804
1805 if let Some(status) = self.check_convergence(
1806 iteration,
1807 cost_before_step,
1808 new_cost,
1809 parameter_norm,
1810 step_norm,
1811 step_result.gradient_norm,
1812 self.config.trust_region_radius,
1813 elapsed,
1814 step_eval.accepted,
1815 ) {
1816 let summary = self.create_summary(
1817 state.initial_cost,
1818 state.current_cost,
1819 iteration + 1,
1820 successful_steps,
1821 unsuccessful_steps,
1822 max_gradient_norm,
1823 final_gradient_norm,
1824 max_parameter_update_norm,
1825 final_parameter_update_norm,
1826 total_cost_reduction,
1827 elapsed,
1828 iteration_stats.clone(),
1829 status.clone(),
1830 );
1831
1832 // Print summary only if debug level is enabled
1833 if tracing::enabled!(tracing::Level::DEBUG) {
1834 debug!("{}", summary);
1835 }
1836
1837 // Compute covariances if enabled
1838 let covariances = if self.config.compute_covariances {
1839 problem.compute_and_set_covariances(
1840 &mut linear_solver,
1841 &mut state.variables,
1842 &state.variable_index_map,
1843 )
1844 } else {
1845 None
1846 };
1847
1848 return Ok(optimizer::SolverResult {
1849 status,
1850 iterations: iteration + 1,
1851 initial_cost: state.initial_cost,
1852 final_cost: state.current_cost,
1853 parameters: state.variables.into_iter().collect(),
1854 elapsed_time: elapsed,
1855 convergence_info: Some(optimizer::ConvergenceInfo {
1856 final_gradient_norm,
1857 final_parameter_update_norm,
1858 cost_evaluations,
1859 jacobian_evaluations,
1860 }),
1861 covariances,
1862 });
1863 }
1864
1865 // Note: Max iterations and timeout checks are now handled inside check_convergence()
1866
1867 iteration += 1;
1868 }
1869 }
1870}
1871
1872impl optimizer::Solver for DogLeg {
1873 type Config = DogLegConfig;
1874 type Error = error::ApexSolverError;
1875
1876 fn new() -> Self {
1877 Self::default()
1878 }
1879
1880 fn optimize(
1881 &mut self,
1882 problem: &problem::Problem,
1883 initial_params: &collections::HashMap<
1884 String,
1885 (manifold::ManifoldType, nalgebra::DVector<f64>),
1886 >,
1887 ) -> Result<
1888 optimizer::SolverResult<collections::HashMap<String, problem::VariableEnum>>,
1889 Self::Error,
1890 > {
1891 self.optimize(problem, initial_params)
1892 }
1893}
1894
1895#[cfg(test)]
1896mod tests {
1897 use super::*;
1898 use crate::{factors, manifold};
1899 use nalgebra;
1900
1901 /// Custom Rosenbrock Factor 1: r1 = 10(x2 - x1²)
1902 /// Demonstrates extensibility - custom factors can be defined outside of factors.rs
1903 #[derive(Debug, Clone)]
1904 struct RosenbrockFactor1;
1905
1906 impl factors::Factor for RosenbrockFactor1 {
1907 fn linearize(
1908 &self,
1909 params: &[nalgebra::DVector<f64>],
1910 compute_jacobian: bool,
1911 ) -> (nalgebra::DVector<f64>, Option<nalgebra::DMatrix<f64>>) {
1912 let x1 = params[0][0];
1913 let x2 = params[1][0];
1914
1915 // Residual: r1 = 10(x2 - x1²)
1916 let residual = nalgebra::dvector![10.0 * (x2 - x1 * x1)];
1917
1918 // Jacobian: ∂r1/∂x1 = -20*x1, ∂r1/∂x2 = 10
1919 let jacobian = if compute_jacobian {
1920 let mut jac = nalgebra::DMatrix::zeros(1, 2);
1921 jac[(0, 0)] = -20.0 * x1;
1922 jac[(0, 1)] = 10.0;
1923 Some(jac)
1924 } else {
1925 None
1926 };
1927
1928 (residual, jacobian)
1929 }
1930
1931 fn get_dimension(&self) -> usize {
1932 1
1933 }
1934 }
1935
1936 /// Custom Rosenbrock Factor 2: r2 = 1 - x1
1937 /// Demonstrates extensibility - custom factors can be defined outside of factors.rs
1938 #[derive(Debug, Clone)]
1939 struct RosenbrockFactor2;
1940
1941 impl factors::Factor for RosenbrockFactor2 {
1942 fn linearize(
1943 &self,
1944 params: &[nalgebra::DVector<f64>],
1945 compute_jacobian: bool,
1946 ) -> (nalgebra::DVector<f64>, Option<nalgebra::DMatrix<f64>>) {
1947 let x1 = params[0][0];
1948
1949 // Residual: r2 = 1 - x1
1950 let residual = nalgebra::dvector![1.0 - x1];
1951
1952 // Jacobian: ∂r2/∂x1 = -1
1953 let jacobian = if compute_jacobian {
1954 Some(nalgebra::DMatrix::from_element(1, 1, -1.0))
1955 } else {
1956 None
1957 };
1958
1959 (residual, jacobian)
1960 }
1961
1962 fn get_dimension(&self) -> usize {
1963 1
1964 }
1965 }
1966
1967 #[test]
1968 fn test_rosenbrock_optimization() -> Result<(), Box<dyn std::error::Error>> {
1969 // Rosenbrock function test:
1970 // Minimize: r1² + r2² where
1971 // r1 = 10(x2 - x1²)
1972 // r2 = 1 - x1
1973 // Starting point: [-1.2, 1.0]
1974 // Expected minimum: [1.0, 1.0]
1975
1976 let mut problem = problem::Problem::new();
1977 let mut initial_values = collections::HashMap::new();
1978
1979 // Add variables using Rn manifold (Euclidean space)
1980 initial_values.insert(
1981 "x1".to_string(),
1982 (manifold::ManifoldType::RN, nalgebra::dvector![-1.2]),
1983 );
1984 initial_values.insert(
1985 "x2".to_string(),
1986 (manifold::ManifoldType::RN, nalgebra::dvector![1.0]),
1987 );
1988
1989 // Add custom factors (demonstrates extensibility!)
1990 problem.add_residual_block(&["x1", "x2"], Box::new(RosenbrockFactor1), None);
1991 problem.add_residual_block(&["x1"], Box::new(RosenbrockFactor2), None);
1992
1993 // Configure Dog Leg optimizer with appropriate trust region
1994 let config = DogLegConfig::new()
1995 .with_max_iterations(100)
1996 .with_cost_tolerance(1e-8)
1997 .with_parameter_tolerance(1e-8)
1998 .with_gradient_tolerance(1e-10)
1999 .with_trust_region_radius(10.0); // Start with larger trust region
2000
2001 let mut solver = DogLeg::with_config(config);
2002 let result = solver.optimize(&problem, &initial_values)?;
2003
2004 // Extract final values
2005 let x1_final = result
2006 .parameters
2007 .get("x1")
2008 .ok_or("x1 not found")?
2009 .to_vector()[0];
2010 let x2_final = result
2011 .parameters
2012 .get("x2")
2013 .ok_or("x2 not found")?
2014 .to_vector()[0];
2015
2016 // Verify convergence to [1.0, 1.0]
2017 assert!(
2018 matches!(
2019 result.status,
2020 optimizer::OptimizationStatus::Converged
2021 | optimizer::OptimizationStatus::CostToleranceReached
2022 | optimizer::OptimizationStatus::ParameterToleranceReached
2023 | optimizer::OptimizationStatus::GradientToleranceReached
2024 ),
2025 "Optimization should converge"
2026 );
2027 assert!(
2028 (x1_final - 1.0).abs() < 1e-4,
2029 "x1 should converge to 1.0, got {}",
2030 x1_final
2031 );
2032 assert!(
2033 (x2_final - 1.0).abs() < 1e-4,
2034 "x2 should converge to 1.0, got {}",
2035 x2_final
2036 );
2037 assert!(
2038 result.final_cost < 1e-6,
2039 "Final cost should be near zero, got {}",
2040 result.final_cost
2041 );
2042 Ok(())
2043 }
2044}