scirs2_integrate/analysis/
advanced.rs

1//! Advanced analysis tools for complex dynamical systems
2//!
3//! This module contains sophisticated analysis methods for dynamical systems including:
4//! - Poincaré section analysis for periodic orbit detection
5//! - Lyapunov exponent calculation for chaos detection
6//! - Fractal dimension analysis for strange attractors
7//! - Recurrence analysis for pattern detection
8//! - Continuation analysis for bifurcation detection
9//! - Monodromy analysis for periodic orbit stability
10
11use crate::analysis::types::*;
12use crate::error::{IntegrateError, IntegrateResult};
13use scirs2_core::ndarray::{Array1, Array2};
14use scirs2_core::numeric::Complex64;
15use scirs2_core::random::Rng;
16use std::collections::HashSet;
17
18/// Poincaré section analysis for periodic orbit detection
19pub struct PoincareAnalyzer {
20    /// Section definition (hyperplane normal)
21    pub section_normal: Array1<f64>,
22    /// Point on the section
23    pub section_point: Array1<f64>,
24    /// Crossing direction (1: positive, -1: negative, 0: both)
25    pub crossing_direction: i8,
26    /// Tolerance for section crossing detection
27    pub crossing_tolerance: f64,
28}
29
30impl PoincareAnalyzer {
31    /// Create a new Poincaré analyzer
32    pub fn new(
33        section_normal: Array1<f64>,
34        section_point: Array1<f64>,
35        crossing_direction: i8,
36    ) -> Self {
37        Self {
38            section_normal,
39            section_point,
40            crossing_direction,
41            crossing_tolerance: 1e-8,
42        }
43    }
44
45    /// Analyze trajectory to find Poincaré map
46    pub fn analyze_trajectory(
47        &self,
48        trajectory: &[Array1<f64>],
49        times: &[f64],
50    ) -> IntegrateResult<PoincareMap> {
51        let mut crossings = Vec::new();
52        let mut crossing_times = Vec::new();
53
54        for i in 1..trajectory.len() {
55            if let Some((crossing_point, crossing_time)) =
56                self.detect_crossing(&trajectory[i - 1], &trajectory[i], times[i - 1], times[i])?
57            {
58                crossings.push(crossing_point);
59                crossing_times.push(crossing_time);
60            }
61        }
62
63        // Compute return map if sufficient crossings
64        let return_map = if crossings.len() > 1 {
65            Some(self.compute_return_map(&crossings)?)
66        } else {
67            None
68        };
69
70        // Detect periodic orbits
71        let periodic_orbits = self.detect_periodic_orbits(&crossings)?;
72
73        Ok(PoincareMap {
74            crossings,
75            crossing_times,
76            return_map,
77            periodic_orbits,
78            section_normal: self.section_normal.clone(),
79            section_point: self.section_point.clone(),
80        })
81    }
82
83    /// Detect crossing of Poincaré section
84    fn detect_crossing(
85        &self,
86        point1: &Array1<f64>,
87        point2: &Array1<f64>,
88        t1: f64,
89        t2: f64,
90    ) -> IntegrateResult<Option<(Array1<f64>, f64)>> {
91        // Calculate distances from section
92        let d1 = self.distance_from_section(point1);
93        let d2 = self.distance_from_section(point2);
94
95        // Check for crossing
96        let crossed = match self.crossing_direction {
97            1 => d1 < 0.0 && d2 > 0.0,  // Positive crossing
98            -1 => d1 > 0.0 && d2 < 0.0, // Negative crossing
99            0 => d1 * d2 < 0.0,         // Any crossing
100            _ => false,
101        };
102
103        if !crossed {
104            return Ok(None);
105        }
106
107        // Interpolate crossing point
108        let alpha = d1.abs() / (d1.abs() + d2.abs());
109        let crossing_point = (1.0 - alpha) * point1 + alpha * point2;
110        let crossing_time = (1.0 - alpha) * t1 + alpha * t2;
111
112        Ok(Some((crossing_point, crossing_time)))
113    }
114
115    /// Calculate distance from point to section
116    fn distance_from_section(&self, point: &Array1<f64>) -> f64 {
117        let relative_pos = point - &self.section_point;
118        relative_pos.dot(&self.section_normal)
119    }
120
121    /// Compute return map from crossings
122    fn compute_return_map(&self, crossings: &[Array1<f64>]) -> IntegrateResult<ReturnMap> {
123        let mut current_points = Vec::new();
124        let mut next_points = Vec::new();
125
126        for i in 0..crossings.len() - 1 {
127            // Project points onto section (remove normal component)
128            let current_projected = self.project_to_section(&crossings[i]);
129            let next_projected = self.project_to_section(&crossings[i + 1]);
130
131            current_points.push(current_projected);
132            next_points.push(next_projected);
133        }
134
135        Ok(ReturnMap {
136            current_points,
137            next_points,
138        })
139    }
140
141    /// Project point onto Poincaré section
142    fn project_to_section(&self, point: &Array1<f64>) -> Array1<f64> {
143        let distance = self.distance_from_section(point);
144        point - distance * &self.section_normal
145    }
146
147    /// Detect periodic orbits from crossings
148    fn detect_periodic_orbits(
149        &self,
150        crossings: &[Array1<f64>],
151    ) -> IntegrateResult<Vec<PeriodicOrbit>> {
152        let mut periodic_orbits = Vec::new();
153        let tolerance = 1e-6;
154
155        // Look for approximate returns to previous crossing points
156        for i in 0..crossings.len() {
157            for j in (i + 2)..crossings.len() {
158                let distance = self.euclidean_distance(&crossings[i], &crossings[j]);
159                if distance < tolerance {
160                    // Found potential periodic orbit
161                    let period_length = j - i;
162                    let representative_point = crossings[i].clone();
163
164                    // Verify periodicity by checking intermediate points
165                    let mut is_periodic = true;
166                    for k in 1..period_length {
167                        if i + k < crossings.len() && j + k < crossings.len() {
168                            let dist =
169                                self.euclidean_distance(&crossings[i + k], &crossings[j + k]);
170                            if dist > tolerance {
171                                is_periodic = false;
172                                break;
173                            }
174                        }
175                    }
176
177                    if is_periodic {
178                        // Calculate approximate period in time
179                        let period = (j - i) as f64; // Simplified period estimate
180
181                        periodic_orbits.push(PeriodicOrbit {
182                            representative_point,
183                            period,
184                            stability: StabilityType::Stable, // Would need proper analysis
185                            floquet_multipliers: vec![],      // Would need computation
186                        });
187                    }
188                }
189            }
190        }
191
192        Ok(periodic_orbits)
193    }
194
195    fn euclidean_distance(&self, a: &Array1<f64>, b: &Array1<f64>) -> f64 {
196        (a - b).iter().map(|&x| x * x).sum::<f64>().sqrt()
197    }
198}
199
200/// Poincaré map data structure
201#[derive(Debug, Clone)]
202pub struct PoincareMap {
203    /// Section crossing points
204    pub crossings: Vec<Array1<f64>>,
205    /// Times of crossings
206    pub crossing_times: Vec<f64>,
207    /// Return map data
208    pub return_map: Option<ReturnMap>,
209    /// Detected periodic orbits
210    pub periodic_orbits: Vec<PeriodicOrbit>,
211    /// Section normal vector
212    pub section_normal: Array1<f64>,
213    /// Point on section
214    pub section_point: Array1<f64>,
215}
216
217/// Return map data
218#[derive(Debug, Clone)]
219pub struct ReturnMap {
220    /// Current points
221    pub current_points: Vec<Array1<f64>>,
222    /// Next points in return map
223    pub next_points: Vec<Array1<f64>>,
224}
225
226/// Lyapunov exponent calculator for chaos detection
227pub struct LyapunovCalculator {
228    /// Number of exponents to calculate
229    pub n_exponents: usize,
230    /// Perturbation magnitude for tangent vectors
231    pub perturbation_magnitude: f64,
232    /// Renormalization interval
233    pub renormalization_interval: usize,
234    /// Integration time step
235    pub dt: f64,
236}
237
238impl LyapunovCalculator {
239    /// Create new Lyapunov calculator
240    pub fn new(nexponents: usize, dt: f64) -> Self {
241        Self {
242            n_exponents: nexponents,
243            perturbation_magnitude: 1e-8,
244            renormalization_interval: 100,
245            dt,
246        }
247    }
248
249    /// Calculate Lyapunov exponents using tangent space evolution
250    pub fn calculate_lyapunov_exponents<F>(
251        &self,
252        system: F,
253        initial_state: &Array1<f64>,
254        total_time: f64,
255    ) -> IntegrateResult<Array1<f64>>
256    where
257        F: Fn(&Array1<f64>) -> Array1<f64> + Send + Sync,
258    {
259        let n_steps = (total_time / self.dt) as usize;
260        let dim = initial_state.len();
261
262        if self.n_exponents > dim {
263            return Err(IntegrateError::ValueError(
264                "Number of exponents cannot exceed system dimension".to_string(),
265            ));
266        }
267
268        // Initialize state and tangent vectors
269        let mut state = initial_state.clone();
270        let mut tangent_vectors = self.initialize_tangent_vectors(dim)?;
271        let mut lyapunov_sums = Array1::zeros(self.n_exponents);
272
273        // Main integration loop
274        for step in 0..n_steps {
275            // Evolve main trajectory
276            let derivative = system(&state);
277            state += &(derivative * self.dt);
278
279            // Evolve tangent vectors
280            for i in 0..self.n_exponents {
281                let jacobian = self.compute_jacobian(&system, &state)?;
282                let tangent_derivative = jacobian.dot(&tangent_vectors.column(i));
283                let old_tangent = tangent_vectors.column(i).to_owned();
284                tangent_vectors
285                    .column_mut(i)
286                    .assign(&(&old_tangent + &(tangent_derivative * self.dt)));
287            }
288
289            // Renormalization to prevent overflow
290            if step % self.renormalization_interval == 0 && step > 0 {
291                let (q, r) = self.qr_decomposition(&tangent_vectors)?;
292                tangent_vectors = q;
293
294                // Add to Lyapunov sum
295                for i in 0..self.n_exponents {
296                    lyapunov_sums[i] += r[[i, i]].abs().ln();
297                }
298            }
299        }
300
301        // Final normalization
302        let lyapunov_exponents = lyapunov_sums / total_time;
303
304        Ok(lyapunov_exponents)
305    }
306
307    /// Initialize orthonormal tangent vectors
308    fn initialize_tangent_vectors(&self, dim: usize) -> IntegrateResult<Array2<f64>> {
309        let mut vectors = Array2::zeros((dim, self.n_exponents));
310
311        // Initialize with random vectors
312        let mut rng = scirs2_core::random::rng();
313        for i in 0..self.n_exponents {
314            for j in 0..dim {
315                vectors[[j, i]] = rng.random::<f64>() - 0.5;
316            }
317        }
318
319        // Gram-Schmidt orthogonalization
320        for i in 0..self.n_exponents {
321            // Orthogonalize against previous vectors
322            for j in 0..i {
323                let projection = vectors.column(i).dot(&vectors.column(j));
324                let col_j = vectors.column(j).to_owned();
325                let mut col_i = vectors.column_mut(i);
326                col_i -= &(projection * &col_j);
327            }
328
329            // Normalize
330            let norm = vectors.column(i).iter().map(|&x| x * x).sum::<f64>().sqrt();
331            if norm > 1e-12 {
332                vectors.column_mut(i).mapv_inplace(|x| x / norm);
333            }
334        }
335
336        Ok(vectors)
337    }
338
339    /// Compute Jacobian matrix using finite differences
340    fn compute_jacobian<F>(&self, system: &F, state: &Array1<f64>) -> IntegrateResult<Array2<f64>>
341    where
342        F: Fn(&Array1<f64>) -> Array1<f64>,
343    {
344        let dim = state.len();
345        let mut jacobian = Array2::zeros((dim, dim));
346        let h = self.perturbation_magnitude;
347
348        for j in 0..dim {
349            let mut state_plus = state.clone();
350            let mut state_minus = state.clone();
351
352            state_plus[j] += h;
353            state_minus[j] -= h;
354
355            let f_plus = system(&state_plus);
356            let f_minus = system(&state_minus);
357
358            for i in 0..dim {
359                jacobian[[i, j]] = (f_plus[i] - f_minus[i]) / (2.0 * h);
360            }
361        }
362
363        Ok(jacobian)
364    }
365
366    /// QR decomposition using Gram-Schmidt
367    fn qr_decomposition(
368        &self,
369        matrix: &Array2<f64>,
370    ) -> IntegrateResult<(Array2<f64>, Array2<f64>)> {
371        let (m, n) = matrix.dim();
372        let mut q = matrix.clone();
373        let mut r = Array2::zeros((n, n));
374
375        for j in 0..n {
376            // Orthogonalize against previous columns
377            for i in 0..j {
378                r[[i, j]] = q.column(j).dot(&q.column(i));
379                let col_i = q.column(i).to_owned();
380                let mut col_j = q.column_mut(j);
381                col_j -= &(r[[i, j]] * &col_i);
382            }
383
384            // Normalize
385            r[[j, j]] = q.column(j).iter().map(|&x| x * x).sum::<f64>().sqrt();
386            if r[[j, j]] > 1e-12 {
387                q.column_mut(j).mapv_inplace(|x| x / r[[j, j]]);
388            }
389        }
390
391        Ok((q, r))
392    }
393
394    /// Calculate largest Lyapunov exponent using Wolf's algorithm
395    pub fn calculate_largest_lyapunov_exponent<F>(
396        &self,
397        system: F,
398        initial_state: &Array1<f64>,
399        total_time: f64,
400        min_separation: f64,
401        max_separation: f64,
402    ) -> IntegrateResult<f64>
403    where
404        F: Fn(&Array1<f64>) -> Array1<f64> + Send + Sync,
405    {
406        let n_steps = (total_time / self.dt) as usize;
407        let dim = initial_state.len();
408
409        // Initialize reference trajectory
410        let mut reference_state = initial_state.clone();
411
412        // Initialize nearby trajectory with small perturbation
413        let mut nearby_state = initial_state.clone();
414        nearby_state[0] += self.perturbation_magnitude;
415
416        let mut lyapunov_sum = 0.0;
417        let mut n_rescales = 0;
418
419        for _step in 0..n_steps {
420            // Evolve both trajectories
421            let ref_derivative = system(&reference_state);
422            let nearby_derivative = system(&nearby_state);
423
424            reference_state += &(ref_derivative * self.dt);
425            nearby_state += &(nearby_derivative * self.dt);
426
427            // Calculate separation
428            let separation_vector = &nearby_state - &reference_state;
429            let separation = separation_vector.iter().map(|&x| x * x).sum::<f64>().sqrt();
430
431            // Check if rescaling is needed
432            if (separation > max_separation || separation < min_separation) && separation > 1e-15 {
433                // Add to Lyapunov sum
434                lyapunov_sum += separation.ln();
435                n_rescales += 1;
436
437                // Rescale the separation vector
438                let scale_factor = self.perturbation_magnitude / separation;
439                nearby_state = &reference_state + &(separation_vector * scale_factor);
440            }
441        }
442
443        if n_rescales > 0 {
444            Ok(lyapunov_sum / total_time)
445        } else {
446            Ok(0.0) // No chaos detected
447        }
448    }
449
450    /// Estimate Lyapunov exponent from time series using delay embedding
451    pub fn estimate_lyapunov_from_timeseries(
452        &self,
453        timeseries: &Array1<f64>,
454        embedding_dimension: usize,
455        delay: usize,
456    ) -> IntegrateResult<f64> {
457        let n = timeseries.len();
458        if n < embedding_dimension * delay + 1 {
459            return Err(IntegrateError::ValueError(
460                "Time series too short for embedding".to_string(),
461            ));
462        }
463
464        // Create delay embedding
465        let n_vectors = n - embedding_dimension * delay;
466        let mut embedded_vectors = Vec::new();
467
468        for i in 0..n_vectors {
469            let mut vector = Array1::zeros(embedding_dimension);
470            for j in 0..embedding_dimension {
471                vector[j] = timeseries[i + j * delay];
472            }
473            embedded_vectors.push(vector);
474        }
475
476        // Calculate nearest neighbor distances and their evolution
477        let mut lyapunov_sum = 0.0;
478        let mut count = 0;
479        let min_time_separation = 2 * delay; // Avoid temporal correlations
480
481        for i in 0..embedded_vectors.len() - 1 {
482            // Find nearest neighbor with sufficient time separation
483            let mut min_distance = f64::INFINITY;
484            let mut nearest_index = None;
485
486            for j in 0..embedded_vectors.len() - 1 {
487                if (j as i32 - i as i32).abs() >= min_time_separation as i32 {
488                    let distance =
489                        self.euclidean_distance_arrays(&embedded_vectors[i], &embedded_vectors[j]);
490                    if distance < min_distance && distance > 1e-12 {
491                        min_distance = distance;
492                        nearest_index = Some(j);
493                    }
494                }
495            }
496
497            if let Some(j) = nearest_index {
498                // Calculate distance after one time step
499                if i + 1 < embedded_vectors.len() && j + 1 < embedded_vectors.len() {
500                    let initial_distance = min_distance;
501                    let final_distance = self.euclidean_distance_arrays(
502                        &embedded_vectors[i + 1],
503                        &embedded_vectors[j + 1],
504                    );
505
506                    if final_distance > 1e-12 && initial_distance > 1e-12 {
507                        lyapunov_sum += (final_distance / initial_distance).ln();
508                        count += 1;
509                    }
510                }
511            }
512        }
513
514        if count > 0 {
515            Ok(lyapunov_sum / (count as f64))
516        } else {
517            Ok(0.0)
518        }
519    }
520
521    /// Helper function for distance calculation between arrays
522    fn euclidean_distance_arrays(&self, a: &Array1<f64>, b: &Array1<f64>) -> f64 {
523        (a - b).iter().map(|&x| x * x).sum::<f64>().sqrt()
524    }
525
526    /// Calculate Lyapunov spectrum with error estimates
527    pub fn calculate_lyapunov_spectrum_with_errors<F>(
528        &self,
529        system: F,
530        initial_state: &Array1<f64>,
531        total_time: f64,
532        n_trials: usize,
533    ) -> IntegrateResult<(Array1<f64>, Array1<f64>)>
534    where
535        F: Fn(&Array1<f64>) -> Array1<f64> + Send + Sync + Clone,
536    {
537        let dim = initial_state.len();
538        let n_exponents = self.n_exponents.min(dim);
539
540        let mut all_exponents = Array2::zeros((n_trials, n_exponents));
541
542        // Calculate Lyapunov exponents multiple times with slightly different initial conditions
543        let mut rng = scirs2_core::random::rng();
544
545        for trial in 0..n_trials {
546            // Add small random perturbation to initial state
547            let mut perturbed_initial = initial_state.clone();
548            for i in 0..dim {
549                perturbed_initial[i] += (rng.random::<f64>() - 0.5) * 1e-6;
550            }
551
552            let exponents =
553                self.calculate_lyapunov_exponents(system.clone(), &perturbed_initial, total_time)?;
554
555            for i in 0..n_exponents {
556                all_exponents[[trial, i]] = exponents[i];
557            }
558        }
559
560        // Calculate mean and standard deviation
561        let mut means = Array1::zeros(n_exponents);
562        let mut std_devs = Array1::zeros(n_exponents);
563
564        for i in 0..n_exponents {
565            let column = all_exponents.column(i);
566            let mean = column.sum() / n_trials as f64;
567            let variance =
568                column.iter().map(|&x| (x - mean).powi(2)).sum::<f64>() / n_trials as f64;
569
570            means[i] = mean;
571            std_devs[i] = variance.sqrt();
572        }
573
574        Ok((means, std_devs))
575    }
576}
577
578/// Fractal dimension analyzer for strange attractors
579pub struct FractalAnalyzer {
580    /// Range of scales to analyze
581    pub scale_range: (f64, f64),
582    /// Number of scale points
583    pub n_scales: usize,
584    /// Box-counting parameters
585    pub box_counting_method: BoxCountingMethod,
586}
587
588/// Box counting methods
589#[derive(Debug, Clone, Copy)]
590pub enum BoxCountingMethod {
591    /// Standard box counting
592    Standard,
593    /// Differential box counting
594    Differential,
595    /// Correlation dimension
596    Correlation,
597}
598
599impl Default for FractalAnalyzer {
600    fn default() -> Self {
601        Self::new()
602    }
603}
604
605impl FractalAnalyzer {
606    /// Create new fractal analyzer
607    pub fn new() -> Self {
608        Self {
609            scale_range: (1e-4, 1e-1),
610            n_scales: 20,
611            box_counting_method: BoxCountingMethod::Standard,
612        }
613    }
614
615    /// Calculate fractal dimension of attractor
616    pub fn calculate_fractal_dimension(
617        &self,
618        attractor_points: &[Array1<f64>],
619    ) -> IntegrateResult<FractalDimension> {
620        match self.box_counting_method {
621            BoxCountingMethod::Standard => self.box_counting_dimension(attractor_points),
622            BoxCountingMethod::Differential => self.differential_box_counting(attractor_points),
623            BoxCountingMethod::Correlation => self.correlation_dimension(attractor_points),
624        }
625    }
626
627    /// Standard box counting dimension
628    fn box_counting_dimension(&self, points: &[Array1<f64>]) -> IntegrateResult<FractalDimension> {
629        if points.is_empty() {
630            return Err(IntegrateError::ValueError(
631                "Cannot analyze empty point set".to_string(),
632            ));
633        }
634
635        let dim = points[0].len();
636
637        // Find bounding box
638        let (min_bounds, max_bounds) = self.find_bounding_box(points);
639        let domain_size = max_bounds
640            .iter()
641            .zip(min_bounds.iter())
642            .map(|(&max, &min)| max - min)
643            .fold(0.0f64, |acc, x| acc.max(x));
644
645        let mut scales = Vec::new();
646        let mut counts = Vec::new();
647
648        // Analyze different scales
649        for i in 0..self.n_scales {
650            let t = i as f64 / (self.n_scales - 1) as f64;
651            let scale = self.scale_range.0 * (self.scale_range.1 / self.scale_range.0).powf(t);
652
653            let box_size = scale * domain_size;
654            let count = self.count_occupied_boxes(points, &min_bounds, box_size, dim)?;
655
656            scales.push(scale);
657            counts.push(count as f64);
658        }
659
660        // Linear regression on log-log plot
661        // For box counting: dimension = -slope of log(count) vs log(scale)
662        let slope = self.calculate_slope_from_log_data(&scales, &counts)?;
663        let dimension = -slope;
664
665        let r_squared = self.calculate_r_squared(&scales, &counts, slope)?;
666
667        Ok(FractalDimension {
668            dimension,
669            method: self.box_counting_method,
670            scales,
671            counts,
672            r_squared,
673        })
674    }
675
676    /// Differential box counting for higher accuracy
677    fn differential_box_counting(
678        &self,
679        points: &[Array1<f64>],
680    ) -> IntegrateResult<FractalDimension> {
681        // Simplified implementation - would need full differential box counting
682        self.box_counting_dimension(points)
683    }
684
685    /// Correlation dimension using Grassberger-Procaccia algorithm
686    fn correlation_dimension(&self, points: &[Array1<f64>]) -> IntegrateResult<FractalDimension> {
687        let n_points = points.len();
688        let mut scales = Vec::new();
689        let mut correlations = Vec::new();
690
691        for i in 0..self.n_scales {
692            let t = i as f64 / (self.n_scales - 1) as f64;
693            let r = self.scale_range.0 * (self.scale_range.1 / self.scale_range.0).powf(t);
694
695            let mut count = 0;
696            for i in 0..n_points {
697                for j in i + 1..n_points {
698                    let distance = self.euclidean_distance(&points[i], &points[j]);
699                    if distance < r {
700                        count += 1;
701                    }
702                }
703            }
704
705            let correlation = 2.0 * count as f64 / (n_points * (n_points - 1)) as f64;
706
707            scales.push(r);
708            correlations.push(correlation);
709        }
710
711        // Filter out zero correlations for log calculation
712        let filtered_data: Vec<(f64, f64)> = scales
713            .iter()
714            .zip(correlations.iter())
715            .filter(|(_, &c)| c > 0.0)
716            .map(|(&s, &c)| (s, c))
717            .collect();
718
719        if filtered_data.len() < 2 {
720            return Err(IntegrateError::ComputationError(
721                "Insufficient data for correlation dimension calculation".to_string(),
722            ));
723        }
724
725        let filtered_scales: Vec<f64> = filtered_data.iter().map(|(s, _)| *s).collect();
726        let filtered_correlations: Vec<f64> = filtered_data.iter().map(|(_, c)| *c).collect();
727
728        let dimension =
729            self.calculate_slope_from_log_data(&filtered_scales, &filtered_correlations)?;
730
731        Ok(FractalDimension {
732            dimension,
733            method: BoxCountingMethod::Correlation,
734            scales,
735            counts: correlations,
736            r_squared: self.calculate_r_squared(
737                &filtered_scales,
738                &filtered_correlations,
739                dimension,
740            )?,
741        })
742    }
743
744    /// Helper functions
745    fn find_bounding_box(&self, points: &[Array1<f64>]) -> (Array1<f64>, Array1<f64>) {
746        let dim = points[0].len();
747        let mut min_bounds = Array1::from_elem(dim, f64::INFINITY);
748        let mut max_bounds = Array1::from_elem(dim, f64::NEG_INFINITY);
749
750        for point in points {
751            for i in 0..dim {
752                min_bounds[i] = min_bounds[i].min(point[i]);
753                max_bounds[i] = max_bounds[i].max(point[i]);
754            }
755        }
756
757        (min_bounds, max_bounds)
758    }
759
760    fn count_occupied_boxes(
761        &self,
762        points: &[Array1<f64>],
763        min_bounds: &Array1<f64>,
764        box_size: f64,
765        dim: usize,
766    ) -> IntegrateResult<usize> {
767        let mut occupied_boxes = HashSet::new();
768
769        for point in points {
770            let mut box_index = Vec::with_capacity(dim);
771            for i in 0..dim {
772                let index = ((point[i] - min_bounds[i]) / box_size).floor() as i64;
773                box_index.push(index);
774            }
775            occupied_boxes.insert(box_index);
776        }
777
778        Ok(occupied_boxes.len())
779    }
780
781    fn calculate_slope_from_log_data(
782        &self,
783        x_data: &[f64],
784        y_data: &[f64],
785    ) -> IntegrateResult<f64> {
786        if x_data.len() != y_data.len() || x_data.len() < 2 {
787            return Err(IntegrateError::ValueError(
788                "Insufficient data for slope calculation".to_string(),
789            ));
790        }
791
792        let n = x_data.len() as f64;
793        let log_x: Vec<f64> = x_data.iter().map(|&x| x.ln()).collect();
794        let log_y: Vec<f64> = y_data.iter().map(|&y| y.ln()).collect();
795
796        let sum_x: f64 = log_x.iter().sum();
797        let sum_y: f64 = log_y.iter().sum();
798        let sum_xy: f64 = log_x.iter().zip(log_y.iter()).map(|(&x, &y)| x * y).sum();
799        let sum_xx: f64 = log_x.iter().map(|&x| x * x).sum();
800
801        let slope = (n * sum_xy - sum_x * sum_y) / (n * sum_xx - sum_x * sum_x);
802
803        Ok(slope)
804    }
805
806    fn calculate_r_squared(
807        &self,
808        x_data: &[f64],
809        y_data: &[f64],
810        slope: f64,
811    ) -> IntegrateResult<f64> {
812        let log_x: Vec<f64> = x_data.iter().map(|&x| x.ln()).collect();
813        let log_y: Vec<f64> = y_data.iter().map(|&y| y.ln()).collect();
814
815        let mean_y = log_y.iter().sum::<f64>() / log_y.len() as f64;
816        let mean_x = log_x.iter().sum::<f64>() / log_x.len() as f64;
817        let intercept = mean_y - slope * mean_x;
818
819        let mut ss_tot = 0.0;
820        let mut ss_res = 0.0;
821
822        for i in 0..log_y.len() {
823            let y_pred = slope * log_x[i] + intercept;
824            ss_res += (log_y[i] - y_pred).powi(2);
825            ss_tot += (log_y[i] - mean_y).powi(2);
826        }
827
828        let r_squared = 1.0 - (ss_res / ss_tot);
829        Ok(r_squared)
830    }
831
832    fn euclidean_distance(&self, a: &Array1<f64>, b: &Array1<f64>) -> f64 {
833        a.iter()
834            .zip(b.iter())
835            .map(|(&x, &y)| (x - y).powi(2))
836            .sum::<f64>()
837            .sqrt()
838    }
839}
840
841/// Fractal dimension result
842#[derive(Debug, Clone)]
843pub struct FractalDimension {
844    /// Calculated dimension
845    pub dimension: f64,
846    /// Method used
847    pub method: BoxCountingMethod,
848    /// Scale values used
849    pub scales: Vec<f64>,
850    /// Count/correlation values
851    pub counts: Vec<f64>,
852    /// Quality of fit (R²)
853    pub r_squared: f64,
854}
855
856/// Recurrence analysis for detecting patterns and periodicities
857pub struct RecurrenceAnalyzer {
858    /// Recurrence threshold
859    pub threshold: f64,
860    /// Embedding dimension for delay coordinate embedding
861    pub embedding_dimension: usize,
862    /// Time delay for embedding
863    pub time_delay: usize,
864    /// Distance metric
865    pub distance_metric: DistanceMetric,
866}
867
868/// Distance metrics for recurrence analysis
869#[derive(Debug, Clone, Copy)]
870pub enum DistanceMetric {
871    /// Euclidean distance
872    Euclidean,
873    /// Maximum (Chebyshev) distance
874    Maximum,
875    /// Manhattan distance
876    Manhattan,
877}
878
879impl RecurrenceAnalyzer {
880    /// Create new recurrence analyzer
881    pub fn new(threshold: f64, embedding_dimension: usize, timedelay: usize) -> Self {
882        Self {
883            threshold,
884            embedding_dimension,
885            time_delay: timedelay,
886            distance_metric: DistanceMetric::Euclidean,
887        }
888    }
889
890    /// Perform recurrence analysis
891    pub fn analyze_recurrence(&self, timeseries: &[f64]) -> IntegrateResult<RecurrenceAnalysis> {
892        // Create delay coordinate embedding
893        let embedded_vectors = self.create_embedding(timeseries)?;
894
895        // Compute recurrence matrix
896        let recurrence_matrix = self.compute_recurrence_matrix(&embedded_vectors)?;
897
898        // Calculate recurrence quantification measures
899        let rqa_measures = self.calculate_rqa_measures(&recurrence_matrix)?;
900
901        Ok(RecurrenceAnalysis {
902            recurrence_matrix,
903            embedded_vectors,
904            rqa_measures,
905            threshold: self.threshold,
906            embedding_dimension: self.embedding_dimension,
907            time_delay: self.time_delay,
908        })
909    }
910
911    /// Create delay coordinate embedding
912    fn create_embedding(&self, timeseries: &[f64]) -> IntegrateResult<Vec<Array1<f64>>> {
913        let n = timeseries.len();
914        let embedded_length = n - (self.embedding_dimension - 1) * self.time_delay;
915
916        if embedded_length <= 0 {
917            return Err(IntegrateError::ValueError(
918                "Time series too short for given embedding parameters".to_string(),
919            ));
920        }
921
922        let mut embedded_vectors = Vec::with_capacity(embedded_length);
923
924        for i in 0..embedded_length {
925            let mut vector = Array1::zeros(self.embedding_dimension);
926            for j in 0..self.embedding_dimension {
927                vector[j] = timeseries[i + j * self.time_delay];
928            }
929            embedded_vectors.push(vector);
930        }
931
932        Ok(embedded_vectors)
933    }
934
935    /// Compute recurrence matrix
936    fn compute_recurrence_matrix(&self, vectors: &[Array1<f64>]) -> IntegrateResult<Array2<bool>> {
937        let n = vectors.len();
938        let mut recurrence_matrix = Array2::from_elem((n, n), false);
939
940        for i in 0..n {
941            for j in 0..n {
942                let distance = self.calculate_distance(&vectors[i], &vectors[j]);
943                recurrence_matrix[[i, j]] = distance <= self.threshold;
944            }
945        }
946
947        Ok(recurrence_matrix)
948    }
949
950    /// Calculate distance between vectors
951    fn calculate_distance(&self, a: &Array1<f64>, b: &Array1<f64>) -> f64 {
952        match self.distance_metric {
953            DistanceMetric::Euclidean => a
954                .iter()
955                .zip(b.iter())
956                .map(|(&x, &y)| (x - y).powi(2))
957                .sum::<f64>()
958                .sqrt(),
959            DistanceMetric::Maximum => a
960                .iter()
961                .zip(b.iter())
962                .map(|(&x, &y)| (x - y).abs())
963                .fold(0.0, f64::max),
964            DistanceMetric::Manhattan => a.iter().zip(b.iter()).map(|(&x, &y)| (x - y).abs()).sum(),
965        }
966    }
967
968    /// Calculate Recurrence Quantification Analysis measures
969    fn calculate_rqa_measures(
970        &self,
971        recurrence_matrix: &Array2<bool>,
972    ) -> IntegrateResult<RQAMeasures> {
973        let n = recurrence_matrix.nrows();
974        let total_points = (n * n) as f64;
975
976        // Recurrence rate
977        let recurrent_points = recurrence_matrix.iter().filter(|&&x| x).count() as f64;
978        let recurrence_rate = recurrent_points / total_points;
979
980        // Determinism (percentage of recurrent points forming diagonal lines)
981        let diagonal_lines = self.find_diagonal_lines(recurrence_matrix, 2)?;
982        let diagonal_points: usize = diagonal_lines.iter().map(|line| line.length).sum();
983        let determinism = diagonal_points as f64 / recurrent_points;
984
985        // Average diagonal line length
986        let avg_diagonal_length = if !diagonal_lines.is_empty() {
987            diagonal_points as f64 / diagonal_lines.len() as f64
988        } else {
989            0.0
990        };
991
992        // Maximum diagonal line length
993        let max_diagonal_length = diagonal_lines
994            .iter()
995            .map(|line| line.length)
996            .max()
997            .unwrap_or(0) as f64;
998
999        // Laminarity (percentage of recurrent points forming vertical lines)
1000        let vertical_lines = self.find_vertical_lines(recurrence_matrix, 2)?;
1001        let vertical_points: usize = vertical_lines.iter().map(|line| line.length).sum();
1002        let laminarity = vertical_points as f64 / recurrent_points;
1003
1004        // Trapping time (average vertical line length)
1005        let trapping_time = if !vertical_lines.is_empty() {
1006            vertical_points as f64 / vertical_lines.len() as f64
1007        } else {
1008            0.0
1009        };
1010
1011        Ok(RQAMeasures {
1012            recurrence_rate,
1013            determinism,
1014            avg_diagonal_length,
1015            max_diagonal_length,
1016            laminarity,
1017            trapping_time,
1018        })
1019    }
1020
1021    /// Find diagonal lines in recurrence matrix
1022    fn find_diagonal_lines(
1023        &self,
1024        matrix: &Array2<bool>,
1025        min_length: usize,
1026    ) -> IntegrateResult<Vec<RecurrentLine>> {
1027        let n = matrix.nrows();
1028        let mut lines = Vec::new();
1029
1030        // Check all diagonals
1031        for k in -(n as i32 - 1)..(n as i32) {
1032            let mut current_length = 0;
1033            let mut start_i = 0;
1034            let mut start_j = 0;
1035
1036            let (start_row, start_col) = if k >= 0 {
1037                (0, k as usize)
1038            } else {
1039                ((-k) as usize, 0)
1040            };
1041
1042            let max_steps = n - start_row.max(start_col);
1043
1044            for step in 0..max_steps {
1045                let i = start_row + step;
1046                let j = start_col + step;
1047
1048                if i < n && j < n && matrix[[i, j]] {
1049                    if current_length == 0 {
1050                        start_i = i;
1051                        start_j = j;
1052                    }
1053                    current_length += 1;
1054                } else {
1055                    if current_length >= min_length {
1056                        lines.push(RecurrentLine {
1057                            start_i,
1058                            start_j,
1059                            length: current_length,
1060                            line_type: LineType::Diagonal,
1061                        });
1062                    }
1063                    current_length = 0;
1064                }
1065            }
1066
1067            // Check end of diagonal
1068            if current_length >= min_length {
1069                lines.push(RecurrentLine {
1070                    start_i,
1071                    start_j,
1072                    length: current_length,
1073                    line_type: LineType::Diagonal,
1074                });
1075            }
1076        }
1077
1078        Ok(lines)
1079    }
1080
1081    /// Find vertical lines in recurrence matrix
1082    fn find_vertical_lines(
1083        &self,
1084        matrix: &Array2<bool>,
1085        min_length: usize,
1086    ) -> IntegrateResult<Vec<RecurrentLine>> {
1087        let n = matrix.nrows();
1088        let mut lines = Vec::new();
1089
1090        for j in 0..n {
1091            let mut current_length = 0;
1092            let mut start_i = 0;
1093
1094            for i in 0..n {
1095                if matrix[[i, j]] {
1096                    if current_length == 0 {
1097                        start_i = i;
1098                    }
1099                    current_length += 1;
1100                } else {
1101                    if current_length >= min_length {
1102                        lines.push(RecurrentLine {
1103                            start_i,
1104                            start_j: j,
1105                            length: current_length,
1106                            line_type: LineType::Vertical,
1107                        });
1108                    }
1109                    current_length = 0;
1110                }
1111            }
1112
1113            // Check end of column
1114            if current_length >= min_length {
1115                lines.push(RecurrentLine {
1116                    start_i,
1117                    start_j: j,
1118                    length: current_length,
1119                    line_type: LineType::Vertical,
1120                });
1121            }
1122        }
1123
1124        Ok(lines)
1125    }
1126}
1127
1128/// Recurrence analysis result
1129#[derive(Debug, Clone)]
1130pub struct RecurrenceAnalysis {
1131    /// Recurrence matrix
1132    pub recurrence_matrix: Array2<bool>,
1133    /// Embedded vectors
1134    pub embedded_vectors: Vec<Array1<f64>>,
1135    /// RQA measures
1136    pub rqa_measures: RQAMeasures,
1137    /// Analysis parameters
1138    pub threshold: f64,
1139    pub embedding_dimension: usize,
1140    pub time_delay: usize,
1141}
1142
1143/// Recurrence Quantification Analysis measures
1144#[derive(Debug, Clone)]
1145pub struct RQAMeasures {
1146    /// Recurrence rate
1147    pub recurrence_rate: f64,
1148    /// Determinism
1149    pub determinism: f64,
1150    /// Average diagonal line length
1151    pub avg_diagonal_length: f64,
1152    /// Maximum diagonal line length
1153    pub max_diagonal_length: f64,
1154    /// Laminarity
1155    pub laminarity: f64,
1156    /// Trapping time
1157    pub trapping_time: f64,
1158}
1159
1160/// Recurrent line structure
1161#[derive(Debug, Clone)]
1162pub struct RecurrentLine {
1163    pub start_i: usize,
1164    pub start_j: usize,
1165    pub length: usize,
1166    pub line_type: LineType,
1167}
1168
1169/// Line types in recurrence plot
1170#[derive(Debug, Clone, Copy)]
1171pub enum LineType {
1172    Diagonal,
1173    Vertical,
1174    Horizontal,
1175}
1176
1177/// Advanced continuation and monodromy analysis for bifurcation detection
1178pub struct ContinuationAnalyzer {
1179    /// Parameter range for continuation
1180    pub param_range: (f64, f64),
1181    /// Number of continuation steps
1182    pub n_steps: usize,
1183    /// Convergence tolerance
1184    pub tol: f64,
1185    /// Maximum Newton iterations
1186    pub max_newton_iter: usize,
1187}
1188
1189impl ContinuationAnalyzer {
1190    /// Create new continuation analyzer
1191    pub fn new(paramrange: (f64, f64), n_steps: usize) -> Self {
1192        Self {
1193            param_range: paramrange,
1194            n_steps,
1195            tol: 1e-8,
1196            max_newton_iter: 50,
1197        }
1198    }
1199
1200    /// Perform parameter continuation to trace bifurcation curves
1201    pub fn trace_bifurcation_curve<F>(
1202        &self,
1203        system: F,
1204        initial_state: &Array1<f64>,
1205    ) -> IntegrateResult<ContinuationResult>
1206    where
1207        F: Fn(&Array1<f64>, f64) -> Array1<f64>,
1208    {
1209        let mut bifurcation_points = Vec::new();
1210        let mut fixed_points = Vec::new();
1211
1212        let (param_start, param_end) = self.param_range;
1213        let step = (param_end - param_start) / self.n_steps as f64;
1214
1215        let mut current_state = initial_state.clone();
1216
1217        for i in 0..=self.n_steps {
1218            let param = param_start + i as f64 * step;
1219
1220            // Find fixed point at current parameter
1221            let fixed_point = self.find_fixed_point(&system, &current_state, param)?;
1222
1223            // Compute stability via numerical Jacobian
1224            let jac = self.numerical_jacobian(&system, &fixed_point, param)?;
1225            let eigenvalues = self.compute_eigenvalues(&jac)?;
1226
1227            // Check for bifurcations
1228            if let Some(bif_type) = self.detect_bifurcation(&eigenvalues) {
1229                bifurcation_points.push(BifurcationPointData {
1230                    parameter: param,
1231                    state: fixed_point.clone(),
1232                    bifurcation_type: bif_type,
1233                });
1234            }
1235
1236            fixed_points.push(FixedPointData {
1237                parameter: param,
1238                state: fixed_point.clone(),
1239                eigenvalues: eigenvalues.clone(),
1240                stability: self.classify_stability(&eigenvalues),
1241            });
1242
1243            current_state = fixed_point;
1244        }
1245
1246        Ok(ContinuationResult {
1247            bifurcation_points,
1248            fixed_points,
1249            parameter_range: self.param_range,
1250        })
1251    }
1252
1253    /// Find fixed point using Newton's method
1254    fn find_fixed_point<F>(
1255        &self,
1256        system: &F,
1257        initial_guess: &Array1<f64>,
1258        parameter: f64,
1259    ) -> IntegrateResult<Array1<f64>>
1260    where
1261        F: Fn(&Array1<f64>, f64) -> Array1<f64>,
1262    {
1263        let mut x = initial_guess.clone();
1264
1265        for _ in 0..self.max_newton_iter {
1266            let f = system(&x, parameter);
1267            let norm_f = f.iter().map(|&v| v * v).sum::<f64>().sqrt();
1268
1269            if norm_f < self.tol {
1270                return Ok(x);
1271            }
1272
1273            let jac = self.numerical_jacobian(system, &x, parameter)?;
1274            let delta_x = self.solve_linear_system(&jac, &f)?;
1275
1276            x = &x - &delta_x;
1277        }
1278
1279        Err(IntegrateError::ConvergenceError(
1280            "Fixed point not found".to_string(),
1281        ))
1282    }
1283
1284    /// Compute numerical Jacobian
1285    fn numerical_jacobian<F>(
1286        &self,
1287        system: &F,
1288        x: &Array1<f64>,
1289        parameter: f64,
1290    ) -> IntegrateResult<Array2<f64>>
1291    where
1292        F: Fn(&Array1<f64>, f64) -> Array1<f64>,
1293    {
1294        let n = x.len();
1295        let mut jac = Array2::zeros((n, n));
1296        let h = 1e-8;
1297
1298        for j in 0..n {
1299            let mut x_plus = x.clone();
1300            let mut x_minus = x.clone();
1301            x_plus[j] += h;
1302            x_minus[j] -= h;
1303
1304            let f_plus = system(&x_plus, parameter);
1305            let f_minus = system(&x_minus, parameter);
1306
1307            for i in 0..n {
1308                jac[[i, j]] = (f_plus[i] - f_minus[i]) / (2.0 * h);
1309            }
1310        }
1311
1312        Ok(jac)
1313    }
1314
1315    /// Solve linear system using Gaussian elimination
1316    fn solve_linear_system(
1317        &self,
1318        a: &Array2<f64>,
1319        b: &Array1<f64>,
1320    ) -> IntegrateResult<Array1<f64>> {
1321        let n = a.nrows();
1322        let mut aug = Array2::zeros((n, n + 1));
1323
1324        for i in 0..n {
1325            for j in 0..n {
1326                aug[[i, j]] = a[[i, j]];
1327            }
1328            aug[[i, n]] = b[i];
1329        }
1330
1331        // Forward elimination
1332        for k in 0..n {
1333            let mut max_row = k;
1334            for i in k + 1..n {
1335                if aug[[i, k]].abs() > aug[[max_row, k]].abs() {
1336                    max_row = i;
1337                }
1338            }
1339
1340            for j in 0..n + 1 {
1341                let temp = aug[[k, j]];
1342                aug[[k, j]] = aug[[max_row, j]];
1343                aug[[max_row, j]] = temp;
1344            }
1345
1346            for i in k + 1..n {
1347                let factor = aug[[i, k]] / aug[[k, k]];
1348                for j in k..n + 1 {
1349                    aug[[i, j]] -= factor * aug[[k, j]];
1350                }
1351            }
1352        }
1353
1354        // Back substitution
1355        let mut x = Array1::zeros(n);
1356        for i in (0..n).rev() {
1357            x[i] = aug[[i, n]];
1358            for j in i + 1..n {
1359                x[i] -= aug[[i, j]] * x[j];
1360            }
1361            x[i] /= aug[[i, i]];
1362        }
1363
1364        Ok(x)
1365    }
1366
1367    /// Compute eigenvalues for 2x2 matrices
1368    fn compute_eigenvalues(&self, matrix: &Array2<f64>) -> IntegrateResult<Vec<Complex64>> {
1369        let n = matrix.nrows();
1370
1371        if n == 2 {
1372            let a = matrix[[0, 0]];
1373            let b = matrix[[0, 1]];
1374            let c = matrix[[1, 0]];
1375            let d = matrix[[1, 1]];
1376
1377            let trace = a + d;
1378            let det = a * d - b * c;
1379            let discriminant = trace * trace - 4.0 * det;
1380
1381            if discriminant >= 0.0 {
1382                let sqrt_disc = discriminant.sqrt();
1383                Ok(vec![
1384                    Complex64::new((trace + sqrt_disc) / 2.0, 0.0),
1385                    Complex64::new((trace - sqrt_disc) / 2.0, 0.0),
1386                ])
1387            } else {
1388                let sqrt_disc = (-discriminant).sqrt();
1389                Ok(vec![
1390                    Complex64::new(trace / 2.0, sqrt_disc / 2.0),
1391                    Complex64::new(trace / 2.0, -sqrt_disc / 2.0),
1392                ])
1393            }
1394        } else {
1395            Err(IntegrateError::InvalidInput(
1396                "Only 2x2 matrices supported".to_string(),
1397            ))
1398        }
1399    }
1400
1401    /// Detect bifurcation types
1402    fn detect_bifurcation(&self, eigenvalues: &[Complex64]) -> Option<BifurcationType> {
1403        for eigenval in eigenvalues {
1404            if eigenval.im == 0.0 && eigenval.re.abs() < 1e-6 {
1405                return Some(BifurcationType::SaddleNode);
1406            }
1407
1408            if eigenval.im != 0.0 && eigenval.re.abs() < 1e-6 {
1409                return Some(BifurcationType::Hopf);
1410            }
1411        }
1412        None
1413    }
1414
1415    /// Classify stability
1416    fn classify_stability(&self, eigenvalues: &[Complex64]) -> StabilityType {
1417        for eigenval in eigenvalues {
1418            if eigenval.re > 1e-12 {
1419                return StabilityType::Unstable;
1420            }
1421            if eigenval.re.abs() < 1e-12 {
1422                return StabilityType::Marginally;
1423            }
1424        }
1425        StabilityType::Stable
1426    }
1427}
1428
1429/// Monodromy matrix analyzer for periodic orbits
1430pub struct MonodromyAnalyzer {
1431    pub period: f64,
1432    pub tol: f64,
1433    pub n_steps: usize,
1434}
1435
1436impl MonodromyAnalyzer {
1437    /// Create new monodromy analyzer
1438    pub fn new(period: f64, nsteps: usize) -> Self {
1439        Self {
1440            period,
1441            tol: 1e-8,
1442            n_steps: nsteps,
1443        }
1444    }
1445
1446    /// Compute monodromy matrix
1447    pub fn compute_monodromy_matrix<F>(
1448        &self,
1449        system: F,
1450        initial_state: &Array1<f64>,
1451    ) -> IntegrateResult<MonodromyResult>
1452    where
1453        F: Fn(&Array1<f64>) -> Array1<f64>,
1454    {
1455        let n = initial_state.len();
1456        let dt = self.period / self.n_steps as f64;
1457
1458        // Integrate fundamental matrix
1459        let mut fundamental_matrix = Array2::eye(n);
1460        let mut current_state = initial_state.clone();
1461
1462        for _ in 0..self.n_steps {
1463            let jac = self.numerical_jacobian(&system, &current_state)?;
1464
1465            // Euler step for fundamental matrix: dΦ/dt = J(t)Φ
1466            fundamental_matrix = &fundamental_matrix + &(jac.dot(&fundamental_matrix) * dt);
1467
1468            // RK4 for state evolution
1469            let k1 = system(&current_state);
1470            let k2 = system(&(&current_state + &(&k1 * (dt / 2.0))));
1471            let k3 = system(&(&current_state + &(&k2 * (dt / 2.0))));
1472            let k4 = system(&(&current_state + &(&k3 * dt)));
1473
1474            current_state = &current_state + &((&k1 + &k2 * 2.0 + &k3 * 2.0 + &k4) * (dt / 6.0));
1475        }
1476
1477        let eigenvalues = self.compute_eigenvalues(&fundamental_matrix)?;
1478        let stability = self.classify_periodic_stability(&eigenvalues);
1479
1480        Ok(MonodromyResult {
1481            monodromy_matrix: fundamental_matrix,
1482            eigenvalues,
1483            stability,
1484            period: self.period,
1485            final_state: current_state,
1486        })
1487    }
1488
1489    /// Compute numerical Jacobian
1490    fn numerical_jacobian<F>(&self, system: &F, x: &Array1<f64>) -> IntegrateResult<Array2<f64>>
1491    where
1492        F: Fn(&Array1<f64>) -> Array1<f64>,
1493    {
1494        let n = x.len();
1495        let mut jac = Array2::zeros((n, n));
1496        let h = 1e-8;
1497
1498        for j in 0..n {
1499            let mut x_plus = x.clone();
1500            let mut x_minus = x.clone();
1501            x_plus[j] += h;
1502            x_minus[j] -= h;
1503
1504            let f_plus = system(&x_plus);
1505            let f_minus = system(&x_minus);
1506
1507            for i in 0..n {
1508                jac[[i, j]] = (f_plus[i] - f_minus[i]) / (2.0 * h);
1509            }
1510        }
1511
1512        Ok(jac)
1513    }
1514
1515    /// Compute eigenvalues
1516    fn compute_eigenvalues(&self, matrix: &Array2<f64>) -> IntegrateResult<Vec<Complex64>> {
1517        let n = matrix.nrows();
1518
1519        if n == 2 {
1520            let a = matrix[[0, 0]];
1521            let b = matrix[[0, 1]];
1522            let c = matrix[[1, 0]];
1523            let d = matrix[[1, 1]];
1524
1525            let trace = a + d;
1526            let det = a * d - b * c;
1527            let discriminant = trace * trace - 4.0 * det;
1528
1529            if discriminant >= 0.0 {
1530                let sqrt_disc = discriminant.sqrt();
1531                Ok(vec![
1532                    Complex64::new((trace + sqrt_disc) / 2.0, 0.0),
1533                    Complex64::new((trace - sqrt_disc) / 2.0, 0.0),
1534                ])
1535            } else {
1536                let sqrt_disc = (-discriminant).sqrt();
1537                Ok(vec![
1538                    Complex64::new(trace / 2.0, sqrt_disc / 2.0),
1539                    Complex64::new(trace / 2.0, -sqrt_disc / 2.0),
1540                ])
1541            }
1542        } else {
1543            Err(IntegrateError::InvalidInput(
1544                "Only 2x2 matrices supported".to_string(),
1545            ))
1546        }
1547    }
1548
1549    /// Classify periodic stability
1550    fn classify_periodic_stability(&self, multipliers: &[Complex64]) -> PeriodicStabilityType {
1551        let max_magnitude = multipliers.iter().map(|m| m.norm()).fold(0.0, f64::max);
1552
1553        if max_magnitude > 1.0 + 1e-6 {
1554            PeriodicStabilityType::Unstable
1555        } else if (max_magnitude - 1.0).abs() < 1e-6 {
1556            PeriodicStabilityType::Marginally
1557        } else {
1558            PeriodicStabilityType::Stable
1559        }
1560    }
1561}
1562
1563/// Continuation analysis result
1564#[derive(Debug, Clone)]
1565pub struct ContinuationResult {
1566    pub bifurcation_points: Vec<BifurcationPointData>,
1567    pub fixed_points: Vec<FixedPointData>,
1568    pub parameter_range: (f64, f64),
1569}
1570
1571/// Fixed point with stability data
1572#[derive(Debug, Clone)]
1573pub struct FixedPointData {
1574    pub parameter: f64,
1575    pub state: Array1<f64>,
1576    pub eigenvalues: Vec<Complex64>,
1577    pub stability: StabilityType,
1578}
1579
1580/// Bifurcation point data
1581#[derive(Debug, Clone)]
1582pub struct BifurcationPointData {
1583    pub parameter: f64,
1584    pub state: Array1<f64>,
1585    pub bifurcation_type: BifurcationType,
1586}
1587
1588/// Monodromy analysis result
1589#[derive(Debug, Clone)]
1590pub struct MonodromyResult {
1591    pub monodromy_matrix: Array2<f64>,
1592    pub eigenvalues: Vec<Complex64>,
1593    pub stability: PeriodicStabilityType,
1594    pub period: f64,
1595    pub final_state: Array1<f64>,
1596}
1597
1598/// Extended bifurcation types
1599#[derive(Debug, Clone, Copy, PartialEq)]
1600pub enum BifurcationType {
1601    SaddleNode,
1602    Hopf,
1603    PeriodDoubling,
1604    Transcritical,
1605    Pitchfork,
1606}
1607
1608/// Periodic orbit stability
1609#[derive(Debug, Clone, Copy, PartialEq)]
1610pub enum PeriodicStabilityType {
1611    Stable,
1612    Unstable,
1613    Marginally,
1614}
1615
1616#[cfg(test)]
1617mod tests {
1618    use super::*;
1619
1620    #[test]
1621    fn test_poincare_analyzer() {
1622        // Test with helical trajectory that crosses z = 0 plane
1623        let mut trajectory = Vec::new();
1624        let mut times = Vec::new();
1625
1626        for i in 0..100 {
1627            let t = i as f64 * 0.1;
1628            let x = t.cos();
1629            let y = t.sin();
1630            let z = 0.1 * t.sin(); // Oscillates around z = 0, creating crossings
1631
1632            trajectory.push(Array1::from_vec(vec![x, y, z]));
1633            times.push(t);
1634        }
1635
1636        // Define Poincaré section as z = 0 plane
1637        let section_normal = Array1::from_vec(vec![0.0, 0.0, 1.0]);
1638        let section_point = Array1::from_vec(vec![0.0, 0.0, 0.0]);
1639
1640        let analyzer = PoincareAnalyzer::new(section_normal, section_point, 1);
1641        let result = analyzer.analyze_trajectory(&trajectory, &times).unwrap();
1642
1643        // Should find crossings for this trajectory
1644        assert!(!result.crossings.is_empty());
1645    }
1646
1647    #[test]
1648    fn test_lyapunov_calculator() {
1649        // Test with simple linear system (should have negative Lyapunov exponent)
1650        let system =
1651            |state: &Array1<f64>| -> Array1<f64> { Array1::from_vec(vec![-state[0], -state[1]]) };
1652
1653        let calculator = LyapunovCalculator::new(2, 0.01);
1654        let initial_state = Array1::from_vec(vec![1.0, 1.0]);
1655
1656        let exponents = calculator
1657            .calculate_lyapunov_exponents(system, &initial_state, 10.0)
1658            .unwrap();
1659
1660        // Both exponents should be negative for stable linear system
1661        assert!(exponents[0] < 0.0);
1662        assert!(exponents[1] < 0.0);
1663    }
1664
1665    #[test]
1666    fn test_fractal_analyzer() {
1667        // Create a simple 2D filled area for testing - should have dimension close to 2
1668        let mut points = Vec::new();
1669        let mut rng = scirs2_core::random::rng();
1670
1671        // Generate points uniformly distributed in a square with some noise
1672        for _i in 0..500 {
1673            let x = rng.random::<f64>() * 2.0 - 1.0; // range [-1, 1]
1674            let y = rng.random::<f64>() * 2.0 - 1.0; // range [-1, 1]
1675            let point = Array1::from_vec(vec![x, y]);
1676            points.push(point);
1677        }
1678
1679        // Optimized analyzer configuration for better performance
1680        let mut analyzer = FractalAnalyzer::new();
1681        analyzer.scale_range = (0.1, 0.5); // Adjusted range for better scale coverage
1682        analyzer.n_scales = 5; // Increased scales for more stable slope calculation
1683        analyzer.box_counting_method = BoxCountingMethod::Standard; // Use standard method
1684
1685        let result = analyzer.calculate_fractal_dimension(&points).unwrap();
1686
1687        // Verify the results are mathematically valid
1688        assert!(result.dimension.is_finite(), "Dimension should be finite");
1689        assert!(
1690            result.dimension > 0.0,
1691            "Dimension should be positive, got: {}",
1692            result.dimension
1693        );
1694        assert!(
1695            result.dimension <= 3.0,
1696            "Dimension should not exceed embedding dimension, got: {}",
1697            result.dimension
1698        );
1699        assert!(
1700            result.dimension >= 1.5 && result.dimension <= 2.5,
1701            "2D filled area should have dimension between 1.5 and 2.5, got: {}",
1702            result.dimension
1703        );
1704        assert!(
1705            result.r_squared >= 0.0 && result.r_squared <= 1.0,
1706            "R-squared should be in [0,1], got: {}",
1707            result.r_squared
1708        );
1709
1710        // Verify that the fractal dimension makes sense for a spiral pattern
1711        println!(
1712            "Fractal dimension: {}, R²: {}",
1713            result.dimension, result.r_squared
1714        );
1715    }
1716
1717    #[test]
1718    fn test_recurrence_analyzer() {
1719        // Test with sinusoidal time series
1720        let timeseries: Vec<f64> = (0..100).map(|i| (i as f64 * 0.1).sin()).collect();
1721
1722        let analyzer = RecurrenceAnalyzer::new(0.1, 3, 1);
1723        let result = analyzer.analyze_recurrence(&timeseries).unwrap();
1724
1725        // Should have reasonable recurrence measures
1726        assert!(result.rqa_measures.recurrence_rate > 0.0);
1727        assert!(result.rqa_measures.recurrence_rate <= 1.0);
1728        assert!(result.rqa_measures.determinism >= 0.0);
1729        assert!(result.rqa_measures.determinism <= 1.0);
1730    }
1731}