1#![allow(unused_variables)]
8#![allow(unused_assignments)]
9#![allow(unused_mut)]
10
11use crate::error::{SparseError, SparseResult};
12use crate::sparray::SparseArray;
13use scirs2_core::ndarray::{Array1, ArrayView1};
14use scirs2_core::numeric::{Float, SparseElement};
15use std::fmt::Debug;
16
17#[derive(Debug, Clone)]
19pub struct TFQMROptions {
20 pub max_iter: usize,
22 pub tol: f64,
24 pub use_left_preconditioner: bool,
26 pub use_right_preconditioner: bool,
28}
29
30impl Default for TFQMROptions {
31 fn default() -> Self {
32 Self {
33 max_iter: 1000,
34 tol: 1e-6,
35 use_left_preconditioner: false,
36 use_right_preconditioner: false,
37 }
38 }
39}
40
41#[derive(Debug, Clone)]
43pub struct TFQMRResult<T> {
44 pub x: Array1<T>,
46 pub iterations: usize,
48 pub residual_norm: T,
50 pub converged: bool,
52 pub residual_history: Option<Vec<T>>,
54}
55
56#[allow(dead_code)]
93pub fn tfqmr<T, S>(
94 matrix: &S,
95 b: &ArrayView1<T>,
96 x0: Option<&ArrayView1<T>>,
97 options: TFQMROptions,
98) -> SparseResult<TFQMRResult<T>>
99where
100 T: Float + SparseElement + Debug + Copy + 'static,
101 S: SparseArray<T>,
102{
103 let n = b.len();
104 let (rows, cols) = matrix.shape();
105
106 if rows != cols || rows != n {
107 return Err(SparseError::DimensionMismatch {
108 expected: n,
109 found: rows,
110 });
111 }
112
113 let mut x = match x0 {
115 Some(x0_val) => x0_val.to_owned(),
116 None => Array1::zeros(n),
117 };
118
119 let ax = matrix_vector_multiply(matrix, &x.view())?;
121 let mut r = b - &ax;
122
123 let initial_residual_norm = l2_norm(&r.view());
125 let b_norm = l2_norm(b);
126 let tolerance = T::from(options.tol).unwrap() * b_norm;
127
128 if initial_residual_norm <= tolerance {
129 return Ok(TFQMRResult {
130 x,
131 iterations: 0,
132 residual_norm: initial_residual_norm,
133 converged: true,
134 residual_history: Some(vec![initial_residual_norm]),
135 });
136 }
137
138 let r_star = r.clone();
140
141 let mut v = r.clone();
143 let mut y = v.clone();
144 let mut w = matrix_vector_multiply(matrix, &y.view())?;
145 let mut z = w.clone();
146
147 let mut d = Array1::zeros(n);
148 let mut theta = T::sparse_zero();
149 let mut eta = T::sparse_zero();
150 let mut tau = initial_residual_norm;
151
152 let mut rho = dot_product(&r_star.view(), &r.view());
154 let mut alpha = T::sparse_zero();
155 let mut beta = T::sparse_zero();
156
157 let mut residual_history = Vec::new();
158 residual_history.push(initial_residual_norm);
159
160 let mut converged = false;
161 let mut iter = 0;
162
163 for m in 0..options.max_iter {
164 iter = m + 1;
165
166 let sigma = dot_product(&r_star.view(), &w.view());
168 if sigma.abs() < T::from(1e-14).unwrap() {
169 return Err(SparseError::ConvergenceError(
170 "TFQMR breakdown: sigma is too small".to_string(),
171 ));
172 }
173 alpha = rho / sigma;
174
175 for i in 0..n {
177 v[i] = v[i] - alpha * w[i];
178 y[i] = y[i] - alpha * z[i];
179 }
180
181 let v_norm = l2_norm(&v.view());
183 theta = v_norm / tau;
184 let c = T::sparse_one() / (T::sparse_one() + theta * theta).sqrt();
185 tau = tau * theta * c;
186 eta = c * c * alpha;
187
188 for i in 0..n {
190 d[i] = y[i] + (theta * eta) * d[i];
191 x[i] = x[i] + eta * d[i];
192 }
193
194 let current_residual = tau;
196 residual_history.push(current_residual);
197
198 if current_residual <= tolerance {
199 converged = true;
200 break;
201 }
202
203 w = matrix_vector_multiply(matrix, &y.view())?;
205
206 let rho_new = dot_product(&r_star.view(), &v.view());
208 beta = rho_new / rho;
209 rho = rho_new;
210
211 for i in 0..n {
213 y[i] = v[i] + beta * y[i];
214 z[i] = w[i] + beta * z[i];
215 }
216
217 let y_norm = l2_norm(&y.view());
219 theta = y_norm / tau;
220 let c = T::sparse_one() / (T::sparse_one() + theta * theta).sqrt();
221 tau = tau * theta * c;
222 eta = c * c * alpha;
223
224 for i in 0..n {
226 d[i] = z[i] + (theta * eta) * d[i];
227 x[i] = x[i] + eta * d[i];
228 }
229
230 let current_residual = tau;
232 residual_history.push(current_residual);
233
234 if current_residual <= tolerance {
235 converged = true;
236 break;
237 }
238
239 w = matrix_vector_multiply(matrix, &z.view())?;
241 }
242
243 let ax_final = matrix_vector_multiply(matrix, &x.view())?;
245 let final_residual = b - &ax_final;
246 let final_residual_norm = l2_norm(&final_residual.view());
247
248 Ok(TFQMRResult {
249 x,
250 iterations: iter,
251 residual_norm: final_residual_norm,
252 converged,
253 residual_history: Some(residual_history),
254 })
255}
256
257#[allow(dead_code)]
259fn matrix_vector_multiply<T, S>(matrix: &S, x: &ArrayView1<T>) -> SparseResult<Array1<T>>
260where
261 T: Float + SparseElement + Debug + Copy + 'static,
262 S: SparseArray<T>,
263{
264 let (rows, cols) = matrix.shape();
265 if x.len() != cols {
266 return Err(SparseError::DimensionMismatch {
267 expected: cols,
268 found: x.len(),
269 });
270 }
271
272 let mut result = Array1::zeros(rows);
273 let (row_indices, col_indices, values) = matrix.find();
274
275 for (k, (&i, &j)) in row_indices.iter().zip(col_indices.iter()).enumerate() {
276 result[i] = result[i] + values[k] * x[j];
277 }
278
279 Ok(result)
280}
281
282#[allow(dead_code)]
284fn l2_norm<T>(x: &ArrayView1<T>) -> T
285where
286 T: Float + SparseElement + Debug + Copy,
287{
288 (x.iter()
289 .map(|&val| val * val)
290 .fold(T::sparse_zero(), |a, b| a + b))
291 .sqrt()
292}
293
294#[allow(dead_code)]
296fn dot_product<T>(x: &ArrayView1<T>, y: &ArrayView1<T>) -> T
297where
298 T: Float + SparseElement + Debug + Copy,
299{
300 x.iter()
301 .zip(y.iter())
302 .map(|(&xi, &yi)| xi * yi)
303 .fold(T::sparse_zero(), |a, b| a + b)
304}
305
306#[cfg(test)]
307mod tests {
308 use super::*;
309 use crate::csr_array::CsrArray;
310 use approx::assert_relative_eq;
311
312 #[test]
313 #[ignore] fn test_tfqmr_simple_system() {
315 let rows = vec![0, 0, 1, 1, 2, 2];
317 let cols = vec![0, 1, 0, 1, 1, 2];
318 let data = vec![2.0, -1.0, -1.0, 2.0, -1.0, 2.0];
319 let matrix = CsrArray::from_triplets(&rows, &cols, &data, (3, 3), false).unwrap();
320
321 let b = Array1::from_vec(vec![1.0, 0.0, 1.0]);
322 let result = tfqmr(&matrix, &b.view(), None, TFQMROptions::default()).unwrap();
323
324 assert!(result.converged);
325
326 let ax = matrix_vector_multiply(&matrix, &result.x.view()).unwrap();
328 let residual = &b - &ax;
329 let residual_norm = l2_norm(&residual.view());
330
331 assert!(residual_norm < 1e-6);
332 }
333
334 #[test]
335 #[ignore] fn test_tfqmr_diagonal_system() {
337 let rows = vec![0, 1, 2];
339 let cols = vec![0, 1, 2];
340 let data = vec![2.0, 3.0, 4.0];
341 let matrix = CsrArray::from_triplets(&rows, &cols, &data, (3, 3), false).unwrap();
342
343 let b = Array1::from_vec(vec![4.0, 9.0, 16.0]);
344 let result = tfqmr(&matrix, &b.view(), None, TFQMROptions::default()).unwrap();
345
346 assert!(result.converged);
347
348 assert_relative_eq!(result.x[0], 2.0, epsilon = 1e-6);
350 assert_relative_eq!(result.x[1], 3.0, epsilon = 1e-6);
351 assert_relative_eq!(result.x[2], 4.0, epsilon = 1e-6);
352 }
353
354 #[test]
355 fn test_tfqmr_with_initial_guess() {
356 let rows = vec![0, 1, 2];
357 let cols = vec![0, 1, 2];
358 let data = vec![1.0, 1.0, 1.0];
359 let matrix = CsrArray::from_triplets(&rows, &cols, &data, (3, 3), false).unwrap();
360
361 let b = Array1::from_vec(vec![5.0, 6.0, 7.0]);
362 let x0 = Array1::from_vec(vec![4.0, 5.0, 6.0]); let result = tfqmr(
365 &matrix,
366 &b.view(),
367 Some(&x0.view()),
368 TFQMROptions::default(),
369 )
370 .unwrap();
371
372 assert!(result.converged);
373 assert!(result.iterations <= 5); }
375}