pub fn gmres_with_truncation<T, F, Tr>(
apply_a: F,
b: &T,
x0: &T,
options: &GmresOptions,
truncate: Tr,
) -> Result<GmresResult<T>>Expand description
Solve A x = b using GMRES with optional truncation after each iteration.
This is an extension of gmres that allows truncating Krylov basis vectors
to control bond dimension growth in tensor network representations.
§Type Parameters
T- A tensor type implementingTensorLikeF- A function that applies the linear operator:F(x) = A xTr- A function that truncates a tensor in-place:Tr(&mut x)
§Arguments
apply_a- Function that applies the linear operator A to a tensorb- Right-hand side tensorx0- Initial guessoptions- Solver optionstruncate- Function that truncates a tensor to control bond dimension
§Note
Truncation is applied after each Gram-Schmidt orthogonalization step and after the final solution update. This helps control the bond dimension growth that would otherwise occur in MPS/MPO representations.
§Examples
Solve 2x = b with a no-op truncation function:
use tensor4all_core::{DynIndex, TensorDynLen, TensorLike, AnyScalar};
use tensor4all_core::krylov::{gmres_with_truncation, GmresOptions};
let i = DynIndex::new_dyn(2);
let b = TensorDynLen::from_dense(vec![i.clone()], vec![4.0, 6.0]).unwrap();
let x0 = TensorDynLen::from_dense(vec![i.clone()], vec![0.0, 0.0]).unwrap();
// Operator A = 2*I (scales input by 2)
let apply_a = |x: &TensorDynLen| x.scale(AnyScalar::new_real(2.0));
// No-op truncation
let truncate = |_x: &mut TensorDynLen| Ok(());
let result = gmres_with_truncation(apply_a, &b, &x0, &GmresOptions::default(), truncate).unwrap();
assert!(result.converged);
// Solution should be [2.0, 3.0]
let expected = TensorDynLen::from_dense(vec![i], vec![2.0, 3.0]).unwrap();
assert!(result.solution.sub(&expected).unwrap().maxabs() < 1e-8);