apex-solver 1.2.1

High-performance nonlinear least squares optimization with Lie group support for SLAM and bundle adjustment
Documentation
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
//! Observer pattern for optimization monitoring.
//!
//! This module provides a clean observer pattern for monitoring optimization progress.
//! Observers can be registered with any optimizer and will be notified at each iteration,
//! enabling real-time visualization, logging, metrics collection, and custom analysis.
//!
//! # Design Philosophy
//!
//! The observer pattern provides complete separation between optimization algorithms
//! and monitoring/visualization logic:
//!
//! - **Decoupling**: Optimization logic is independent of how progress is monitored
//! - **Extensibility**: Easy to add new observers (Rerun, CSV, metrics, dashboards)
//! - **Composability**: Multiple observers can run simultaneously
//! - **Zero overhead**: When no observers are registered, notification is a no-op
//!
//! # Architecture
//!
//! ```text
//! ┌─────────────────┐
//! │   Optimizer     │
//! │  (LM/GN/DogLeg) │
//! └────────┬────────┘
//!          │ observers.notify(values, iteration)
//!          ├──────────────┬──────────────┬──────────────┐
//!          ▼              ▼              ▼              ▼
//!    ┌──────────┐  ┌──────────┐  ┌──────────┐  ┌──────────┐
//!    │  Rerun   │  │   CSV    │  │ Metrics  │  │  Custom  │
//!    │ Observer │  │ Observer │  │ Observer │  │ Observer │
//!    └──────────┘  └──────────┘  └──────────┘  └──────────┘
//! ```
//!
//! # Examples
//!
//! ## Single Observer
//!
//! ```no_run
//! use apex_solver::{LevenbergMarquardt, LevenbergMarquardtConfig};
//! use apex_solver::observers::OptObserver;
//! # use apex_solver::core::problem::Problem;
//! # use std::collections::HashMap;
//!
//! # fn main() -> Result<(), Box<dyn std::error::Error>> {
//! # let problem = Problem::new();
//! # let initial_values = HashMap::new();
//!
//! let config = LevenbergMarquardtConfig::new().with_max_iterations(100);
//! let mut solver = LevenbergMarquardt::with_config(config);
//!
//! #[cfg(feature = "visualization")]
//! {
//!     use apex_solver::observers::RerunObserver;
//!     let rerun_observer = RerunObserver::new(true)?;
//!     solver.add_observer(rerun_observer);
//! }
//!
//! let result = solver.optimize(&problem, &initial_values)?;
//! # Ok(())
//! # }
//! ```
//!
//! ## Multiple Observers
//!
//! ```no_run
//! # use apex_solver::{LevenbergMarquardt, LevenbergMarquardtConfig};
//! # use apex_solver::core::problem::{Problem, VariableEnum};
//! # use apex_solver::observers::OptObserver;
//! # use std::collections::HashMap;
//!
//! // Custom observer that logs to CSV
//! struct CsvObserver {
//!     file: std::fs::File,
//! }
//!
//! impl OptObserver for CsvObserver {
//!     fn on_step(&self, _values: &HashMap<String, VariableEnum>, iteration: usize) {
//!         // Write iteration data to CSV
//!         // ... implementation ...
//!     }
//! }
//!
//! # fn main() -> Result<(), Box<dyn std::error::Error>> {
//! # let problem = Problem::new();
//! # let initial_values = HashMap::new();
//! let mut solver = LevenbergMarquardt::new();
//!
//! // Add Rerun visualization
//! #[cfg(feature = "visualization")]
//! {
//!     use apex_solver::observers::RerunObserver;
//!     solver.add_observer(RerunObserver::new(true)?);
//! }
//!
//! // Add CSV logging
//! // solver.add_observer(CsvObserver { file: ... });
//!
//! let result = solver.optimize(&problem, &initial_values)?;
//! # Ok(())
//! # }
//! ```
//!
//! ## Custom Observer
//!
//! ```no_run
//! use apex_solver::observers::OptObserver;
//! use apex_solver::core::problem::VariableEnum;
//! use std::collections::HashMap;
//!
//! struct MetricsObserver {
//!     max_variables_seen: std::cell::RefCell<usize>,
//! }
//!
//! impl OptObserver for MetricsObserver {
//!     fn on_step(&self, values: &HashMap<String, VariableEnum>, iteration: usize) {
//!         let count = values.len();
//!         let mut max = self.max_variables_seen.borrow_mut();
//!         *max = (*max).max(count);
//!     }
//! }
//! ```

