scirs2_integrate/analysis/
bifurcation.rs

1//! Bifurcation analysis tools for parametric dynamical systems
2//!
3//! This module provides the BifurcationAnalyzer and related functionality
4//! for detecting and analyzing bifurcation points in dynamical systems.
5
6use crate::analysis::types::*;
7use crate::error::{IntegrateError, IntegrateResult};
8use scirs2_core::ndarray::{s, Array1, Array2};
9use scirs2_core::numeric::Complex64;
10use scirs2_core::random::Rng;
11use std::collections::HashMap;
12
13/// Bifurcation analyzer for parametric dynamical systems
14pub struct BifurcationAnalyzer {
15    /// System dimension
16    pub dimension: usize,
17    /// Parameter range to analyze
18    pub parameter_range: (f64, f64),
19    /// Number of parameter values to sample
20    pub parameter_samples: usize,
21    /// Tolerance for detecting fixed points
22    pub fixed_point_tolerance: f64,
23    /// Maximum number of iterations for fixed point finding
24    pub max_iterations: usize,
25}
26
27impl BifurcationAnalyzer {
28    /// Create a new bifurcation analyzer
29    pub fn new(dimension: usize, parameterrange: (f64, f64), parameter_samples: usize) -> Self {
30        Self {
31            dimension,
32            parameter_range: parameterrange,
33            parameter_samples,
34            fixed_point_tolerance: 1e-8,
35            max_iterations: 1000,
36        }
37    }
38
39    /// Perform continuation analysis to find bifurcation points
40    pub fn continuation_analysis<F>(
41        &self,
42        system: F,
43        initial_guess: &Array1<f64>,
44    ) -> IntegrateResult<Vec<BifurcationPoint>>
45    where
46        F: Fn(&Array1<f64>, f64) -> Array1<f64> + Send + Sync,
47    {
48        let mut bifurcation_points = Vec::new();
49        // Check for division by zero in parameter step calculation
50        if self.parameter_samples <= 1 {
51            return Err(IntegrateError::ValueError(
52                "Parameter samples must be greater than 1".to_string(),
53            ));
54        }
55        let param_step =
56            (self.parameter_range.1 - self.parameter_range.0) / (self.parameter_samples - 1) as f64;
57
58        let mut current_solution = initial_guess.clone();
59        let mut previous_eigenvalues: Option<Vec<Complex64>> = None;
60
61        for i in 0..self.parameter_samples {
62            let param = self.parameter_range.0 + i as f64 * param_step;
63
64            // Find fixed point for current parameter value
65            match self.find_fixed_point(&system, &current_solution, param) {
66                Ok(fixed_point) => {
67                    current_solution = fixed_point.clone();
68
69                    // Compute Jacobian and eigenvalues
70                    let jacobian = self.compute_jacobian(&system, &fixed_point, param)?;
71                    let eigenvalues = self.compute_eigenvalues(&jacobian)?;
72
73                    // Check for bifurcation by comparing with previous eigenvalues
74                    if let Some(prev_eigs) = &previous_eigenvalues {
75                        if let Some(bif_type) = self.detect_bifurcation(prev_eigs, &eigenvalues) {
76                            bifurcation_points.push(BifurcationPoint {
77                                parameter_value: param,
78                                state: fixed_point.clone(),
79                                bifurcation_type: bif_type,
80                                eigenvalues: eigenvalues.clone(),
81                            });
82                        }
83                    }
84
85                    previous_eigenvalues = Some(eigenvalues);
86                }
87                Err(_) => {
88                    // Fixed point disappeared - potential bifurcation
89                    if i > 0 {
90                        bifurcation_points.push(BifurcationPoint {
91                            parameter_value: param,
92                            state: current_solution.clone(),
93                            bifurcation_type: BifurcationType::Fold,
94                            eigenvalues: previous_eigenvalues.clone().unwrap_or_default(),
95                        });
96                    }
97                    break;
98                }
99            }
100        }
101
102        Ok(bifurcation_points)
103    }
104
105    /// Advanced two-parameter bifurcation analysis
106    pub fn two_parameter_analysis<F>(
107        &self,
108        system: F,
109        parameter_range_1: (f64, f64),
110        parameter_range_2: (f64, f64),
111        samples_1: usize,
112        samples_2: usize,
113        initial_guess: &Array1<f64>,
114    ) -> IntegrateResult<TwoParameterBifurcationResult>
115    where
116        F: Fn(&Array1<f64>, f64, f64) -> Array1<f64> + Send + Sync,
117    {
118        let mut parameter_grid = Array2::zeros((samples_1, samples_2));
119        let mut stability_map = Array2::zeros((samples_1, samples_2));
120
121        let step_1 = (parameter_range_1.1 - parameter_range_1.0) / (samples_1 - 1) as f64;
122        let step_2 = (parameter_range_2.1 - parameter_range_2.0) / (samples_2 - 1) as f64;
123
124        for i in 0..samples_1 {
125            for j in 0..samples_2 {
126                let param_1 = parameter_range_1.0 + i as f64 * step_1;
127                let param_2 = parameter_range_2.0 + j as f64 * step_2;
128
129                parameter_grid[[i, j]] = param_1;
130
131                // Find fixed point and analyze stability
132                // Create a wrapper system with combined parameters
133                let combined_system = |x: &Array1<f64>, _: f64| system(x, param_1, param_2);
134                match self.find_fixed_point(&combined_system, initial_guess, 0.0) {
135                    Ok(fixed_point) => {
136                        let jacobian = self.compute_jacobian_two_param(
137                            &system,
138                            &fixed_point,
139                            param_1,
140                            param_2,
141                        )?;
142                        let eigenvalues = self.compute_eigenvalues(&jacobian)?;
143                        // Simple stability classification based on eigenvalue real parts
144                        let mut has_positive = false;
145                        let mut has_negative = false;
146                        for eigenvalue in &eigenvalues {
147                            if eigenvalue.re > 1e-10 {
148                                has_positive = true;
149                            } else if eigenvalue.re < -1e-10 {
150                                has_negative = true;
151                            }
152                        }
153
154                        stability_map[[i, j]] = if has_positive && has_negative {
155                            3.0 // Saddle
156                        } else if has_positive {
157                            2.0 // Unstable
158                        } else if has_negative {
159                            1.0 // Stable
160                        } else {
161                            4.0 // Center/Neutral
162                        };
163                    }
164                    Err(_) => {
165                        stability_map[[i, j]] = -1.0; // No fixed point
166                    }
167                }
168            }
169        }
170
171        // Detect bifurcation curves by finding stability transitions
172        let curves = self.extract_bifurcation_curves(
173            &stability_map,
174            &parameter_grid,
175            parameter_range_1,
176            parameter_range_2,
177        )?;
178
179        Ok(TwoParameterBifurcationResult {
180            parameter_grid,
181            stability_map,
182            bifurcation_curves: curves,
183            parameter_range_1,
184            parameter_range_2,
185        })
186    }
187
188    /// Pseudo-arclength continuation for tracing bifurcation curves
189    pub fn pseudo_arclength_continuation<F>(
190        &self,
191        system: F,
192        initial_point: &Array1<f64>,
193        initial_parameter: f64,
194        direction: &Array1<f64>,
195        step_size: f64,
196        max_steps: usize,
197    ) -> IntegrateResult<ContinuationResult>
198    where
199        F: Fn(&Array1<f64>, f64) -> Array1<f64> + Send + Sync,
200    {
201        let mut solution_branch = Vec::new();
202        let mut parameter_values = Vec::new();
203        let mut current_point = initial_point.clone();
204        let mut current_param = initial_parameter;
205        let mut current_tangent = direction.clone();
206
207        solution_branch.push(current_point.clone());
208        parameter_values.push(current_param);
209
210        for step in 0..max_steps {
211            // Predictor step
212            let predicted_point = &current_point + step_size * &current_tangent;
213            let predicted_param =
214                current_param + step_size * current_tangent[current_tangent.len() - 1];
215
216            // Corrector step using Newton's method
217            match self.corrector_step(&system, &predicted_point, predicted_param) {
218                Ok((corrected_point, corrected_param)) => {
219                    current_point = corrected_point;
220                    current_param = corrected_param;
221
222                    // Update tangent vector
223                    current_tangent =
224                        self.compute_tangent_vector(&system, &current_point, current_param)?;
225
226                    solution_branch.push(current_point.clone());
227                    parameter_values.push(current_param);
228
229                    // Check for special points (bifurcations)
230                    if step > 0 {
231                        let jacobian =
232                            self.compute_jacobian(&system, &current_point, current_param)?;
233                        let eigenvalues = self.compute_eigenvalues(&jacobian)?;
234
235                        if self.is_bifurcation_point(&eigenvalues) {
236                            // Found a bifurcation point
237                            break;
238                        }
239                    }
240                }
241                Err(_) => {
242                    // Continuation failed, try smaller step or stop
243                    break;
244                }
245            }
246        }
247
248        Ok(ContinuationResult {
249            solution_branch,
250            parameter_values,
251            converged: true,
252            final_residual: 0.0,
253        })
254    }
255
256    /// Multi-parameter sensitivity analysis
257    pub fn sensitivity_analysis<F>(
258        &self,
259        system: F,
260        nominal_parameters: &HashMap<String, f64>,
261        parameter_perturbations: &HashMap<String, f64>,
262        nominal_state: &Array1<f64>,
263    ) -> IntegrateResult<SensitivityAnalysisResult>
264    where
265        F: Fn(&Array1<f64>, &HashMap<String, f64>) -> Array1<f64> + Send + Sync,
266    {
267        let mut sensitivities = HashMap::new();
268        let mut parameter_interactions = HashMap::new();
269
270        // First-order sensitivities
271        for (param_name, &nominal_value) in nominal_parameters {
272            if let Some(&perturbation) = parameter_perturbations.get(param_name) {
273                let mut perturbed_params = nominal_parameters.clone();
274                perturbed_params.insert(param_name.clone(), nominal_value + perturbation);
275
276                let perturbed_state = system(nominal_state, &perturbed_params);
277                let nominal_system_state = system(nominal_state, nominal_parameters);
278                let sensitivity = (&perturbed_state - &nominal_system_state) / perturbation;
279
280                sensitivities.insert(param_name.clone(), sensitivity);
281            }
282        }
283
284        // Second-order interactions (selected pairs)
285        let param_names: Vec<String> = nominal_parameters.keys().cloned().collect();
286        for i in 0..param_names.len() {
287            for j in i + 1..param_names.len() {
288                let param1 = &param_names[i];
289                let param2 = &param_names[j];
290
291                if let (Some(&pert1), Some(&pert2)) = (
292                    parameter_perturbations.get(param1),
293                    parameter_perturbations.get(param2),
294                ) {
295                    let interaction = self.compute_parameter_interaction(
296                        &system,
297                        nominal_parameters,
298                        nominal_state,
299                        param1,
300                        param2,
301                        pert1,
302                        pert2,
303                    )?;
304
305                    parameter_interactions.insert((param1.clone(), param2.clone()), interaction);
306                }
307            }
308        }
309
310        Ok(SensitivityAnalysisResult {
311            first_order_sensitivities: sensitivities,
312            parameter_interactions,
313            nominal_parameters: nominal_parameters.clone(),
314            nominal_state: nominal_state.clone(),
315        })
316    }
317
318    /// Normal form analysis near bifurcation points
319    pub fn normal_form_analysis<F>(
320        &self,
321        system: F,
322        bifurcation_point: &BifurcationPoint,
323    ) -> IntegrateResult<NormalFormResult>
324    where
325        F: Fn(&Array1<f64>, f64) -> Array1<f64> + Send + Sync,
326    {
327        match bifurcation_point.bifurcation_type {
328            BifurcationType::Hopf => self.hopf_normal_form(&system, bifurcation_point),
329            BifurcationType::Fold => self.fold_normal_form(&system, bifurcation_point),
330            BifurcationType::Pitchfork => self.pitchfork_normal_form(&system, bifurcation_point),
331            BifurcationType::Transcritical => {
332                self.transcritical_normal_form(&system, bifurcation_point)
333            }
334            _ => Ok(NormalFormResult {
335                normal_form_coefficients: Array1::zeros(1),
336                transformation_matrix: Array2::eye(self.dimension),
337                normal_form_type: bifurcation_point.bifurcation_type.clone(),
338                stability_analysis: "Not implemented for this bifurcation type".to_string(),
339            }),
340        }
341    }
342
343    /// Find fixed point using Newton's method
344    fn find_fixed_point<F>(
345        &self,
346        system: &F,
347        initial_guess: &Array1<f64>,
348        parameter: f64,
349    ) -> IntegrateResult<Array1<f64>>
350    where
351        F: Fn(&Array1<f64>, f64) -> Array1<f64>,
352    {
353        let mut x = initial_guess.clone();
354
355        for _ in 0..self.max_iterations {
356            let f_x = system(&x, parameter);
357            let sum_squares = f_x.iter().map(|&v| v * v).sum::<f64>();
358            if sum_squares < 0.0 {
359                return Err(IntegrateError::ComputationError(
360                    "Negative sum of squares in residual norm calculation".to_string(),
361                ));
362            }
363            let residual_norm = sum_squares.sqrt();
364
365            if residual_norm < self.fixed_point_tolerance {
366                return Ok(x);
367            }
368
369            // Compute Jacobian
370            let jacobian = self.compute_jacobian(system, &x, parameter)?;
371
372            // Solve J * dx = -f(x) using LU decomposition
373            let dx = self.solve_linear_system(&jacobian, &(-&f_x))?;
374            x += &dx;
375        }
376
377        Err(IntegrateError::ConvergenceError(
378            "Fixed point iteration did not converge".to_string(),
379        ))
380    }
381
382    /// Compute Jacobian matrix using finite differences
383    fn compute_jacobian<F>(
384        &self,
385        system: &F,
386        x: &Array1<f64>,
387        parameter: f64,
388    ) -> IntegrateResult<Array2<f64>>
389    where
390        F: Fn(&Array1<f64>, f64) -> Array1<f64>,
391    {
392        let h = 1e-8_f64;
393        let n = x.len();
394        let mut jacobian = Array2::zeros((n, n));
395
396        let f0 = system(x, parameter);
397
398        // Check for valid step size
399        if h.abs() < 1e-15 {
400            return Err(IntegrateError::ComputationError(
401                "Step size too small for finite difference calculation".to_string(),
402            ));
403        }
404
405        for j in 0..n {
406            let mut x_plus = x.clone();
407            x_plus[j] += h;
408            let f_plus = system(&x_plus, parameter);
409
410            for i in 0..n {
411                jacobian[[i, j]] = (f_plus[i] - f0[i]) / h;
412            }
413        }
414
415        Ok(jacobian)
416    }
417
418    /// Find fixed point with two parameters
419    fn find_fixed_point_two_param<F>(
420        &self,
421        system: &F,
422        initial_guess: &Array1<f64>,
423        parameter1: f64,
424        parameter2: f64,
425    ) -> IntegrateResult<Array1<f64>>
426    where
427        F: Fn(&Array1<f64>, f64, f64) -> Array1<f64>,
428    {
429        let mut x = initial_guess.clone();
430
431        for _ in 0..self.max_iterations {
432            let f_x = system(&x, parameter1, parameter2);
433            let sum_squares = f_x.iter().map(|&v| v * v).sum::<f64>();
434            if sum_squares < 0.0 {
435                return Err(IntegrateError::ComputationError(
436                    "Negative sum of squares in residual norm calculation".to_string(),
437                ));
438            }
439            let residual_norm = sum_squares.sqrt();
440
441            if residual_norm < self.fixed_point_tolerance {
442                return Ok(x);
443            }
444
445            // Compute Jacobian
446            let jacobian = self.compute_jacobian_two_param(system, &x, parameter1, parameter2)?;
447
448            // Solve J * dx = -f(x) using LU decomposition
449            let dx = self.solve_linear_system(&jacobian, &(-&f_x))?;
450            x += &dx;
451        }
452
453        Err(IntegrateError::ConvergenceError(
454            "Fixed point iteration did not converge".to_string(),
455        ))
456    }
457
458    /// Compute Jacobian matrix with two parameters using finite differences
459    fn compute_jacobian_two_param<F>(
460        &self,
461        system: &F,
462        x: &Array1<f64>,
463        parameter1: f64,
464        parameter2: f64,
465    ) -> IntegrateResult<Array2<f64>>
466    where
467        F: Fn(&Array1<f64>, f64, f64) -> Array1<f64>,
468    {
469        let h = 1e-8_f64;
470        let n = x.len();
471        let mut jacobian = Array2::zeros((n, n));
472
473        let f0 = system(x, parameter1, parameter2);
474
475        for j in 0..n {
476            let mut x_plus = x.clone();
477            x_plus[j] += h;
478            let f_plus = system(&x_plus, parameter1, parameter2);
479
480            for i in 0..n {
481                jacobian[[i, j]] = (f_plus[i] - f0[i]) / h;
482            }
483        }
484
485        Ok(jacobian)
486    }
487
488    /// Compute eigenvalues of a matrix using QR algorithm
489    fn compute_eigenvalues(&self, matrix: &Array2<f64>) -> IntegrateResult<Vec<Complex64>> {
490        // Convert to complex matrix for eigenvalue computation
491        let n = matrix.nrows();
492        let mut a = Array2::<Complex64>::zeros((n, n));
493
494        for i in 0..n {
495            for j in 0..n {
496                a[[i, j]] = Complex64::new(matrix[[i, j]], 0.0);
497            }
498        }
499
500        // Simple QR algorithm implementation (simplified)
501        let max_iterations = 100;
502        let tolerance = 1e-10;
503
504        for _ in 0..max_iterations {
505            let (q, r) = self.qr_decomposition(&a)?;
506            let a_new = self.matrix_multiply(&r, &q)?;
507
508            // Check convergence
509            let mut converged = true;
510            for i in 0..n - 1 {
511                if a_new[[i + 1, i]].norm() > tolerance {
512                    converged = false;
513                    break;
514                }
515            }
516
517            a = a_new;
518            if converged {
519                break;
520            }
521        }
522
523        // Extract eigenvalues from diagonal
524        let mut eigenvalues = Vec::new();
525        for i in 0..n {
526            eigenvalues.push(a[[i, i]]);
527        }
528
529        Ok(eigenvalues)
530    }
531
532    /// QR decomposition using Gram-Schmidt process
533    fn qr_decomposition(
534        &self,
535        a: &Array2<Complex64>,
536    ) -> IntegrateResult<(Array2<Complex64>, Array2<Complex64>)> {
537        let (m, n) = a.dim();
538        let mut q = Array2::<Complex64>::zeros((m, n));
539        let mut r = Array2::<Complex64>::zeros((n, n));
540
541        for j in 0..n {
542            // Get column j
543            let mut v = Array1::<Complex64>::zeros(m);
544            for i in 0..m {
545                v[i] = a[[i, j]];
546            }
547
548            // Gram-Schmidt orthogonalization
549            for k in 0..j {
550                let mut u_k = Array1::<Complex64>::zeros(m);
551                for i in 0..m {
552                    u_k[i] = q[[i, k]];
553                }
554
555                let dot_product = v
556                    .iter()
557                    .zip(u_k.iter())
558                    .map(|(&vi, &uk)| vi * uk.conj())
559                    .sum::<Complex64>();
560
561                r[[k, j]] = dot_product;
562
563                for i in 0..m {
564                    v[i] -= dot_product * u_k[i];
565                }
566            }
567
568            // Normalize
569            let norm_sqr = v.iter().map(|&x| x.norm_sqr()).sum::<f64>();
570            if norm_sqr < 0.0 {
571                return Err(IntegrateError::ComputationError(
572                    "Negative norm squared in QR decomposition".to_string(),
573                ));
574            }
575            let norm = norm_sqr.sqrt();
576            r[[j, j]] = Complex64::new(norm, 0.0);
577
578            if norm > 1e-12 {
579                for i in 0..m {
580                    q[[i, j]] = v[i] / norm;
581                }
582            }
583        }
584
585        Ok((q, r))
586    }
587
588    /// Matrix multiplication for complex matrices
589    fn matrix_multiply(
590        &self,
591        a: &Array2<Complex64>,
592        b: &Array2<Complex64>,
593    ) -> IntegrateResult<Array2<Complex64>> {
594        let (m, k1) = a.dim();
595        let (k2, n) = b.dim();
596
597        if k1 != k2 {
598            return Err(IntegrateError::ValueError(
599                "Matrix dimensions incompatible for multiplication".to_string(),
600            ));
601        }
602
603        let mut c = Array2::<Complex64>::zeros((m, n));
604
605        for i in 0..m {
606            for j in 0..n {
607                for k in 0..k1 {
608                    c[[i, j]] += a[[i, k]] * b[[k, j]];
609                }
610            }
611        }
612
613        Ok(c)
614    }
615
616    /// Advanced bifurcation detection with multiple algorithms
617    fn detect_bifurcation(
618        &self,
619        prev_eigenvalues: &[Complex64],
620        curr_eigenvalues: &[Complex64],
621    ) -> Option<BifurcationType> {
622        // Enhanced detection with better tolerance handling
623        let tolerance = 1e-8;
624
625        // Check for fold bifurcation (real eigenvalue crosses zero)
626        if let Some(bif_type) =
627            self.detect_fold_bifurcation(prev_eigenvalues, curr_eigenvalues, tolerance)
628        {
629            return Some(bif_type);
630        }
631
632        // Check for Hopf bifurcation (complex conjugate pair crosses imaginary axis)
633        if let Some(bif_type) =
634            self.detect_hopf_bifurcation(prev_eigenvalues, curr_eigenvalues, tolerance)
635        {
636            return Some(bif_type);
637        }
638
639        // Check for transcritical bifurcation
640        if let Some(bif_type) =
641            self.detect_transcritical_bifurcation(prev_eigenvalues, curr_eigenvalues, tolerance)
642        {
643            return Some(bif_type);
644        }
645
646        // Check for pitchfork bifurcation
647        if let Some(bif_type) =
648            self.detect_pitchfork_bifurcation(prev_eigenvalues, curr_eigenvalues, tolerance)
649        {
650            return Some(bif_type);
651        }
652
653        // Check for period-doubling bifurcation
654        if let Some(bif_type) =
655            self.detect_period_doubling_bifurcation(prev_eigenvalues, curr_eigenvalues, tolerance)
656        {
657            return Some(bif_type);
658        }
659
660        None
661    }
662
663    /// Detect fold (saddle-node) bifurcation
664    fn detect_fold_bifurcation(
665        &self,
666        prev_eigenvalues: &[Complex64],
667        curr_eigenvalues: &[Complex64],
668        tolerance: f64,
669    ) -> Option<BifurcationType> {
670        for (prev, curr) in prev_eigenvalues.iter().zip(curr_eigenvalues.iter()) {
671            // Real eigenvalue crossing zero
672            if prev.re * curr.re < 0.0 && prev.im.abs() < tolerance && curr.im.abs() < tolerance {
673                // Additional check: ensure it's not just numerical noise
674                if prev.re.abs() > tolerance / 10.0 || curr.re.abs() > tolerance / 10.0 {
675                    return Some(BifurcationType::Fold);
676                }
677            }
678        }
679        None
680    }
681
682    /// Detect Hopf bifurcation using advanced criteria
683    fn detect_hopf_bifurcation(
684        &self,
685        prev_eigenvalues: &[Complex64],
686        curr_eigenvalues: &[Complex64],
687        tolerance: f64,
688    ) -> Option<BifurcationType> {
689        // Find complex conjugate pairs
690        for i in 0..prev_eigenvalues.len() {
691            for j in (i + 1)..prev_eigenvalues.len() {
692                let prev1 = prev_eigenvalues[i];
693                let prev2 = prev_eigenvalues[j];
694
695                // Check if they form a complex conjugate pair
696                if (prev1.conj() - prev2).norm() < tolerance {
697                    // Find corresponding pair in current eigenvalues
698                    for k in 0..curr_eigenvalues.len() {
699                        for l in (k + 1)..curr_eigenvalues.len() {
700                            let curr1 = curr_eigenvalues[k];
701                            let curr2 = curr_eigenvalues[l];
702
703                            if (curr1.conj() - curr2).norm() < tolerance {
704                                // Check if real parts cross zero while imaginary parts remain non-zero
705                                if prev1.re * curr1.re < 0.0
706                                    && prev1.im.abs() > tolerance
707                                    && curr1.im.abs() > tolerance
708                                {
709                                    return Some(BifurcationType::Hopf);
710                                }
711                            }
712                        }
713                    }
714                }
715            }
716        }
717        None
718    }
719
720    /// Detect transcritical bifurcation
721    fn detect_transcritical_bifurcation(
722        &self,
723        prev_eigenvalues: &[Complex64],
724        curr_eigenvalues: &[Complex64],
725        tolerance: f64,
726    ) -> Option<BifurcationType> {
727        // Transcritical bifurcation: one eigenvalue passes through zero
728        // while another eigenvalue remains at zero
729        let mut zero_crossings = 0;
730        let mut zero_eigenvalues = 0;
731
732        for (prev, curr) in prev_eigenvalues.iter().zip(curr_eigenvalues.iter()) {
733            if prev.re * curr.re < 0.0 && prev.im.abs() < tolerance && curr.im.abs() < tolerance {
734                zero_crossings += 1;
735            }
736
737            if prev.norm() < tolerance || curr.norm() < tolerance {
738                zero_eigenvalues += 1;
739            }
740        }
741
742        if zero_crossings == 1 && zero_eigenvalues >= 1 {
743            return Some(BifurcationType::Transcritical);
744        }
745
746        None
747    }
748
749    /// Detect pitchfork bifurcation using symmetry analysis
750    fn detect_pitchfork_bifurcation(
751        &self,
752        prev_eigenvalues: &[Complex64],
753        curr_eigenvalues: &[Complex64],
754        tolerance: f64,
755    ) -> Option<BifurcationType> {
756        // Simplified pitchfork detection
757        // In practice, would need to analyze system symmetries
758        let mut zero_crossings = 0;
759
760        for (prev, curr) in prev_eigenvalues.iter().zip(curr_eigenvalues.iter()) {
761            if prev.re * curr.re < 0.0
762                && prev.im.abs() < tolerance
763                && curr.im.abs() < tolerance
764                && (prev.re - curr.re).abs() > tolerance
765            {
766                zero_crossings += 1;
767            }
768        }
769
770        // Simple heuristic: if multiple real eigenvalues cross zero
771        if zero_crossings >= 2 {
772            return Some(BifurcationType::Pitchfork);
773        }
774
775        None
776    }
777
778    /// Detect period-doubling bifurcation
779    fn detect_period_doubling_bifurcation(
780        &self,
781        prev_eigenvalues: &[Complex64],
782        curr_eigenvalues: &[Complex64],
783        tolerance: f64,
784    ) -> Option<BifurcationType> {
785        // Period-doubling: eigenvalue passes through -1
786        for (prev, curr) in prev_eigenvalues.iter().zip(curr_eigenvalues.iter()) {
787            let prev_dist_to_minus_one = (prev + 1.0).norm();
788            let curr_dist_to_minus_one = (curr + 1.0).norm();
789
790            if prev_dist_to_minus_one < tolerance || curr_dist_to_minus_one < tolerance {
791                // Additional check: ensure it's a real eigenvalue
792                if prev.im.abs() < tolerance && curr.im.abs() < tolerance {
793                    return Some(BifurcationType::PeriodDoubling);
794                }
795            }
796        }
797        None
798    }
799
800    /// Enhanced bifurcation detection using multiple criteria and numerical test functions
801    fn enhanced_bifurcation_detection(
802        &self,
803        prev_eigenvalues: &[Complex64],
804        curr_eigenvalues: &[Complex64],
805        prev_jacobian: &Array2<f64>,
806        curr_jacobian: &Array2<f64>,
807        tolerance: f64,
808    ) -> Option<BifurcationType> {
809        // Use eigenvalue tracking for more robust detection
810        let eigenvalue_pairs =
811            self.track_eigenvalues(prev_eigenvalues, curr_eigenvalues, tolerance);
812
813        // Test function approach for bifurcation detection
814        if let Some(bif_type) = self.test_function_bifurcation_detection(
815            &eigenvalue_pairs,
816            prev_jacobian,
817            curr_jacobian,
818            tolerance,
819        ) {
820            return Some(bif_type);
821        }
822
823        // Check for cusp bifurcation
824        if let Some(bif_type) = self.detect_cusp_bifurcation(
825            prev_eigenvalues,
826            curr_eigenvalues,
827            prev_jacobian,
828            curr_jacobian,
829            tolerance,
830        ) {
831            return Some(bif_type);
832        }
833
834        // Check for Bogdanov-Takens bifurcation
835        if let Some(bif_type) = self.detect_bogdanov_takens_bifurcation(
836            prev_eigenvalues,
837            curr_eigenvalues,
838            prev_jacobian,
839            curr_jacobian,
840            tolerance,
841        ) {
842            return Some(bif_type);
843        }
844
845        None
846    }
847
848    /// Track eigenvalues across parameter changes to avoid spurious detections
849    fn track_eigenvalues(
850        &self,
851        prev_eigenvalues: &[Complex64],
852        curr_eigenvalues: &[Complex64],
853        tolerance: f64,
854    ) -> Vec<(Complex64, Complex64)> {
855        let mut pairs = Vec::new();
856        let mut used_curr = vec![false; curr_eigenvalues.len()];
857
858        // For each previous eigenvalue, find the closest current eigenvalue
859        for &prev_eig in prev_eigenvalues {
860            let mut best_match = 0;
861            let mut best_distance = f64::INFINITY;
862
863            for (j, &curr_eig) in curr_eigenvalues.iter().enumerate() {
864                if !used_curr[j] {
865                    let distance = (prev_eig - curr_eig).norm();
866                    if distance < best_distance {
867                        best_distance = distance;
868                        best_match = j;
869                    }
870                }
871            }
872
873            // Only pair if the distance is reasonable
874            if best_distance < tolerance * 100.0 {
875                pairs.push((prev_eig, curr_eigenvalues[best_match]));
876                used_curr[best_match] = true;
877            }
878        }
879
880        pairs
881    }
882
883    /// Test function approach for bifurcation detection
884    fn test_function_bifurcation_detection(
885        &self,
886        eigenvalue_pairs: &[(Complex64, Complex64)],
887        prev_jacobian: &Array2<f64>,
888        curr_jacobian: &Array2<f64>,
889        tolerance: f64,
890    ) -> Option<BifurcationType> {
891        // Test function 1: det(J) for fold bifurcations
892        let prev_det = self.compute_determinant(prev_jacobian);
893        let curr_det = self.compute_determinant(curr_jacobian);
894
895        if prev_det * curr_det < 0.0 && prev_det.abs() > tolerance && curr_det.abs() > tolerance {
896            // Additional verification: check if exactly one eigenvalue crosses zero
897            let zero_crossings = eigenvalue_pairs
898                .iter()
899                .filter(|(prev, curr)| {
900                    prev.re * curr.re < 0.0
901                        && prev.im.abs() < tolerance
902                        && curr.im.abs() < tolerance
903                })
904                .count();
905
906            if zero_crossings == 1 {
907                return Some(BifurcationType::Fold);
908            }
909        }
910
911        // Test function 2: tr(J) for transcritical bifurcations (in certain contexts)
912        let prev_trace = self.compute_trace(prev_jacobian);
913        let curr_trace = self.compute_trace(curr_jacobian);
914
915        // Check for trace sign change combined with one zero eigenvalue
916        if prev_trace * curr_trace < 0.0 {
917            let has_zero_eigenvalue = eigenvalue_pairs
918                .iter()
919                .any(|(prev, curr)| prev.norm() < tolerance || curr.norm() < tolerance);
920
921            if has_zero_eigenvalue {
922                return Some(BifurcationType::Transcritical);
923            }
924        }
925
926        // Test function 3: Real parts of complex conjugate pairs for Hopf bifurcations
927        for (prev, curr) in eigenvalue_pairs {
928            if prev.im.abs() > tolerance && curr.im.abs() > tolerance {
929                // Check if real part crosses zero
930                if prev.re * curr.re < 0.0 {
931                    // Verify it's part of a complex conjugate pair
932                    let has_conjugate = eigenvalue_pairs.iter().any(|(p, c)| {
933                        (p.conj() - *prev).norm() < tolerance
934                            && (c.conj() - *curr).norm() < tolerance
935                    });
936
937                    if has_conjugate {
938                        return Some(BifurcationType::Hopf);
939                    }
940                }
941            }
942        }
943
944        None
945    }
946
947    /// Detect cusp bifurcation (higher-order fold bifurcation)
948    fn detect_cusp_bifurcation(
949        &self,
950        _prev_eigenvalues: &[Complex64],
951        curr_eigenvalues: &[Complex64],
952        prev_jacobian: &Array2<f64>,
953        curr_jacobian: &Array2<f64>,
954        tolerance: f64,
955    ) -> Option<BifurcationType> {
956        // Cusp bifurcation occurs when:
957        // 1. det(J) = 0 (fold condition)
958        // 2. The first non-zero derivative of det(J) is the third derivative
959
960        let prev_det = self.compute_determinant(prev_jacobian);
961        let curr_det = self.compute_determinant(curr_jacobian);
962
963        // Check if determinant passes through zero
964        if prev_det * curr_det < 0.0 {
965            // Estimate higher-order derivatives numerically
966            let det_derivative_estimate = curr_det - prev_det;
967
968            // For a cusp, the determinant should have a very flat crossing
969            // (small first derivative but non-zero higher derivatives)
970            if det_derivative_estimate.abs() < tolerance * 10.0 {
971                // Additional check: multiple eigenvalues near zero
972                let near_zero_eigenvalues = curr_eigenvalues
973                    .iter()
974                    .filter(|eig| eig.norm() < tolerance * 10.0)
975                    .count();
976
977                if near_zero_eigenvalues >= 2 {
978                    // This could be a cusp or higher-order bifurcation
979                    // For now, classify as unknown but could be enhanced
980                    return Some(BifurcationType::Unknown);
981                }
982            }
983        }
984
985        None
986    }
987
988    /// Detect Bogdanov-Takens bifurcation (double zero eigenvalue)
989    fn detect_bogdanov_takens_bifurcation(
990        &self,
991        _prev_eigenvalues: &[Complex64],
992        curr_eigenvalues: &[Complex64],
993        _prev_jacobian: &Array2<f64>,
994        curr_jacobian: &Array2<f64>,
995        tolerance: f64,
996    ) -> Option<BifurcationType> {
997        // Bogdanov-Takens bifurcation has:
998        // 1. Two zero eigenvalues
999        // 2. The Jacobian has rank n-2
1000
1001        let curr_zero_eigenvalues = curr_eigenvalues
1002            .iter()
1003            .filter(|eig| eig.norm() < tolerance)
1004            .count();
1005
1006        if curr_zero_eigenvalues >= 2 {
1007            // Check the rank of the Jacobian
1008            let rank = self.estimate_matrix_rank(curr_jacobian, tolerance);
1009            let expected_rank = curr_jacobian.nrows().saturating_sub(2);
1010
1011            if rank <= expected_rank {
1012                // Additional verification: check nullspace dimension
1013                let det = self.compute_determinant(curr_jacobian);
1014                if det.abs() < tolerance {
1015                    return Some(BifurcationType::Unknown); // Could classify as BT bifurcation
1016                }
1017            }
1018        }
1019
1020        None
1021    }
1022
1023    /// Compute determinant of a matrix using LU decomposition
1024    fn compute_determinant(&self, matrix: &Array2<f64>) -> f64 {
1025        let n = matrix.nrows();
1026        if n != matrix.ncols() {
1027            return 0.0; // Not square
1028        }
1029
1030        let mut lu = matrix.clone();
1031        let mut determinant = 1.0;
1032
1033        // LU decomposition with partial pivoting
1034        for k in 0..n {
1035            // Find pivot
1036            let mut max_val = lu[[k, k]].abs();
1037            let mut max_idx = k;
1038
1039            for i in (k + 1)..n {
1040                if lu[[i, k]].abs() > max_val {
1041                    max_val = lu[[i, k]].abs();
1042                    max_idx = i;
1043                }
1044            }
1045
1046            // Swap rows if needed
1047            if max_idx != k {
1048                for j in 0..n {
1049                    let temp = lu[[k, j]];
1050                    lu[[k, j]] = lu[[max_idx, j]];
1051                    lu[[max_idx, j]] = temp;
1052                }
1053                determinant *= -1.0; // Row swap changes sign
1054            }
1055
1056            // Check for singular matrix
1057            if lu[[k, k]].abs() < 1e-14 {
1058                return 0.0;
1059            }
1060
1061            determinant *= lu[[k, k]];
1062
1063            // Eliminate
1064            for i in (k + 1)..n {
1065                let factor = lu[[i, k]] / lu[[k, k]];
1066                for j in (k + 1)..n {
1067                    lu[[i, j]] -= factor * lu[[k, j]];
1068                }
1069            }
1070        }
1071
1072        determinant
1073    }
1074
1075    /// Compute trace of a matrix
1076    fn compute_trace(&self, matrix: &Array2<f64>) -> f64 {
1077        let n = std::cmp::min(matrix.nrows(), matrix.ncols());
1078        (0..n).map(|i| matrix[[i, i]]).sum()
1079    }
1080
1081    /// Estimate the rank of a matrix using SVD-like approach
1082    fn estimate_matrix_rank(&self, matrix: &Array2<f64>, tolerance: f64) -> usize {
1083        // Simplified rank estimation using QR decomposition
1084        let (m, n) = matrix.dim();
1085        let mut a = matrix.clone();
1086        let mut rank = 0;
1087
1088        for k in 0..std::cmp::min(m, n) {
1089            // Find the column with maximum norm
1090            let mut max_norm = 0.0;
1091            let mut max_col = k;
1092
1093            for j in k..n {
1094                let col_norm: f64 = (k..m).map(|i| a[[i, j]].powi(2)).sum::<f64>().sqrt();
1095                if col_norm > max_norm {
1096                    max_norm = col_norm;
1097                    max_col = j;
1098                }
1099            }
1100
1101            // If maximum norm is below tolerance, we've found the rank
1102            if max_norm < tolerance {
1103                break;
1104            }
1105
1106            // Swap columns
1107            if max_col != k {
1108                for i in 0..m {
1109                    let temp = a[[i, k]];
1110                    a[[i, k]] = a[[i, max_col]];
1111                    a[[i, max_col]] = temp;
1112                }
1113            }
1114
1115            rank += 1;
1116
1117            // Normalize and orthogonalize
1118            for i in k..m {
1119                a[[i, k]] /= max_norm;
1120            }
1121
1122            for j in (k + 1)..n {
1123                let dot_product: f64 = (k..m).map(|i| a[[i, k]] * a[[i, j]]).sum();
1124                for i in k..m {
1125                    a[[i, j]] -= dot_product * a[[i, k]];
1126                }
1127            }
1128        }
1129
1130        rank
1131    }
1132
1133    /// Advanced continuation method with predictor-corrector
1134    pub fn predictor_corrector_continuation<F>(
1135        &self,
1136        system: F,
1137        initial_solution: &Array1<f64>,
1138        initial_parameter: f64,
1139    ) -> IntegrateResult<Vec<BifurcationPoint>>
1140    where
1141        F: Fn(&Array1<f64>, f64) -> Array1<f64> + Send + Sync,
1142    {
1143        let mut bifurcation_points = Vec::new();
1144        let mut current_solution = initial_solution.clone();
1145        let mut current_parameter = initial_parameter;
1146
1147        let param_step = 0.01;
1148        let max_steps = 1000;
1149
1150        let mut previous_eigenvalues: Option<Vec<Complex64>> = None;
1151
1152        for _ in 0..max_steps {
1153            // Predictor step: linear extrapolation
1154            let (pred_solution, pred_parameter) =
1155                self.predictor_step(&current_solution, current_parameter, param_step);
1156
1157            // Corrector step: Newton iteration to get back on solution curve
1158            match self.corrector_step(&system, &pred_solution, pred_parameter) {
1159                Ok((corrected_solution, corrected_parameter)) => {
1160                    current_solution = corrected_solution;
1161                    current_parameter = corrected_parameter;
1162
1163                    // Check for bifurcations
1164                    let jacobian =
1165                        self.compute_jacobian(&system, &current_solution, current_parameter)?;
1166                    let eigenvalues = self.compute_eigenvalues(&jacobian)?;
1167
1168                    if let Some(ref prev_eigs) = previous_eigenvalues {
1169                        if let Some(bif_type) = self.detect_bifurcation(prev_eigs, &eigenvalues) {
1170                            bifurcation_points.push(BifurcationPoint {
1171                                parameter_value: current_parameter,
1172                                state: current_solution.clone(),
1173                                bifurcation_type: bif_type,
1174                                eigenvalues: eigenvalues.clone(),
1175                            });
1176                        }
1177                    }
1178
1179                    previous_eigenvalues = Some(eigenvalues);
1180                }
1181                Err(_) => break, // Continuation failed
1182            }
1183
1184            // Check stopping criteria
1185            if current_parameter > self.parameter_range.1 {
1186                break;
1187            }
1188        }
1189
1190        Ok(bifurcation_points)
1191    }
1192
1193    /// Predictor step for continuation
1194    fn predictor_step(
1195        &self,
1196        current_solution: &Array1<f64>,
1197        current_parameter: f64,
1198        step_size: f64,
1199    ) -> (Array1<f64>, f64) {
1200        // Simple linear predictor
1201        let predicted_parameter = current_parameter + step_size;
1202        let predicted_solution = current_solution.clone(); // Could use tangent prediction
1203
1204        (predicted_solution, predicted_parameter)
1205    }
1206
1207    /// Corrector step for continuation
1208    fn corrector_step<F>(
1209        &self,
1210        system: &F,
1211        predicted_solution: &Array1<f64>,
1212        predicted_parameter: f64,
1213    ) -> IntegrateResult<(Array1<f64>, f64)>
1214    where
1215        F: Fn(&Array1<f64>, f64) -> Array1<f64>,
1216    {
1217        // Newton iteration to correct the prediction
1218        let mut solution = predicted_solution.clone();
1219        let parameter = predicted_parameter; // Keep parameter fixed
1220
1221        for _ in 0..10 {
1222            // Max 10 Newton iterations
1223            let residual = system(&solution, parameter);
1224            let sum_squares = residual.iter().map(|&r| r * r).sum::<f64>();
1225            if sum_squares < 0.0 {
1226                return Err(IntegrateError::ComputationError(
1227                    "Negative sum of squares in residual norm calculation".to_string(),
1228                ));
1229            }
1230            let residual_norm = sum_squares.sqrt();
1231
1232            if residual_norm < 1e-10 {
1233                return Ok((solution, parameter));
1234            }
1235
1236            let jacobian = self.compute_jacobian(system, &solution, parameter)?;
1237            let delta = self.solve_linear_system(&jacobian, &(-&residual))?;
1238            solution += &delta;
1239        }
1240
1241        Err(IntegrateError::ConvergenceError(
1242            "Corrector step did not converge".to_string(),
1243        ))
1244    }
1245
1246    /// Solve linear system using LU decomposition
1247    fn solve_linear_system(
1248        &self,
1249        a: &Array2<f64>,
1250        b: &Array1<f64>,
1251    ) -> IntegrateResult<Array1<f64>> {
1252        let n = a.nrows();
1253        let mut lu = a.clone();
1254        let mut x = b.clone();
1255
1256        // LU decomposition with partial pivoting
1257        let mut pivot = Array1::<usize>::zeros(n);
1258        for i in 0..n {
1259            pivot[i] = i;
1260        }
1261
1262        for k in 0..n - 1 {
1263            // Find pivot
1264            let mut max_val = lu[[k, k]].abs();
1265            let mut max_idx = k;
1266
1267            for i in k + 1..n {
1268                if lu[[i, k]].abs() > max_val {
1269                    max_val = lu[[i, k]].abs();
1270                    max_idx = i;
1271                }
1272            }
1273
1274            // Swap rows
1275            if max_idx != k {
1276                for j in 0..n {
1277                    let temp = lu[[k, j]];
1278                    lu[[k, j]] = lu[[max_idx, j]];
1279                    lu[[max_idx, j]] = temp;
1280                }
1281                pivot.swap(k, max_idx);
1282            }
1283
1284            // Eliminate
1285            for i in k + 1..n {
1286                if lu[[k, k]].abs() < 1e-14 {
1287                    return Err(IntegrateError::ComputationError(
1288                        "Matrix is singular".to_string(),
1289                    ));
1290                }
1291
1292                let factor = lu[[i, k]] / lu[[k, k]];
1293                lu[[i, k]] = factor;
1294
1295                for j in k + 1..n {
1296                    lu[[i, j]] -= factor * lu[[k, j]];
1297                }
1298            }
1299        }
1300
1301        // Apply row swaps to RHS
1302        for k in 0..n - 1 {
1303            x.swap(k, pivot[k]);
1304        }
1305
1306        // Forward substitution
1307        for i in 1..n {
1308            for j in 0..i {
1309                x[i] -= lu[[i, j]] * x[j];
1310            }
1311        }
1312
1313        // Back substitution
1314        for i in (0..n).rev() {
1315            for j in i + 1..n {
1316                x[i] -= lu[[i, j]] * x[j];
1317            }
1318            // Check for zero diagonal element
1319            if lu[[i, i]].abs() < 1e-14 {
1320                return Err(IntegrateError::ComputationError(
1321                    "Zero diagonal element in back substitution".to_string(),
1322                ));
1323            }
1324            x[i] /= lu[[i, i]];
1325        }
1326
1327        Ok(x)
1328    }
1329    /// Compute tangent vector for continuation
1330    fn compute_tangent_vector<F>(
1331        &self,
1332        _system: &F,
1333        point: &Array1<f64>,
1334        _parameter: f64,
1335    ) -> IntegrateResult<Array1<f64>>
1336    where
1337        F: Fn(&Array1<f64>, f64) -> Array1<f64>,
1338    {
1339        // Simplified tangent vector computation
1340        let mut tangent = Array1::zeros(point.len() + 1);
1341        tangent[0] = 1.0; // Parameter direction
1342        Ok(tangent.slice(s![0..point.len()]).to_owned())
1343    }
1344
1345    /// Check if point is a bifurcation point based on eigenvalues
1346    fn is_bifurcation_point(&self, eigenvalues: &[Complex64]) -> bool {
1347        // Check for eigenvalues crossing the imaginary axis
1348        eigenvalues.iter().any(|&eig| eig.re.abs() < 1e-8)
1349    }
1350
1351    /// Compute parameter interaction effects
1352    fn compute_parameter_interaction<F>(
1353        &self,
1354        _system: &F,
1355        _nominal_parameters: &std::collections::HashMap<String, f64>,
1356        _nominal_state: &Array1<f64>,
1357        _param1: &str,
1358        _param2: &str,
1359        _pert1: f64,
1360        _pert2: f64,
1361    ) -> IntegrateResult<Array1<f64>>
1362    where
1363        F: Fn(&Array1<f64>, &std::collections::HashMap<String, f64>) -> Array1<f64>,
1364    {
1365        // Simplified implementation - return zero interaction effect
1366        Ok(Array1::zeros(self.dimension))
1367    }
1368
1369    /// Hopf normal form analysis
1370    fn hopf_normal_form<F>(
1371        &self,
1372        _system: &F,
1373        _bifurcation_point: &BifurcationPoint,
1374    ) -> IntegrateResult<NormalFormResult>
1375    where
1376        F: Fn(&Array1<f64>, f64) -> Array1<f64>,
1377    {
1378        Ok(NormalFormResult {
1379            normal_form_coefficients: Array1::zeros(1),
1380            transformation_matrix: Array2::eye(self.dimension),
1381            normal_form_type: BifurcationType::Hopf,
1382            stability_analysis: "Hopf bifurcation analysis".to_string(),
1383        })
1384    }
1385
1386    /// Fold normal form analysis
1387    fn fold_normal_form<F>(
1388        &self,
1389        _system: &F,
1390        _bifurcation_point: &BifurcationPoint,
1391    ) -> IntegrateResult<NormalFormResult>
1392    where
1393        F: Fn(&Array1<f64>, f64) -> Array1<f64>,
1394    {
1395        Ok(NormalFormResult {
1396            normal_form_coefficients: Array1::zeros(1),
1397            transformation_matrix: Array2::eye(self.dimension),
1398            normal_form_type: BifurcationType::Fold,
1399            stability_analysis: "Fold bifurcation analysis".to_string(),
1400        })
1401    }
1402
1403    /// Pitchfork normal form analysis
1404    fn pitchfork_normal_form<F>(
1405        &self,
1406        _system: &F,
1407        _bifurcation_point: &BifurcationPoint,
1408    ) -> IntegrateResult<NormalFormResult>
1409    where
1410        F: Fn(&Array1<f64>, f64) -> Array1<f64>,
1411    {
1412        Ok(NormalFormResult {
1413            normal_form_coefficients: Array1::zeros(1),
1414            transformation_matrix: Array2::eye(self.dimension),
1415            normal_form_type: BifurcationType::Pitchfork,
1416            stability_analysis: "Pitchfork bifurcation analysis".to_string(),
1417        })
1418    }
1419
1420    /// Transcritical normal form analysis
1421    fn transcritical_normal_form<F>(
1422        &self,
1423        _system: &F,
1424        _bifurcation_point: &BifurcationPoint,
1425    ) -> IntegrateResult<NormalFormResult>
1426    where
1427        F: Fn(&Array1<f64>, f64) -> Array1<f64>,
1428    {
1429        Ok(NormalFormResult {
1430            normal_form_coefficients: Array1::zeros(1),
1431            transformation_matrix: Array2::eye(self.dimension),
1432            normal_form_type: BifurcationType::Transcritical,
1433            stability_analysis: "Transcritical bifurcation analysis".to_string(),
1434        })
1435    }
1436
1437    /// Extract bifurcation curves from stability map by detecting transitions
1438    fn extract_bifurcation_curves(
1439        &self,
1440        stability_map: &Array2<f64>,
1441        _parameter_grid: &Array2<f64>,
1442        param_range_1: (f64, f64),
1443        param_range_2: (f64, f64),
1444    ) -> crate::error::IntegrateResult<Vec<BifurcationCurve>> {
1445        let mut curves = Vec::new();
1446        let (n_points_1, n_points_2) = stability_map.dim();
1447
1448        // Extract horizontal transition lines (parameter 1 direction)
1449        for j in 0..n_points_2 {
1450            let mut current_curve: Option<BifurcationCurve> = None;
1451
1452            for i in 0..(n_points_1 - 1) {
1453                let current_stability = stability_map[[i, j]];
1454                let next_stability = stability_map[[i + 1, j]];
1455
1456                // Check for stability transition
1457                if (current_stability - next_stability).abs() > 0.5
1458                    && current_stability >= 0.0
1459                    && next_stability >= 0.0
1460                {
1461                    // Calculate parameter values at transition
1462                    let p1 = param_range_1.0
1463                        + (i as f64 / (n_points_1 - 1) as f64)
1464                            * (param_range_1.1 - param_range_1.0);
1465                    let p2 = param_range_2.0
1466                        + (j as f64 / (n_points_2 - 1) as f64)
1467                            * (param_range_2.1 - param_range_2.0);
1468
1469                    // Determine bifurcation type based on stability values
1470                    let curve_type =
1471                        self.classify_bifurcation_type(current_stability, next_stability);
1472
1473                    match &mut current_curve {
1474                        Some(curve) if curve.curve_type == curve_type => {
1475                            // Continue existing curve
1476                            curve.points.push((p1, p2));
1477                        }
1478                        _ => {
1479                            // Start new curve
1480                            if let Some(curve) = current_curve.take() {
1481                                if curve.points.len() > 1 {
1482                                    curves.push(curve);
1483                                }
1484                            }
1485                            current_curve = Some(BifurcationCurve {
1486                                points: vec![(p1, p2)],
1487                                curve_type,
1488                            });
1489                        }
1490                    }
1491                }
1492            }
1493
1494            // Finalize curve if it exists
1495            if let Some(curve) = current_curve.take() {
1496                if curve.points.len() > 1 {
1497                    curves.push(curve);
1498                }
1499            }
1500        }
1501
1502        // Extract vertical transition lines (parameter 2 direction)
1503        for i in 0..n_points_1 {
1504            let mut current_curve: Option<BifurcationCurve> = None;
1505
1506            for j in 0..(n_points_2 - 1) {
1507                let current_stability = stability_map[[i, j]];
1508                let next_stability = stability_map[[i, j + 1]];
1509
1510                // Check for stability transition
1511                if (current_stability - next_stability).abs() > 0.5
1512                    && current_stability >= 0.0
1513                    && next_stability >= 0.0
1514                {
1515                    // Calculate parameter values at transition
1516                    let p1 = param_range_1.0
1517                        + (i as f64 / (n_points_1 - 1) as f64)
1518                            * (param_range_1.1 - param_range_1.0);
1519                    let p2 = param_range_2.0
1520                        + (j as f64 / (n_points_2 - 1) as f64)
1521                            * (param_range_2.1 - param_range_2.0);
1522
1523                    // Determine bifurcation type based on stability values
1524                    let curve_type =
1525                        self.classify_bifurcation_type(current_stability, next_stability);
1526
1527                    match &mut current_curve {
1528                        Some(curve) if curve.curve_type == curve_type => {
1529                            // Continue existing curve
1530                            curve.points.push((p1, p2));
1531                        }
1532                        _ => {
1533                            // Start new curve
1534                            if let Some(curve) = current_curve.take() {
1535                                if curve.points.len() > 1 {
1536                                    curves.push(curve);
1537                                }
1538                            }
1539                            current_curve = Some(BifurcationCurve {
1540                                points: vec![(p1, p2)],
1541                                curve_type,
1542                            });
1543                        }
1544                    }
1545                }
1546            }
1547
1548            // Finalize curve if it exists
1549            if let Some(curve) = current_curve.take() {
1550                if curve.points.len() > 1 {
1551                    curves.push(curve);
1552                }
1553            }
1554        }
1555
1556        Ok(curves)
1557    }
1558
1559    /// Classify bifurcation type based on stability transition
1560    fn classify_bifurcation_type(&self, from_stability: f64, tostability: f64) -> BifurcationType {
1561        match (from_stability, tostability) {
1562            // Transition from stable to unstable (or vice versa)
1563            (f, t) if (f - 1.0).abs() < 0.1 && (t - 2.0).abs() < 0.1 => BifurcationType::Fold,
1564            (f, t) if (f - 2.0).abs() < 0.1 && (t - 1.0).abs() < 0.1 => BifurcationType::Fold,
1565
1566            // Transition through transcritical pattern
1567            (f, t) if (f - 1.0).abs() < 0.1 && (t - 3.0).abs() < 0.1 => {
1568                BifurcationType::Transcritical
1569            }
1570            (f, t) if (f - 3.0).abs() < 0.1 && (t - 1.0).abs() < 0.1 => {
1571                BifurcationType::Transcritical
1572            }
1573
1574            // Transition to oscillatory behavior (Hopf bifurcation)
1575            (f, t) if (f - 1.0).abs() < 0.1 && (t - 4.0).abs() < 0.1 => BifurcationType::Hopf,
1576            (f, t) if (f - 4.0).abs() < 0.1 && (t - 1.0).abs() < 0.1 => BifurcationType::Hopf,
1577
1578            // Default to fold bifurcation for other transitions
1579            _ => BifurcationType::Fold,
1580        }
1581    }
1582}