scirs2_sparse/linalg/
preconditioners.rs1use crate::csr::CsrMatrix;
4use crate::error::{SparseError, SparseResult};
5use crate::linalg::interface::LinearOperator;
6use scirs2_core::numeric::{Float, NumAssign, SparseElement};
7use std::fmt::Debug;
8use std::iter::Sum;
9
10pub struct JacobiPreconditioner<F> {
15 inv_diagonal: Vec<F>,
16 size: usize,
17}
18
19impl<F: Float + SparseElement> JacobiPreconditioner<F> {
20 pub fn new(matrix: &CsrMatrix<F>) -> SparseResult<Self>
22 where
23 F: Debug,
24 {
25 let n = matrix.rows();
26 let mut inv_diagonal = vec![F::sparse_zero(); n];
27
28 for (i, diag_inv) in inv_diagonal.iter_mut().enumerate().take(n) {
30 let diag_val = matrix.get(i, i);
31 if diag_val.abs() < F::epsilon() {
32 return Err(SparseError::ValueError(format!(
33 "Zero diagonal element at position {i}"
34 )));
35 }
36 *diag_inv = F::sparse_one() / diag_val;
37 }
38
39 Ok(Self {
40 inv_diagonal,
41 size: n,
42 })
43 }
44
45 pub fn from_diagonal(diagonal: Vec<F>) -> SparseResult<Self> {
47 let size = diagonal.len();
48 let mut inv_diagonal = vec![F::sparse_zero(); size];
49
50 for (i, &d) in diagonal.iter().enumerate() {
51 if d.abs() < F::epsilon() {
52 return Err(SparseError::ValueError(format!(
53 "Zero _diagonal element at position {i}"
54 )));
55 }
56 inv_diagonal[i] = F::sparse_one() / d;
57 }
58
59 Ok(Self { inv_diagonal, size })
60 }
61}
62
63impl<F: Float + NumAssign + SparseElement> LinearOperator<F> for JacobiPreconditioner<F> {
64 fn shape(&self) -> (usize, usize) {
65 (self.size, self.size)
66 }
67
68 fn matvec(&self, x: &[F]) -> SparseResult<Vec<F>> {
69 if x.len() != self.size {
70 return Err(SparseError::DimensionMismatch {
71 expected: self.size,
72 found: x.len(),
73 });
74 }
75
76 Ok(x.iter()
77 .zip(&self.inv_diagonal)
78 .map(|(&xi, &di)| xi * di)
79 .collect())
80 }
81
82 fn has_adjoint(&self) -> bool {
83 true
84 }
85
86 fn rmatvec(&self, x: &[F]) -> SparseResult<Vec<F>> {
87 self.matvec(x)
89 }
90}
91
92pub struct SSORPreconditioner<F> {
97 matrix: CsrMatrix<F>,
98 omega: F,
99 diagonal: Vec<F>,
100}
101
102impl<F: Float + NumAssign + Sum + Debug + SparseElement + 'static> SSORPreconditioner<F> {
103 pub fn new(matrix: CsrMatrix<F>, omega: F) -> SparseResult<Self> {
105 if omega <= F::sparse_zero()
106 || omega >= F::from(2.0).expect("Failed to convert constant to float")
107 {
108 return Err(SparseError::ValueError(
109 "Relaxation parameter omega must be in (0, 2)".to_string(),
110 ));
111 }
112
113 let n = matrix.rows();
114 let mut diagonal = vec![F::sparse_zero(); n];
115
116 for (i, diag) in diagonal.iter_mut().enumerate().take(n) {
118 *diag = matrix.get(i, i);
119 if diag.abs() < F::epsilon() {
120 return Err(SparseError::ValueError(format!(
121 "Zero diagonal element at position {i}"
122 )));
123 }
124 }
125
126 Ok(Self {
127 matrix,
128 omega,
129 diagonal,
130 })
131 }
132}
133
134impl<F: Float + NumAssign + Sum + Debug + SparseElement + 'static> LinearOperator<F>
135 for SSORPreconditioner<F>
136{
137 fn shape(&self) -> (usize, usize) {
138 self.matrix.shape()
139 }
140
141 fn matvec(&self, x: &[F]) -> SparseResult<Vec<F>> {
142 let n = self.matrix.rows();
143 if x.len() != n {
144 return Err(SparseError::DimensionMismatch {
145 expected: n,
146 found: x.len(),
147 });
148 }
149
150 let mut y = vec![F::sparse_zero(); n];
152 for i in 0..n {
153 let mut sum = x[i];
154 let row_range = self.matrix.row_range(i);
155 let row_indices = &self.matrix.indices[row_range.clone()];
156 let row_data = &self.matrix.data[row_range];
157
158 for (idx, &j) in row_indices.iter().enumerate() {
159 if j < i {
160 sum -= self.omega * row_data[idx] * y[j];
161 }
162 }
163 y[i] = sum / self.diagonal[i];
164 }
165
166 let mut z = vec![F::sparse_zero(); n];
168 for i in 0..n {
169 z[i] = y[i] * self.diagonal[i];
170 }
171
172 let mut w = vec![F::sparse_zero(); n];
174 for i in (0..n).rev() {
175 let mut sum = z[i];
176 let row_range = self.matrix.row_range(i);
177 let row_indices = &self.matrix.indices[row_range.clone()];
178 let row_data = &self.matrix.data[row_range];
179
180 for (idx, &j) in row_indices.iter().enumerate() {
181 if j > i {
182 sum -= self.omega * row_data[idx] * w[j];
183 }
184 }
185 w[i] = sum / self.diagonal[i];
186 }
187
188 Ok(w)
189 }
190
191 fn has_adjoint(&self) -> bool {
192 false }
194}
195
196pub struct ILU0Preconditioner<F> {
201 l_data: Vec<F>,
202 u_data: Vec<F>,
203 l_indices: Vec<usize>,
204 u_indices: Vec<usize>,
205 l_indptr: Vec<usize>,
206 u_indptr: Vec<usize>,
207 n: usize,
208}
209
210impl<F: Float + NumAssign + Sum + Debug + SparseElement + 'static> ILU0Preconditioner<F> {
211 pub fn new(matrix: &CsrMatrix<F>) -> SparseResult<Self> {
213 let n = matrix.rows();
214
215 let mut data = matrix.data.clone();
217 let indices = matrix.indices.clone();
218 let indptr = matrix.indptr.clone();
219
220 for k in 0..n {
222 let k_diag_idx = find_diagonal_index(&indices, &indptr, k)?;
223 let k_diag = data[k_diag_idx];
224
225 if k_diag.abs() < F::epsilon() {
226 return Err(SparseError::ValueError(format!(
227 "Zero pivot at position {k}"
228 )));
229 }
230
231 for i in (k + 1)..n {
233 let row_start = indptr[i];
234 let row_end = indptr[i + 1];
235
236 let mut k_pos = None;
238 for (&col, j) in indices[row_start..row_end].iter().zip(row_start..row_end) {
239 if col == k {
240 k_pos = Some(j);
241 break;
242 }
243 if col > k {
244 break;
245 }
246 }
247
248 if let Some(ki_idx) = k_pos {
249 let mult = data[ki_idx] / k_diag;
251 data[ki_idx] = mult;
252
253 let k_row_start = indptr[k];
255 let k_row_end = indptr[k + 1];
256
257 for kj_idx in k_row_start..k_row_end {
258 let j = indices[kj_idx];
259 if j <= k {
260 continue;
261 }
262
263 for ij_idx in row_start..row_end {
265 if indices[ij_idx] == j {
266 data[ij_idx] = data[ij_idx] - mult * data[kj_idx];
267 break;
268 }
269 }
270 }
271 }
272 }
273 }
274
275 let mut l_data = Vec::new();
277 let mut u_data = Vec::new();
278 let mut l_indices = Vec::new();
279 let mut u_indices = Vec::new();
280 let mut l_indptr = vec![0];
281 let mut u_indptr = vec![0];
282
283 for i in 0..n {
284 let row_start = indptr[i];
285 let row_end = indptr[i + 1];
286
287 for j in row_start..row_end {
288 let col = indices[j];
289 let val = data[j];
290
291 match col.cmp(&i) {
292 std::cmp::Ordering::Less => {
293 l_indices.push(col);
295 l_data.push(val);
296 }
297 std::cmp::Ordering::Equal => {
298 u_indices.push(col);
300 u_data.push(val);
301 }
302 std::cmp::Ordering::Greater => {
303 u_indices.push(col);
305 u_data.push(val);
306 }
307 }
308 }
309
310 l_indptr.push(l_indices.len());
311 u_indptr.push(u_indices.len());
312 }
313
314 Ok(Self {
315 l_data,
316 u_data,
317 l_indices,
318 u_indices,
319 l_indptr,
320 u_indptr,
321 n,
322 })
323 }
324}
325
326impl<F: Float + NumAssign + Sum + Debug + SparseElement + 'static> LinearOperator<F>
327 for ILU0Preconditioner<F>
328{
329 fn shape(&self) -> (usize, usize) {
330 (self.n, self.n)
331 }
332
333 fn matvec(&self, x: &[F]) -> SparseResult<Vec<F>> {
334 if x.len() != self.n {
335 return Err(SparseError::DimensionMismatch {
336 expected: self.n,
337 found: x.len(),
338 });
339 }
340
341 let mut y = vec![F::sparse_zero(); self.n];
343 for i in 0..self.n {
344 y[i] = x[i];
345 let row_start = self.l_indptr[i];
346 let row_end = self.l_indptr[i + 1];
347
348 for j in row_start..row_end {
349 let col = self.l_indices[j];
350 y[i] = y[i] - self.l_data[j] * y[col];
351 }
352 }
353
354 let mut z = vec![F::sparse_zero(); self.n];
356 for i in (0..self.n).rev() {
357 z[i] = y[i];
358 let row_start = self.u_indptr[i];
359 let row_end = self.u_indptr[i + 1];
360
361 let mut diag_val = F::sparse_one();
362 for j in row_start..row_end {
363 let col = self.u_indices[j];
364 match col.cmp(&i) {
365 std::cmp::Ordering::Equal => diag_val = self.u_data[j],
366 std::cmp::Ordering::Greater => z[i] = z[i] - self.u_data[j] * z[col],
367 std::cmp::Ordering::Less => {}
368 }
369 }
370
371 z[i] /= diag_val;
372 }
373
374 Ok(z)
375 }
376
377 fn has_adjoint(&self) -> bool {
378 false }
380}
381
382#[allow(dead_code)]
384fn find_diagonal_index(indices: &[usize], indptr: &[usize], row: usize) -> SparseResult<usize> {
385 let row_start = indptr[row];
386 let row_end = indptr[row + 1];
387
388 for (idx, &col) in indices[row_start..row_end]
389 .iter()
390 .enumerate()
391 .map(|(i, col)| (i + row_start, col))
392 {
393 if col == row {
394 return Ok(idx);
395 }
396 }
397
398 Err(SparseError::ValueError(format!(
399 "Missing diagonal element at row {row}"
400 )))
401}
402
403#[cfg(test)]
404mod tests {
405 use super::*;
406
407 #[test]
408 fn test_jacobi_preconditioner() {
409 let diagonal = vec![2.0, 3.0, 4.0];
411 let precond = JacobiPreconditioner::from_diagonal(diagonal).expect("Operation failed");
412
413 let x = vec![2.0, 6.0, 12.0];
415 let result = precond.matvec(&x).expect("Operation failed");
416
417 assert!((result[0] - 1.0).abs() < 1e-10);
419 assert!((result[1] - 2.0).abs() < 1e-10);
420 assert!((result[2] - 3.0).abs() < 1e-10);
421 }
422}