Skip to main content

apex_solver/core/
problem.rs

1//! Optimization problem definition and sparse Jacobian computation.
2//!
3//! The `Problem` struct is the central component that defines a factor graph optimization problem.
4//! It manages residual blocks (constraints), variables, and the construction of sparse Jacobian
5//! matrices for efficient nonlinear least squares optimization.
6//!
7//! # Factor Graph Representation
8//!
9//! The optimization problem is represented as a bipartite factor graph:
10//!
11//! ```text
12//! Variables:  x0 --- x1 --- x2 --- x3
13//!              |      |      |      |
14//! Factors:    f0     f1     f2     f3 (constraints/measurements)
15//! ```
16//!
17//! Each factor connects one or more variables and contributes a residual (error) term to the
18//! overall cost function:
19//!
20//! ```text
21//! minimize Σ_i ρ(||r_i(x)||²)
22//! ```
23//!
24//! where `r_i(x)` is the residual for factor i, and ρ is an optional robust loss function.
25//!
26//! # Key Responsibilities
27//!
28//! 1. **Residual Block Management**: Add/remove factors and track their structure
29//! 2. **Variable Management**: Initialize variables with manifold types and constraints
30//! 3. **Sparsity Pattern**: Build symbolic structure for efficient sparse linear algebra
31//! 4. **Linearization**: Compute residuals and Jacobians in parallel
32//! 5. **Covariance**: Extract per-variable uncertainty estimates after optimization
33//!
34//! # Sparse Jacobian Structure
35//!
36//! The Jacobian matrix `J = ∂r/∂x` is sparse because each factor only depends on a small
37//! subset of variables. For example, a between factor connecting x0 and x1 contributes
38//! a 3×6 block (SE2) or 6×12 block (SE3) to the Jacobian, leaving the rest as zeros.
39//!
40//! The Problem pre-computes the sparsity pattern once, then efficiently fills in the
41//! numerical values during each iteration.
42//!
43//! # Mixed Manifold Support
44//!
45//! The Problem supports mixed manifold types in a single optimization problem via
46//! [`VariableEnum`]. This allows:
47//! - SE2 and SE3 poses in the same graph
48//! - SO3 rotations with R³ landmarks
49//! - Any combination of supported manifolds
50//!
51//! # Example: Building a Problem
52//!
53//! ```
54//! use apex_solver::core::problem::Problem;
55//! use apex_solver::factors::{BetweenFactor, PriorFactor};
56//! use apex_solver::core::loss_functions::HuberLoss;
57//! use apex_solver::manifold::ManifoldType;
58//! use nalgebra::{DVector, dvector};
59//! use std::collections::HashMap;
60//! use apex_solver::manifold::se2::SE2;
61//! use apex_solver::JacobianMode;
62//! # use apex_solver::core::CoreResult;
63//! # fn example() -> CoreResult<()> {
64//!
65//! let mut problem = Problem::new(JacobianMode::Sparse);
66//!
67//! // Add prior factor to anchor the first pose
68//! let prior = Box::new(PriorFactor {
69//!     data: dvector![0.0, 0.0, 0.0],
70//! });
71//! problem.add_residual_block(&["x0"], prior, None);
72//!
73//! // Add between factor with robust loss
74//! let between = Box::new(BetweenFactor::new(SE2::from_xy_angle(1.0, 0.0, 0.1)));
75//! let loss: Option<Box<dyn apex_solver::core::loss_functions::LossFunction + Send>> =
76//!     Some(Box::new(HuberLoss::new(1.0)?));
77//! problem.add_residual_block(&["x0", "x1"], between, loss);
78//!
79//! // Initialize variables
80//! let mut initial_values = HashMap::new();
81//! initial_values.insert("x0".to_string(), (ManifoldType::SE2, dvector![0.0, 0.0, 0.0]));
82//! initial_values.insert("x1".to_string(), (ManifoldType::SE2, dvector![0.9, 0.1, 0.12]));
83//!
84//! let variables = problem.initialize_variables(&initial_values);
85//! # Ok(())
86//! # }
87//! # example().unwrap();
88//! ```
89
90use std::{
91    collections::{HashMap, HashSet},
92    fs::File,
93    io::{Error, Write},
94    sync::{Arc, Mutex},
95};
96
97use faer::{Col, Mat, MatRef, sparse::SparseColMat};
98use nalgebra::DVector;
99use rayon::prelude::*;
100use tracing::warn;
101
102use crate::error::ErrorLogging;
103use crate::{
104    core::CoreResult,
105    core::{
106        CoreError, corrector::Corrector, loss_functions::LossFunction,
107        residual_block::ResidualBlock, variable::Variable,
108    },
109    factors::Factor,
110    linalg::{JacobianMode, LinearSolver, SparseMode, extract_variable_covariances},
111};
112use apex_manifolds::{ManifoldType, rn, se2, se3, se23, sgal3, sim3, so2, so3};
113
114// Re-export SymbolicStructure from the assembly module for backward compatibility
115pub use crate::linearizer::cpu::sparse::SymbolicStructure;
116
117/// Enum to handle mixed manifold variable types
118#[derive(Clone)]
119pub enum VariableEnum {
120    Rn(Variable<rn::Rn>),
121    SE2(Variable<se2::SE2>),
122    SE3(Variable<se3::SE3>),
123    SE23(Variable<se23::SE23>),
124    SGal3(Variable<sgal3::SGal3>),
125    Sim3(Variable<sim3::Sim3>),
126    SO2(Variable<so2::SO2>),
127    SO3(Variable<so3::SO3>),
128}
129
130impl VariableEnum {
131    /// Get the manifold type for this variable
132    pub fn manifold_type(&self) -> ManifoldType {
133        match self {
134            VariableEnum::Rn(_) => ManifoldType::RN,
135            VariableEnum::SE2(_) => ManifoldType::SE2,
136            VariableEnum::SE3(_) => ManifoldType::SE3,
137            VariableEnum::SE23(_) => ManifoldType::SE23,
138            VariableEnum::SGal3(_) => ManifoldType::SGal3,
139            VariableEnum::Sim3(_) => ManifoldType::Sim3,
140            VariableEnum::SO2(_) => ManifoldType::SO2,
141            VariableEnum::SO3(_) => ManifoldType::SO3,
142        }
143    }
144
145    /// Get the tangent space size for this variable
146    pub fn get_size(&self) -> usize {
147        match self {
148            VariableEnum::Rn(var) => var.get_size(),
149            VariableEnum::SE2(var) => var.get_size(),
150            VariableEnum::SE3(var) => var.get_size(),
151            VariableEnum::SE23(var) => var.get_size(),
152            VariableEnum::SGal3(var) => var.get_size(),
153            VariableEnum::Sim3(var) => var.get_size(),
154            VariableEnum::SO2(var) => var.get_size(),
155            VariableEnum::SO3(var) => var.get_size(),
156        }
157    }
158
159    /// Convert to DVector for use with Factor trait
160    pub fn to_vector(&self) -> DVector<f64> {
161        match self {
162            VariableEnum::Rn(var) => var.value.clone().into(),
163            VariableEnum::SE2(var) => var.value.clone().into(),
164            VariableEnum::SE3(var) => var.value.clone().into(),
165            VariableEnum::SE23(var) => var.value.clone().into(),
166            VariableEnum::SGal3(var) => var.value.clone().into(),
167            VariableEnum::Sim3(var) => var.value.clone().into(),
168            VariableEnum::SO2(var) => var.value.clone().into(),
169            VariableEnum::SO3(var) => var.value.clone().into(),
170        }
171    }
172
173    /// Apply a tangent space step to update this variable.
174    ///
175    /// This method applies a manifold plus operation: x_new = x ⊞ δx
176    /// where δx is a tangent vector. It supports all manifold types.
177    ///
178    /// # Arguments
179    /// * `step_slice` - View into the full step vector for this variable's DOF
180    ///
181    /// # Implementation Notes
182    /// Uses explicit clone instead of unsafe memory copy (`IntoNalgebra`) for small vectors.
183    /// This is safe and performant for typical manifold dimensions (1-6 DOF).
184    ///
185    pub fn apply_tangent_step(&mut self, step_slice: MatRef<f64>) {
186        match self {
187            VariableEnum::SE3(var) => {
188                let mut step_data = DVector::from_fn(6, |i, _| step_slice[(i, 0)]);
189                for &fixed_idx in &var.fixed_indices {
190                    if fixed_idx < 6 {
191                        step_data[fixed_idx] = 0.0;
192                    }
193                }
194                let tangent = se3::SE3Tangent::from(step_data);
195                let new_value = var.plus(&tangent);
196                var.set_value(new_value);
197            }
198            VariableEnum::SE2(var) => {
199                let mut step_data = DVector::from_fn(3, |i, _| step_slice[(i, 0)]);
200                for &fixed_idx in &var.fixed_indices {
201                    if fixed_idx < 3 {
202                        step_data[fixed_idx] = 0.0;
203                    }
204                }
205                let tangent = se2::SE2Tangent::from(step_data);
206                let new_value = var.plus(&tangent);
207                var.set_value(new_value);
208            }
209            VariableEnum::SO3(var) => {
210                let mut step_data = DVector::from_fn(3, |i, _| step_slice[(i, 0)]);
211                for &fixed_idx in &var.fixed_indices {
212                    if fixed_idx < 3 {
213                        step_data[fixed_idx] = 0.0;
214                    }
215                }
216                let tangent = so3::SO3Tangent::from(step_data);
217                let new_value = var.plus(&tangent);
218                var.set_value(new_value);
219            }
220            VariableEnum::SO2(var) => {
221                let mut step_data = DVector::from_fn(1, |i, _| step_slice[(i, 0)]);
222                for &fixed_idx in &var.fixed_indices {
223                    if fixed_idx < 1 {
224                        step_data[fixed_idx] = 0.0;
225                    }
226                }
227                let tangent = so2::SO2Tangent::from(step_data);
228                let new_value = var.plus(&tangent);
229                var.set_value(new_value);
230            }
231            VariableEnum::SE23(var) => {
232                // SE_2(3) has 9 DOF in tangent space
233                let mut step_data = DVector::from_fn(9, |i, _| step_slice[(i, 0)]);
234
235                // Enforce fixed indices: zero out step components for fixed DOF
236                for &fixed_idx in &var.fixed_indices {
237                    if fixed_idx < 9 {
238                        step_data[fixed_idx] = 0.0;
239                    }
240                }
241
242                let tangent = se23::SE23Tangent::from(step_data);
243                let new_value = var.plus(&tangent);
244                var.set_value(new_value);
245            }
246            VariableEnum::SGal3(var) => {
247                // SGal(3) has 10 DOF in tangent space
248                let mut step_data = DVector::from_fn(10, |i, _| step_slice[(i, 0)]);
249
250                // Enforce fixed indices: zero out step components for fixed DOF
251                for &fixed_idx in &var.fixed_indices {
252                    if fixed_idx < 10 {
253                        step_data[fixed_idx] = 0.0;
254                    }
255                }
256
257                let tangent = sgal3::SGal3Tangent::from(step_data);
258                let new_value = var.plus(&tangent);
259                var.set_value(new_value);
260            }
261            VariableEnum::Sim3(var) => {
262                // Sim(3) has 7 DOF in tangent space
263                let mut step_data = DVector::from_fn(7, |i, _| step_slice[(i, 0)]);
264
265                // Enforce fixed indices: zero out step components for fixed DOF
266                for &fixed_idx in &var.fixed_indices {
267                    if fixed_idx < 7 {
268                        step_data[fixed_idx] = 0.0;
269                    }
270                }
271
272                let tangent = sim3::Sim3Tangent::from(step_data);
273                let new_value = var.plus(&tangent);
274                var.set_value(new_value);
275            }
276            VariableEnum::Rn(var) => {
277                let size = var.get_size();
278                let mut step_data = DVector::from_fn(size, |i, _| step_slice[(i, 0)]);
279                for &fixed_idx in &var.fixed_indices {
280                    if fixed_idx < size {
281                        step_data[fixed_idx] = 0.0;
282                    }
283                }
284                let tangent = rn::RnTangent::new(step_data);
285                let new_value = var.plus(&tangent);
286                var.set_value(new_value);
287            }
288        }
289    }
290
291    /// Get the covariance matrix for this variable (if computed).
292    ///
293    /// Returns `None` if covariance has not been computed.
294    ///
295    /// # Returns
296    /// Reference to the covariance matrix in tangent space
297    pub fn get_covariance(&self) -> Option<&Mat<f64>> {
298        match self {
299            VariableEnum::Rn(var) => var.get_covariance(),
300            VariableEnum::SE2(var) => var.get_covariance(),
301            VariableEnum::SE3(var) => var.get_covariance(),
302            VariableEnum::SE23(var) => var.get_covariance(),
303            VariableEnum::SGal3(var) => var.get_covariance(),
304            VariableEnum::Sim3(var) => var.get_covariance(),
305            VariableEnum::SO2(var) => var.get_covariance(),
306            VariableEnum::SO3(var) => var.get_covariance(),
307        }
308    }
309
310    /// Set the covariance matrix for this variable.
311    ///
312    /// The covariance matrix should be square with dimension equal to
313    /// the tangent space dimension of this variable.
314    ///
315    /// # Arguments
316    /// * `cov` - Covariance matrix in tangent space
317    pub fn set_covariance(&mut self, cov: Mat<f64>) {
318        match self {
319            VariableEnum::Rn(var) => var.set_covariance(cov),
320            VariableEnum::SE2(var) => var.set_covariance(cov),
321            VariableEnum::SE3(var) => var.set_covariance(cov),
322            VariableEnum::SE23(var) => var.set_covariance(cov),
323            VariableEnum::SGal3(var) => var.set_covariance(cov),
324            VariableEnum::Sim3(var) => var.set_covariance(cov),
325            VariableEnum::SO2(var) => var.set_covariance(cov),
326            VariableEnum::SO3(var) => var.set_covariance(cov),
327        }
328    }
329
330    /// Clear the covariance matrix for this variable.
331    pub fn clear_covariance(&mut self) {
332        match self {
333            VariableEnum::Rn(var) => var.clear_covariance(),
334            VariableEnum::SE2(var) => var.clear_covariance(),
335            VariableEnum::SE3(var) => var.clear_covariance(),
336            VariableEnum::SE23(var) => var.clear_covariance(),
337            VariableEnum::SGal3(var) => var.clear_covariance(),
338            VariableEnum::Sim3(var) => var.clear_covariance(),
339            VariableEnum::SO2(var) => var.clear_covariance(),
340            VariableEnum::SO3(var) => var.clear_covariance(),
341        }
342    }
343
344    /// Get the bounds for this variable.
345    ///
346    /// Returns a reference to the bounds map where keys are indices and values are (lower, upper) pairs.
347    pub fn get_bounds(&self) -> &HashMap<usize, (f64, f64)> {
348        match self {
349            VariableEnum::Rn(var) => &var.bounds,
350            VariableEnum::SE2(var) => &var.bounds,
351            VariableEnum::SE3(var) => &var.bounds,
352            VariableEnum::SE23(var) => &var.bounds,
353            VariableEnum::SGal3(var) => &var.bounds,
354            VariableEnum::Sim3(var) => &var.bounds,
355            VariableEnum::SO2(var) => &var.bounds,
356            VariableEnum::SO3(var) => &var.bounds,
357        }
358    }
359
360    /// Get the fixed indices for this variable.
361    ///
362    /// Returns a reference to the set of indices that should remain fixed during optimization.
363    pub fn get_fixed_indices(&self) -> &HashSet<usize> {
364        match self {
365            VariableEnum::Rn(var) => &var.fixed_indices,
366            VariableEnum::SE2(var) => &var.fixed_indices,
367            VariableEnum::SE3(var) => &var.fixed_indices,
368            VariableEnum::SE23(var) => &var.fixed_indices,
369            VariableEnum::SGal3(var) => &var.fixed_indices,
370            VariableEnum::Sim3(var) => &var.fixed_indices,
371            VariableEnum::SO2(var) => &var.fixed_indices,
372            VariableEnum::SO3(var) => &var.fixed_indices,
373        }
374    }
375
376    /// Set the value of this variable from a vector representation.
377    ///
378    /// This is used to update the variable after applying constraints (bounds and fixed indices).
379    pub fn set_from_vector(&mut self, vec: &DVector<f64>) {
380        match self {
381            VariableEnum::Rn(var) => {
382                var.set_value(rn::Rn::new(vec.clone()));
383            }
384            VariableEnum::SE2(var) => {
385                let new_se2: se2::SE2 = vec.clone().into();
386                var.set_value(new_se2);
387            }
388            VariableEnum::SE3(var) => {
389                let new_se3: se3::SE3 = vec.clone().into();
390                var.set_value(new_se3);
391            }
392            VariableEnum::SE23(var) => {
393                let new_se23: se23::SE23 = vec.clone().into();
394                var.set_value(new_se23);
395            }
396            VariableEnum::SGal3(var) => {
397                let new_sgal3: sgal3::SGal3 = vec.clone().into();
398                var.set_value(new_sgal3);
399            }
400            VariableEnum::Sim3(var) => {
401                let new_sim3: sim3::Sim3 = vec.clone().into();
402                var.set_value(new_sim3);
403            }
404            VariableEnum::SO2(var) => {
405                let new_so2: so2::SO2 = vec.clone().into();
406                var.set_value(new_so2);
407            }
408            VariableEnum::SO3(var) => {
409                let new_so3: so3::SO3 = vec.clone().into();
410                var.set_value(new_so3);
411            }
412        }
413    }
414}
415
416/// The optimization problem definition for factor graph optimization.
417///
418/// Manages residual blocks (factors/constraints), variables, and the sparse Jacobian structure.
419/// Supports mixed manifold types (SE2, SE3, SO2, SO3, Rn) in a single problem and provides
420/// efficient parallel residual/Jacobian computation.
421///
422/// # Architecture
423///
424/// The Problem acts as a container and coordinator:
425/// - Stores all residual blocks (factors with optional loss functions)
426/// - Tracks the global structure (which variables connect to which factors)
427/// - Builds and maintains the sparse Jacobian pattern
428/// - Provides parallel residual/Jacobian evaluation using rayon
429/// - Manages variable constraints (fixed indices, bounds)
430///
431/// # Workflow
432///
433/// 1. **Construction**: Create a new Problem with `Problem::new(JacobianMode::Sparse)`
434/// 2. **Add Factors**: Use `add_residual_block()` to add constraints
435/// 3. **Initialize Variables**: Use `initialize_variables()` with initial values
436/// 4. **Build Sparsity**: Use `build_symbolic_structure()` once before optimization
437/// 5. **Linearize**: Call `compute_residual_and_jacobian_sparse()` each iteration
438/// 6. **Extract Covariance**: Use `compute_and_set_covariances()` after convergence
439///
440/// # Example
441///
442/// ```
443/// use apex_solver::core::problem::Problem;
444/// use apex_solver::factors::BetweenFactor;
445/// use apex_solver::manifold::ManifoldType;
446/// use apex_solver::manifold::se2::SE2;
447/// use apex_solver::JacobianMode;
448/// use nalgebra::dvector;
449/// use std::collections::HashMap;
450///
451/// let mut problem = Problem::new(JacobianMode::Sparse);
452///
453/// // Add a between factor
454/// let factor = Box::new(BetweenFactor::new(SE2::from_xy_angle(1.0, 0.0, 0.1)));
455/// problem.add_residual_block(&["x0", "x1"], factor, None);
456///
457/// // Initialize variables
458/// let mut initial = HashMap::new();
459/// initial.insert("x0".to_string(), (ManifoldType::SE2, dvector![0.0, 0.0, 0.0]));
460/// initial.insert("x1".to_string(), (ManifoldType::SE2, dvector![1.0, 0.0, 0.1]));
461///
462/// let variables = problem.initialize_variables(&initial);
463/// assert_eq!(variables.len(), 2);
464/// ```
465pub struct Problem {
466    /// Total dimension of the stacked residual vector (sum of all residual block dimensions)
467    pub(crate) total_residual_dimension: usize,
468
469    /// Controls which Jacobian assembly strategy is used (sparse or dense).
470    ///
471    /// Set via `Problem::new(JacobianMode::Sparse)` (sparse, default) or
472    /// `Problem::new(JacobianMode::Dense)` (dense).
473    /// The optimizer reads this field to select the assembly path.
474    pub(crate) jacobian_mode: JacobianMode,
475
476    /// Counter for assigning unique IDs to residual blocks
477    residual_id_count: usize,
478
479    /// Map from residual block ID to ResidualBlock instance
480    residual_blocks: HashMap<usize, ResidualBlock>,
481
482    /// Variables with fixed indices (e.g., fix first pose's x,y coordinates)
483    /// Maps variable name -> set of indices to fix
484    pub(crate) fixed_variable_indexes: HashMap<String, HashSet<usize>>,
485
486    /// Variable bounds (box constraints on individual DOF)
487    /// Maps variable name -> (index -> (lower_bound, upper_bound))
488    pub(crate) variable_bounds: HashMap<String, HashMap<usize, (f64, f64)>>,
489}
490impl Default for Problem {
491    fn default() -> Self {
492        Self::new(JacobianMode::Sparse)
493    }
494}
495
496impl Problem {
497    /// Create a new empty optimization problem with the given Jacobian assembly mode.
498    ///
499    /// # Arguments
500    ///
501    /// * `jacobian_mode` - Assembly strategy: [`JacobianMode::Sparse`] (default, large problems)
502    ///   or [`JacobianMode::Dense`] (small-to-medium problems < ~500 DOF)
503    ///
504    /// # Example
505    ///
506    /// ```
507    /// use apex_solver::core::problem::Problem;
508    /// use apex_solver::JacobianMode;
509    ///
510    /// // Sparse (default for large-scale problems)
511    /// let sparse_problem = Problem::new(JacobianMode::Sparse);
512    /// assert_eq!(sparse_problem.num_residual_blocks(), 0);
513    ///
514    /// // Dense (for small-to-medium problems)
515    /// let dense_problem = Problem::new(JacobianMode::Dense);
516    /// assert_eq!(dense_problem.num_residual_blocks(), 0);
517    /// ```
518    pub fn new(jacobian_mode: JacobianMode) -> Self {
519        Self {
520            total_residual_dimension: 0,
521            jacobian_mode,
522            residual_id_count: 0,
523            residual_blocks: HashMap::new(),
524            fixed_variable_indexes: HashMap::new(),
525            variable_bounds: HashMap::new(),
526        }
527    }
528
529    /// Add a residual block (factor with optional loss function) to the problem.
530    ///
531    /// This is the primary method for building the factor graph. Each call adds one constraint
532    /// connecting one or more variables.
533    ///
534    /// # Arguments
535    ///
536    /// * `variable_key_size_list` - Names of the variables this factor connects (order matters)
537    /// * `factor` - The factor implementation that computes residuals and Jacobians
538    /// * `loss_func` - Optional robust loss function for outlier rejection
539    ///
540    /// # Returns
541    ///
542    /// The unique ID assigned to this residual block
543    ///
544    /// # Example
545    ///
546    /// ```
547    /// use apex_solver::core::problem::Problem;
548    /// use apex_solver::factors::{BetweenFactor, PriorFactor};
549    /// use apex_solver::core::loss_functions::HuberLoss;
550    /// use apex_solver::JacobianMode;
551    /// use nalgebra::dvector;
552    /// use apex_solver::manifold::se2::SE2;
553    /// # use apex_solver::core::CoreResult;
554    /// # fn example() -> CoreResult<()> {
555    ///
556    /// let mut problem = Problem::new(JacobianMode::Sparse);
557    ///
558    /// // Add prior factor (unary constraint)
559    /// let prior = Box::new(PriorFactor { data: dvector![0.0, 0.0, 0.0] });
560    /// let id1 = problem.add_residual_block(&["x0"], prior, None);
561    ///
562    /// // Add between factor with robust loss (binary constraint)
563    /// let between = Box::new(BetweenFactor::new(SE2::from_xy_angle(1.0, 0.0, 0.1)));
564    /// let loss: Option<Box<dyn apex_solver::core::loss_functions::LossFunction + Send>> =
565    ///     Some(Box::new(HuberLoss::new(1.0)?));
566    /// let id2 = problem.add_residual_block(&["x0", "x1"], between, loss);
567    ///
568    /// assert_eq!(id1, 0);
569    /// assert_eq!(id2, 1);
570    /// assert_eq!(problem.num_residual_blocks(), 2);
571    /// # Ok(())
572    /// # }
573    /// # example().unwrap();
574    /// ```
575    pub fn add_residual_block(
576        &mut self,
577        variable_key_size_list: &[&str],
578        factor: Box<dyn Factor + Send>,
579        loss_func: Option<Box<dyn LossFunction + Send>>,
580    ) -> usize {
581        let new_residual_dimension = factor.get_dimension();
582        self.residual_blocks.insert(
583            self.residual_id_count,
584            ResidualBlock::new(
585                self.residual_id_count,
586                self.total_residual_dimension,
587                variable_key_size_list,
588                factor,
589                loss_func,
590            ),
591        );
592        let block_id = self.residual_id_count;
593        self.residual_id_count += 1;
594
595        self.total_residual_dimension += new_residual_dimension;
596
597        block_id
598    }
599
600    pub fn remove_residual_block(&mut self, block_id: usize) -> Option<ResidualBlock> {
601        if let Some(residual_block) = self.residual_blocks.remove(&block_id) {
602            self.total_residual_dimension -= residual_block.factor.get_dimension();
603            Some(residual_block)
604        } else {
605            None
606        }
607    }
608
609    pub fn fix_variable(&mut self, var_to_fix: &str, idx: usize) {
610        if let Some(var_mut) = self.fixed_variable_indexes.get_mut(var_to_fix) {
611            var_mut.insert(idx);
612        } else {
613            self.fixed_variable_indexes
614                .insert(var_to_fix.to_owned(), HashSet::from([idx]));
615        }
616    }
617
618    pub fn unfix_variable(&mut self, var_to_unfix: &str) {
619        self.fixed_variable_indexes.remove(var_to_unfix);
620    }
621
622    pub fn set_variable_bounds(
623        &mut self,
624        var_to_bound: &str,
625        idx: usize,
626        lower_bound: f64,
627        upper_bound: f64,
628    ) {
629        if lower_bound > upper_bound {
630            warn!("lower bound is larger than upper bound");
631        } else if let Some(var_mut) = self.variable_bounds.get_mut(var_to_bound) {
632            var_mut.insert(idx, (lower_bound, upper_bound));
633        } else {
634            self.variable_bounds.insert(
635                var_to_bound.to_owned(),
636                HashMap::from([(idx, (lower_bound, upper_bound))]),
637            );
638        }
639    }
640
641    pub fn remove_variable_bounds(&mut self, var_to_unbound: &str) {
642        self.variable_bounds.remove(var_to_unbound);
643    }
644
645    /// Initialize variables from initial values with manifold types.
646    ///
647    /// Converts raw initial values into typed `Variable<M>` instances wrapped in `VariableEnum`.
648    /// This method also applies any fixed indices or bounds that were set via `fix_variable()`
649    /// or `set_variable_bounds()`.
650    ///
651    /// # Arguments
652    ///
653    /// * `initial_values` - Map from variable name to (manifold type, initial value vector)
654    ///
655    /// # Returns
656    ///
657    /// Map from variable name to `VariableEnum` (typed variables ready for optimization)
658    ///
659    /// # Manifold Formats
660    ///
661    /// - **SE2**: `[x, y, theta]` (3 elements)
662    /// - **SE3**: `[tx, ty, tz, qw, qx, qy, qz]` (7 elements)
663    /// - **SO2**: `[theta]` (1 element)
664    /// - **SO3**: `[qw, qx, qy, qz]` (4 elements)
665    /// - **Rn**: `[x1, x2, ..., xn]` (n elements)
666    ///
667    /// # Example
668    ///
669    /// ```
670    /// use apex_solver::core::problem::Problem;
671    /// use apex_solver::manifold::ManifoldType;
672    /// use apex_solver::JacobianMode;
673    /// use nalgebra::dvector;
674    /// use std::collections::HashMap;
675    ///
676    /// let problem = Problem::new(JacobianMode::Sparse);
677    ///
678    /// let mut initial = HashMap::new();
679    /// initial.insert("pose0".to_string(), (ManifoldType::SE2, dvector![0.0, 0.0, 0.0]));
680    /// initial.insert("pose1".to_string(), (ManifoldType::SE2, dvector![1.0, 0.0, 0.1]));
681    /// initial.insert("landmark".to_string(), (ManifoldType::RN, dvector![5.0, 3.0]));
682    ///
683    /// let variables = problem.initialize_variables(&initial);
684    /// assert_eq!(variables.len(), 3);
685    /// ```
686    pub fn initialize_variables(
687        &self,
688        initial_values: &HashMap<String, (ManifoldType, DVector<f64>)>,
689    ) -> HashMap<String, VariableEnum> {
690        let variables: HashMap<String, VariableEnum> = initial_values
691            .iter()
692            .map(|(k, v)| {
693                let variable_enum = match v.0 {
694                    ManifoldType::SO2 => {
695                        assert_eq!(
696                            v.1.len(),
697                            1,
698                            "Variable '{}': SO2 expects 1 element, got {}",
699                            k,
700                            v.1.len()
701                        );
702                        let mut var = Variable::new(so2::SO2::from(v.1.clone()));
703                        if let Some(indexes) = self.fixed_variable_indexes.get(k) {
704                            var.fixed_indices = indexes.clone();
705                        }
706                        if let Some(bounds) = self.variable_bounds.get(k) {
707                            var.bounds = bounds.clone();
708                        }
709                        VariableEnum::SO2(var)
710                    }
711                    ManifoldType::SO3 => {
712                        assert_eq!(
713                            v.1.len(),
714                            4,
715                            "Variable '{}': SO3 expects 4 elements, got {}",
716                            k,
717                            v.1.len()
718                        );
719                        let mut var = Variable::new(so3::SO3::from(v.1.clone()));
720                        if let Some(indexes) = self.fixed_variable_indexes.get(k) {
721                            var.fixed_indices = indexes.clone();
722                        }
723                        if let Some(bounds) = self.variable_bounds.get(k) {
724                            var.bounds = bounds.clone();
725                        }
726                        VariableEnum::SO3(var)
727                    }
728                    ManifoldType::SE2 => {
729                        assert_eq!(
730                            v.1.len(),
731                            3,
732                            "Variable '{}': SE2 expects 3 elements, got {}",
733                            k,
734                            v.1.len()
735                        );
736                        let mut var = Variable::new(se2::SE2::from(v.1.clone()));
737                        if let Some(indexes) = self.fixed_variable_indexes.get(k) {
738                            var.fixed_indices = indexes.clone();
739                        }
740                        if let Some(bounds) = self.variable_bounds.get(k) {
741                            var.bounds = bounds.clone();
742                        }
743                        VariableEnum::SE2(var)
744                    }
745                    ManifoldType::SE3 => {
746                        assert_eq!(
747                            v.1.len(),
748                            7,
749                            "Variable '{}': SE3 expects 7 elements, got {}",
750                            k,
751                            v.1.len()
752                        );
753                        let mut var = Variable::new(se3::SE3::from(v.1.clone()));
754                        if let Some(indexes) = self.fixed_variable_indexes.get(k) {
755                            var.fixed_indices = indexes.clone();
756                        }
757                        if let Some(bounds) = self.variable_bounds.get(k) {
758                            var.bounds = bounds.clone();
759                        }
760                        VariableEnum::SE3(var)
761                    }
762                    ManifoldType::RN => {
763                        let mut var = Variable::new(rn::Rn::from(v.1.clone()));
764                        if let Some(indexes) = self.fixed_variable_indexes.get(k) {
765                            var.fixed_indices = indexes.clone();
766                        }
767                        if let Some(bounds) = self.variable_bounds.get(k) {
768                            var.bounds = bounds.clone();
769                        }
770                        VariableEnum::Rn(var)
771                    }
772                    ManifoldType::SE23 => {
773                        let mut var = Variable::new(se23::SE23::from(v.1.clone()));
774                        if let Some(indexes) = self.fixed_variable_indexes.get(k) {
775                            var.fixed_indices = indexes.clone();
776                        }
777                        if let Some(bounds) = self.variable_bounds.get(k) {
778                            var.bounds = bounds.clone();
779                        }
780                        VariableEnum::SE23(var)
781                    }
782                    ManifoldType::SGal3 => {
783                        let mut var = Variable::new(sgal3::SGal3::from(v.1.clone()));
784                        if let Some(indexes) = self.fixed_variable_indexes.get(k) {
785                            var.fixed_indices = indexes.clone();
786                        }
787                        if let Some(bounds) = self.variable_bounds.get(k) {
788                            var.bounds = bounds.clone();
789                        }
790                        VariableEnum::SGal3(var)
791                    }
792                    ManifoldType::Sim3 => {
793                        let mut var = Variable::new(sim3::Sim3::from(v.1.clone()));
794                        if let Some(indexes) = self.fixed_variable_indexes.get(k) {
795                            var.fixed_indices = indexes.clone();
796                        }
797                        if let Some(bounds) = self.variable_bounds.get(k) {
798                            var.bounds = bounds.clone();
799                        }
800                        VariableEnum::Sim3(var)
801                    }
802                };
803
804                (k.to_owned(), variable_enum)
805            })
806            .collect();
807        variables
808    }
809
810    /// Get the number of residual blocks
811    pub fn num_residual_blocks(&self) -> usize {
812        self.residual_blocks.len()
813    }
814
815    /// Get a reference to the residual blocks map (crate-internal for assembly modules)
816    pub(crate) fn residual_blocks(&self) -> &HashMap<usize, ResidualBlock> {
817        &self.residual_blocks
818    }
819
820    /// Compute only the residual vector for the current variable values.
821    ///
822    /// This is an optimized version that skips Jacobian computation when only the cost
823    /// function value is needed (e.g., during initialization or step evaluation).
824    ///
825    /// # Arguments
826    ///
827    /// * `variables` - Current variable values (from `initialize_variables()` or updated)
828    ///
829    /// # Returns
830    ///
831    /// Residual vector as N×1 column matrix (N = total residual dimension)
832    ///
833    /// # Performance
834    ///
835    /// Approximately **2x faster** than `compute_residual_and_jacobian_sparse()` since it:
836    /// - Skips Jacobian computation for each residual block
837    /// - Avoids Jacobian matrix assembly and storage
838    /// - Only parallelizes residual evaluation
839    ///
840    /// # When to Use
841    ///
842    /// - **Initial cost computation**: When setting up optimization state
843    /// - **Step evaluation**: When computing new cost after applying parameter updates
844    /// - **Cost-only queries**: When you don't need gradients
845    ///
846    /// Use `compute_residual_and_jacobian_sparse()` when you need both residual and Jacobian
847    /// (e.g., in the main optimization iteration loop for linearization).
848    ///
849    /// # Example
850    ///
851    /// ```no_run
852    /// # use apex_solver::core::problem::Problem;
853    /// # use apex_solver::JacobianMode;
854    /// # use std::collections::HashMap;
855    /// # fn example() -> Result<(), Box<dyn std::error::Error>> {
856    /// # let problem = Problem::new(JacobianMode::Sparse);
857    /// # let variables = HashMap::new();
858    /// // Initial cost evaluation (no Jacobian needed)
859    /// let residual = problem.compute_residual_sparse(&variables)?;
860    /// let initial_cost = residual.norm_l2() * residual.norm_l2();
861    /// # Ok(())
862    /// # }
863    /// ```
864    pub fn compute_residual_sparse(
865        &self,
866        variables: &HashMap<String, VariableEnum>,
867    ) -> CoreResult<Mat<f64>> {
868        let total_residual = Arc::new(Mutex::new(Col::<f64>::zeros(self.total_residual_dimension)));
869
870        // Compute residuals in parallel (no Jacobian needed)
871        let result: Result<Vec<()>, CoreError> = self
872            .residual_blocks
873            .par_iter()
874            .map(|(_, residual_block)| {
875                self.compute_residual_block(residual_block, variables, &total_residual)
876            })
877            .collect();
878
879        result?;
880
881        let total_residual = Arc::try_unwrap(total_residual)
882            .map_err(|_| {
883                CoreError::ParallelComputation(
884                    "Failed to unwrap Arc for total residual".to_string(),
885                )
886                .log()
887            })?
888            .into_inner()
889            .map_err(|e| {
890                CoreError::ParallelComputation(
891                    "Failed to extract mutex inner value for total residual".to_string(),
892                )
893                .log_with_source(e)
894            })?;
895
896        // Convert faer Col to Mat (column vector as n×1 matrix)
897        let residual_faer = total_residual.as_ref().as_mat().to_owned();
898        Ok(residual_faer)
899    }
900
901    /// Compute residual vector and sparse Jacobian matrix for the current variable values.
902    ///
903    /// This is the core linearization method called during each optimization iteration. It:
904    /// 1. Evaluates all residual blocks in parallel using rayon
905    /// 2. Assembles the full residual vector
906    /// 3. Constructs the sparse Jacobian matrix using the precomputed symbolic structure
907    ///
908    /// # Arguments
909    ///
910    /// * `variables` - Current variable values (from `initialize_variables()` or updated)
911    /// * `variable_index_sparce_matrix` - Map from variable name to starting column in Jacobian
912    /// * `symbolic_structure` - Precomputed sparsity pattern (from `build_symbolic_structure()`)
913    ///
914    /// # Returns
915    ///
916    /// Tuple `(residual, jacobian)` where:
917    /// - `residual`: N×1 column matrix (total residual dimension)
918    /// - `jacobian`: N×M sparse matrix (N = residual dim, M = total DOF)
919    ///
920    /// # Performance
921    ///
922    /// This method is highly optimized:
923    /// - **Parallel evaluation**: Each residual block is evaluated independently using rayon
924    /// - **Sparse storage**: Only non-zero Jacobian entries are stored and computed
925    /// - **Memory efficient**: Preallocated sparse structure avoids dynamic allocations
926    ///
927    /// Typically accounts for 40-60% of total optimization time (including sparse matrix ops).
928    ///
929    /// # When to Use
930    ///
931    /// Use this method in the main optimization loop when you need both residual and Jacobian
932    /// for linearization. For cost-only evaluation, use `compute_residual_sparse()` instead.
933    ///
934    /// # Example
935    ///
936    /// ```
937    /// # use apex_solver::core::problem::Problem;
938    /// # use std::collections::HashMap;
939    /// // Inside optimizer loop, compute both residual and Jacobian for linearization
940    /// // let (residual, jacobian) = problem.compute_residual_and_jacobian_sparse(
941    /// //     &variables,
942    /// //     &variable_index_map,
943    /// //     &symbolic_structure,
944    /// // )?;
945    /// //
946    /// // Use for linear system: J^T J dx = -J^T r
947    /// ```
948    /// Compute residuals and sparse Jacobian.
949    ///
950    /// Delegates to [`assembly::sparse::assemble_sparse()`](super::assembly::sparse::assemble_sparse).
951    pub fn compute_residual_and_jacobian_sparse(
952        &self,
953        variables: &HashMap<String, VariableEnum>,
954        variable_index_sparce_matrix: &HashMap<String, usize>,
955        symbolic_structure: &SymbolicStructure,
956    ) -> CoreResult<(Mat<f64>, SparseColMat<usize, f64>)> {
957        Ok(crate::linearizer::cpu::sparse::assemble_sparse(
958            self,
959            variables,
960            variable_index_sparce_matrix,
961            symbolic_structure,
962        )?)
963    }
964
965    /// Compute residuals and dense Jacobian.
966    ///
967    /// Delegates to [`assembly::dense::assemble_dense()`](super::assembly::dense::assemble_dense).
968    pub fn compute_residual_and_jacobian_dense(
969        &self,
970        variables: &HashMap<String, VariableEnum>,
971        variable_index_map: &HashMap<String, usize>,
972        total_dof: usize,
973    ) -> CoreResult<(Mat<f64>, Mat<f64>)> {
974        Ok(crate::linearizer::cpu::dense::assemble_dense(
975            self,
976            variables,
977            variable_index_map,
978            total_dof,
979        )?)
980    }
981
982    /// Compute only the residual for a single residual block (no Jacobian).
983    ///
984    /// Helper method for `compute_residual_sparse()`.
985    fn compute_residual_block(
986        &self,
987        residual_block: &ResidualBlock,
988        variables: &HashMap<String, VariableEnum>,
989        total_residual: &Arc<Mutex<Col<f64>>>,
990    ) -> CoreResult<()> {
991        let mut param_vectors: Vec<DVector<f64>> = Vec::new();
992
993        for var_key in &residual_block.variable_key_list {
994            if let Some(variable) = variables.get(var_key) {
995                param_vectors.push(variable.to_vector());
996            }
997        }
998
999        // Compute only residual (linearize still computes Jacobian internally,
1000        // but we don't extract/store it)
1001        let (mut res, _) = residual_block.factor.linearize(&param_vectors, false);
1002
1003        // Apply loss function if present (critical for robust optimization)
1004        if let Some(loss_func) = &residual_block.loss_func {
1005            let squared_norm = res.dot(&res);
1006            let corrector = Corrector::new(loss_func.as_ref(), squared_norm);
1007            corrector.correct_residuals(&mut res);
1008        }
1009
1010        let mut total_residual = total_residual.lock().map_err(|e| {
1011            CoreError::ParallelComputation("Failed to acquire lock on total residual".to_string())
1012                .log_with_source(e)
1013        })?;
1014
1015        // Copy residual values from nalgebra DVector to faer Col
1016        let start_idx = residual_block.residual_row_start_idx;
1017        let dim = residual_block.factor.get_dimension();
1018        let mut total_residual_mut = total_residual.as_mut();
1019        for i in 0..dim {
1020            total_residual_mut[start_idx + i] = res[i];
1021        }
1022
1023        Ok(())
1024    }
1025
1026    /// Log residual vector to a text file
1027    pub fn log_residual_to_file(
1028        &self,
1029        residual: &nalgebra::DVector<f64>,
1030        filename: &str,
1031    ) -> Result<(), Error> {
1032        let mut file = File::create(filename)?;
1033        writeln!(file, "# Residual vector - {} elements", residual.len())?;
1034        for (i, &value) in residual.iter().enumerate() {
1035            writeln!(file, "{}: {:.12}", i, value)?;
1036        }
1037        Ok(())
1038    }
1039
1040    /// Log sparse Jacobian matrix to a text file
1041    pub fn log_sparse_jacobian_to_file(
1042        &self,
1043        jacobian: &SparseColMat<usize, f64>,
1044        filename: &str,
1045    ) -> Result<(), Error> {
1046        let mut file = File::create(filename)?;
1047        writeln!(
1048            file,
1049            "# Sparse Jacobian matrix - {} x {} ({} non-zeros)",
1050            jacobian.nrows(),
1051            jacobian.ncols(),
1052            jacobian.compute_nnz()
1053        )?;
1054        writeln!(file, "# Matrix saved as dimensions and non-zero count only")?;
1055        writeln!(file, "# For detailed access, convert to dense matrix first")?;
1056        Ok(())
1057    }
1058
1059    /// Log variables to a text file
1060    pub fn log_variables_to_file(
1061        &self,
1062        variables: &HashMap<String, VariableEnum>,
1063        filename: &str,
1064    ) -> Result<(), Error> {
1065        let mut file = File::create(filename)?;
1066        writeln!(file, "# Variables - {} total", variables.len())?;
1067        writeln!(file, "# Format: variable_name: [values...]")?;
1068
1069        let mut sorted_vars: Vec<_> = variables.keys().collect();
1070        sorted_vars.sort();
1071
1072        for var_name in sorted_vars {
1073            let var_vector = variables[var_name].to_vector();
1074            write!(file, "{}: [", var_name)?;
1075            for (i, &value) in var_vector.iter().enumerate() {
1076                write!(file, "{:.12}", value)?;
1077                if i < var_vector.len() - 1 {
1078                    write!(file, ", ")?;
1079                }
1080            }
1081            writeln!(file, "]")?;
1082        }
1083        Ok(())
1084    }
1085
1086    /// Compute per-variable covariances and set them in Variable objects
1087    ///
1088    /// This method computes the full covariance matrix by inverting the Hessian
1089    /// from the linear solver, then extracts per-variable covariance blocks and
1090    /// stores them in the corresponding Variable objects.
1091    ///
1092    /// # Arguments
1093    /// * `linear_solver` - Mutable reference to the linear solver containing the cached Hessian
1094    /// * `variables` - Mutable map of variables where covariances will be stored
1095    /// * `variable_index_map` - Map from variable names to their starting column indices
1096    ///
1097    /// # Returns
1098    /// `Some(HashMap)` containing per-variable covariance matrices if successful, `None` otherwise
1099    ///
1100    pub fn compute_and_set_covariances(
1101        &self,
1102        linear_solver: &mut Box<dyn LinearSolver<SparseMode>>,
1103        variables: &mut HashMap<String, VariableEnum>,
1104        variable_index_map: &HashMap<String, usize>,
1105    ) -> Option<HashMap<String, Mat<f64>>> {
1106        // Compute the full covariance matrix (H^{-1}) using the linear solver
1107        linear_solver.compute_covariance_matrix()?;
1108        let full_cov = linear_solver.get_covariance_matrix()?.clone();
1109
1110        // Extract per-variable covariance blocks from the full matrix
1111        let per_var_covariances =
1112            extract_variable_covariances(&full_cov, variables, variable_index_map);
1113
1114        // Set covariances in Variable objects for easy access
1115        for (var_name, cov) in &per_var_covariances {
1116            if let Some(var) = variables.get_mut(var_name) {
1117                var.set_covariance(cov.clone());
1118            }
1119        }
1120
1121        Some(per_var_covariances)
1122    }
1123
1124    /// Compute and set covariances using a generic linear solver.
1125    ///
1126    /// This is the generic version of [`compute_and_set_covariances`] that works
1127    /// with any assembly mode (sparse or dense).
1128    pub fn compute_and_set_covariances_generic<M: crate::linalg::LinearizationMode>(
1129        &self,
1130        linear_solver: &mut dyn crate::linalg::LinearSolver<M>,
1131        variables: &mut HashMap<String, VariableEnum>,
1132        variable_index_map: &HashMap<String, usize>,
1133    ) -> Option<HashMap<String, Mat<f64>>> {
1134        linear_solver.compute_covariance_matrix()?;
1135        let full_cov = linear_solver.get_covariance_matrix()?.clone();
1136
1137        let per_var_covariances =
1138            extract_variable_covariances(&full_cov, variables, variable_index_map);
1139
1140        for (var_name, cov) in &per_var_covariances {
1141            if let Some(var) = variables.get_mut(var_name) {
1142                var.set_covariance(cov.clone());
1143            }
1144        }
1145
1146        Some(per_var_covariances)
1147    }
1148}
1149
1150#[cfg(test)]
1151mod tests {
1152    use super::*;
1153    use crate::core::loss_functions::HuberLoss;
1154    use crate::factors::{BetweenFactor, PriorFactor};
1155    use apex_manifolds::{ManifoldType, se2::SE2, se3::SE3};
1156    use nalgebra::{Quaternion, Vector3, dvector};
1157    use std::collections::HashMap;
1158
1159    type TestResult = Result<(), Box<dyn std::error::Error>>;
1160    type TestProblemResult = Result<
1161        (
1162            Problem,
1163            HashMap<String, (ManifoldType, nalgebra::DVector<f64>)>,
1164        ),
1165        Box<dyn std::error::Error>,
1166    >;
1167
1168    /// Create a test SE2 dataset with 10 vertices in a loop
1169    fn create_se2_test_problem() -> TestProblemResult {
1170        let mut problem = Problem::new(JacobianMode::Sparse);
1171        let mut initial_values = HashMap::new();
1172
1173        // Create 10 SE2 poses in a rough circle pattern
1174        let poses = vec![
1175            (0.0, 0.0, 0.0),    // x0: origin
1176            (1.0, 0.0, 0.1),    // x1: move right
1177            (1.5, 1.0, 0.5),    // x2: move up-right
1178            (1.0, 2.0, 1.0),    // x3: move up
1179            (0.0, 2.5, 1.5),    // x4: move up-left
1180            (-1.0, 2.0, 2.0),   // x5: move left
1181            (-1.5, 1.0, 2.5),   // x6: move down-left
1182            (-1.0, 0.0, 3.0),   // x7: move down
1183            (-0.5, -0.5, -2.8), // x8: move down-right
1184            (0.5, -0.5, -2.3),  // x9: back towards origin
1185        ];
1186
1187        // Add vertices using [x, y, theta] ordering
1188        for (i, (x, y, theta)) in poses.iter().enumerate() {
1189            let var_name = format!("x{}", i);
1190            let se2_data = dvector![*x, *y, *theta];
1191            initial_values.insert(var_name, (ManifoldType::SE2, se2_data));
1192        }
1193
1194        // Add chain of between factors
1195        for i in 0..9 {
1196            let from_pose = poses[i];
1197            let to_pose = poses[i + 1];
1198
1199            // Compute relative transformation
1200            let dx = to_pose.0 - from_pose.0;
1201            let dy = to_pose.1 - from_pose.1;
1202            let dtheta = to_pose.2 - from_pose.2;
1203
1204            let between_factor = BetweenFactor::new(SE2::from_xy_angle(dx, dy, dtheta));
1205            problem.add_residual_block(
1206                &[&format!("x{}", i), &format!("x{}", i + 1)],
1207                Box::new(between_factor),
1208                Some(Box::new(HuberLoss::new(1.0)?)),
1209            );
1210        }
1211
1212        // Add loop closure from x9 back to x0
1213        let dx = poses[0].0 - poses[9].0;
1214        let dy = poses[0].1 - poses[9].1;
1215        let dtheta = poses[0].2 - poses[9].2;
1216
1217        let loop_closure = BetweenFactor::new(SE2::from_xy_angle(dx, dy, dtheta));
1218        problem.add_residual_block(
1219            &["x9", "x0"],
1220            Box::new(loop_closure),
1221            Some(Box::new(HuberLoss::new(1.0)?)),
1222        );
1223
1224        // Add prior factor for x0
1225        let prior_factor = PriorFactor {
1226            data: dvector![0.0, 0.0, 0.0],
1227        };
1228        problem.add_residual_block(&["x0"], Box::new(prior_factor), None);
1229
1230        Ok((problem, initial_values))
1231    }
1232
1233    /// Create a test SE3 dataset with 8 vertices in a 3D pattern
1234    fn create_se3_test_problem() -> TestProblemResult {
1235        let mut problem = Problem::new(JacobianMode::Sparse);
1236        let mut initial_values = HashMap::new();
1237
1238        // Create 8 SE3 poses in a rough 3D cube pattern
1239        let poses = [
1240            // Bottom face of cube
1241            (0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0),   // x0: origin
1242            (1.0, 0.0, 0.0, 0.0, 0.0, 0.1, 0.995), // x1: +X
1243            (1.0, 1.0, 0.0, 0.0, 0.0, 0.2, 0.98),  // x2: +X+Y
1244            (0.0, 1.0, 0.0, 0.0, 0.0, 0.3, 0.955), // x3: +Y
1245            // Top face of cube
1246            (0.0, 0.0, 1.0, 0.1, 0.0, 0.0, 0.995), // x4: +Z
1247            (1.0, 0.0, 1.0, 0.1, 0.0, 0.1, 0.99),  // x5: +X+Z
1248            (1.0, 1.0, 1.0, 0.1, 0.0, 0.2, 0.975), // x6: +X+Y+Z
1249            (0.0, 1.0, 1.0, 0.1, 0.0, 0.3, 0.95),  // x7: +Y+Z
1250        ];
1251
1252        // Add vertices using [tx, ty, tz, qw, qx, qy, qz] ordering
1253        for (i, (tx, ty, tz, qx, qy, qz, qw)) in poses.iter().enumerate() {
1254            let var_name = format!("x{}", i);
1255            let se3_data = dvector![*tx, *ty, *tz, *qw, *qx, *qy, *qz];
1256            initial_values.insert(var_name, (ManifoldType::SE3, se3_data));
1257        }
1258
1259        // Add between factors connecting the cube edges
1260        let edges = vec![
1261            (0, 1),
1262            (1, 2),
1263            (2, 3),
1264            (3, 0), // Bottom face
1265            (4, 5),
1266            (5, 6),
1267            (6, 7),
1268            (7, 4), // Top face
1269            (0, 4),
1270            (1, 5),
1271            (2, 6),
1272            (3, 7), // Vertical edges
1273        ];
1274
1275        for (from_idx, to_idx) in edges {
1276            let from_pose = poses[from_idx];
1277            let to_pose = poses[to_idx];
1278
1279            // Create a simple relative transformation (simplified for testing)
1280            let relative_se3 = SE3::from_translation_quaternion(
1281                Vector3::new(
1282                    to_pose.0 - from_pose.0, // dx
1283                    to_pose.1 - from_pose.1, // dy
1284                    to_pose.2 - from_pose.2, // dz
1285                ),
1286                Quaternion::new(1.0, 0.0, 0.0, 0.0), // identity quaternion
1287            );
1288
1289            let between_factor = BetweenFactor::new(relative_se3);
1290            problem.add_residual_block(
1291                &[&format!("x{}", from_idx), &format!("x{}", to_idx)],
1292                Box::new(between_factor),
1293                Some(Box::new(HuberLoss::new(1.0)?)),
1294            );
1295        }
1296
1297        // Add prior factor for x0
1298        let prior_factor = PriorFactor {
1299            data: dvector![0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0],
1300        };
1301        problem.add_residual_block(&["x0"], Box::new(prior_factor), None);
1302
1303        Ok((problem, initial_values))
1304    }
1305
1306    #[test]
1307    fn test_problem_construction_se2() -> TestResult {
1308        let (problem, initial_values) = create_se2_test_problem()?;
1309
1310        // Test basic problem properties
1311        assert_eq!(problem.num_residual_blocks(), 11); // 9 between + 1 loop closure + 1 prior
1312        assert_eq!(problem.total_residual_dimension, 33); // 11 * 3
1313        assert_eq!(initial_values.len(), 10);
1314
1315        Ok(())
1316    }
1317
1318    #[test]
1319    fn test_problem_construction_se3() -> TestResult {
1320        let (problem, initial_values) = create_se3_test_problem()?;
1321
1322        // Test basic problem properties
1323        assert_eq!(problem.num_residual_blocks(), 13); // 12 between + 1 prior
1324        assert_eq!(problem.total_residual_dimension, 79); // 12 * 6 + 1 * 7 (SE3 between factors are 6-dim, prior factor is 7-dim)
1325        assert_eq!(initial_values.len(), 8);
1326
1327        Ok(())
1328    }
1329
1330    #[test]
1331    fn test_variable_initialization_se2() -> TestResult {
1332        let (problem, initial_values) = create_se2_test_problem()?;
1333
1334        // Test variable initialization
1335        let variables = problem.initialize_variables(&initial_values);
1336        assert_eq!(variables.len(), 10);
1337
1338        // Test variable sizes
1339        for (name, var) in &variables {
1340            assert_eq!(
1341                var.get_size(),
1342                3,
1343                "SE2 variable {} should have size 3",
1344                name
1345            );
1346        }
1347
1348        // Test conversion to DVector
1349        for (name, var) in &variables {
1350            let vec = var.to_vector();
1351            assert_eq!(
1352                vec.len(),
1353                3,
1354                "SE2 variable {} vector should have length 3",
1355                name
1356            );
1357        }
1358
1359        Ok(())
1360    }
1361
1362    #[test]
1363    fn test_variable_initialization_se3() -> TestResult {
1364        let (problem, initial_values) = create_se3_test_problem()?;
1365
1366        // Test variable initialization
1367        let variables = problem.initialize_variables(&initial_values);
1368        assert_eq!(variables.len(), 8);
1369
1370        // Test variable sizes
1371        for (name, var) in &variables {
1372            assert_eq!(
1373                var.get_size(),
1374                6,
1375                "SE3 variable {} should have size 6 (DOF)",
1376                name
1377            );
1378        }
1379
1380        // Test conversion to DVector
1381        for (name, var) in &variables {
1382            let vec = var.to_vector();
1383            assert_eq!(
1384                vec.len(),
1385                7,
1386                "SE3 variable {} vector should have length 7",
1387                name
1388            );
1389        }
1390
1391        Ok(())
1392    }
1393
1394    #[test]
1395    fn test_column_mapping_se2() -> TestResult {
1396        let (problem, initial_values) = create_se2_test_problem()?;
1397        let variables = problem.initialize_variables(&initial_values);
1398
1399        // Create column mapping for variables
1400        let mut variable_index_sparce_matrix = HashMap::new();
1401        let mut col_offset = 0;
1402        let mut sorted_vars: Vec<_> = variables.keys().collect();
1403        sorted_vars.sort(); // Ensure consistent ordering
1404
1405        for var_name in sorted_vars {
1406            variable_index_sparce_matrix.insert(var_name.clone(), col_offset);
1407            col_offset += variables[var_name].get_size();
1408        }
1409
1410        // Test total degrees of freedom
1411        let total_dof: usize = variables.values().map(|v| v.get_size()).sum();
1412        assert_eq!(total_dof, 30); // 10 variables * 3 DOF each
1413        assert_eq!(col_offset, 30);
1414
1415        // Test each variable has correct column mapping
1416        for (var_name, &col_idx) in &variable_index_sparce_matrix {
1417            assert!(
1418                col_idx < total_dof,
1419                "Column index {} for {} should be < {}",
1420                col_idx,
1421                var_name,
1422                total_dof
1423            );
1424        }
1425
1426        Ok(())
1427    }
1428
1429    #[test]
1430    fn test_symbolic_structure_se2() -> TestResult {
1431        let (problem, initial_values) = create_se2_test_problem()?;
1432        let variables = problem.initialize_variables(&initial_values);
1433
1434        // Create column mapping
1435        let mut variable_index_sparce_matrix = HashMap::new();
1436        let mut col_offset = 0;
1437        let mut sorted_vars: Vec<_> = variables.keys().collect();
1438        sorted_vars.sort();
1439
1440        for var_name in sorted_vars {
1441            variable_index_sparce_matrix.insert(var_name.clone(), col_offset);
1442            col_offset += variables[var_name].get_size();
1443        }
1444
1445        // Build symbolic structure
1446        let symbolic_structure = crate::linearizer::cpu::sparse::build_symbolic_structure(
1447            &problem,
1448            &variables,
1449            &variable_index_sparce_matrix,
1450            col_offset,
1451        )?;
1452
1453        // Test symbolic structure dimensions
1454        assert_eq!(
1455            symbolic_structure.pattern.nrows(),
1456            problem.total_residual_dimension
1457        );
1458        assert_eq!(symbolic_structure.pattern.ncols(), 30); // total DOF
1459
1460        Ok(())
1461    }
1462
1463    #[test]
1464    fn test_residual_jacobian_computation_se2() -> TestResult {
1465        let (problem, initial_values) = create_se2_test_problem()?;
1466        let variables = problem.initialize_variables(&initial_values);
1467
1468        // Create column mapping
1469        let mut variable_index_sparce_matrix = HashMap::new();
1470        let mut col_offset = 0;
1471        let mut sorted_vars: Vec<_> = variables.keys().collect();
1472        sorted_vars.sort();
1473
1474        for var_name in sorted_vars {
1475            variable_index_sparce_matrix.insert(var_name.clone(), col_offset);
1476            col_offset += variables[var_name].get_size();
1477        }
1478
1479        // Test sparse computation
1480        let symbolic_structure = crate::linearizer::cpu::sparse::build_symbolic_structure(
1481            &problem,
1482            &variables,
1483            &variable_index_sparce_matrix,
1484            col_offset,
1485        )?;
1486        let (residual_sparse, jacobian_sparse) = problem.compute_residual_and_jacobian_sparse(
1487            &variables,
1488            &variable_index_sparce_matrix,
1489            &symbolic_structure,
1490        )?;
1491
1492        // Test sparse dimensions
1493        assert_eq!(residual_sparse.nrows(), problem.total_residual_dimension);
1494        assert_eq!(jacobian_sparse.nrows(), problem.total_residual_dimension);
1495        assert_eq!(jacobian_sparse.ncols(), 30);
1496
1497        Ok(())
1498    }
1499
1500    #[test]
1501    fn test_residual_jacobian_computation_se3() -> TestResult {
1502        let (problem, initial_values) = create_se3_test_problem()?;
1503        let variables = problem.initialize_variables(&initial_values);
1504
1505        // Create column mapping
1506        let mut variable_index_sparce_matrix = HashMap::new();
1507        let mut col_offset = 0;
1508        let mut sorted_vars: Vec<_> = variables.keys().collect();
1509        sorted_vars.sort();
1510
1511        for var_name in sorted_vars {
1512            variable_index_sparce_matrix.insert(var_name.clone(), col_offset);
1513            col_offset += variables[var_name].get_size();
1514        }
1515
1516        // Test sparse computation
1517        let symbolic_structure = crate::linearizer::cpu::sparse::build_symbolic_structure(
1518            &problem,
1519            &variables,
1520            &variable_index_sparce_matrix,
1521            col_offset,
1522        )?;
1523        let (residual_sparse, jacobian_sparse) = problem.compute_residual_and_jacobian_sparse(
1524            &variables,
1525            &variable_index_sparce_matrix,
1526            &symbolic_structure,
1527        )?;
1528
1529        // Test sparse dimensions match
1530        assert_eq!(residual_sparse.nrows(), problem.total_residual_dimension);
1531        assert_eq!(jacobian_sparse.nrows(), problem.total_residual_dimension);
1532        assert_eq!(jacobian_sparse.ncols(), 48); // 8 variables * 6 DOF each
1533
1534        Ok(())
1535    }
1536
1537    #[test]
1538    fn test_residual_block_operations() -> TestResult {
1539        let mut problem = Problem::new(JacobianMode::Sparse);
1540
1541        // Test adding residual blocks
1542        let block_id1 = problem.add_residual_block(
1543            &["x0", "x1"],
1544            Box::new(BetweenFactor::new(SE2::from_xy_angle(1.0, 0.0, 0.1))),
1545            Some(Box::new(HuberLoss::new(1.0)?)),
1546        );
1547
1548        let block_id2 = problem.add_residual_block(
1549            &["x0"],
1550            Box::new(PriorFactor {
1551                data: dvector![0.0, 0.0, 0.0],
1552            }),
1553            None,
1554        );
1555
1556        assert_eq!(problem.num_residual_blocks(), 2);
1557        assert_eq!(problem.total_residual_dimension, 6); // 3 + 3
1558        assert_eq!(block_id1, 0);
1559        assert_eq!(block_id2, 1);
1560
1561        // Test removing residual blocks
1562        let removed_block = problem.remove_residual_block(block_id1);
1563        assert!(removed_block.is_some());
1564        assert_eq!(problem.num_residual_blocks(), 1);
1565        assert_eq!(problem.total_residual_dimension, 3); // Only prior factor remains
1566
1567        // Test removing non-existent block
1568        let non_existent = problem.remove_residual_block(999);
1569        assert!(non_existent.is_none());
1570
1571        Ok(())
1572    }
1573
1574    #[test]
1575    fn test_variable_constraints() -> TestResult {
1576        let mut problem = Problem::new(JacobianMode::Sparse);
1577
1578        // Test fixing variables
1579        problem.fix_variable("x0", 0);
1580        problem.fix_variable("x0", 1);
1581        problem.fix_variable("x1", 2);
1582
1583        assert!(problem.fixed_variable_indexes.contains_key("x0"));
1584        assert!(problem.fixed_variable_indexes.contains_key("x1"));
1585        assert_eq!(problem.fixed_variable_indexes["x0"].len(), 2);
1586        assert_eq!(problem.fixed_variable_indexes["x1"].len(), 1);
1587
1588        // Test unfixing variables
1589        problem.unfix_variable("x0");
1590        assert!(!problem.fixed_variable_indexes.contains_key("x0"));
1591        assert!(problem.fixed_variable_indexes.contains_key("x1"));
1592
1593        // Test variable bounds
1594        problem.set_variable_bounds("x2", 0, -1.0, 1.0);
1595        problem.set_variable_bounds("x2", 1, -2.0, 2.0);
1596        problem.set_variable_bounds("x3", 0, 0.0, 5.0);
1597
1598        assert!(problem.variable_bounds.contains_key("x2"));
1599        assert!(problem.variable_bounds.contains_key("x3"));
1600        assert_eq!(problem.variable_bounds["x2"].len(), 2);
1601        assert_eq!(problem.variable_bounds["x3"].len(), 1);
1602
1603        // Test removing bounds
1604        problem.remove_variable_bounds("x2");
1605        assert!(!problem.variable_bounds.contains_key("x2"));
1606        assert!(problem.variable_bounds.contains_key("x3"));
1607
1608        Ok(())
1609    }
1610
1611    // -------------------------------------------------------------------------
1612    // New tests for previously uncovered code paths
1613    // -------------------------------------------------------------------------
1614
1615    #[test]
1616    fn test_problem_default_equals_new_sparse() {
1617        let default_problem = Problem::default();
1618        let new_problem = Problem::new(JacobianMode::Sparse);
1619        assert_eq!(default_problem.jacobian_mode, new_problem.jacobian_mode);
1620        assert_eq!(
1621            default_problem.num_residual_blocks(),
1622            new_problem.num_residual_blocks()
1623        );
1624    }
1625
1626    /// Test that SO2 variables can be initialized correctly.
1627    #[test]
1628    fn test_initialize_so2_variable() -> TestResult {
1629        let problem = Problem::new(JacobianMode::Sparse);
1630        let mut initial = HashMap::new();
1631        initial.insert("angle".to_string(), (ManifoldType::SO2, dvector![0.5]));
1632        let variables = problem.initialize_variables(&initial);
1633        assert_eq!(variables.len(), 1);
1634        let var = variables.get("angle").ok_or("angle not found")?;
1635        assert_eq!(var.manifold_type(), ManifoldType::SO2);
1636        assert_eq!(var.get_size(), 1);
1637        assert_eq!(var.to_vector().len(), 1);
1638        Ok(())
1639    }
1640
1641    /// Test that SO3 variables can be initialized correctly.
1642    #[test]
1643    fn test_initialize_so3_variable() -> TestResult {
1644        let problem = Problem::new(JacobianMode::Sparse);
1645        let mut initial = HashMap::new();
1646        // SO3: [qw, qx, qy, qz]
1647        initial.insert(
1648            "rot".to_string(),
1649            (ManifoldType::SO3, dvector![1.0, 0.0, 0.0, 0.0]),
1650        );
1651        let variables = problem.initialize_variables(&initial);
1652        assert_eq!(variables.len(), 1);
1653        let var = variables.get("rot").ok_or("rot not found")?;
1654        assert_eq!(var.manifold_type(), ManifoldType::SO3);
1655        assert_eq!(var.get_size(), 3); // SO3 tangent space is 3-dimensional
1656        Ok(())
1657    }
1658
1659    /// Test that RN variables can be initialized correctly (arbitrary dimension).
1660    #[test]
1661    fn test_initialize_rn_variable() -> TestResult {
1662        let problem = Problem::new(JacobianMode::Sparse);
1663        let mut initial = HashMap::new();
1664        initial.insert(
1665            "pt".to_string(),
1666            (ManifoldType::RN, dvector![1.0, 2.0, 3.0]),
1667        );
1668        let variables = problem.initialize_variables(&initial);
1669        let var = variables.get("pt").ok_or("pt not found")?;
1670        assert_eq!(var.manifold_type(), ManifoldType::RN);
1671        assert_eq!(var.get_size(), 3);
1672        let vec = var.to_vector();
1673        assert!((vec[0] - 1.0).abs() < 1e-12);
1674        assert!((vec[1] - 2.0).abs() < 1e-12);
1675        assert!((vec[2] - 3.0).abs() < 1e-12);
1676        Ok(())
1677    }
1678
1679    /// Test VariableEnum covariance lifecycle: set → get → clear.
1680    #[test]
1681    fn test_variable_enum_covariance_lifecycle() -> TestResult {
1682        use faer::Mat;
1683        let problem = Problem::new(JacobianMode::Sparse);
1684        let mut initial = HashMap::new();
1685        initial.insert(
1686            "x0".to_string(),
1687            (ManifoldType::SE2, dvector![0.0, 0.0, 0.0]),
1688        );
1689        let mut variables = problem.initialize_variables(&initial);
1690        let var = variables.get_mut("x0").ok_or("x0 not found")?;
1691
1692        // Before setting: no covariance
1693        assert!(var.get_covariance().is_none());
1694
1695        // Set a 3×3 identity covariance
1696        let cov = Mat::identity(3, 3);
1697        var.set_covariance(cov);
1698        let retrieved = var.get_covariance().ok_or("covariance not set")?;
1699        assert_eq!(retrieved.nrows(), 3);
1700        assert!((retrieved[(0, 0)] - 1.0).abs() < 1e-12);
1701
1702        // Clear covariance
1703        var.clear_covariance();
1704        assert!(var.get_covariance().is_none());
1705        Ok(())
1706    }
1707
1708    /// Test VariableEnum::get_bounds() propagation via initialize_variables.
1709    #[test]
1710    fn test_variable_enum_bounds_propagation() -> TestResult {
1711        let mut problem = Problem::new(JacobianMode::Sparse);
1712        problem.set_variable_bounds("x0", 0, -1.0, 1.0);
1713        problem.set_variable_bounds("x0", 1, -2.0, 2.0);
1714        let mut initial = HashMap::new();
1715        initial.insert(
1716            "x0".to_string(),
1717            (ManifoldType::SE2, dvector![0.0, 0.0, 0.0]),
1718        );
1719        let variables = problem.initialize_variables(&initial);
1720        let var = variables.get("x0").ok_or("x0 not found")?;
1721        let bounds = var.get_bounds();
1722        assert_eq!(bounds.len(), 2);
1723        let (lo, hi) = bounds[&0];
1724        assert!((lo + 1.0).abs() < 1e-12);
1725        assert!((hi - 1.0).abs() < 1e-12);
1726        Ok(())
1727    }
1728
1729    /// Test VariableEnum::get_fixed_indices() propagation via initialize_variables.
1730    #[test]
1731    fn test_variable_enum_fixed_indices_propagation() -> TestResult {
1732        let mut problem = Problem::new(JacobianMode::Sparse);
1733        problem.fix_variable("x0", 0);
1734        problem.fix_variable("x0", 2);
1735        let mut initial = HashMap::new();
1736        initial.insert(
1737            "x0".to_string(),
1738            (ManifoldType::SE2, dvector![0.0, 0.0, 0.0]),
1739        );
1740        let variables = problem.initialize_variables(&initial);
1741        let var = variables.get("x0").ok_or("x0 not found")?;
1742        let fixed = var.get_fixed_indices();
1743        assert_eq!(fixed.len(), 2);
1744        assert!(fixed.contains(&0));
1745        assert!(fixed.contains(&2));
1746        Ok(())
1747    }
1748
1749    /// Test VariableEnum::set_from_vector() for all five manifold variants.
1750    #[test]
1751    fn test_variable_enum_set_from_vector_all_variants() -> TestResult {
1752        let problem = Problem::new(JacobianMode::Sparse);
1753        let mut initial = HashMap::new();
1754        initial.insert(
1755            "se2".to_string(),
1756            (ManifoldType::SE2, dvector![0.0, 0.0, 0.0]),
1757        );
1758        initial.insert(
1759            "se3".to_string(),
1760            (
1761                ManifoldType::SE3,
1762                dvector![0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0],
1763            ),
1764        );
1765        initial.insert("so2".to_string(), (ManifoldType::SO2, dvector![0.0]));
1766        initial.insert(
1767            "so3".to_string(),
1768            (ManifoldType::SO3, dvector![1.0, 0.0, 0.0, 0.0]),
1769        );
1770        initial.insert("rn".to_string(), (ManifoldType::RN, dvector![1.0, 2.0]));
1771        let mut variables = problem.initialize_variables(&initial);
1772
1773        // SE2: update to [1, 2, 0.5]
1774        let new_se2 = dvector![1.0, 2.0, 0.5];
1775        variables
1776            .get_mut("se2")
1777            .ok_or("se2 not found")?
1778            .set_from_vector(&new_se2);
1779        let got = variables.get("se2").ok_or("se2 not found")?.to_vector();
1780        assert!((got[0] - 1.0).abs() < 1e-10);
1781
1782        // SE3: update translation part
1783        let new_se3 = dvector![1.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0];
1784        variables
1785            .get_mut("se3")
1786            .ok_or("se3 not found")?
1787            .set_from_vector(&new_se3);
1788        let got = variables.get("se3").ok_or("se3 not found")?.to_vector();
1789        assert!((got[0] - 1.0).abs() < 1e-10);
1790
1791        // SO2
1792        let new_so2 = dvector![0.5];
1793        variables
1794            .get_mut("so2")
1795            .ok_or("so2 not found")?
1796            .set_from_vector(&new_so2);
1797        assert_eq!(variables.get("so2").ok_or("so2 not found")?.get_size(), 1);
1798
1799        // SO3
1800        let new_so3 = dvector![1.0, 0.0, 0.0, 0.0];
1801        variables
1802            .get_mut("so3")
1803            .ok_or("so3 not found")?
1804            .set_from_vector(&new_so3);
1805        assert_eq!(variables.get("so3").ok_or("so3 not found")?.get_size(), 3);
1806
1807        // RN
1808        let new_rn = dvector![5.0, 6.0];
1809        variables
1810            .get_mut("rn")
1811            .ok_or("rn not found")?
1812            .set_from_vector(&new_rn);
1813        let got = variables.get("rn").ok_or("rn not found")?.to_vector();
1814        assert!((got[0] - 5.0).abs() < 1e-10);
1815        assert!((got[1] - 6.0).abs() < 1e-10);
1816        Ok(())
1817    }
1818
1819    /// Test apply_tangent_step() zeros out fixed indices for SE2.
1820    #[test]
1821    fn test_apply_tangent_step_se2_fixed_index() -> TestResult {
1822        use faer::Mat;
1823        let mut problem = Problem::new(JacobianMode::Sparse);
1824        // Fix SE2 index 0 (x translation)
1825        problem.fix_variable("x0", 0);
1826        let mut initial = HashMap::new();
1827        initial.insert(
1828            "x0".to_string(),
1829            (ManifoldType::SE2, dvector![0.0, 0.0, 0.0]),
1830        );
1831        let mut variables = problem.initialize_variables(&initial);
1832
1833        // Step vector: [1, 2, 3] — but index 0 should be zeroed
1834        let mut step = Mat::zeros(3, 1);
1835        step[(0, 0)] = 1.0;
1836        step[(1, 0)] = 2.0;
1837        step[(2, 0)] = 3.0;
1838
1839        let var = variables.get_mut("x0").ok_or("x0 not found")?;
1840        var.apply_tangent_step(step.as_ref());
1841        // After update, get_size() should still be 3
1842        assert_eq!(var.get_size(), 3);
1843        Ok(())
1844    }
1845
1846    /// Test apply_tangent_step() zeros out fixed indices for SE3.
1847    #[test]
1848    fn test_apply_tangent_step_se3_fixed_index() -> TestResult {
1849        use faer::Mat;
1850        let mut problem = Problem::new(JacobianMode::Sparse);
1851        problem.fix_variable("p", 0);
1852        let mut initial = HashMap::new();
1853        initial.insert(
1854            "p".to_string(),
1855            (
1856                ManifoldType::SE3,
1857                dvector![0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0],
1858            ),
1859        );
1860        let mut variables = problem.initialize_variables(&initial);
1861
1862        let mut step = Mat::zeros(6, 1);
1863        for i in 0..6 {
1864            step[(i, 0)] = (i + 1) as f64;
1865        }
1866
1867        let var = variables.get_mut("p").ok_or("p not found")?;
1868        var.apply_tangent_step(step.as_ref());
1869        assert_eq!(var.get_size(), 6);
1870        Ok(())
1871    }
1872
1873    /// Test apply_tangent_step for SO2 and SO3 and RN.
1874    #[test]
1875    fn test_apply_tangent_step_remaining_manifolds() -> TestResult {
1876        use faer::Mat;
1877        let problem = Problem::new(JacobianMode::Sparse);
1878        let mut initial = HashMap::new();
1879        initial.insert("so2".to_string(), (ManifoldType::SO2, dvector![0.0]));
1880        initial.insert(
1881            "so3".to_string(),
1882            (ManifoldType::SO3, dvector![1.0, 0.0, 0.0, 0.0]),
1883        );
1884        initial.insert("rn".to_string(), (ManifoldType::RN, dvector![0.0, 0.0]));
1885        let mut variables = problem.initialize_variables(&initial);
1886
1887        // SO2: 1-DOF step
1888        let mut step_so2 = Mat::zeros(1, 1);
1889        step_so2[(0, 0)] = 0.1;
1890        variables
1891            .get_mut("so2")
1892            .ok_or("so2 not found")?
1893            .apply_tangent_step(step_so2.as_ref());
1894
1895        // SO3: 3-DOF step
1896        let mut step_so3 = Mat::zeros(3, 1);
1897        step_so3[(0, 0)] = 0.01;
1898        variables
1899            .get_mut("so3")
1900            .ok_or("so3 not found")?
1901            .apply_tangent_step(step_so3.as_ref());
1902
1903        // RN: 2-DOF step
1904        let mut step_rn = Mat::zeros(2, 1);
1905        step_rn[(0, 0)] = 1.0;
1906        step_rn[(1, 0)] = 2.0;
1907        variables
1908            .get_mut("rn")
1909            .ok_or("rn not found")?
1910            .apply_tangent_step(step_rn.as_ref());
1911        let vec = variables.get("rn").ok_or("rn not found")?.to_vector();
1912        assert!((vec[0] - 1.0).abs() < 1e-10);
1913        assert!((vec[1] - 2.0).abs() < 1e-10);
1914        Ok(())
1915    }
1916
1917    /// Test compute_residual_sparse returns a non-negative squared norm.
1918    #[test]
1919    fn test_compute_residual_sparse_smoke() -> TestResult {
1920        let (problem, initial_values) = create_se2_test_problem()?;
1921        let variables = problem.initialize_variables(&initial_values);
1922        let residual = problem.compute_residual_sparse(&variables)?;
1923        // Squared norm must be non-negative
1924        let norm_sq: f64 = (0..residual.nrows())
1925            .map(|i| residual[(i, 0)].powi(2))
1926            .sum();
1927        assert!(norm_sq >= 0.0);
1928        assert_eq!(residual.nrows(), problem.total_residual_dimension);
1929        Ok(())
1930    }
1931
1932    /// Test log_residual_to_file writes a file without error.
1933    #[test]
1934    fn test_log_residual_to_file() -> TestResult {
1935        let problem = Problem::new(JacobianMode::Sparse);
1936        let residual = nalgebra::dvector![1.0, 2.0, 3.0];
1937        let path = std::env::temp_dir().join("apex_test_residual.txt");
1938        let path_str = path.to_str().ok_or("temp path is not valid UTF-8")?;
1939        problem.log_residual_to_file(&residual, path_str)?;
1940        assert!(path.exists());
1941        Ok(())
1942    }
1943
1944    /// Test log_variables_to_file writes a file without error.
1945    #[test]
1946    fn test_log_variables_to_file() -> TestResult {
1947        let problem = Problem::new(JacobianMode::Sparse);
1948        let mut initial = HashMap::new();
1949        initial.insert(
1950            "x0".to_string(),
1951            (ManifoldType::SE2, dvector![1.0, 2.0, 0.3]),
1952        );
1953        let variables = problem.initialize_variables(&initial);
1954        let path = std::env::temp_dir().join("apex_test_variables.txt");
1955        let path_str = path.to_str().ok_or("temp path is not valid UTF-8")?;
1956        problem.log_variables_to_file(&variables, path_str)?;
1957        assert!(path.exists());
1958        Ok(())
1959    }
1960
1961    /// Test log_sparse_jacobian_to_file writes a file without error.
1962    #[test]
1963    fn test_log_sparse_jacobian_to_file() -> TestResult {
1964        use faer::sparse::SparseColMat;
1965        let problem = Problem::new(JacobianMode::Sparse);
1966        let triplets = vec![faer::sparse::Triplet::new(0usize, 0usize, 1.0f64)];
1967        let jacobian =
1968            SparseColMat::try_new_from_triplets(1, 1, &triplets).map_err(|e| format!("{e:?}"))?;
1969        let path = std::env::temp_dir().join("apex_test_jacobian.txt");
1970        let path_str = path.to_str().ok_or("temp path is not valid UTF-8")?;
1971        problem.log_sparse_jacobian_to_file(&jacobian, path_str)?;
1972        assert!(path.exists());
1973        Ok(())
1974    }
1975
1976    /// Test set_variable_bounds with lower > upper logs a warning but does not update bounds.
1977    #[test]
1978    fn test_set_variable_bounds_invalid_order() {
1979        let mut problem = Problem::new(JacobianMode::Sparse);
1980        // lower > upper — should warn but NOT insert into the bounds map
1981        problem.set_variable_bounds("x0", 0, 5.0, 1.0);
1982        // Since lower > upper, the bounds map should NOT be updated
1983        assert!(!problem.variable_bounds.contains_key("x0"));
1984    }
1985
1986    /// Test VariableEnum::manifold_type() returns the correct type for each variant.
1987    #[test]
1988    fn test_variable_enum_manifold_type() {
1989        let problem = Problem::new(JacobianMode::Sparse);
1990        let mut initial = HashMap::new();
1991        initial.insert(
1992            "se2".to_string(),
1993            (ManifoldType::SE2, dvector![0.0, 0.0, 0.0]),
1994        );
1995        initial.insert(
1996            "se3".to_string(),
1997            (
1998                ManifoldType::SE3,
1999                dvector![0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0],
2000            ),
2001        );
2002        initial.insert("so2".to_string(), (ManifoldType::SO2, dvector![0.0]));
2003        initial.insert(
2004            "so3".to_string(),
2005            (ManifoldType::SO3, dvector![1.0, 0.0, 0.0, 0.0]),
2006        );
2007        initial.insert("rn".to_string(), (ManifoldType::RN, dvector![0.0]));
2008        let variables = problem.initialize_variables(&initial);
2009
2010        assert_eq!(variables["se2"].manifold_type(), ManifoldType::SE2);
2011        assert_eq!(variables["se3"].manifold_type(), ManifoldType::SE3);
2012        assert_eq!(variables["so2"].manifold_type(), ManifoldType::SO2);
2013        assert_eq!(variables["so3"].manifold_type(), ManifoldType::SO3);
2014        assert_eq!(variables["rn"].manifold_type(), ManifoldType::RN);
2015    }
2016}