// Visualization-specific submodules (feature-gated)
#[cfg(feature = "visualization")]
pub mod conversions;
#[cfg(feature = "visualization")]
pub mod visualization;

// Re-export RerunObserver when visualization is enabled
#[cfg(feature = "visualization")]
pub use visualization::{RerunObserver, VisualizationConfig, VisualizationMode};

// Re-export conversion traits for ergonomic use
#[cfg(feature = "visualization")]
pub use conversions::{
    CollectRerunArrows2D, CollectRerunArrows3D, CollectRerunPoints2D, CollectRerunPoints3D,
    ToRerunArrows2D, ToRerunArrows3D, ToRerunPoints2D, ToRerunPoints3D, ToRerunTransform3D,
    ToRerunTransform3DFrom2D, ToRerunVec2D, ToRerunVec3D,
};

use crate::core::problem::VariableEnum;
use faer::Mat;
use faer::sparse;
use std::collections::HashMap;
use thiserror::Error;
use tracing::error;

/// Observer-specific error types for apex-solver
#[derive(Debug, Clone, Error)]
pub enum ObserverError {
    /// Failed to initialize Rerun recording stream
    #[error("Failed to initialize Rerun recording stream: {0}")]
    RerunInitialization(String),

    /// Failed to spawn Rerun viewer process
    #[error("Failed to spawn Rerun viewer: {0}")]
    ViewerSpawnFailed(String),

    /// Failed to save recording to file
    #[error("Failed to save recording to file '{path}': {reason}")]
    RecordingSaveFailed { path: String, reason: String },

    /// Failed to log data to Rerun
    #[error("Failed to log data to Rerun at '{entity_path}': {reason}")]
    LoggingFailed { entity_path: String, reason: String },

    /// Failed to convert matrix to visualization format
    #[error("Failed to convert matrix to image: {0}")]
    MatrixVisualizationFailed(String),

    /// Failed to convert tensor data
    #[error("Failed to create tensor data: {0}")]
    TensorConversionFailed(String),

    /// Recording stream is in invalid state
    #[error("Recording stream is in invalid state: {0}")]
    InvalidState(String),

    /// Mutex was poisoned (thread panicked while holding lock)
    #[error("Mutex poisoned in {context}: {reason}")]
    MutexPoisoned { context: String, reason: String },
}

impl ObserverError {
    /// Log the error with tracing::error and return self for chaining
    ///
    /// This method allows for a consistent error logging pattern throughout
    /// the observers module, ensuring all errors are properly recorded.
    ///
    /// # Example
    /// ```
    /// # use apex_solver::observers::ObserverError;
    /// # fn operation() -> Result<(), ObserverError> { Ok(()) }
    /// # fn example() -> Result<(), ObserverError> {
    /// operation()
    ///     .map_err(|e| e.log())?;
    /// # Ok(())
    /// # }
    /// ```
    #[must_use]
    pub fn log(self) -> Self {
        error!("{}", self);
        self
    }

    /// Log the error with the original source error from a third-party library
    ///
    /// This method logs both the ObserverError and the underlying error
    /// from external libraries (e.g., Rerun's errors). This provides full
    /// debugging context when errors occur in third-party code.
    ///
    /// # Arguments
    /// * `source_error` - The original error from the third-party library (must implement Debug)
    ///
    /// # Example
    /// ```no_run
    /// # use apex_solver::observers::ObserverError;
    /// # fn rec_log() -> Result<(), std::io::Error> { Ok(()) }
    /// # fn example() -> Result<(), ObserverError> {
    /// rec_log()
    ///     .map_err(|e| {
    ///         ObserverError::LoggingFailed {
    ///             entity_path: "world/points".to_string(),
    ///             reason: format!("{}", e)
    ///         }
    ///         .log_with_source(e)
    ///     })?;
    /// # Ok(())
    /// # }
    /// ```
    #[must_use]
    pub fn log_with_source<E: std::fmt::Debug>(self, source_error: E) -> Self {
        error!("{} | Source: {:?}", self, source_error);
        self
    }
}

/// Result type for observer operations
pub type ObserverResult<T> = Result<T, ObserverError>;

