1use crate::preconditioner::{IdentityPreconditioner, Preconditioner};
11use crate::sparse::SparseMatrix;
12use crate::Scalar;
13use faer::{ComplexField, Conjugate, Entity, SimpleEntity};
14use numra_core::LinalgError;
15
16pub struct IterativeOptions<S: Scalar> {
18 pub tol: S,
20 pub max_iter: usize,
22 pub restart: usize,
24}
25
26impl<S: Scalar> Default for IterativeOptions<S> {
27 fn default() -> Self {
28 Self {
29 tol: S::from_f64(1e-8),
30 max_iter: 1000,
31 restart: 30,
32 }
33 }
34}
35
36impl<S: Scalar> IterativeOptions<S> {
37 pub fn tol(mut self, tol: S) -> Self {
39 self.tol = tol;
40 self
41 }
42
43 pub fn max_iter(mut self, max_iter: usize) -> Self {
45 self.max_iter = max_iter;
46 self
47 }
48
49 pub fn restart(mut self, restart: usize) -> Self {
51 self.restart = restart;
52 self
53 }
54}
55
56pub struct IterativeResult<S: Scalar> {
58 pub x: Vec<S>,
60 pub residual_norm: S,
62 pub iterations: usize,
64 pub converged: bool,
66}
67
68fn dot<S: Scalar>(a: &[S], b: &[S]) -> S {
73 a.iter()
74 .zip(b.iter())
75 .map(|(&ai, &bi)| ai * bi)
76 .fold(S::ZERO, |acc, x| acc + x)
77}
78
79fn norm<S: Scalar>(v: &[S]) -> S {
80 S::from_f64(dot(v, v).to_f64().sqrt())
81}
82
83fn axpy<S: Scalar>(a: S, x: &[S], y: &mut [S]) {
85 for (yi, &xi) in y.iter_mut().zip(x.iter()) {
86 *yi += a * xi;
87 }
88}
89
90pub fn cg<S: Scalar + SimpleEntity + Conjugate<Canonical = S> + ComplexField + Entity>(
100 a: &SparseMatrix<S>,
101 b: &[S],
102 x0: Option<&[S]>,
103 opts: &IterativeOptions<S>,
104) -> Result<IterativeResult<S>, LinalgError> {
105 pcg(a, b, &IdentityPreconditioner, x0, opts)
106}
107
108pub fn pcg<S: Scalar + SimpleEntity + Conjugate<Canonical = S> + ComplexField + Entity>(
113 a: &SparseMatrix<S>,
114 b: &[S],
115 precond: &dyn Preconditioner<S>,
116 x0: Option<&[S]>,
117 opts: &IterativeOptions<S>,
118) -> Result<IterativeResult<S>, LinalgError> {
119 let n = b.len();
120 if a.nrows() != n || a.ncols() != n {
121 return Err(LinalgError::DimensionMismatch {
122 expected: (n, n),
123 actual: (a.nrows(), a.ncols()),
124 });
125 }
126
127 let mut x = match x0 {
128 Some(x0) => x0.to_vec(),
129 None => vec![S::ZERO; n],
130 };
131
132 let ax = a.mul_vec(&x)?;
134 let mut r: Vec<S> = b.iter().zip(ax.iter()).map(|(&bi, &ai)| bi - ai).collect();
135
136 let b_norm = norm(b);
137 if b_norm.to_f64() < S::EPSILON.to_f64() {
138 return Ok(IterativeResult {
139 x: vec![S::ZERO; n],
140 residual_norm: S::ZERO,
141 iterations: 0,
142 converged: true,
143 });
144 }
145
146 let mut z = vec![S::ZERO; n];
148 precond.apply(&r, &mut z);
149
150 let mut p = z.clone();
151 let mut rz = dot(&r, &z);
152
153 for iter in 0..opts.max_iter {
154 let r_norm = norm(&r);
155 if r_norm.to_f64() / b_norm.to_f64() < opts.tol.to_f64() {
156 return Ok(IterativeResult {
157 x,
158 residual_norm: r_norm,
159 iterations: iter,
160 converged: true,
161 });
162 }
163
164 let ap = a.mul_vec(&p)?;
166 let pap = dot(&p, &ap);
167
168 if pap.to_f64().abs() < S::EPSILON.to_f64() * 100.0 {
169 return Ok(IterativeResult {
170 x,
171 residual_norm: r_norm,
172 iterations: iter,
173 converged: false,
174 });
175 }
176
177 let alpha = S::from_f64(rz.to_f64() / pap.to_f64());
178
179 axpy(alpha, &p, &mut x);
181 axpy(S::ZERO - alpha, &ap, &mut r);
183
184 precond.apply(&r, &mut z);
186 let rz_new = dot(&r, &z);
187
188 let beta = S::from_f64(rz_new.to_f64() / rz.to_f64());
189 for i in 0..n {
191 p[i] = z[i] + beta * p[i];
192 }
193 rz = rz_new;
194 }
195
196 let r_norm = norm(&r);
197 Ok(IterativeResult {
198 x,
199 residual_norm: r_norm,
200 iterations: opts.max_iter,
201 converged: r_norm.to_f64() / b_norm.to_f64() < opts.tol.to_f64(),
202 })
203}
204
205pub fn gmres<S: Scalar + SimpleEntity + Conjugate<Canonical = S> + ComplexField + Entity>(
215 a: &SparseMatrix<S>,
216 b: &[S],
217 x0: Option<&[S]>,
218 opts: &IterativeOptions<S>,
219) -> Result<IterativeResult<S>, LinalgError> {
220 let n = b.len();
221 if a.nrows() != n || a.ncols() != n {
222 return Err(LinalgError::DimensionMismatch {
223 expected: (n, n),
224 actual: (a.nrows(), a.ncols()),
225 });
226 }
227
228 let b_norm = norm(b);
229 if b_norm.to_f64() < S::EPSILON.to_f64() {
230 return Ok(IterativeResult {
231 x: vec![S::ZERO; n],
232 residual_norm: S::ZERO,
233 iterations: 0,
234 converged: true,
235 });
236 }
237
238 let mut x = match x0 {
239 Some(x0) => x0.to_vec(),
240 None => vec![S::ZERO; n],
241 };
242
243 let m = opts.restart.min(n); let mut total_iter = 0;
245
246 while total_iter < opts.max_iter {
247 let ax = a.mul_vec(&x)?;
249 let r: Vec<S> = b.iter().zip(ax.iter()).map(|(&bi, &ai)| bi - ai).collect();
250 let r_norm = norm(&r);
251
252 if r_norm.to_f64() / b_norm.to_f64() < opts.tol.to_f64() {
253 return Ok(IterativeResult {
254 x,
255 residual_norm: r_norm,
256 iterations: total_iter,
257 converged: true,
258 });
259 }
260
261 let mut v: Vec<Vec<S>> = Vec::with_capacity(m + 1); let inv_r_norm = S::from_f64(1.0 / r_norm.to_f64());
264 v.push(r.iter().map(|&ri| ri * inv_r_norm).collect());
265
266 let mut h = vec![vec![S::ZERO; m + 1]; m];
268
269 let mut cs = vec![S::ZERO; m]; let mut sn = vec![S::ZERO; m]; let mut g = vec![S::ZERO; m + 1]; g[0] = r_norm;
274
275 let mut k = 0; for j in 0..m {
277 if total_iter >= opts.max_iter {
278 break;
279 }
280
281 let w_tmp = a.mul_vec(&v[j])?;
283 let mut w = w_tmp;
284
285 for i in 0..=j {
287 h[j][i] = dot(&w, &v[i]);
288 axpy(S::ZERO - h[j][i], &v[i], &mut w);
289 }
290 h[j][j + 1] = norm(&w);
291
292 if h[j][j + 1].to_f64().abs() < S::EPSILON.to_f64() * 100.0 {
294 k = j + 1;
295 apply_givens_to_col(&mut h[j], &cs, &sn, j);
297 let (c, s) = givens_rotation(h[j][j], h[j][j + 1]);
299 cs[j] = c;
300 sn[j] = s;
301 apply_givens_pair(&mut h[j], j, c, s);
302 let gj = g[j];
303 let gjp1 = g[j + 1];
304 g[j] = c * gj + s * gjp1;
305 g[j + 1] = S::ZERO - s * gj + c * gjp1;
306 total_iter += 1;
307 break;
308 }
309
310 let inv_hjj1 = S::from_f64(1.0 / h[j][j + 1].to_f64());
312 v.push(w.iter().map(|&wi| wi * inv_hjj1).collect());
313
314 apply_givens_to_col(&mut h[j], &cs, &sn, j);
316
317 let (c, s) = givens_rotation(h[j][j], h[j][j + 1]);
319 cs[j] = c;
320 sn[j] = s;
321 apply_givens_pair(&mut h[j], j, c, s);
322
323 let gj = g[j];
325 let gjp1 = g[j + 1];
326 g[j] = c * gj + s * gjp1;
327 g[j + 1] = S::ZERO - s * gj + c * gjp1;
328
329 total_iter += 1;
330 k = j + 1;
331
332 if g[j + 1].abs().to_f64() / b_norm.to_f64() < opts.tol.to_f64() {
334 break;
335 }
336 }
337
338 if k == 0 {
339 break;
340 }
341
342 let mut y = vec![S::ZERO; k];
344 for i in (0..k).rev() {
345 let mut sum = g[i];
346 for j in (i + 1)..k {
347 sum -= h[j][i] * y[j];
348 }
349 if h[i][i].to_f64().abs() < S::EPSILON.to_f64() * 100.0 {
350 break;
351 }
352 y[i] = S::from_f64(sum.to_f64() / h[i][i].to_f64());
353 }
354
355 for j in 0..k {
357 axpy(y[j], &v[j], &mut x);
358 }
359
360 let ax = a.mul_vec(&x)?;
362 let r_final: Vec<S> = b.iter().zip(ax.iter()).map(|(&bi, &ai)| bi - ai).collect();
363 let r_final_norm = norm(&r_final);
364 if r_final_norm.to_f64() / b_norm.to_f64() < opts.tol.to_f64() {
365 return Ok(IterativeResult {
366 x,
367 residual_norm: r_final_norm,
368 iterations: total_iter,
369 converged: true,
370 });
371 }
372 }
373
374 let ax = a.mul_vec(&x)?;
376 let r: Vec<S> = b.iter().zip(ax.iter()).map(|(&bi, &ai)| bi - ai).collect();
377 let r_norm = norm(&r);
378 Ok(IterativeResult {
379 x,
380 residual_norm: r_norm,
381 iterations: total_iter,
382 converged: r_norm.to_f64() / b_norm.to_f64() < opts.tol.to_f64(),
383 })
384}
385
386fn givens_rotation<S: Scalar>(a: S, b: S) -> (S, S) {
388 let af = a.to_f64();
389 let bf = b.to_f64();
390 if bf.abs() < f64::EPSILON * 100.0 {
391 return (S::ONE, S::ZERO);
392 }
393 if af.abs() < f64::EPSILON * 100.0 {
394 return (S::ZERO, S::ONE);
395 }
396 let r = (af * af + bf * bf).sqrt();
397 (S::from_f64(af / r), S::from_f64(bf / r))
398}
399
400fn apply_givens_to_col<S: Scalar>(col: &mut [S], cs: &[S], sn: &[S], j: usize) {
402 for i in 0..j {
403 let a = col[i];
404 let b = col[i + 1];
405 col[i] = cs[i] * a + sn[i] * b;
406 col[i + 1] = (S::ZERO - sn[i]) * a + cs[i] * b;
407 }
408}
409
410fn apply_givens_pair<S: Scalar>(col: &mut [S], i: usize, c: S, s: S) {
412 let a = col[i];
413 let b = col[i + 1];
414 col[i] = c * a + s * b;
415 col[i + 1] = (S::ZERO - s) * a + c * b;
416}
417
418pub fn bicgstab<S: Scalar + SimpleEntity + Conjugate<Canonical = S> + ComplexField + Entity>(
427 a: &SparseMatrix<S>,
428 b: &[S],
429 x0: Option<&[S]>,
430 opts: &IterativeOptions<S>,
431) -> Result<IterativeResult<S>, LinalgError> {
432 let n = b.len();
433 if a.nrows() != n || a.ncols() != n {
434 return Err(LinalgError::DimensionMismatch {
435 expected: (n, n),
436 actual: (a.nrows(), a.ncols()),
437 });
438 }
439
440 let b_norm = norm(b);
441 if b_norm.to_f64() < S::EPSILON.to_f64() {
442 return Ok(IterativeResult {
443 x: vec![S::ZERO; n],
444 residual_norm: S::ZERO,
445 iterations: 0,
446 converged: true,
447 });
448 }
449
450 let mut x = match x0 {
451 Some(x0) => x0.to_vec(),
452 None => vec![S::ZERO; n],
453 };
454
455 let ax = a.mul_vec(&x)?;
457 let mut r: Vec<S> = b.iter().zip(ax.iter()).map(|(&bi, &ai)| bi - ai).collect();
458
459 let r0_hat = r.clone(); let mut rho = S::ONE;
461 let mut alpha = S::ONE;
462 let mut omega = S::ONE;
463
464 let mut p = vec![S::ZERO; n];
465 let mut v = vec![S::ZERO; n];
466
467 for iter in 0..opts.max_iter {
468 let r_norm = norm(&r);
469 if r_norm.to_f64() / b_norm.to_f64() < opts.tol.to_f64() {
470 return Ok(IterativeResult {
471 x,
472 residual_norm: r_norm,
473 iterations: iter,
474 converged: true,
475 });
476 }
477
478 let rho_new = dot(&r0_hat, &r);
479 if rho_new.to_f64().abs() < S::EPSILON.to_f64() * 100.0 {
480 return Ok(IterativeResult {
482 x,
483 residual_norm: r_norm,
484 iterations: iter,
485 converged: false,
486 });
487 }
488
489 let beta =
490 S::from_f64((rho_new.to_f64() / rho.to_f64()) * (alpha.to_f64() / omega.to_f64()));
491
492 for i in 0..n {
494 p[i] = r[i] + beta * (p[i] - omega * v[i]);
495 }
496
497 v = a.mul_vec(&p)?;
499
500 let r0v = dot(&r0_hat, &v);
501 if r0v.to_f64().abs() < S::EPSILON.to_f64() * 100.0 {
502 return Ok(IterativeResult {
503 x,
504 residual_norm: r_norm,
505 iterations: iter,
506 converged: false,
507 });
508 }
509 alpha = S::from_f64(rho_new.to_f64() / r0v.to_f64());
510
511 let s: Vec<S> = r
513 .iter()
514 .zip(v.iter())
515 .map(|(&ri, &vi)| ri - alpha * vi)
516 .collect();
517
518 let s_norm = norm(&s);
519 if s_norm.to_f64() / b_norm.to_f64() < opts.tol.to_f64() {
520 axpy(alpha, &p, &mut x);
522 return Ok(IterativeResult {
523 x,
524 residual_norm: s_norm,
525 iterations: iter + 1,
526 converged: true,
527 });
528 }
529
530 let t = a.mul_vec(&s)?;
532
533 let tt = dot(&t, &t);
534 if tt.to_f64().abs() < S::EPSILON.to_f64() * 100.0 {
535 axpy(alpha, &p, &mut x);
536 r = s;
537 } else {
538 omega = S::from_f64(dot(&t, &s).to_f64() / tt.to_f64());
539
540 axpy(alpha, &p, &mut x);
542 axpy(omega, &s, &mut x);
543
544 for i in 0..n {
546 r[i] = s[i] - omega * t[i];
547 }
548 }
549
550 rho = rho_new;
551 }
552
553 let r_norm = norm(&r);
554 Ok(IterativeResult {
555 x,
556 residual_norm: r_norm,
557 iterations: opts.max_iter,
558 converged: r_norm.to_f64() / b_norm.to_f64() < opts.tol.to_f64(),
559 })
560}
561
562pub fn minres<S: Scalar + SimpleEntity + Conjugate<Canonical = S> + ComplexField + Entity>(
575 a: &SparseMatrix<S>,
576 b: &[S],
577 x0: Option<&[S]>,
578 opts: &IterativeOptions<S>,
579) -> Result<IterativeResult<S>, LinalgError> {
580 let n = b.len();
581 if a.nrows() != n || a.ncols() != n {
582 return Err(LinalgError::DimensionMismatch {
583 expected: (n, n),
584 actual: (a.nrows(), a.ncols()),
585 });
586 }
587
588 let b_norm = norm(b);
589 if b_norm.to_f64() < S::EPSILON.to_f64() {
590 return Ok(IterativeResult {
591 x: vec![S::ZERO; n],
592 residual_norm: S::ZERO,
593 iterations: 0,
594 converged: true,
595 });
596 }
597
598 let mut x = match x0 {
599 Some(x0) => x0.to_vec(),
600 None => vec![S::ZERO; n],
601 };
602
603 let ax = a.mul_vec(&x)?;
605 let r: Vec<S> = b.iter().zip(ax.iter()).map(|(&bi, &ai)| bi - ai).collect();
606
607 let mut beta_curr = norm(&r);
608 if beta_curr.to_f64() / b_norm.to_f64() < opts.tol.to_f64() {
609 return Ok(IterativeResult {
610 x,
611 residual_norm: beta_curr,
612 iterations: 0,
613 converged: true,
614 });
615 }
616
617 let mut v_prev = vec![S::ZERO; n];
619 let mut v_curr: Vec<S> = r
620 .iter()
621 .map(|&ri| S::from_f64(ri.to_f64() / beta_curr.to_f64()))
622 .collect();
623
624 let mut c_old = S::ONE;
626 let mut s_old = S::ZERO;
627 let mut c_curr = S::ONE;
628 let mut s_curr = S::ZERO;
629
630 let mut w_curr = vec![S::ZERO; n];
632 let mut w_prev = vec![S::ZERO; n];
633
634 let mut phi_bar = beta_curr;
636
637 for iter in 0..opts.max_iter {
638 let av = a.mul_vec(&v_curr)?;
640 let alpha = dot(&v_curr, &av);
641
642 let mut v_next: Vec<S> = Vec::with_capacity(n);
644 for i in 0..n {
645 v_next.push(av[i] - alpha * v_curr[i] - beta_curr * v_prev[i]);
646 }
647 let beta_next = norm(&v_next);
648
649 if beta_next.to_f64().abs() > S::EPSILON.to_f64() * 100.0 {
651 let inv_bn = S::from_f64(1.0 / beta_next.to_f64());
652 for vi in &mut v_next {
653 *vi *= inv_bn;
654 }
655 }
656
657 let epsilon = s_old * beta_curr;
663 let delta_hat = c_old * beta_curr;
664
665 let delta_tilde = c_curr * delta_hat + s_curr * alpha;
667 let gamma_bar = (S::ZERO - s_curr) * delta_hat + c_curr * alpha;
668
669 let gb = gamma_bar.to_f64();
671 let bn = beta_next.to_f64();
672 let gamma = (gb * gb + bn * bn).sqrt();
673
674 let (c_new, s_new) = if gamma.abs() < f64::EPSILON * 100.0 {
675 (S::ONE, S::ZERO)
676 } else {
677 (S::from_f64(gb / gamma), S::from_f64(bn / gamma))
678 };
679 let gamma_s = S::from_f64(gamma);
680
681 let phi = c_new * phi_bar;
683 let phi_bar_new = (S::ZERO - s_new) * phi_bar;
684
685 if gamma.abs() < f64::EPSILON * 100.0 {
687 return Ok(IterativeResult {
689 x,
690 residual_norm: phi_bar.abs(),
691 iterations: iter + 1,
692 converged: false,
693 });
694 }
695
696 let mut w_next = Vec::with_capacity(n);
697 for i in 0..n {
698 w_next.push((v_curr[i] - delta_tilde * w_curr[i] - epsilon * w_prev[i]) / gamma_s);
699 }
700
701 axpy(phi, &w_next, &mut x);
703
704 phi_bar = phi_bar_new;
706 v_prev = v_curr;
707 v_curr = v_next;
708 w_prev = w_curr;
709 w_curr = w_next;
710 beta_curr = beta_next;
711 c_old = c_curr;
712 s_old = s_curr;
713 c_curr = c_new;
714 s_curr = s_new;
715
716 if phi_bar.abs().to_f64() / b_norm.to_f64() < opts.tol.to_f64() {
718 return Ok(IterativeResult {
719 x,
720 residual_norm: phi_bar.abs(),
721 iterations: iter + 1,
722 converged: true,
723 });
724 }
725 }
726
727 Ok(IterativeResult {
728 x,
729 residual_norm: phi_bar.abs(),
730 iterations: opts.max_iter,
731 converged: phi_bar.abs().to_f64() / b_norm.to_f64() < opts.tol.to_f64(),
732 })
733}
734
735#[cfg(test)]
740mod tests {
741 use super::*;
742 use crate::preconditioner::{Ilu0, Jacobi, Ssor};
743
744 fn laplacian_1d(n: usize) -> SparseMatrix<f64> {
746 let mut triplets = Vec::new();
747 for i in 0..n {
748 triplets.push((i, i, 2.0));
749 if i > 0 {
750 triplets.push((i, i - 1, -1.0));
751 }
752 if i < n - 1 {
753 triplets.push((i, i + 1, -1.0));
754 }
755 }
756 SparseMatrix::from_triplets(n, n, &triplets).unwrap()
757 }
758
759 fn convection_diffusion(n: usize, eps: f64) -> SparseMatrix<f64> {
761 let mut triplets = Vec::new();
762 for i in 0..n {
763 triplets.push((i, i, 2.0));
764 if i > 0 {
765 triplets.push((i, i - 1, -1.0 - eps));
766 }
767 if i < n - 1 {
768 triplets.push((i, i + 1, -1.0 + eps));
769 }
770 }
771 SparseMatrix::from_triplets(n, n, &triplets).unwrap()
772 }
773
774 #[test]
775 fn test_cg_laplacian() {
776 let n = 50;
777 let a = laplacian_1d(n);
778 let x_exact: Vec<f64> = (0..n).map(|i| (i + 1) as f64).collect();
779 let b = a.mul_vec(&x_exact).unwrap();
780
781 let opts = IterativeOptions::default().tol(1e-10);
782 let result = cg(&a, &b, None, &opts).unwrap();
783
784 assert!(
785 result.converged,
786 "CG did not converge after {} iters",
787 result.iterations
788 );
789 for i in 0..n {
790 assert!(
791 (result.x[i] - x_exact[i]).abs() < 1e-6,
792 "CG error at {}: {} vs {}",
793 i,
794 result.x[i],
795 x_exact[i]
796 );
797 }
798 }
799
800 #[test]
801 fn test_cg_with_initial_guess() {
802 let n = 20;
803 let a = laplacian_1d(n);
804 let x_exact: Vec<f64> = (0..n).map(|i| (i + 1) as f64).collect();
805 let b = a.mul_vec(&x_exact).unwrap();
806
807 let x0: Vec<f64> = (0..n).map(|i| i as f64 + 0.9).collect();
809 let opts = IterativeOptions::default().tol(1e-10);
810 let result = cg(&a, &b, Some(&x0), &opts).unwrap();
811
812 assert!(result.converged);
813 }
814
815 #[test]
816 fn test_pcg_jacobi() {
817 let n = 50;
818 let a = laplacian_1d(n);
819 let x_exact: Vec<f64> = (0..n).map(|i| (i + 1) as f64).collect();
820 let b = a.mul_vec(&x_exact).unwrap();
821
822 let jac = Jacobi::new(&a);
823 let opts = IterativeOptions::default().tol(1e-10);
824
825 let result_plain = cg(&a, &b, None, &opts).unwrap();
826 let result_pcg = pcg(&a, &b, &jac, None, &opts).unwrap();
827
828 assert!(result_pcg.converged);
829 assert!(
832 result_pcg.iterations <= result_plain.iterations + 5,
833 "PCG: {} vs CG: {}",
834 result_pcg.iterations,
835 result_plain.iterations
836 );
837 }
838
839 #[test]
840 fn test_pcg_ilu0() {
841 let n = 50;
842 let a = laplacian_1d(n);
843 let x_exact: Vec<f64> = (0..n).map(|i| (i + 1) as f64).collect();
844 let b = a.mul_vec(&x_exact).unwrap();
845
846 let ilu = Ilu0::new(&a).unwrap();
847 let opts = IterativeOptions::default().tol(1e-10);
848
849 let result_plain = cg(&a, &b, None, &opts).unwrap();
850 let result_pcg = pcg(&a, &b, &ilu, None, &opts).unwrap();
851
852 assert!(result_pcg.converged);
853 assert!(
855 result_pcg.iterations <= 2,
856 "ILU(0) PCG should converge immediately on tridiagonal: {} iters",
857 result_pcg.iterations
858 );
859 assert!(
860 result_pcg.iterations < result_plain.iterations,
861 "ILU(0) should beat unpreconditioned CG"
862 );
863 }
864
865 #[test]
866 fn test_gmres_nonsymmetric() {
867 let n = 30;
868 let a = convection_diffusion(n, 0.1);
869 let x_exact: Vec<f64> = (0..n).map(|i| (i as f64).sin()).collect();
870 let b = a.mul_vec(&x_exact).unwrap();
871
872 let opts = IterativeOptions::default().tol(1e-10).restart(30);
873 let result = gmres(&a, &b, None, &opts).unwrap();
874
875 assert!(
876 result.converged,
877 "GMRES did not converge after {} iterations",
878 result.iterations
879 );
880 for i in 0..n {
881 assert!(
882 (result.x[i] - x_exact[i]).abs() < 1e-6,
883 "GMRES error at {}: {} vs {}",
884 i,
885 result.x[i],
886 x_exact[i]
887 );
888 }
889 }
890
891 #[test]
892 fn test_gmres_spd_matches_cg() {
893 let n = 20;
895 let a = laplacian_1d(n);
896 let b: Vec<f64> = (0..n).map(|i| (i + 1) as f64).collect();
897
898 let opts = IterativeOptions::default().tol(1e-10);
899 let cg_result = cg(&a, &b, None, &opts).unwrap();
900 let gmres_result = gmres(&a, &b, None, &opts).unwrap();
901
902 assert!(cg_result.converged);
903 assert!(gmres_result.converged);
904 for i in 0..n {
905 assert!(
906 (cg_result.x[i] - gmres_result.x[i]).abs() < 1e-6,
907 "CG vs GMRES mismatch at {}: {} vs {}",
908 i,
909 cg_result.x[i],
910 gmres_result.x[i]
911 );
912 }
913 }
914
915 #[test]
916 fn test_bicgstab_nonsymmetric() {
917 let n = 30;
918 let a = convection_diffusion(n, 0.2);
919 let x_exact: Vec<f64> = (0..n).map(|i| (i as f64).cos()).collect();
920 let b = a.mul_vec(&x_exact).unwrap();
921
922 let opts = IterativeOptions::default().tol(1e-10);
923 let result = bicgstab(&a, &b, None, &opts).unwrap();
924
925 assert!(
926 result.converged,
927 "BiCGSTAB did not converge after {} iters",
928 result.iterations
929 );
930 for i in 0..n {
931 assert!(
932 (result.x[i] - x_exact[i]).abs() < 1e-6,
933 "BiCGSTAB error at {}: {} vs {}",
934 i,
935 result.x[i],
936 x_exact[i]
937 );
938 }
939 }
940
941 #[test]
942 fn test_bicgstab_spd_matches_cg() {
943 let n = 20;
944 let a = laplacian_1d(n);
945 let b: Vec<f64> = (0..n).map(|i| (i + 1) as f64).collect();
946
947 let opts = IterativeOptions::default().tol(1e-10);
948 let cg_result = cg(&a, &b, None, &opts).unwrap();
949 let bicg_result = bicgstab(&a, &b, None, &opts).unwrap();
950
951 assert!(cg_result.converged);
952 assert!(bicg_result.converged);
953 for i in 0..n {
954 assert!(
955 (cg_result.x[i] - bicg_result.x[i]).abs() < 1e-6,
956 "CG vs BiCGSTAB mismatch at {}: {} vs {}",
957 i,
958 cg_result.x[i],
959 bicg_result.x[i]
960 );
961 }
962 }
963
964 #[test]
965 fn test_minres_symmetric_indefinite() {
966 let n = 5;
968 let mut triplets = Vec::new();
969 for i in 0..n {
970 let sign = if i % 2 == 0 { -1.0 } else { 2.0 };
971 triplets.push((i, i, sign * 3.0));
972 if i > 0 {
973 triplets.push((i, i - 1, 1.0));
974 }
975 if i < n - 1 {
976 triplets.push((i, i + 1, 1.0));
977 }
978 }
979 let a = SparseMatrix::from_triplets(n, n, &triplets).unwrap();
980 let x_exact = vec![1.0, 2.0, 3.0, 4.0, 5.0];
981 let b = a.mul_vec(&x_exact).unwrap();
982
983 let opts = IterativeOptions::default().tol(1e-10);
984 let result = minres(&a, &b, None, &opts).unwrap();
985
986 assert!(
987 result.converged,
988 "MINRES did not converge after {} iters",
989 result.iterations
990 );
991 for i in 0..n {
992 assert!(
993 (result.x[i] - x_exact[i]).abs() < 1e-4,
994 "MINRES error at {}: {} vs {}",
995 i,
996 result.x[i],
997 x_exact[i]
998 );
999 }
1000 }
1001
1002 #[test]
1003 fn test_minres_spd() {
1004 let n = 20;
1005 let a = laplacian_1d(n);
1006 let x_exact: Vec<f64> = (0..n).map(|i| (i + 1) as f64).collect();
1007 let b = a.mul_vec(&x_exact).unwrap();
1008
1009 let opts = IterativeOptions::default().tol(1e-10);
1010 let result = minres(&a, &b, None, &opts).unwrap();
1011
1012 assert!(result.converged, "MINRES did not converge");
1013 for i in 0..n {
1014 assert!(
1015 (result.x[i] - x_exact[i]).abs() < 1e-6,
1016 "MINRES error at {}: {} vs {}",
1017 i,
1018 result.x[i],
1019 x_exact[i]
1020 );
1021 }
1022 }
1023
1024 #[test]
1025 fn test_cg_zero_rhs() {
1026 let n = 10;
1027 let a = laplacian_1d(n);
1028 let b = vec![0.0; n];
1029 let opts = IterativeOptions::default();
1030 let result = cg(&a, &b, None, &opts).unwrap();
1031 assert!(result.converged);
1032 assert_eq!(result.iterations, 0);
1033 }
1034
1035 #[test]
1036 fn test_dimension_mismatch() {
1037 let a = laplacian_1d(5);
1038 let b = vec![1.0; 3]; let opts = IterativeOptions::default();
1040 assert!(cg(&a, &b, None, &opts).is_err());
1041 assert!(gmres(&a, &b, None, &opts).is_err());
1042 assert!(bicgstab(&a, &b, None, &opts).is_err());
1043 assert!(minres(&a, &b, None, &opts).is_err());
1044 }
1045
1046 #[test]
1047 fn test_pcg_ssor() {
1048 let n = 20;
1049 let a = laplacian_1d(n);
1050 let x_exact: Vec<f64> = (0..n).map(|i| (i + 1) as f64).collect();
1051 let b = a.mul_vec(&x_exact).unwrap();
1052
1053 let ssor = Ssor::new(&a, 1.0);
1055 let opts = IterativeOptions::default().tol(1e-8);
1056 let result = pcg(&a, &b, &ssor, None, &opts).unwrap();
1057
1058 assert!(
1059 result.converged,
1060 "SSOR-PCG did not converge after {} iters",
1061 result.iterations
1062 );
1063 for i in 0..n {
1064 assert!(
1065 (result.x[i] - x_exact[i]).abs() < 1e-4,
1066 "SSOR-PCG error at {}: {} vs {}",
1067 i,
1068 result.x[i],
1069 x_exact[i]
1070 );
1071 }
1072 }
1073
1074 #[test]
1075 fn test_iterative_vs_direct() {
1076 let n = 100;
1078 let a = laplacian_1d(n);
1079 let b: Vec<f64> = (0..n).map(|i| (i as f64 * 0.1).sin()).collect();
1080
1081 let lu = crate::SparseLU::new(&a).unwrap();
1083 let x_direct = lu.solve(&b).unwrap();
1084
1085 let opts = IterativeOptions::default().tol(1e-12);
1087 let x_cg = cg(&a, &b, None, &opts).unwrap();
1088 assert!(x_cg.converged);
1089
1090 for i in 0..n {
1091 assert!(
1092 (x_cg.x[i] - x_direct[i]).abs() < 1e-8,
1093 "CG vs direct mismatch at {}: {} vs {}",
1094 i,
1095 x_cg.x[i],
1096 x_direct[i]
1097 );
1098 }
1099 }
1100
1101 #[test]
1102 fn test_cg_f32() {
1103 let n = 10;
1104 let mut triplets: Vec<(usize, usize, f32)> = Vec::new();
1105 for i in 0..n {
1106 triplets.push((i, i, 2.0f32));
1107 if i > 0 {
1108 triplets.push((i, i - 1, -1.0f32));
1109 }
1110 if i < n - 1 {
1111 triplets.push((i, i + 1, -1.0f32));
1112 }
1113 }
1114 let a = SparseMatrix::from_triplets(n, n, &triplets).unwrap();
1115 let b: Vec<f32> = (0..n).map(|i| (i + 1) as f32).collect();
1116
1117 let opts = IterativeOptions::default().tol(1e-5f32);
1118 let result = cg(&a, &b, None, &opts).unwrap();
1119 assert!(result.converged);
1120 }
1121
1122 #[test]
1123 fn test_gmres_restart() {
1124 let n = 30;
1126 let a = convection_diffusion(n, 0.1);
1127 let x_exact: Vec<f64> = (0..n).map(|i| (i as f64).sin()).collect();
1128 let b = a.mul_vec(&x_exact).unwrap();
1129
1130 let opts = IterativeOptions::default().tol(1e-10).restart(5); let result = gmres(&a, &b, None, &opts).unwrap();
1132
1133 assert!(result.converged, "Restarted GMRES did not converge");
1134 for i in 0..n {
1135 assert!(
1136 (result.x[i] - x_exact[i]).abs() < 1e-5,
1137 "Restarted GMRES error at {}",
1138 i
1139 );
1140 }
1141 }
1142}