Skip to main content

crossinterpolate2

Function crossinterpolate2 

Source
pub fn crossinterpolate2<T, F, B>(
    f: F,
    batched_f: Option<B>,
    local_dims: Vec<usize>,
    initial_pivots: Vec<MultiIndex>,
    options: TCI2Options,
) -> Result<(TensorCI2<T>, Vec<usize>, Vec<f64>)>
where T: Scalar + TTScalar + Default + MatrixLuciScalar, DenseFaerLuKernel: PivotKernel<T>, LazyBlockRookKernel: PivotKernel<T>, F: Fn(&MultiIndex) -> T, B: Fn(&[MultiIndex]) -> Vec<T>,
Expand description

Approximate a function as a tensor train using the TCI2 algorithm.

This is the main entry point for tensor cross interpolation. It performs alternating two-site sweeps with global pivot search until convergence, then cleans up with a final one-site sweep.

§Arguments

  • f – Function to interpolate. Takes a multi-index &Vec<usize> where each element is in 0..local_dims[i] (0-indexed) and returns a scalar.
  • batched_f – Optional batch evaluation function for efficiency. Takes &[Vec<usize>] and returns Vec<T>. Pass None to use element-wise evaluation only.
  • local_dims – Number of values each index can take. Must have at least 2 elements (TCI requires at least 2 sites).
  • initial_pivots – Starting multi-indices for the algorithm. At least one pivot must have a non-zero function value. Choose pivots where |f| is large for best convergence. If empty, defaults to the all-zeros index.
  • options – Algorithm configuration; see TCI2Options.

§Returns

A tuple (tci, ranks, errors):

  • tci: TensorCI2<T> – The interpolation state. Call to_tensor_train to get a TensorTrain.
  • ranks: Vec<usize> – Bond dimension after each half-sweep.
  • errors: Vec<f64> – Normalized error estimate after each half-sweep. The last value should be below tolerance if the algorithm converged.

§Errors

Returns TCIError::DimensionMismatch if local_dims has fewer than 2 elements, or TCIError::InvalidPivot if all initial pivots evaluate to zero.

§Examples

Basic usage: interpolate f(i, j) = i + j + 1 on a 4x4 grid.

use tensor4all_tensorci::{crossinterpolate2, TCI2Options};
use tensor4all_simplett::AbstractTensorTrain;

let f = |idx: &Vec<usize>| (idx[0] + idx[1] + 1) as f64;
let local_dims = vec![4, 4];
let initial_pivots = vec![vec![3, 3]]; // pick where |f| is large

let (tci, _ranks, errors) =
    crossinterpolate2::<f64, _, fn(&[Vec<usize>]) -> Vec<f64>>(
        f,
        None,
        local_dims,
        initial_pivots,
        TCI2Options {
            tolerance: 1e-10,
            seed: Some(42),
            ..TCI2Options::default()
        },
    )
    .unwrap();

// Verify convergence
assert!(*errors.last().unwrap() < 1e-10);

// Evaluate through the tensor train
let tt = tci.to_tensor_train().unwrap();
let val = tt.evaluate(&[2, 3]).unwrap();
assert!((val - 6.0).abs() < 1e-10); // f(2,3) = 2+3+1 = 6

// Bond dimensions are available
assert!(!tci.link_dims().is_empty());