/// Observer trait for monitoring optimization progress.
///
/// Implement this trait to create custom observers that are notified at each
/// optimization iteration. Observers receive the current variable values and
/// iteration number, enabling real-time monitoring, visualization, logging,
/// or custom analysis.
///
/// # Design Notes
///
/// - Observers should be lightweight and non-blocking
/// - Errors in observers should not crash optimization (handle internally)
/// - For expensive operations (file I/O, network), consider buffering
/// - Observers receive immutable references (cannot modify optimization state)
///
/// # Thread Safety
///
/// Observers must be `Send` to support parallel optimization in the future.
/// Use interior mutability (`RefCell`, `Mutex`) if you need to mutate state.
pub trait OptObserver: Send {
    /// Called after each optimization iteration.
    ///
    /// # Arguments
    ///
    /// * `values` - Current variable values (manifold states)
    /// * `iteration` - Current iteration number (0 = initial values, 1+ = after steps)
    ///
    /// # Implementation Guidelines
    ///
    /// - Keep this method fast to avoid slowing optimization
    /// - Handle errors internally (log warnings, don't panic)
    /// - Don't mutate `values` (you receive `&HashMap`)
    /// - Consider buffering expensive operations
    ///
    /// # Examples
    ///
    /// ```no_run
    /// use apex_solver::observers::OptObserver;
    /// use apex_solver::core::problem::VariableEnum;
    /// use std::collections::HashMap;
    ///
    /// struct SimpleLogger;
    ///
    /// impl OptObserver for SimpleLogger {
    ///     fn on_step(&self, values: &HashMap<String, VariableEnum>, iteration: usize) {
    ///         // Track optimization progress
    ///     }
    /// }
    /// ```
    fn on_step(&self, values: &HashMap<String, VariableEnum>, iteration: usize);

    /// Set iteration metrics for visualization and monitoring.
    ///
    /// This method is called before `on_step` to provide optimization metrics
    /// such as cost, gradient norm, damping parameter, etc. Observers can use
    /// this data for visualization, logging, or analysis.
    ///
    /// # Arguments
    ///
    /// * `cost` - Current cost function value
    /// * `gradient_norm` - L2 norm of the gradient vector
    /// * `damping` - Damping parameter (for Levenberg-Marquardt, may be None for other solvers)
    /// * `step_norm` - L2 norm of the parameter update step
    /// * `step_quality` - Step quality metric (e.g., rho for trust region methods)
    ///
    /// # Default Implementation
    ///
    /// The default implementation does nothing, allowing simple observers to ignore metrics.
    fn set_iteration_metrics(
        &self,
        _cost: f64,
        _gradient_norm: f64,
        _damping: Option<f64>,
        _step_norm: f64,
        _step_quality: Option<f64>,
    ) {
        // Default implementation does nothing
    }

    /// Set matrix data for advanced visualization.
    ///
    /// This method provides access to the Hessian matrix and gradient vector
    /// for observers that want to visualize matrix structure or perform
    /// advanced analysis.
    ///
    /// # Arguments
    ///
    /// * `hessian` - Sparse Hessian matrix (J^T * J)
    /// * `gradient` - Gradient vector (J^T * r)
    ///
    /// # Default Implementation
    ///
    /// The default implementation does nothing, allowing simple observers to ignore matrices.
    fn set_matrix_data(
        &self,
        _hessian: Option<sparse::SparseColMat<usize, f64>>,
        _gradient: Option<Mat<f64>>,
    ) {
        // Default implementation does nothing
    }

    /// Called when optimization completes.
    ///
    /// This method is called once at the end of optimization, after all iterations
    /// are complete. Use this for final visualization, cleanup, or summary logging.
    ///
    /// # Arguments
    ///
    /// * `values` - Final optimized variable values
    /// * `iterations` - Total number of iterations performed
    ///
    /// # Default Implementation
    ///
    /// The default implementation does nothing, allowing simple observers to ignore completion.
    ///
    /// # Examples
    ///
    /// ```no_run
    /// use apex_solver::observers::OptObserver;
    /// use apex_solver::core::problem::VariableEnum;
    /// use std::collections::HashMap;
    ///
    /// struct FinalStateLogger;
    ///
    /// impl OptObserver for FinalStateLogger {
    ///     fn on_step(&self, _values: &HashMap<String, VariableEnum>, _iteration: usize) {}
    ///
    ///     fn on_optimization_complete(&self, values: &HashMap<String, VariableEnum>, iterations: usize) {
    ///         println!("Optimization completed after {} iterations with {} variables",
    ///                  iterations, values.len());
    ///     }
    /// }
    /// ```
    fn on_optimization_complete(
        &self,
        _values: &HashMap<String, VariableEnum>,
        _iterations: usize,
    ) {
        // Default implementation does nothing
    }
}

