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