1use crate::dae::types::DAEIndex;
13use crate::dae::utils::{compute_constraint_jacobian, is_singular_matrix};
14use crate::error::{IntegrateError, IntegrateResult};
15use scirs2_core::ndarray::{Array1, Array2, ArrayView1, ArrayView2, ScalarOperand};
16use scirs2_core::numeric::{Float, FromPrimitive};
17use std::fmt::{Debug, Display, LowerExp};
18
19#[derive(Debug, Clone)]
21pub struct DAEStructure<
22 F: Float
23 + FromPrimitive
24 + Debug
25 + ScalarOperand
26 + std::ops::AddAssign
27 + std::ops::SubAssign
28 + std::ops::MulAssign
29 + std::ops::DivAssign
30 + Display
31 + LowerExp
32 + std::iter::Sum
33 + crate::common::IntegrateFloat,
34> {
35 pub n_differential: usize,
37
38 pub n_algebraic: usize,
40
41 pub n_diff_eqs: usize,
43
44 pub n_alg_eqs: usize,
46
47 pub index: DAEIndex,
49
50 pub diff_index: usize,
52
53 pub constraint_jacobian: Option<Array2<F>>,
55
56 pub derivative_jacobian: Option<Array2<F>>,
58
59 pub incidence_matrix: Option<Array2<bool>>,
61
62 pub variable_dependencies: Option<Vec<Vec<usize>>>,
64
65 pub equation_dependencies: Option<Vec<Vec<usize>>>,
67}
68
69impl<
70 F: Float
71 + FromPrimitive
72 + Debug
73 + ScalarOperand
74 + std::ops::AddAssign
75 + std::ops::SubAssign
76 + std::ops::MulAssign
77 + std::ops::DivAssign
78 + Display
79 + LowerExp
80 + std::iter::Sum
81 + crate::common::IntegrateFloat,
82 > Default for DAEStructure<F>
83{
84 fn default() -> Self {
85 DAEStructure {
86 n_differential: 0,
87 n_algebraic: 0,
88 n_diff_eqs: 0,
89 n_alg_eqs: 0,
90 index: DAEIndex::Index1,
91 diff_index: 1,
92 constraint_jacobian: None,
93 derivative_jacobian: None,
94 incidence_matrix: None,
95 variable_dependencies: None,
96 equation_dependencies: None,
97 }
98 }
99}
100
101impl<
102 F: Float
103 + FromPrimitive
104 + Debug
105 + ScalarOperand
106 + std::ops::AddAssign
107 + std::ops::SubAssign
108 + std::ops::MulAssign
109 + std::ops::DivAssign
110 + Display
111 + LowerExp
112 + std::iter::Sum
113 + crate::common::IntegrateFloat,
114 > DAEStructure<F>
115{
116 pub fn new_semi_explicit(n_differential: usize, nalgebraic: usize) -> Self {
118 DAEStructure {
119 n_differential,
120 n_algebraic: nalgebraic,
121 n_diff_eqs: n_differential,
122 n_alg_eqs: nalgebraic,
123 index: DAEIndex::Index1, diff_index: 1, constraint_jacobian: None,
126 derivative_jacobian: None,
127 incidence_matrix: None,
128 variable_dependencies: None,
129 equation_dependencies: None,
130 }
131 }
132
133 pub fn new_fully_implicit(n_equations: usize, nvariables: usize) -> Self {
135 DAEStructure {
136 n_differential: nvariables, n_algebraic: 0,
138 n_diff_eqs: n_equations,
139 n_alg_eqs: 0,
140 index: DAEIndex::Index1, diff_index: 1, constraint_jacobian: None,
143 derivative_jacobian: None,
144 incidence_matrix: None,
145 variable_dependencies: None,
146 equation_dependencies: None,
147 }
148 }
149
150 pub fn compute_index<FFunc, GFunc>(
152 &mut self,
153 t: F,
154 x: ArrayView1<F>,
155 y: ArrayView1<F>,
156 f: &FFunc,
157 g: &GFunc,
158 ) -> IntegrateResult<DAEIndex>
159 where
160 FFunc: Fn(F, ArrayView1<F>, ArrayView1<F>) -> Array1<F>,
161 GFunc: Fn(F, ArrayView1<F>, ArrayView1<F>) -> Array1<F>,
162 {
163 let x_slice: Vec<F> = x.to_vec();
167 let y_slice: Vec<F> = y.to_vec();
168 let g_y = compute_constraint_jacobian(
169 &|t, x, y| g(t, ArrayView1::from(x), ArrayView1::from(y)).to_vec(),
170 t,
171 &x_slice,
172 &y_slice,
173 );
174 self.constraint_jacobian = Some(g_y.clone());
175
176 let singular = is_singular_matrix(g_y.view());
178
179 if !singular {
180 self.index = DAEIndex::Index1;
182 self.diff_index = 1;
183 return Ok(DAEIndex::Index1);
184 }
185
186 let _f_x = compute_jacobian_for_variables(f, t, x, y, 0, self.n_differential)?;
191 let f_y =
192 compute_jacobian_for_variables(f, t, x, y, self.n_differential, self.n_algebraic)?;
193 let g_x = compute_jacobian_for_variables(g, t, x, y, 0, self.n_differential)?;
194
195 let mut index2_matrix = Array2::<F>::zeros((self.n_alg_eqs, self.n_algebraic));
197 let g_x_f_y = g_x.dot(&f_y);
198
199 for i in 0..self.n_alg_eqs {
200 for j in 0..self.n_algebraic {
201 index2_matrix[[i, j]] = g_y[[i, j]];
202 if j < g_x_f_y.dim().1 {
203 index2_matrix[[i, j]] += g_x_f_y[[i, j]];
204 }
205 }
206 }
207
208 let index2_singular = is_singular_matrix(index2_matrix.view());
210
211 if !index2_singular {
212 self.index = DAEIndex::Index2;
214 self.diff_index = 2;
215 return Ok(DAEIndex::Index2);
216 }
217
218 self.index = DAEIndex::Index3;
224 self.diff_index = 3;
225
226 Ok(DAEIndex::Index3)
227 }
228}
229
230pub struct PantelidesReducer<
235 F: Float
236 + FromPrimitive
237 + Debug
238 + ScalarOperand
239 + std::ops::AddAssign
240 + std::ops::SubAssign
241 + std::ops::MulAssign
242 + std::ops::DivAssign
243 + Display
244 + LowerExp
245 + std::iter::Sum
246 + crate::common::IntegrateFloat,
247> {
248 pub structure: DAEStructure<F>,
250
251 pub max_diff_steps: usize,
253
254 pub tol: F,
256}
257
258impl<
259 F: Float
260 + FromPrimitive
261 + Debug
262 + ScalarOperand
263 + std::ops::AddAssign
264 + std::ops::SubAssign
265 + std::ops::MulAssign
266 + std::ops::DivAssign
267 + Display
268 + LowerExp
269 + std::iter::Sum
270 + crate::common::IntegrateFloat,
271 > PantelidesReducer<F>
272{
273 pub fn new(structure: DAEStructure<F>) -> Self {
275 PantelidesReducer {
276 structure,
277 max_diff_steps: 5, tol: F::from_f64(1e-10).unwrap(),
279 }
280 }
281
282 pub fn initialize_incidence_matrix<FFunc, GFunc>(
286 &mut self,
287 t: F,
288 x: ArrayView1<F>,
289 y: ArrayView1<F>,
290 f: &FFunc,
291 g: &GFunc,
292 ) -> IntegrateResult<()>
293 where
294 FFunc: Fn(F, ArrayView1<F>, ArrayView1<F>) -> Array1<F>,
295 GFunc: Fn(F, ArrayView1<F>, ArrayView1<F>) -> Array1<F>,
296 {
297 let n_diff = self.structure.n_differential;
298 let n_alg = self.structure.n_algebraic;
299 let n_diff_eqs = self.structure.n_diff_eqs;
300 let n_alg_eqs = self.structure.n_alg_eqs;
301 let n_vars = n_diff + n_alg;
302 let n_eqs = n_diff_eqs + n_alg_eqs;
303
304 let mut incidence = Array2::<bool>::from_elem((n_eqs, n_vars), false);
306
307 let f_x = compute_jacobian_for_variables(f, t, x, y, 0, n_diff)?;
309 let f_y = compute_jacobian_for_variables(f, t, x, y, n_diff, n_alg)?;
310 let g_x = compute_jacobian_for_variables(g, t, x, y, 0, n_diff)?;
311 let g_y = compute_jacobian_for_variables(g, t, x, y, n_diff, n_alg)?;
312
313 for i in 0..n_diff_eqs {
317 for j in 0..n_diff {
319 if f_x[[i, j]].abs() > self.tol {
320 incidence[[i, j]] = true;
321 }
322 }
323
324 for j in 0..n_alg {
326 if j < f_y.dim().1 && f_y[[i, j]].abs() > self.tol {
327 incidence[[i, n_diff + j]] = true;
328 }
329 }
330 }
331
332 for i in 0..n_alg_eqs {
334 for j in 0..n_diff {
336 if j < g_x.dim().1 && g_x[[i, j]].abs() > self.tol {
337 incidence[[n_diff_eqs + i, j]] = true;
338 }
339 }
340
341 for j in 0..n_alg {
343 if j < g_y.dim().1 && g_y[[i, j]].abs() > self.tol {
344 incidence[[n_diff_eqs + i, n_diff + j]] = true;
345 }
346 }
347 }
348
349 self.structure.incidence_matrix = Some(incidence);
350
351 let mut var_deps = Vec::with_capacity(n_vars);
353 let mut eq_deps = Vec::with_capacity(n_eqs);
354
355 for j in 0..n_vars {
357 let mut deps = Vec::new();
358 for i in 0..n_eqs {
359 if let Some(incidence) = &self.structure.incidence_matrix {
360 if incidence[[i, j]] {
361 deps.push(i);
362 }
363 }
364 }
365 var_deps.push(deps);
366 }
367
368 for i in 0..n_eqs {
370 let mut deps = Vec::new();
371 for j in 0..n_vars {
372 if let Some(incidence) = &self.structure.incidence_matrix {
373 if incidence[[i, j]] {
374 deps.push(j);
375 }
376 }
377 }
378 eq_deps.push(deps);
379 }
380
381 self.structure.variable_dependencies = Some(var_deps);
382 self.structure.equation_dependencies = Some(eq_deps);
383
384 Ok(())
385 }
386
387 pub fn reduce_index<FFunc, GFunc>(
392 &mut self,
393 t: F,
394 x: ArrayView1<F>,
395 y: ArrayView1<F>,
396 f: &FFunc,
397 g: &GFunc,
398 ) -> IntegrateResult<()>
399 where
400 FFunc: Fn(F, ArrayView1<F>, ArrayView1<F>) -> Array1<F>,
401 GFunc: Fn(F, ArrayView1<F>, ArrayView1<F>) -> Array1<F>,
402 {
403 self.structure.compute_index(t, x, y, f, g)?;
405
406 if self.structure.index == DAEIndex::Index1 {
408 return Ok(());
409 }
410
411 if self.structure.incidence_matrix.is_none() {
413 self.initialize_incidence_matrix(t, x, y, f, g)?;
414 }
415
416 let max_iterations = 10;
420 let mut current_iteration = 0;
421
422 while self.structure.index != DAEIndex::Index1 && current_iteration < max_iterations {
423 current_iteration += 1;
424
425 let singular_subsets = self.find_singular_subsets()?;
427
428 if singular_subsets.is_empty() {
429 break;
432 }
433
434 let differentiated_equations =
436 self.differentiate_singular_equations(&singular_subsets)?;
437
438 self.update_structure_after_differentiation(&differentiated_equations)?;
440
441 self.structure.compute_index(t, x, y, f, g)?;
443
444 if self.structure.diff_index > 0 {
446 self.structure.diff_index -= 1;
447 }
448
449 self.structure.index = match self.structure.diff_index {
451 0 | 1 => DAEIndex::Index1,
452 2 => DAEIndex::Index2,
453 3 => DAEIndex::Index3,
454 _ => DAEIndex::HigherIndex,
455 };
456 }
457
458 if self.structure.index != DAEIndex::Index1 {
459 return Err(IntegrateError::ConvergenceError(format!(
460 "Failed to reduce DAE to index-1 after {max_iterations} iterations"
461 )));
462 }
463
464 Ok(())
465 }
466
467 fn find_singular_subsets(&self) -> IntegrateResult<Vec<Vec<usize>>> {
469 let mut singular_subsets = Vec::new();
470
471 if let Some(ref incidence) = self.structure.incidence_matrix {
472 let n_eqs = incidence.nrows();
473 let n_vars = incidence.ncols();
474
475 for subset_size in 2..=std::cmp::min(n_eqs, 6) {
480 let equation_combinations = generate_combinations(n_eqs, subset_size);
482
483 for eq_subset in equation_combinations {
484 let mut involved_vars = std::collections::HashSet::new();
486
487 for &eq_idx in &eq_subset {
488 for var_idx in 0..n_vars {
489 if incidence[[eq_idx, var_idx]] {
490 involved_vars.insert(var_idx);
491 }
492 }
493 }
494
495 if eq_subset.len() > involved_vars.len() {
497 singular_subsets.push(eq_subset);
498 }
499 }
500 }
501 }
502
503 Ok(singular_subsets)
504 }
505
506 fn differentiate_singular_equations(
508 &self,
509 singular_subsets: &[Vec<usize>],
510 ) -> IntegrateResult<Vec<usize>> {
511 let mut equations_to_differentiate = Vec::new();
512
513 for subset in singular_subsets {
515 if !subset.is_empty() {
516 equations_to_differentiate.push(subset[0]);
519 }
520 }
521
522 Ok(equations_to_differentiate)
523 }
524
525 fn update_structure_after_differentiation(
527 &mut self,
528 differentiated_equations: &[usize],
529 ) -> IntegrateResult<()> {
530 if differentiated_equations.is_empty() {
531 return Ok(());
532 }
533
534 let num_new_variables = differentiated_equations.len();
536 let old_n_diff = self.structure.n_differential;
537 let old_n_alg = self.structure.n_algebraic;
538 let old_n_vars = old_n_diff + old_n_alg;
539 let old_n_eqs = self.structure.n_diff_eqs + self.structure.n_alg_eqs;
540
541 self.structure.n_differential += num_new_variables;
543 self.structure.n_diff_eqs += num_new_variables;
544
545 let new_n_vars = self.structure.n_differential + self.structure.n_algebraic;
546 let new_n_eqs = self.structure.n_diff_eqs + self.structure.n_alg_eqs;
547
548 if let Some(ref old_incidence) = self.structure.incidence_matrix.clone() {
550 let mut new_incidence = Array2::<bool>::from_elem((new_n_eqs, new_n_vars), false);
551
552 for i in 0..old_n_eqs {
554 for j in 0..old_n_vars {
555 new_incidence[[i, j]] = old_incidence[[i, j]];
556 }
557 }
558
559 for (new_eq_idx, &orig_eq_idx) in differentiated_equations.iter().enumerate() {
561 let new_eq_row = old_n_eqs + new_eq_idx;
562
563 let new_var_col = old_n_vars + new_eq_idx;
565 new_incidence[[new_eq_row, new_var_col]] = true;
566
567 for j in 0..old_n_vars {
569 if orig_eq_idx < old_n_eqs && old_incidence[[orig_eq_idx, j]] {
570 new_incidence[[new_eq_row, j]] = true;
571 }
572 }
573
574 for j in 0..old_n_diff {
577 if orig_eq_idx < old_n_eqs && old_incidence[[orig_eq_idx, j]] {
578 new_incidence[[new_eq_row, j]] = true;
580 }
581 }
582 }
583
584 self.structure.incidence_matrix = Some(new_incidence);
585 }
586
587 self.update_dependencies_after_structure_change()?;
589
590 if self.structure.diff_index > 1 {
592 self.structure.diff_index -= 1;
593 }
594
595 self.structure.constraint_jacobian = None;
597 self.structure.derivative_jacobian = None;
598
599 Ok(())
600 }
601
602 fn update_dependencies_after_structure_change(&mut self) -> IntegrateResult<()> {
604 if let Some(ref incidence) = self.structure.incidence_matrix {
605 let n_eqs = incidence.nrows();
606 let n_vars = incidence.ncols();
607
608 let mut var_deps = Vec::new();
610 for j in 0..n_vars {
611 let mut deps = Vec::new();
612 for i in 0..n_eqs {
613 if incidence[[i, j]] {
614 deps.push(i);
615 }
616 }
617 var_deps.push(deps);
618 }
619
620 let mut eq_deps = Vec::new();
622 for i in 0..n_eqs {
623 let mut deps = Vec::new();
624 for j in 0..n_vars {
625 if incidence[[i, j]] {
626 deps.push(j);
627 }
628 }
629 eq_deps.push(deps);
630 }
631
632 self.structure.variable_dependencies = Some(var_deps);
633 self.structure.equation_dependencies = Some(eq_deps);
634 }
635
636 Ok(())
637 }
638}
639
640#[derive(Debug, Clone)]
642struct DummySelection<
643 F: Float
644 + FromPrimitive
645 + Debug
646 + ScalarOperand
647 + std::ops::AddAssign
648 + std::ops::SubAssign
649 + std::ops::MulAssign
650 + std::ops::DivAssign
651 + Display
652 + LowerExp
653 + std::iter::Sum
654 + crate::common::IntegrateFloat,
655> {
656 dummy_vars: Vec<usize>,
658 dummy_eqs: Vec<usize>,
660 q: Array2<F>,
662 r: Array2<F>,
664}
665
666pub struct DummyDerivativeReducer<
671 F: Float
672 + FromPrimitive
673 + Debug
674 + ScalarOperand
675 + std::ops::AddAssign
676 + std::ops::SubAssign
677 + std::ops::MulAssign
678 + std::ops::DivAssign
679 + Display
680 + LowerExp
681 + std::iter::Sum
682 + crate::common::IntegrateFloat,
683> {
684 pub structure: DAEStructure<F>,
686
687 pub dummy_variables: Vec<usize>,
689
690 pub dummy_equations: Vec<usize>,
692}
693
694impl<
695 F: Float
696 + FromPrimitive
697 + Debug
698 + ScalarOperand
699 + std::ops::AddAssign
700 + std::ops::SubAssign
701 + std::ops::MulAssign
702 + std::ops::DivAssign
703 + Display
704 + LowerExp
705 + std::iter::Sum
706 + crate::common::IntegrateFloat,
707 > DummyDerivativeReducer<F>
708{
709 pub fn new(structure: DAEStructure<F>) -> Self {
711 DummyDerivativeReducer {
712 structure,
713 dummy_variables: Vec::new(),
714 dummy_equations: Vec::new(),
715 }
716 }
717
718 pub fn reduce_index<FFunc, GFunc>(
723 &mut self,
724 t: F,
725 x: ArrayView1<F>,
726 y: ArrayView1<F>,
727 f: &FFunc,
728 g: &GFunc,
729 ) -> IntegrateResult<()>
730 where
731 FFunc: Fn(F, ArrayView1<F>, ArrayView1<F>) -> Array1<F>,
732 GFunc: Fn(F, ArrayView1<F>, ArrayView1<F>) -> Array1<F>,
733 {
734 self.structure.compute_index(t, x, y, f, g)?;
736
737 if self.structure.index == DAEIndex::Index1 {
739 return Ok(());
740 }
741
742 let constraint_derivative = self.differentiate_constraints(t, x, y, g)?;
744
745 let extended_jacobian =
747 self.build_extended_jacobian(t, x, y, f, g, &constraint_derivative)?;
748
749 let dummy_selection = self.select_dummy_variables(&extended_jacobian)?;
751
752 self.dummy_variables = dummy_selection.dummy_vars;
754 self.dummy_equations = dummy_selection.dummy_eqs;
755
756 self.update_structure_with_dummies()?;
758
759 let final_index = self.verify_reduced_index(t, x, y, f, g)?;
761
762 if final_index != DAEIndex::Index1 {
763 return Err(IntegrateError::ComputationError(format!(
764 "Dummy derivative method failed to reduce to index-1. Final index: {final_index:?}"
765 )));
766 }
767
768 self.structure.index = DAEIndex::Index1;
769 self.structure.diff_index = 1;
770
771 Ok(())
772 }
773
774 fn differentiate_constraints<GFunc>(
776 &self,
777 t: F,
778 x: ArrayView1<F>,
779 y: ArrayView1<F>,
780 g: &GFunc,
781 ) -> IntegrateResult<Array1<F>>
782 where
783 GFunc: Fn(F, ArrayView1<F>, ArrayView1<F>) -> Array1<F>,
784 {
785 let h = F::from_f64(1e-8).unwrap();
786
787 let g0 = g(t, x, y);
789
790 let g_plus_t = g(t + h, x, y);
792
793 let dg_dt = (&g_plus_t - &g0) / h;
795
796 Ok(dg_dt)
797 }
798
799 fn build_extended_jacobian<FFunc, GFunc>(
801 &self,
802 t: F,
803 x: ArrayView1<F>,
804 y: ArrayView1<F>,
805 f: &FFunc,
806 g: &GFunc,
807 constraint_derivative: &Array1<F>,
808 ) -> IntegrateResult<Array2<F>>
809 where
810 FFunc: Fn(F, ArrayView1<F>, ArrayView1<F>) -> Array1<F>,
811 GFunc: Fn(F, ArrayView1<F>, ArrayView1<F>) -> Array1<F>,
812 {
813 let n_diff = self.structure.n_differential;
814 let n_alg = self.structure.n_algebraic;
815 let n_constraints = constraint_derivative.len();
816
817 let extended_size = n_diff + n_alg + n_constraints;
819 let mut extended_jacobian = Array2::<F>::zeros((extended_size, extended_size));
820
821 let f_x = compute_jacobian_for_variables(f, t, x, y, 0, n_diff)?;
823 let f_y = compute_jacobian_for_variables(f, t, x, y, n_diff, n_alg)?;
824 let g_x = compute_jacobian_for_variables(g, t, x, y, 0, n_diff)?;
825 let g_y = compute_jacobian_for_variables(g, t, x, y, n_diff, n_alg)?;
826
827 for i in 0..n_diff {
835 for j in 0..n_diff {
836 if j < f_x.dim().1 {
837 extended_jacobian[[i, j]] = f_x[[i, j]];
838 }
839 }
840 for j in 0..n_alg {
841 if j < f_y.dim().1 {
842 extended_jacobian[[i, n_diff + j]] = f_y[[i, j]];
843 }
844 }
845 }
846
847 for i in 0..n_alg {
849 for j in 0..n_diff {
850 if j < g_x.dim().1 && i < g_x.dim().0 {
851 extended_jacobian[[n_diff + i, j]] = g_x[[i, j]];
852 }
853 }
854 for j in 0..n_alg {
855 if j < g_y.dim().1 && i < g_y.dim().0 {
856 extended_jacobian[[n_diff + i, n_diff + j]] = g_y[[i, j]];
857 }
858 }
859 }
860
861 for i in 0..n_constraints {
864 for j in 0..n_diff {
866 if j < g_x.dim().1 && i < g_x.dim().0 {
867 extended_jacobian[[n_diff + n_alg + i, j]] = g_x[[i, j]];
868 }
869 }
870 for j in 0..n_alg {
871 if j < g_y.dim().1 && i < g_y.dim().0 {
872 extended_jacobian[[n_diff + n_alg + i, n_diff + j]] = g_y[[i, j]];
873 }
874 }
875 if i < extended_size - n_diff - n_alg {
877 extended_jacobian[[n_diff + n_alg + i, n_diff + n_alg + i]] = F::one();
878 }
879 }
880
881 Ok(extended_jacobian)
882 }
883
884 fn select_dummy_variables(
886 &self,
887 extended_jacobian: &Array2<F>,
888 ) -> IntegrateResult<DummySelection<F>> {
889 let n_diff = self.structure.n_differential;
890 let n_alg = self.structure.n_algebraic;
891 let total_vars = n_diff + n_alg;
892
893 let (q, r, pivots) = self.qr_decomposition_with_pivoting(extended_jacobian)?;
895
896 let rank = Self::compute_matrix_rank(&r)?;
898
899 let n_dummy = total_vars.saturating_sub(rank);
901
902 let mut dummy_vars = Vec::new();
904 let mut dummy_eqs = Vec::new();
905
906 for i in rank..total_vars.min(pivots.len()) {
908 if pivots[i] < total_vars {
909 dummy_vars.push(pivots[i]);
910 dummy_eqs.push(n_diff + n_alg + dummy_eqs.len());
911 }
912 }
913
914 while dummy_vars.len() < n_dummy && dummy_vars.len() < total_vars {
916 for i in 0..total_vars {
917 if !dummy_vars.contains(&i) && dummy_vars.len() < n_dummy {
918 dummy_vars.push(i);
919 dummy_eqs.push(n_diff + n_alg + dummy_eqs.len());
920 }
921 }
922 }
923
924 Ok(DummySelection {
925 dummy_vars,
926 dummy_eqs,
927 q,
928 r,
929 })
930 }
931
932 fn qr_decomposition_with_pivoting(
934 &self,
935 matrix: &Array2<F>,
936 ) -> IntegrateResult<(Array2<F>, Array2<F>, Vec<usize>)> {
937 let (m, n) = matrix.dim();
938 let mut a = matrix.clone();
939 let q = Array2::<F>::eye(m);
940 let mut pivots: Vec<usize> = (0..n).collect();
941
942 for k in 0..std::cmp::min(m, n) {
944 let mut max_norm = F::zero();
946 let mut max_col = k;
947
948 for j in k..n {
949 let col_norm = (k..m)
950 .map(|i| a[[i, j]].powi(2))
951 .fold(F::zero(), |acc, x| acc + x)
952 .sqrt();
953 if col_norm > max_norm {
954 max_norm = col_norm;
955 max_col = j;
956 }
957 }
958
959 if max_col != k {
961 for i in 0..m {
962 let temp = a[[i, k]];
963 a[[i, k]] = a[[i, max_col]];
964 a[[i, max_col]] = temp;
965 }
966 pivots.swap(k, max_col);
967 }
968
969 if max_norm > F::from_f64(1e-12).unwrap() {
971 for i in k..m {
973 a[[i, k]] /= max_norm;
974 }
975
976 for j in (k + 1)..n {
978 let dot_product = (k..m)
979 .map(|i| a[[i, k]] * a[[i, j]])
980 .fold(F::zero(), |acc, x| acc + x);
981 for i in k..m {
982 a[[i, j]] = a[[i, j]] - dot_product * a[[i, k]];
983 }
984 }
985 }
986 }
987
988 let r = a.clone();
990
991 Ok((q, r, pivots))
992 }
993
994 fn compute_matrix_rank(r: &Array2<F>) -> IntegrateResult<usize> {
996 let tolerance = F::from_f64(1e-10).unwrap();
997 let min_dim = std::cmp::min(r.nrows(), r.ncols());
998
999 let mut rank = 0;
1000 for i in 0..min_dim {
1001 if r[[i, i]].abs() > tolerance {
1002 rank += 1;
1003 }
1004 }
1005
1006 Ok(rank)
1007 }
1008
1009 fn update_structure_with_dummies(&mut self) -> IntegrateResult<()> {
1011 let n_dummy = self.dummy_variables.len();
1012
1013 self.structure.n_differential += n_dummy;
1015 self.structure.n_diff_eqs += n_dummy;
1016
1017 self.structure.constraint_jacobian = None;
1019 self.structure.derivative_jacobian = None;
1020 self.structure.incidence_matrix = None;
1021
1022 Ok(())
1023 }
1024
1025 fn verify_reduced_index<FFunc, GFunc>(
1027 &self,
1028 t: F,
1029 x: ArrayView1<F>,
1030 y: ArrayView1<F>,
1031 _f: &FFunc,
1032 g: &GFunc,
1033 ) -> IntegrateResult<DAEIndex>
1034 where
1035 FFunc: Fn(F, ArrayView1<F>, ArrayView1<F>) -> Array1<F>,
1036 GFunc: Fn(F, ArrayView1<F>, ArrayView1<F>) -> Array1<F>,
1037 {
1038 let x_slice: Vec<F> = x.to_vec();
1040 let y_slice: Vec<F> = y.to_vec();
1041
1042 let g_y = compute_constraint_jacobian(
1043 &|t, x, y| g(t, ArrayView1::from(x), ArrayView1::from(y)).to_vec(),
1044 t,
1045 &x_slice,
1046 &y_slice,
1047 );
1048
1049 let singular = is_singular_matrix(g_y.view());
1051
1052 if !singular {
1053 Ok(DAEIndex::Index1)
1054 } else {
1055 Ok(DAEIndex::Index2) }
1057 }
1058}
1059
1060pub struct ProjectionMethod<
1065 F: Float
1066 + FromPrimitive
1067 + Debug
1068 + ScalarOperand
1069 + std::ops::AddAssign
1070 + std::ops::SubAssign
1071 + std::ops::MulAssign
1072 + std::ops::DivAssign
1073 + Display
1074 + LowerExp
1075 + std::iter::Sum
1076 + crate::common::IntegrateFloat,
1077> {
1078 pub structure: DAEStructure<F>,
1080
1081 pub project_after_step: bool,
1083
1084 pub consistent_initialization: bool,
1086
1087 pub constraint_tol: F,
1089}
1090
1091impl<
1092 F: Float
1093 + FromPrimitive
1094 + Debug
1095 + ScalarOperand
1096 + std::ops::AddAssign
1097 + std::ops::SubAssign
1098 + std::ops::MulAssign
1099 + std::ops::DivAssign
1100 + Display
1101 + LowerExp
1102 + std::iter::Sum
1103 + crate::common::IntegrateFloat,
1104 > ProjectionMethod<F>
1105{
1106 pub fn new(structure: DAEStructure<F>) -> Self {
1108 ProjectionMethod {
1109 structure,
1110 project_after_step: true,
1111 consistent_initialization: true,
1112 constraint_tol: F::from_f64(1e-8).unwrap(),
1113 }
1114 }
1115
1116 pub fn project_solution<GFunc>(
1120 &self,
1121 t: F,
1122 x: &mut Array1<F>,
1123 y: &mut Array1<F>,
1124 g: &GFunc,
1125 ) -> IntegrateResult<()>
1126 where
1127 GFunc: Fn(F, ArrayView1<F>, ArrayView1<F>) -> Array1<F>,
1128 {
1129 let g_val = g(t, x.view(), y.view());
1131
1132 let g_norm = g_val
1134 .iter()
1135 .fold(F::zero(), |acc, &val| acc + val.powi(2))
1136 .sqrt();
1137
1138 if g_norm <= self.constraint_tol {
1139 return Ok(());
1141 }
1142
1143 let x_slice: Vec<F> = x.to_vec();
1146 let y_slice: Vec<F> = y.to_vec();
1147 let g_y = compute_constraint_jacobian(
1148 &|t, x, y| g(t, ArrayView1::from(x), ArrayView1::from(y)).to_vec(),
1149 t,
1150 &x_slice,
1151 &y_slice,
1152 );
1153
1154 let neg_g = g_val.mapv(|val| -val);
1156
1157 let (delta_y, residual) = solve_constrained_least_squares(g_y.view(), neg_g.view())?;
1159
1160 *y = y.clone() + delta_y;
1162
1163 if residual > self.constraint_tol {
1165 return Err(IntegrateError::ComputationError(format!(
1166 "Failed to project solution onto constraint manifold. Residual: {residual}"
1167 )));
1168 }
1169
1170 Ok(())
1171 }
1172
1173 pub fn make_consistent<GFunc>(
1177 &self,
1178 t: F,
1179 x: &mut Array1<F>,
1180 y: &mut Array1<F>,
1181 g: &GFunc,
1182 ) -> IntegrateResult<()>
1183 where
1184 GFunc: Fn(F, ArrayView1<F>, ArrayView1<F>) -> Array1<F>,
1185 {
1186 let max_iter = 10;
1189 let mut iter = 0;
1190
1191 while iter < max_iter {
1192 let g_val = g(t, x.view(), y.view());
1194
1195 let g_norm = g_val
1197 .iter()
1198 .fold(F::zero(), |acc, &val| acc + val.powi(2))
1199 .sqrt();
1200
1201 if g_norm <= self.constraint_tol {
1202 return Ok(());
1204 }
1205
1206 self.project_solution(t, x, y, g)?;
1208
1209 iter += 1;
1210 }
1211
1212 let g_val = g(t, x.view(), y.view());
1214 let g_norm = g_val
1215 .iter()
1216 .fold(F::zero(), |acc, &val| acc + val.powi(2))
1217 .sqrt();
1218
1219 if g_norm > self.constraint_tol {
1220 return Err(IntegrateError::ComputationError(format!(
1221 "Failed to find consistent initial conditions after {max_iter} iterations. Residual: {g_norm}"
1222 )));
1223 }
1224
1225 Ok(())
1226 }
1227}
1228
1229#[allow(dead_code)]
1233fn compute_jacobian_for_variables<F, Func>(
1234 f: &Func,
1235 t: F,
1236 x: ArrayView1<F>,
1237 y: ArrayView1<F>,
1238 start_idx: usize,
1239 n_vars: usize,
1240) -> IntegrateResult<Array2<F>>
1241where
1242 F: Float
1243 + FromPrimitive
1244 + Debug
1245 + ScalarOperand
1246 + std::ops::AddAssign
1247 + std::ops::SubAssign
1248 + std::ops::MulAssign
1249 + std::ops::DivAssign
1250 + Display
1251 + LowerExp
1252 + std::iter::Sum
1253 + crate::common::IntegrateFloat,
1254 Func: Fn(F, ArrayView1<F>, ArrayView1<F>) -> Array1<F>,
1255{
1256 let f_val = f(t, x, y);
1258 let n_eqs = f_val.len();
1259
1260 let mut jac = Array2::<F>::zeros((n_eqs, n_vars));
1262
1263 let epsilon = F::from_f64(1e-8).unwrap();
1265
1266 if start_idx < x.len() {
1267 let end_idx = (start_idx + n_vars).min(x.len());
1269 let n_x_vars = end_idx - start_idx;
1270
1271 let mut x_perturbed = x.to_owned();
1272
1273 for j in 0..n_x_vars {
1274 let _idx = start_idx + j;
1275 let h = epsilon.max(x[_idx].abs() * epsilon);
1276
1277 x_perturbed[_idx] = x[_idx] + h;
1278 let f_plus = f(t, x_perturbed.view(), y);
1279 x_perturbed[_idx] = x[_idx];
1280
1281 let df = (f_plus - f_val.view()) / h;
1282
1283 for i in 0..n_eqs {
1284 jac[[i, j]] = df[i];
1285 }
1286 }
1287 }
1288
1289 if start_idx + n_vars > x.len() {
1290 let n_x_vars = x.len().saturating_sub(start_idx);
1292 let n_y_vars = n_vars - n_x_vars;
1293 let y_start = start_idx.saturating_sub(x.len());
1294 let y_end = (y_start + n_y_vars).min(y.len());
1295
1296 let mut y_perturbed = y.to_owned();
1297
1298 for j in 0..y_end.saturating_sub(y_start) {
1299 let _idx = y_start + j;
1300 let h = epsilon.max(y[_idx].abs() * epsilon);
1301
1302 y_perturbed[_idx] = y[_idx] + h;
1303 let f_plus = f(t, x, y_perturbed.view());
1304 y_perturbed[_idx] = y[_idx];
1305
1306 let df = (f_plus - f_val.view()) / h;
1307
1308 for i in 0..n_eqs {
1309 jac[[i, n_x_vars + j]] = df[i];
1310 }
1311 }
1312 }
1313
1314 Ok(jac)
1315}
1316
1317#[allow(dead_code)]
1321fn solve_constrained_least_squares<F>(
1322 j: ArrayView2<F>,
1323 b: ArrayView1<F>,
1324) -> IntegrateResult<(Array1<F>, F)>
1325where
1326 F: Float
1327 + FromPrimitive
1328 + Debug
1329 + ScalarOperand
1330 + std::ops::AddAssign
1331 + std::ops::SubAssign
1332 + std::ops::MulAssign
1333 + std::ops::DivAssign
1334 + Display
1335 + LowerExp
1336 + std::iter::Sum
1337 + crate::common::IntegrateFloat,
1338{
1339 use crate::dae::utils::linear_solvers::solve_linear_system;
1340
1341 let j_jt = j.dot(&j.t());
1347
1348 let lambda = match solve_linear_system(&j_jt.view(), &b.view()) {
1350 Ok(sol) => sol,
1351 Err(e) => {
1352 return Err(IntegrateError::ComputationError(format!(
1353 "Failed to solve least squares system: {e}"
1354 )))
1355 }
1356 };
1357
1358 let delta_y = j.t().dot(&lambda);
1360
1361 let residual = (&j.dot(&delta_y) - &b)
1363 .iter()
1364 .fold(F::zero(), |acc, &val| acc + val.powi(2))
1365 .sqrt();
1366
1367 Ok((delta_y, residual))
1368}
1369
1370#[allow(dead_code)]
1372fn generate_combinations(n: usize, k: usize) -> Vec<Vec<usize>> {
1373 if k == 0 || k > n {
1374 return vec![];
1375 }
1376
1377 let mut combinations = Vec::new();
1378 let mut current = Vec::new();
1379
1380 generate_combinations_recursive(0, n, k, &mut current, &mut combinations);
1381
1382 combinations
1383}
1384
1385#[allow(dead_code)]
1387fn generate_combinations_recursive(
1388 start: usize,
1389 n: usize,
1390 k: usize,
1391 current: &mut Vec<usize>,
1392 combinations: &mut Vec<Vec<usize>>,
1393) {
1394 if current.len() == k {
1395 combinations.push(current.clone());
1396 return;
1397 }
1398
1399 for i in start..n {
1400 if current.len() + (n - i) >= k {
1401 current.push(i);
1402 generate_combinations_recursive(i + 1, n, k, current, combinations);
1403 current.pop();
1404 }
1405 }
1406}