/// Collection of observers for optimization monitoring.
///
/// This struct manages a vector of observers and provides a convenient
/// `notify()` method to call all observers at once. Optimizers use this
/// internally to manage their observers.
///
/// # Usage
///
/// Typically you don't create this directly - use the `add_observer()` method
/// on optimizers. However, you can use it for custom optimization algorithms:
///
/// ```no_run
/// use apex_solver::observers::{OptObserver, OptObserverVec};
/// use apex_solver::core::problem::VariableEnum;
/// use std::collections::HashMap;
///
/// struct MyOptimizer {
///     observers: OptObserverVec,
///     // ... other fields ...
/// }
///
/// impl MyOptimizer {
///     fn step(&mut self, values: &HashMap<String, VariableEnum>, iteration: usize) {
///         // ... optimization logic ...
///
///         // Notify all observers
///         self.observers.notify(values, iteration);
///     }
/// }
/// ```
#[derive(Default)]
pub struct OptObserverVec {
    observers: Vec<Box<dyn OptObserver>>,
}

impl OptObserverVec {
    /// Create a new empty observer collection.
    pub fn new() -> Self {
        Self {
            observers: Vec::new(),
        }
    }

    /// Add an observer to the collection.
    ///
    /// The observer will be called at each optimization iteration in the order
    /// it was added.
    ///
    /// # Arguments
    ///
    /// * `observer` - Any type implementing `OptObserver`
    ///
    /// # Examples
    ///
    /// ```no_run
    /// use apex_solver::observers::{OptObserver, OptObserverVec};
    /// use apex_solver::core::problem::VariableEnum;
    /// use std::collections::HashMap;
    ///
    /// struct MyObserver;
    /// impl OptObserver for MyObserver {
    ///     fn on_step(&self, _values: &HashMap<String, VariableEnum>, _iteration: usize) {
    ///         // Handle optimization step
    ///     }
    /// }
    ///
    /// let mut observers = OptObserverVec::new();
    /// observers.add(MyObserver);
    /// ```
    pub fn add(&mut self, observer: impl OptObserver + 'static) {
        self.observers.push(Box::new(observer));
    }

    /// Set iteration metrics for all observers.
    ///
    /// Calls `set_iteration_metrics()` on each registered observer. This should
    /// be called before `notify()` to provide optimization metrics.
    ///
    /// # Arguments
    ///
    /// * `cost` - Current cost function value
    /// * `gradient_norm` - L2 norm of the gradient vector
    /// * `damping` - Damping parameter (may be None)
    /// * `step_norm` - L2 norm of the parameter update step
    /// * `step_quality` - Step quality metric (may be None)
    #[inline]
    pub fn set_iteration_metrics(
        &self,
        cost: f64,
        gradient_norm: f64,
        damping: Option<f64>,
        step_norm: f64,
        step_quality: Option<f64>,
    ) {
        for observer in &self.observers {
            observer.set_iteration_metrics(cost, gradient_norm, damping, step_norm, step_quality);
        }
    }

    /// Set matrix data for all observers.
    ///
    /// Calls `set_matrix_data()` on each registered observer. This should
    /// be called before `notify()` to provide matrix data for visualization.
    ///
    /// # Arguments
    ///
    /// * `hessian` - Sparse Hessian matrix
    /// * `gradient` - Gradient vector
    #[inline]
    pub fn set_matrix_data(
        &self,
        hessian: Option<sparse::SparseColMat<usize, f64>>,
        gradient: Option<Mat<f64>>,
    ) {
        for observer in &self.observers {
            observer.set_matrix_data(hessian.clone(), gradient.clone());
        }
    }

    /// Notify all observers with current optimization state.
    ///
    /// Calls `on_step()` on each registered observer in order. If no observers
    /// are registered, this is a no-op with zero overhead.
    ///
    /// # Arguments
    ///
    /// * `values` - Current variable values
    /// * `iteration` - Current iteration number
    ///
    /// # Examples
    ///
    /// ```no_run
    /// use apex_solver::observers::OptObserverVec;
    /// use std::collections::HashMap;
    ///
    /// let observers = OptObserverVec::new();
    /// let values = HashMap::new();
    ///
    /// // Notify all observers (safe even if empty)
    /// observers.notify(&values, 0);
    /// ```
    #[inline]
    pub fn notify(&self, values: &HashMap<String, VariableEnum>, iteration: usize) {
        for observer in &self.observers {
            observer.on_step(values, iteration);
        }
    }

