solvers/preconditioners/
amg.rs

1//! Algebraic Multigrid (AMG) Preconditioner
2//!
3//! This module implements an algebraic multigrid preconditioner inspired by
4//! hypre's BoomerAMG, designed for better parallel scalability across CPU cores.
5//!
6//! ## Features
7//!
8//! - **Parallel coarsening**: Classical Ruge-Stüben (RS) and PMIS algorithms
9//! - **Interpolation**: Standard and extended interpolation operators
10//! - **Smoothers**: Jacobi (fully parallel) and symmetric Gauss-Seidel
11//! - **V-cycle**: Standard V(ν₁, ν₂) cycling with configurable pre/post smoothing
12//!
13//! ## Scalability
14//!
15//! The AMG preconditioner scales better than ILU across multiple cores because:
16//! - Coarsening can be parallelized (PMIS is inherently parallel)
17//! - Jacobi smoothing is embarrassingly parallel
18//! - Each level's operations can be parallelized independently
19//!
20//! ## Usage
21//!
22//! ```ignore
23//! use solvers::{AmgPreconditioner, AmgConfig, CsrMatrix};
24//!
25//! let config = AmgConfig::default();
26//! let precond = AmgPreconditioner::from_csr(&matrix, config);
27//!
28//! // Use with GMRES
29//! let z = precond.apply(&residual);
30//! ```
31
32#[cfg(any(feature = "native", feature = "wasm"))]
33use crate::parallel::{parallel_enumerate_map, parallel_map_indexed};
34use crate::sparse::CsrMatrix;
35use crate::traits::{ComplexField, Preconditioner};
36use ndarray::Array1;
37use num_traits::FromPrimitive;
38use std::collections::HashMap;
39
40/// Coarsening algorithm for AMG hierarchy construction
41#[derive(Debug, Clone, Copy, PartialEq, Eq, Default)]
42pub enum AmgCoarsening {
43    /// Classical Ruge-Stüben coarsening
44    /// Good quality but limited parallelism in the selection phase
45    #[default]
46    RugeStuben,
47
48    /// Parallel Modified Independent Set (PMIS)
49    /// Better parallel scalability, may produce slightly larger coarse grids
50    Pmis,
51
52    /// Hybrid MIS (HMIS) - PMIS in first pass, then RS cleanup
53    /// Balance between quality and parallelism
54    Hmis,
55}
56
57/// Interpolation operator type
58#[derive(Debug, Clone, Copy, PartialEq, Eq, Default)]
59pub enum AmgInterpolation {
60    /// Standard interpolation - direct interpolation from coarse neighbors
61    #[default]
62    Standard,
63
64    /// Extended interpolation - includes indirect (distance-2) connections
65    /// Better for some problem types but more expensive
66    Extended,
67
68    /// Direct interpolation - simplest, only immediate strong connections
69    /// Fastest but may have poor convergence for hard problems
70    Direct,
71}
72
73/// Smoother type for AMG relaxation
74#[derive(Debug, Clone, Copy, PartialEq, Eq, Default)]
75pub enum AmgSmoother {
76    /// Jacobi relaxation - fully parallel, requires damping (ω ≈ 0.6-0.8)
77    #[default]
78    Jacobi,
79
80    /// l1-Jacobi - Jacobi with l1 norm scaling, more robust
81    L1Jacobi,
82
83    /// Symmetric Gauss-Seidel - forward then backward sweep
84    /// Better convergence but limited parallelism
85    SymmetricGaussSeidel,
86
87    /// Chebyshev polynomial smoother - fully parallel, no damping needed
88    Chebyshev,
89}
90
91/// AMG cycle type
92#[derive(Debug, Clone, Copy, PartialEq, Eq, Default)]
93pub enum AmgCycle {
94    /// V-cycle: one visit to each level
95    #[default]
96    VCycle,
97
98    /// W-cycle: two visits to coarser levels (more expensive)
99    WCycle,
100
101    /// F-cycle: hybrid between V and W
102    FCycle,
103}
104
105/// Configuration for AMG preconditioner
106#[derive(Debug, Clone)]
107pub struct AmgConfig {
108    /// Coarsening algorithm
109    pub coarsening: AmgCoarsening,
110
111    /// Interpolation operator type
112    pub interpolation: AmgInterpolation,
113
114    /// Smoother for pre- and post-relaxation
115    pub smoother: AmgSmoother,
116
117    /// Cycle type (V, W, or F)
118    pub cycle: AmgCycle,
119
120    /// Strong connection threshold (default: 0.25)
121    /// Connections with |a_ij| >= θ * max_k |a_ik| are considered strong
122    pub strong_threshold: f64,
123
124    /// Maximum number of levels in the hierarchy
125    pub max_levels: usize,
126
127    /// Coarsest level size - switch to direct solve below this
128    pub coarse_size: usize,
129
130    /// Number of pre-smoothing sweeps (ν₁)
131    pub num_pre_smooth: usize,
132
133    /// Number of post-smoothing sweeps (ν₂)
134    pub num_post_smooth: usize,
135
136    /// Jacobi damping parameter (ω)
137    pub jacobi_weight: f64,
138
139    /// Truncation factor for interpolation (drop small weights)
140    pub trunc_factor: f64,
141
142    /// Maximum interpolation stencil size per row
143    pub max_interp_elements: usize,
144
145    /// Enable aggressive coarsening on first few levels
146    pub aggressive_coarsening_levels: usize,
147}
148
149impl Default for AmgConfig {
150    fn default() -> Self {
151        Self {
152            coarsening: AmgCoarsening::default(),
153            interpolation: AmgInterpolation::default(),
154            smoother: AmgSmoother::default(),
155            cycle: AmgCycle::default(),
156            strong_threshold: 0.25,
157            max_levels: 25,
158            coarse_size: 50,
159            num_pre_smooth: 1,
160            num_post_smooth: 1,
161            jacobi_weight: 0.6667, // 2/3 is optimal for Poisson
162            trunc_factor: 0.0,
163            max_interp_elements: 4,
164            aggressive_coarsening_levels: 0,
165        }
166    }
167}
168
169impl AmgConfig {
170    /// Configuration optimized for BEM systems
171    ///
172    /// BEM matrices are typically denser and less sparse than FEM,
173    /// requiring adjusted thresholds.
174    pub fn for_bem() -> Self {
175        Self {
176            strong_threshold: 0.5,           // Higher for denser BEM matrices
177            coarsening: AmgCoarsening::Pmis, // Better parallel scalability
178            smoother: AmgSmoother::L1Jacobi, // More robust for BEM
179            max_interp_elements: 6,
180            ..Default::default()
181        }
182    }
183
184    /// Configuration optimized for FEM systems
185    pub fn for_fem() -> Self {
186        Self {
187            strong_threshold: 0.25,
188            coarsening: AmgCoarsening::RugeStuben,
189            smoother: AmgSmoother::SymmetricGaussSeidel,
190            ..Default::default()
191        }
192    }
193
194    /// Configuration optimized for maximum parallel scalability
195    pub fn for_parallel() -> Self {
196        Self {
197            coarsening: AmgCoarsening::Pmis,
198            smoother: AmgSmoother::Jacobi,
199            jacobi_weight: 0.8,
200            num_pre_smooth: 2,
201            num_post_smooth: 2,
202            ..Default::default()
203        }
204    }
205
206    /// Configuration for difficult/ill-conditioned problems
207    pub fn for_difficult_problems() -> Self {
208        Self {
209            coarsening: AmgCoarsening::RugeStuben,
210            interpolation: AmgInterpolation::Extended,
211            smoother: AmgSmoother::SymmetricGaussSeidel,
212            strong_threshold: 0.25,
213            max_interp_elements: 8,
214            num_pre_smooth: 2,
215            num_post_smooth: 2,
216            ..Default::default()
217        }
218    }
219}
220
221/// Point classification in coarsening
222#[derive(Debug, Clone, Copy, PartialEq, Eq)]
223enum PointType {
224    /// Undecided
225    Undecided,
226    /// Coarse point (C-point)
227    Coarse,
228    /// Fine point (F-point)
229    Fine,
230}
231
232/// Single level in the AMG hierarchy
233#[derive(Debug, Clone)]
234struct AmgLevel<T: ComplexField> {
235    /// System matrix A at this level (CSR format)
236    matrix: CsrMatrix<T>,
237
238    /// Prolongation operator P: coarse -> fine
239    prolongation: Option<CsrMatrix<T>>,
240
241    /// Restriction operator R: fine -> coarse (typically R = P^T)
242    restriction: Option<CsrMatrix<T>>,
243
244    /// Inverse diagonal for Jacobi smoothing
245    diag_inv: Array1<T>,
246
247    /// Coarse-to-fine mapping
248    coarse_to_fine: Vec<usize>,
249
250    /// Number of DOFs at this level
251    num_dofs: usize,
252}
253
254/// Algebraic Multigrid Preconditioner
255///
256/// Implements a classical AMG V-cycle preconditioner with configurable
257/// coarsening, interpolation, and smoothing strategies.
258#[derive(Debug, Clone)]
259pub struct AmgPreconditioner<T: ComplexField> {
260    /// AMG hierarchy (finest to coarsest)
261    levels: Vec<AmgLevel<T>>,
262
263    /// Configuration
264    config: AmgConfig,
265
266    /// Statistics
267    setup_time_ms: f64,
268    grid_complexity: f64,
269    operator_complexity: f64,
270}
271
272impl<T: ComplexField> AmgPreconditioner<T> {
273    /// Create AMG preconditioner from a CSR matrix
274    pub fn from_csr(matrix: &CsrMatrix<T>, config: AmgConfig) -> Self {
275        let start = std::time::Instant::now();
276
277        let mut levels = Vec::new();
278        let mut current_matrix = matrix.clone();
279
280        // Extract diagonal for first level
281        let diag_inv = Self::compute_diag_inv(&current_matrix);
282
283        levels.push(AmgLevel {
284            matrix: current_matrix.clone(),
285            prolongation: None,
286            restriction: None,
287            diag_inv,
288            coarse_to_fine: Vec::new(),
289            num_dofs: current_matrix.num_rows,
290        });
291
292        // Build hierarchy
293        for _level_idx in 0..config.max_levels - 1 {
294            let n = current_matrix.num_rows;
295            if n <= config.coarse_size {
296                break;
297            }
298
299            // Compute strength matrix
300            let strong_connections =
301                Self::compute_strength_matrix(&current_matrix, config.strong_threshold);
302
303            // Coarsening
304            let (point_types, coarse_to_fine) = match config.coarsening {
305                AmgCoarsening::RugeStuben => {
306                    Self::coarsen_ruge_stuben(&current_matrix, &strong_connections)
307                }
308                AmgCoarsening::Pmis | AmgCoarsening::Hmis => {
309                    Self::coarsen_pmis(&current_matrix, &strong_connections)
310                }
311            };
312
313            let num_coarse = coarse_to_fine.len();
314            if num_coarse == 0 || num_coarse >= n {
315                // Can't coarsen further
316                break;
317            }
318
319            // Build interpolation operator
320            let prolongation = Self::build_interpolation(
321                &current_matrix,
322                &strong_connections,
323                &point_types,
324                &coarse_to_fine,
325                &config,
326            );
327
328            // Restriction is transpose of prolongation
329            let restriction = Self::transpose_csr(&prolongation);
330
331            // Galerkin coarse grid: A_c = R * A * P
332            let coarse_matrix =
333                Self::galerkin_product(&restriction, &current_matrix, &prolongation);
334
335            // Extract diagonal for new level
336            let coarse_diag_inv = Self::compute_diag_inv(&coarse_matrix);
337
338            // Update level with P and R
339            if let Some(last) = levels.last_mut() {
340                last.prolongation = Some(prolongation);
341                last.restriction = Some(restriction);
342                last.coarse_to_fine = coarse_to_fine;
343            }
344
345            // Add coarse level
346            levels.push(AmgLevel {
347                matrix: coarse_matrix.clone(),
348                prolongation: None,
349                restriction: None,
350                diag_inv: coarse_diag_inv,
351                coarse_to_fine: Vec::new(),
352                num_dofs: num_coarse,
353            });
354
355            current_matrix = coarse_matrix;
356        }
357
358        let setup_time_ms = start.elapsed().as_secs_f64() * 1000.0;
359
360        // Compute complexities
361        let (grid_complexity, operator_complexity) = Self::compute_complexities(&levels);
362
363        Self {
364            levels,
365            config,
366            setup_time_ms,
367            grid_complexity,
368            operator_complexity,
369        }
370    }
371
372    /// Get number of levels in hierarchy
373    pub fn num_levels(&self) -> usize {
374        self.levels.len()
375    }
376
377    /// Get setup time in milliseconds
378    pub fn setup_time_ms(&self) -> f64 {
379        self.setup_time_ms
380    }
381
382    /// Get grid complexity (sum of DOFs / fine DOFs)
383    pub fn grid_complexity(&self) -> f64 {
384        self.grid_complexity
385    }
386
387    /// Get operator complexity (sum of nnz / fine nnz)
388    pub fn operator_complexity(&self) -> f64 {
389        self.operator_complexity
390    }
391
392    /// Get configuration
393    pub fn config(&self) -> &AmgConfig {
394        &self.config
395    }
396
397    /// Compute inverse diagonal for Jacobi smoothing
398    fn compute_diag_inv(matrix: &CsrMatrix<T>) -> Array1<T> {
399        let n = matrix.num_rows;
400        let mut diag_inv = Array1::from_elem(n, T::one());
401
402        for i in 0..n {
403            let diag = matrix.get(i, i);
404            let tol = T::Real::from_f64(1e-15).unwrap();
405            if diag.norm() > tol {
406                diag_inv[i] = diag.inv();
407            }
408        }
409
410        diag_inv
411    }
412
413    /// Compute strength of connection matrix
414    ///
415    /// Entry (i,j) is strong if |a_ij| >= θ * max_k!=i |a_ik|
416    fn compute_strength_matrix(matrix: &CsrMatrix<T>, theta: f64) -> Vec<Vec<usize>> {
417        let n = matrix.num_rows;
418
419        #[cfg(any(feature = "native", feature = "wasm"))]
420        {
421            parallel_map_indexed(n, |i| {
422                // Find max off-diagonal magnitude in row i
423                let mut max_off_diag = T::Real::from_f64(0.0).unwrap();
424                for (j, val) in matrix.row_entries(i) {
425                    if i != j {
426                        let norm = val.norm();
427                        if norm > max_off_diag {
428                            max_off_diag = norm;
429                        }
430                    }
431                }
432
433                let threshold = T::Real::from_f64(theta).unwrap() * max_off_diag;
434
435                // Collect strong connections
436                let mut row_strong = Vec::new();
437                for (j, val) in matrix.row_entries(i) {
438                    if i != j && val.norm() >= threshold {
439                        row_strong.push(j);
440                    }
441                }
442                row_strong
443            })
444        }
445
446        #[cfg(not(any(feature = "native", feature = "wasm")))]
447        {
448            let mut strong: Vec<Vec<usize>> = vec![Vec::new(); n];
449            for (i, row_strong) in strong.iter_mut().enumerate().take(n) {
450                // Find max off-diagonal magnitude in row i
451                let mut max_off_diag = T::Real::from_f64(0.0).unwrap();
452                for (j, val) in matrix.row_entries(i) {
453                    if i != j {
454                        let norm = val.norm();
455                        if norm > max_off_diag {
456                            max_off_diag = norm;
457                        }
458                    }
459                }
460
461                let threshold = T::Real::from_f64(theta).unwrap() * max_off_diag;
462
463                // Collect strong connections
464                for (j, val) in matrix.row_entries(i) {
465                    if i != j && val.norm() >= threshold {
466                        row_strong.push(j);
467                    }
468                }
469            }
470            strong
471        }
472    }
473
474    /// Classical Ruge-Stüben coarsening
475    fn coarsen_ruge_stuben(
476        matrix: &CsrMatrix<T>,
477        strong: &[Vec<usize>],
478    ) -> (Vec<PointType>, Vec<usize>) {
479        let n = matrix.num_rows;
480        let mut point_types = vec![PointType::Undecided; n];
481
482        // Compute influence measure λ_i = |S_i^T| (how many points strongly depend on i)
483        let mut lambda: Vec<usize> = vec![0; n];
484        for row in strong.iter().take(n) {
485            for &j in row {
486                lambda[j] += 1;
487            }
488        }
489
490        // Build priority queue (we use a simple approach: process by decreasing lambda)
491        let mut order: Vec<usize> = (0..n).collect();
492        order.sort_by(|&a, &b| lambda[b].cmp(&lambda[a]));
493
494        // First pass: select C-points
495        for &i in &order {
496            if point_types[i] != PointType::Undecided {
497                continue;
498            }
499
500            // Make i a C-point
501            point_types[i] = PointType::Coarse;
502
503            // All points that strongly depend on i become F-points
504            for j in 0..n {
505                if point_types[j] == PointType::Undecided && strong[j].contains(&i) {
506                    point_types[j] = PointType::Fine;
507                    // Update lambda for neighbors
508                    for &k in &strong[j] {
509                        if point_types[k] == PointType::Undecided {
510                            lambda[k] = lambda[k].saturating_sub(1);
511                        }
512                    }
513                }
514            }
515        }
516
517        // Ensure all remaining undecided become fine
518        for pt in &mut point_types {
519            if *pt == PointType::Undecided {
520                *pt = PointType::Fine;
521            }
522        }
523
524        // Build coarse-to-fine mapping
525        let coarse_to_fine: Vec<usize> = (0..n)
526            .filter(|&i| point_types[i] == PointType::Coarse)
527            .collect();
528
529        (point_types, coarse_to_fine)
530    }
531
532    /// Parallel Modified Independent Set (PMIS) coarsening
533    fn coarsen_pmis(matrix: &CsrMatrix<T>, strong: &[Vec<usize>]) -> (Vec<PointType>, Vec<usize>) {
534        let n = matrix.num_rows;
535        let mut point_types = vec![PointType::Undecided; n];
536
537        // Compute weights based on number of strong connections
538        let weights: Vec<f64> = (0..n)
539            .map(|i| {
540                // Weight = |S_i| + random tie-breaker
541                strong[i].len() as f64 + (i as f64 * 0.0001) % 0.001
542            })
543            .collect();
544
545        // Iterative independent set selection
546        let mut changed = true;
547        let mut iteration = 0;
548        const MAX_ITERATIONS: usize = 100;
549
550        while changed && iteration < MAX_ITERATIONS {
551            changed = false;
552            iteration += 1;
553
554            #[cfg(any(feature = "native", feature = "wasm"))]
555            {
556                // Parallel pass: determine new C-points and F-points
557                let updates: Vec<(usize, PointType)> =
558                    parallel_enumerate_map(&point_types, |i, pt| {
559                        if *pt != PointType::Undecided {
560                            return (i, *pt);
561                        }
562
563                        // Check if i has maximum weight among undecided strong neighbors
564                        let mut is_max = true;
565                        for &j in &strong[i] {
566                            if point_types[j] == PointType::Undecided && weights[j] > weights[i] {
567                                is_max = false;
568                                break;
569                            }
570                        }
571
572                        // Check if any strong neighbor is already C
573                        let has_c_neighbor = strong[i]
574                            .iter()
575                            .any(|&j| point_types[j] == PointType::Coarse);
576
577                        if has_c_neighbor {
578                            (i, PointType::Fine)
579                        } else if is_max {
580                            (i, PointType::Coarse)
581                        } else {
582                            (i, PointType::Undecided)
583                        }
584                    });
585
586                for (i, new_type) in updates {
587                    if point_types[i] != new_type {
588                        point_types[i] = new_type;
589                        changed = true;
590                    }
591                }
592            }
593
594            #[cfg(not(any(feature = "native", feature = "wasm")))]
595            {
596                // Sequential fallback
597                let old_types = point_types.clone();
598                for i in 0..n {
599                    if old_types[i] != PointType::Undecided {
600                        continue;
601                    }
602
603                    // Check if i has maximum weight among undecided strong neighbors
604                    let mut is_max = true;
605                    for &j in &strong[i] {
606                        if old_types[j] == PointType::Undecided && weights[j] > weights[i] {
607                            is_max = false;
608                            break;
609                        }
610                    }
611
612                    // Check if any strong neighbor is already C
613                    let has_c_neighbor =
614                        strong[i].iter().any(|&j| old_types[j] == PointType::Coarse);
615
616                    if has_c_neighbor {
617                        point_types[i] = PointType::Fine;
618                        changed = true;
619                    } else if is_max {
620                        point_types[i] = PointType::Coarse;
621                        changed = true;
622                    }
623                }
624            }
625        }
626
627        // Any remaining undecided become coarse
628        for pt in &mut point_types {
629            if *pt == PointType::Undecided {
630                *pt = PointType::Coarse;
631            }
632        }
633
634        // Build coarse-to-fine mapping
635        let coarse_to_fine: Vec<usize> = (0..n)
636            .filter(|&i| point_types[i] == PointType::Coarse)
637            .collect();
638
639        (point_types, coarse_to_fine)
640    }
641
642    /// Build interpolation operator P
643    fn build_interpolation(
644        matrix: &CsrMatrix<T>,
645        strong: &[Vec<usize>],
646        point_types: &[PointType],
647        coarse_to_fine: &[usize],
648        config: &AmgConfig,
649    ) -> CsrMatrix<T> {
650        let n_fine = matrix.num_rows;
651        let n_coarse = coarse_to_fine.len();
652
653        // Build fine-to-coarse mapping
654        let mut fine_to_coarse = vec![usize::MAX; n_fine];
655        for (coarse_idx, &fine_idx) in coarse_to_fine.iter().enumerate() {
656            fine_to_coarse[fine_idx] = coarse_idx;
657        }
658
659        // Build P row by row
660        let mut triplets: Vec<(usize, usize, T)> = Vec::new();
661
662        for i in 0..n_fine {
663            match point_types[i] {
664                PointType::Coarse => {
665                    // C-point: identity mapping P_ij = 1 if j = coarse_index(i)
666                    let coarse_idx = fine_to_coarse[i];
667                    triplets.push((i, coarse_idx, T::one()));
668                }
669                PointType::Fine => {
670                    // F-point: interpolate from strong C-neighbors
671                    let a_ii = matrix.get(i, i);
672
673                    // Collect strong C-neighbors
674                    let c_neighbors: Vec<usize> = strong[i]
675                        .iter()
676                        .copied()
677                        .filter(|&j| point_types[j] == PointType::Coarse)
678                        .collect();
679
680                    if c_neighbors.is_empty() {
681                        continue;
682                    }
683
684                    // Standard interpolation weights
685                    match config.interpolation {
686                        AmgInterpolation::Direct | AmgInterpolation::Standard => {
687                            let mut weights: Vec<(usize, T)> = Vec::new();
688                            let mut sum_weights = T::zero();
689
690                            for &j in &c_neighbors {
691                                let a_ij = matrix.get(i, j);
692                                let tol = T::Real::from_f64(1e-15).unwrap();
693                                if a_ii.norm() > tol {
694                                    let w = T::zero() - a_ij * a_ii.inv();
695                                    weights.push((fine_to_coarse[j], w));
696                                    sum_weights += w;
697                                }
698                            }
699
700                            // Add weak connections contribution (standard interpolation)
701                            if config.interpolation == AmgInterpolation::Standard {
702                                let mut weak_sum = T::zero();
703                                for (j, val) in matrix.row_entries(i) {
704                                    if j != i && !c_neighbors.contains(&j) {
705                                        weak_sum += val;
706                                    }
707                                }
708
709                                let tol = T::Real::from_f64(1e-15).unwrap();
710                                if sum_weights.norm() > tol && weak_sum.norm() > tol {
711                                    let scale = T::one() + weak_sum * (a_ii * sum_weights).inv();
712                                    for (_, w) in &mut weights {
713                                        *w *= scale;
714                                    }
715                                }
716                            }
717
718                            // Truncate small weights if configured
719                            if config.trunc_factor > 0.0 {
720                                let max_w = weights.iter().map(|(_, w)| w.norm()).fold(
721                                    T::Real::from_f64(0.0).unwrap(),
722                                    |a, b| {
723                                        if a > b { a } else { b }
724                                    },
725                                );
726                                let threshold =
727                                    T::Real::from_f64(config.trunc_factor).unwrap() * max_w;
728                                weights.retain(|(_, w)| w.norm() >= threshold);
729
730                                if weights.len() > config.max_interp_elements {
731                                    weights.sort_by(|a, b| {
732                                        b.1.norm().partial_cmp(&a.1.norm()).unwrap()
733                                    });
734                                    weights.truncate(config.max_interp_elements);
735                                }
736                            }
737
738                            for (coarse_idx, w) in weights {
739                                triplets.push((i, coarse_idx, w));
740                            }
741                        }
742                        AmgInterpolation::Extended => {
743                            let mut weights: Vec<(usize, T)> = Vec::new();
744
745                            // Direct C-neighbors
746                            for &j in &c_neighbors {
747                                let a_ij = matrix.get(i, j);
748                                let tol = T::Real::from_f64(1e-15).unwrap();
749                                if a_ii.norm() > tol {
750                                    let w = T::zero() - a_ij * a_ii.inv();
751                                    weights.push((fine_to_coarse[j], w));
752                                }
753                            }
754
755                            // F-neighbors contribute through their C-neighbors
756                            let f_neighbors: Vec<usize> = strong[i]
757                                .iter()
758                                .copied()
759                                .filter(|&j| point_types[j] == PointType::Fine)
760                                .collect();
761
762                            for &k in &f_neighbors {
763                                let a_ik = matrix.get(i, k);
764                                let a_kk = matrix.get(k, k);
765
766                                let tol = T::Real::from_f64(1e-15).unwrap();
767                                if a_kk.norm() < tol {
768                                    continue;
769                                }
770
771                                for &j in &strong[k] {
772                                    if point_types[j] == PointType::Coarse {
773                                        let a_kj = matrix.get(k, j);
774                                        let w = T::zero() - a_ik * a_kj * (a_ii * a_kk).inv();
775
776                                        let coarse_j = fine_to_coarse[j];
777                                        if let Some((_, existing)) =
778                                            weights.iter_mut().find(|(idx, _)| *idx == coarse_j)
779                                        {
780                                            *existing += w;
781                                        } else {
782                                            weights.push((coarse_j, w));
783                                        }
784                                    }
785                                }
786                            }
787
788                            if weights.len() > config.max_interp_elements {
789                                weights
790                                    .sort_by(|a, b| b.1.norm().partial_cmp(&a.1.norm()).unwrap());
791                                weights.truncate(config.max_interp_elements);
792                            }
793
794                            for (coarse_idx, w) in weights {
795                                triplets.push((i, coarse_idx, w));
796                            }
797                        }
798                    }
799                }
800                PointType::Undecided => {}
801            }
802        }
803
804        CsrMatrix::from_triplets(n_fine, n_coarse, triplets)
805    }
806
807    /// Transpose a CSR matrix
808    fn transpose_csr(matrix: &CsrMatrix<T>) -> CsrMatrix<T> {
809        let m = matrix.num_rows;
810        let n = matrix.num_cols;
811
812        let mut triplets: Vec<(usize, usize, T)> = Vec::new();
813        for i in 0..m {
814            for (j, val) in matrix.row_entries(i) {
815                triplets.push((j, i, val));
816            }
817        }
818
819        CsrMatrix::from_triplets(n, m, triplets)
820    }
821
822    /// Compute Galerkin coarse grid operator: A_c = R * A * P
823    fn galerkin_product(r: &CsrMatrix<T>, a: &CsrMatrix<T>, p: &CsrMatrix<T>) -> CsrMatrix<T> {
824        let ap = Self::sparse_matmul(a, p);
825        Self::sparse_matmul(r, &ap)
826    }
827
828    /// Sparse matrix multiplication
829    fn sparse_matmul(a: &CsrMatrix<T>, b: &CsrMatrix<T>) -> CsrMatrix<T> {
830        assert_eq!(a.num_cols, b.num_rows, "Matrix dimension mismatch");
831
832        let m = a.num_rows;
833        let n = b.num_cols;
834        let mut triplets: Vec<(usize, usize, T)> = Vec::new();
835
836        for i in 0..m {
837            let mut row_acc: HashMap<usize, T> = HashMap::new();
838
839            for (k, a_ik) in a.row_entries(i) {
840                for (j, b_kj) in b.row_entries(k) {
841                    *row_acc.entry(j).or_insert(T::zero()) += a_ik * b_kj;
842                }
843            }
844
845            let tol = T::Real::from_f64(1e-15).unwrap();
846            for (j, val) in row_acc {
847                if val.norm() > tol {
848                    triplets.push((i, j, val));
849                }
850            }
851        }
852
853        CsrMatrix::from_triplets(m, n, triplets)
854    }
855
856    /// Compute grid and operator complexities
857    fn compute_complexities(levels: &[AmgLevel<T>]) -> (f64, f64) {
858        if levels.is_empty() {
859            return (1.0, 1.0);
860        }
861
862        let fine_dofs = levels[0].num_dofs as f64;
863        let fine_nnz = levels[0].matrix.nnz() as f64;
864
865        let total_dofs: f64 = levels.iter().map(|l| l.num_dofs as f64).sum();
866        let total_nnz: f64 = levels.iter().map(|l| l.matrix.nnz() as f64).sum();
867
868        let grid_complexity = total_dofs / fine_dofs;
869        let operator_complexity = total_nnz / fine_nnz;
870
871        (grid_complexity, operator_complexity)
872    }
873
874    /// Apply Jacobi smoothing: x = x + ω * D^{-1} * (b - A*x)
875    fn smooth_jacobi(
876        matrix: &CsrMatrix<T>,
877        diag_inv: &Array1<T>,
878        x: &mut Array1<T>,
879        b: &Array1<T>,
880        omega: f64,
881        num_sweeps: usize,
882    ) {
883        let omega = T::from_real(T::Real::from_f64(omega).unwrap());
884        let n = x.len();
885
886        for _ in 0..num_sweeps {
887            let r = b - &matrix.matvec(x);
888
889            #[cfg(any(feature = "native", feature = "wasm"))]
890            {
891                let updates: Vec<T> = parallel_map_indexed(n, |i| omega * diag_inv[i] * r[i]);
892                for (i, delta) in updates.into_iter().enumerate() {
893                    x[i] += delta;
894                }
895            }
896
897            #[cfg(not(any(feature = "native", feature = "wasm")))]
898            {
899                for i in 0..n {
900                    x[i] += omega * diag_inv[i] * r[i];
901                }
902            }
903        }
904    }
905
906    /// Apply l1-Jacobi smoothing
907    fn smooth_l1_jacobi(
908        matrix: &CsrMatrix<T>,
909        x: &mut Array1<T>,
910        b: &Array1<T>,
911        num_sweeps: usize,
912    ) {
913        let n = x.len();
914
915        let l1_diag: Vec<T::Real> = (0..n)
916            .map(|i| {
917                let mut sum = T::Real::from_f64(0.0).unwrap();
918                for (_, val) in matrix.row_entries(i) {
919                    sum += val.norm();
920                }
921                let tol = T::Real::from_f64(1e-15).unwrap();
922                if sum > tol {
923                    sum
924                } else {
925                    T::Real::from_f64(1.0).unwrap()
926                }
927            })
928            .collect();
929
930        for _ in 0..num_sweeps {
931            let r = b - &matrix.matvec(x);
932
933            #[cfg(any(feature = "native", feature = "wasm"))]
934            {
935                let updates: Vec<T> =
936                    parallel_map_indexed(n, |i| r[i] * T::from_real(l1_diag[i]).inv());
937                for (i, delta) in updates.into_iter().enumerate() {
938                    x[i] += delta;
939                }
940            }
941
942            #[cfg(not(any(feature = "native", feature = "wasm")))]
943            {
944                for i in 0..n {
945                    x[i] += r[i] * T::from_real(l1_diag[i]).inv();
946                }
947            }
948        }
949    }
950
951    /// Apply symmetric Gauss-Seidel smoothing
952    fn smooth_sym_gauss_seidel(
953        matrix: &CsrMatrix<T>,
954        x: &mut Array1<T>,
955        b: &Array1<T>,
956        num_sweeps: usize,
957    ) {
958        let n = x.len();
959        let tol = T::Real::from_f64(1e-15).unwrap();
960
961        for _ in 0..num_sweeps {
962            // Forward sweep
963            for i in 0..n {
964                let mut sum = b[i];
965                let mut diag = T::one();
966
967                for (j, val) in matrix.row_entries(i) {
968                    if j == i {
969                        diag = val;
970                    } else {
971                        sum -= val * x[j];
972                    }
973                }
974
975                if diag.norm() > tol {
976                    x[i] = sum * diag.inv();
977                }
978            }
979
980            // Backward sweep
981            for i in (0..n).rev() {
982                let mut sum = b[i];
983                let mut diag = T::one();
984
985                for (j, val) in matrix.row_entries(i) {
986                    if j == i {
987                        diag = val;
988                    } else {
989                        sum -= val * x[j];
990                    }
991                }
992
993                if diag.norm() > tol {
994                    x[i] = sum * diag.inv();
995                }
996            }
997        }
998    }
999
1000    /// Apply V-cycle
1001    fn v_cycle(&self, level: usize, x: &mut Array1<T>, b: &Array1<T>) {
1002        let lvl = &self.levels[level];
1003
1004        // Coarsest level: direct solve (or many smoothing iterations)
1005        if level == self.levels.len() - 1 || lvl.prolongation.is_none() {
1006            match self.config.smoother {
1007                AmgSmoother::Jacobi | AmgSmoother::Chebyshev => {
1008                    Self::smooth_jacobi(
1009                        &lvl.matrix,
1010                        &lvl.diag_inv,
1011                        x,
1012                        b,
1013                        self.config.jacobi_weight,
1014                        20,
1015                    );
1016                }
1017                AmgSmoother::L1Jacobi => {
1018                    Self::smooth_l1_jacobi(&lvl.matrix, x, b, 20);
1019                }
1020                AmgSmoother::SymmetricGaussSeidel => {
1021                    Self::smooth_sym_gauss_seidel(&lvl.matrix, x, b, 10);
1022                }
1023            }
1024            return;
1025        }
1026
1027        // Pre-smoothing
1028        match self.config.smoother {
1029            AmgSmoother::Jacobi | AmgSmoother::Chebyshev => {
1030                Self::smooth_jacobi(
1031                    &lvl.matrix,
1032                    &lvl.diag_inv,
1033                    x,
1034                    b,
1035                    self.config.jacobi_weight,
1036                    self.config.num_pre_smooth,
1037                );
1038            }
1039            AmgSmoother::L1Jacobi => {
1040                Self::smooth_l1_jacobi(&lvl.matrix, x, b, self.config.num_pre_smooth);
1041            }
1042            AmgSmoother::SymmetricGaussSeidel => {
1043                Self::smooth_sym_gauss_seidel(&lvl.matrix, x, b, self.config.num_pre_smooth);
1044            }
1045        }
1046
1047        // Compute residual: r = b - A*x
1048        let r = b - &lvl.matrix.matvec(x);
1049
1050        // Restrict residual to coarse grid: r_c = R * r
1051        let r_coarse = lvl.restriction.as_ref().unwrap().matvec(&r);
1052
1053        // Initialize coarse correction
1054        let n_coarse = self.levels[level + 1].num_dofs;
1055        let mut e_coarse = Array1::from_elem(n_coarse, T::zero());
1056
1057        // Recursive call
1058        self.v_cycle(level + 1, &mut e_coarse, &r_coarse);
1059
1060        // Prolongate correction: e = P * e_c
1061        let e = lvl.prolongation.as_ref().unwrap().matvec(&e_coarse);
1062
1063        // Apply correction: x = x + e
1064        *x = x.clone() + e;
1065
1066        // Post-smoothing
1067        match self.config.smoother {
1068            AmgSmoother::Jacobi | AmgSmoother::Chebyshev => {
1069                Self::smooth_jacobi(
1070                    &lvl.matrix,
1071                    &lvl.diag_inv,
1072                    x,
1073                    b,
1074                    self.config.jacobi_weight,
1075                    self.config.num_post_smooth,
1076                );
1077            }
1078            AmgSmoother::L1Jacobi => {
1079                Self::smooth_l1_jacobi(&lvl.matrix, x, b, self.config.num_post_smooth);
1080            }
1081            AmgSmoother::SymmetricGaussSeidel => {
1082                Self::smooth_sym_gauss_seidel(&lvl.matrix, x, b, self.config.num_post_smooth);
1083            }
1084        }
1085    }
1086}
1087
1088impl<T: ComplexField> Preconditioner<T> for AmgPreconditioner<T> {
1089    fn apply(&self, r: &Array1<T>) -> Array1<T> {
1090        if self.levels.is_empty() {
1091            return r.clone();
1092        }
1093
1094        let n = self.levels[0].num_dofs;
1095        if r.len() != n {
1096            return r.clone();
1097        }
1098
1099        let mut z = Array1::from_elem(n, T::zero());
1100
1101        match self.config.cycle {
1102            AmgCycle::VCycle => {
1103                self.v_cycle(0, &mut z, r);
1104            }
1105            AmgCycle::WCycle => {
1106                self.v_cycle(0, &mut z, r);
1107                self.v_cycle(0, &mut z, r);
1108            }
1109            AmgCycle::FCycle => {
1110                self.v_cycle(0, &mut z, r);
1111                let residual = r - &self.levels[0].matrix.matvec(&z);
1112                let mut correction = Array1::from_elem(n, T::zero());
1113                self.v_cycle(0, &mut correction, &residual);
1114                z = z + correction;
1115            }
1116        }
1117
1118        z
1119    }
1120}
1121
1122/// Diagnostic information about AMG setup
1123#[derive(Debug, Clone)]
1124pub struct AmgDiagnostics {
1125    /// Number of levels
1126    pub num_levels: usize,
1127    /// Grid complexity
1128    pub grid_complexity: f64,
1129    /// Operator complexity
1130    pub operator_complexity: f64,
1131    /// Setup time in milliseconds
1132    pub setup_time_ms: f64,
1133    /// DOFs per level
1134    pub level_dofs: Vec<usize>,
1135    /// NNZ per level
1136    pub level_nnz: Vec<usize>,
1137}
1138
1139impl<T: ComplexField> AmgPreconditioner<T> {
1140    /// Get diagnostic information
1141    pub fn diagnostics(&self) -> AmgDiagnostics {
1142        AmgDiagnostics {
1143            num_levels: self.levels.len(),
1144            grid_complexity: self.grid_complexity,
1145            operator_complexity: self.operator_complexity,
1146            setup_time_ms: self.setup_time_ms,
1147            level_dofs: self.levels.iter().map(|l| l.num_dofs).collect(),
1148            level_nnz: self.levels.iter().map(|l| l.matrix.nnz()).collect(),
1149        }
1150    }
1151}
1152
1153#[cfg(test)]
1154mod tests {
1155    use super::*;
1156    use num_complex::Complex64;
1157
1158    /// Create a simple 1D Laplacian matrix for testing
1159    fn create_1d_laplacian(n: usize) -> CsrMatrix<Complex64> {
1160        let mut triplets: Vec<(usize, usize, Complex64)> = Vec::new();
1161
1162        for i in 0..n {
1163            triplets.push((i, i, Complex64::new(2.0, 0.0)));
1164            if i > 0 {
1165                triplets.push((i, i - 1, Complex64::new(-1.0, 0.0)));
1166            }
1167            if i < n - 1 {
1168                triplets.push((i, i + 1, Complex64::new(-1.0, 0.0)));
1169            }
1170        }
1171
1172        CsrMatrix::from_triplets(n, n, triplets)
1173    }
1174
1175    #[test]
1176    fn test_amg_creation() {
1177        let matrix = create_1d_laplacian(100);
1178        let config = AmgConfig::default();
1179
1180        let amg = AmgPreconditioner::from_csr(&matrix, config);
1181
1182        assert!(amg.num_levels() >= 2);
1183        assert!(amg.grid_complexity() >= 1.0);
1184        assert!(amg.operator_complexity() >= 1.0);
1185    }
1186
1187    #[test]
1188    fn test_amg_apply() {
1189        let matrix = create_1d_laplacian(50);
1190        let config = AmgConfig::default();
1191        let amg = AmgPreconditioner::from_csr(&matrix, config);
1192
1193        let r = Array1::from_vec((0..50).map(|i| Complex64::new(i as f64, 0.0)).collect());
1194
1195        let z = amg.apply(&r);
1196
1197        assert_eq!(z.len(), r.len());
1198
1199        let diff: f64 = (&z - &r).iter().map(|x| x.norm()).sum();
1200        assert!(diff > 1e-10, "Preconditioner should modify the vector");
1201    }
1202
1203    #[test]
1204    fn test_amg_pmis_coarsening() {
1205        let matrix = create_1d_laplacian(100);
1206        let config = AmgConfig {
1207            coarsening: AmgCoarsening::Pmis,
1208            ..Default::default()
1209        };
1210
1211        let amg = AmgPreconditioner::from_csr(&matrix, config);
1212        assert!(amg.num_levels() >= 2);
1213    }
1214
1215    #[test]
1216    fn test_amg_different_smoothers() {
1217        let matrix = create_1d_laplacian(50);
1218        let r = Array1::from_vec((0..50).map(|i| Complex64::new(i as f64, 0.0)).collect());
1219
1220        for smoother in [
1221            AmgSmoother::Jacobi,
1222            AmgSmoother::L1Jacobi,
1223            AmgSmoother::SymmetricGaussSeidel,
1224        ] {
1225            let config = AmgConfig {
1226                smoother,
1227                ..Default::default()
1228            };
1229            let amg = AmgPreconditioner::from_csr(&matrix, config);
1230
1231            let z = amg.apply(&r);
1232            assert_eq!(z.len(), r.len());
1233        }
1234    }
1235
1236    #[test]
1237    fn test_amg_reduces_residual() {
1238        let n = 64;
1239        let matrix = create_1d_laplacian(n);
1240        let config = AmgConfig::default();
1241        let amg = AmgPreconditioner::from_csr(&matrix, config);
1242
1243        let b = Array1::from_vec(
1244            (0..n)
1245                .map(|i| Complex64::new((i as f64).sin(), 0.0))
1246                .collect(),
1247        );
1248
1249        let mut x = Array1::from_elem(n, Complex64::new(0.0, 0.0));
1250
1251        let r0 = &b - &matrix.matvec(&x);
1252        let norm_r0: f64 = r0.iter().map(|v| v.norm_sqr()).sum::<f64>().sqrt();
1253
1254        for _ in 0..10 {
1255            let r = &b - &matrix.matvec(&x);
1256            let z = amg.apply(&r);
1257            x = x + z;
1258        }
1259
1260        let rf = &b - &matrix.matvec(&x);
1261        let norm_rf: f64 = rf.iter().map(|v| v.norm_sqr()).sum::<f64>().sqrt();
1262
1263        assert!(
1264            norm_rf < norm_r0 * 0.1,
1265            "AMG should significantly reduce residual: {} -> {}",
1266            norm_r0,
1267            norm_rf
1268        );
1269    }
1270
1271    #[test]
1272    fn test_diagnostics() {
1273        let matrix = create_1d_laplacian(100);
1274        let amg = AmgPreconditioner::from_csr(&matrix, AmgConfig::default());
1275
1276        let diag = amg.diagnostics();
1277
1278        assert!(diag.num_levels >= 2);
1279        assert_eq!(diag.level_dofs.len(), diag.num_levels);
1280        assert_eq!(diag.level_nnz.len(), diag.num_levels);
1281        assert!(diag.grid_complexity >= 1.0);
1282        assert!(diag.setup_time_ms >= 0.0);
1283    }
1284}