pub fn restart_gmres_with_truncation<T, F, Tr>(
apply_a: F,
b: &T,
x0: Option<&T>,
options: &RestartGmresOptions,
truncate: Tr,
) -> Result<RestartGmresResult<T>>Expand description
Solve A x = b using restarted GMRES with truncation.
This wraps gmres_with_truncation with an outer loop that recomputes the true residual
at each restart. This is particularly useful for MPS/MPO computations where truncation
can cause the inner GMRES residual to be inaccurate.
§Algorithm
for outer_iter in 0..max_outer_iters:
r = b - A*x0 // Compute true residual
r = truncate(r)
if ||r|| / ||b|| < rtol:
return x0 // Converged
x' = gmres_with_truncation(A, r, 0, inner_options, truncate)
x0 = truncate(x0 + x')§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 guess (if None, starts from zero)options- Solver optionstruncate- Function that truncates a tensor to control bond dimension
§Returns
A RestartGmresResult containing the solution and convergence information.
§Examples
Solve 5x = b with no truncation:
use tensor4all_core::{DynIndex, TensorDynLen, TensorLike, AnyScalar};
use tensor4all_core::krylov::{restart_gmres_with_truncation, RestartGmresOptions};
let i = DynIndex::new_dyn(3);
let b = TensorDynLen::from_dense(vec![i.clone()], vec![5.0, 10.0, 15.0]).unwrap();
let apply_a = |x: &TensorDynLen| x.scale(AnyScalar::new_real(5.0));
let truncate = |_x: &mut TensorDynLen| Ok(());
let result = restart_gmres_with_truncation(
apply_a, &b, None, &RestartGmresOptions::default(), truncate,
).unwrap();
assert!(result.converged);
let expected = TensorDynLen::from_dense(vec![i], vec![1.0, 2.0, 3.0]).unwrap();
assert!(result.solution.sub(&expected).unwrap().maxabs() < 1e-8);