    /// Notify all observers that optimization is complete.
    ///
    /// Calls `on_optimization_complete()` on each registered observer. This should
    /// be called once at the end of optimization, after all iterations are done.
    ///
    /// # Arguments
    ///
    /// * `values` - Final optimized variable values
    /// * `iterations` - Total number of iterations performed
    ///
    /// # Examples
    ///
    /// ```no_run
    /// use apex_solver::observers::OptObserverVec;
    /// use std::collections::HashMap;
    ///
    /// let observers = OptObserverVec::new();
    /// let values = HashMap::new();
    ///
    /// // Notify all observers that optimization is complete
    /// observers.notify_complete(&values, 50);
    /// ```
    #[inline]
    pub fn notify_complete(&self, values: &HashMap<String, VariableEnum>, iterations: usize) {
        for observer in &self.observers {
            observer.on_optimization_complete(values, iterations);
        }
    }

    /// Check if any observers are registered.
    ///
    /// Useful for conditional logic or debugging.
    #[inline]
    pub fn is_empty(&self) -> bool {
        self.observers.is_empty()
    }

    /// Get the number of registered observers.
    #[inline]
    pub fn len(&self) -> usize {
        self.observers.len()
    }
}

#[cfg(test)]
mod tests {
    use super::*;
    use std::sync::{Arc, Mutex};

    #[derive(Clone)]
    struct TestObserver {
        calls: Arc<Mutex<Vec<usize>>>,
    }

    impl OptObserver for TestObserver {
        fn on_step(&self, _values: &HashMap<String, VariableEnum>, iteration: usize) {
            // In test code, we log and ignore mutex poisoning errors since they indicate test bugs
            if let Ok(mut guard) = self.calls.lock().map_err(|e| {
                ObserverError::MutexPoisoned {
                    context: "TestObserver::on_step".to_string(),
                    reason: e.to_string(),
                }
                .log()
            }) {
                guard.push(iteration);
            }
        }
    }

    #[test]
    fn test_empty_observers() {
        let observers = OptObserverVec::new();
        assert!(observers.is_empty());
        assert_eq!(observers.len(), 0);

        // Should not panic with no observers
        observers.notify(&HashMap::new(), 0);
    }

    #[test]
    fn test_single_observer() -> Result<(), ObserverError> {
        let calls = Arc::new(Mutex::new(Vec::new()));
        let observer = TestObserver {
            calls: calls.clone(),
        };

        let mut observers = OptObserverVec::new();
        observers.add(observer);

        assert_eq!(observers.len(), 1);

        observers.notify(&HashMap::new(), 0);
        observers.notify(&HashMap::new(), 1);
        observers.notify(&HashMap::new(), 2);

        let guard = calls.lock().map_err(|e| {
            ObserverError::MutexPoisoned {
                context: "test_single_observer".to_string(),
                reason: e.to_string(),
            }
            .log()
        })?;
        assert_eq!(*guard, vec![0, 1, 2]);
        Ok(())
    }

    #[test]
    fn test_multiple_observers() -> Result<(), ObserverError> {
        let calls1 = Arc::new(Mutex::new(Vec::new()));
        let calls2 = Arc::new(Mutex::new(Vec::new()));

        let observer1 = TestObserver {
            calls: calls1.clone(),
        };
        let observer2 = TestObserver {
            calls: calls2.clone(),
        };

        let mut observers = OptObserverVec::new();
        observers.add(observer1);
        observers.add(observer2);

        assert_eq!(observers.len(), 2);

        observers.notify(&HashMap::new(), 5);

        let guard1 = calls1.lock().map_err(|e| {
            ObserverError::MutexPoisoned {
                context: "test_multiple_observers (calls1)".to_string(),
                reason: e.to_string(),
            }
            .log()
        })?;
        assert_eq!(*guard1, vec![5]);

        let guard2 = calls2.lock().map_err(|e| {
            ObserverError::MutexPoisoned {
                context: "test_multiple_observers (calls2)".to_string(),
                reason: e.to_string(),
            }
            .log()
        })?;
        assert_eq!(*guard2, vec![5]);
        Ok(())
    }
}