Skip to main content

crossinterpolate1

Function crossinterpolate1 

Source
pub fn crossinterpolate1<T, F>(
    f: F,
    local_dims: Vec<usize>,
    first_pivot: MultiIndex,
    options: TCI1Options,
) -> Result<(TensorCI1<T>, Vec<usize>, Vec<f64>)>
where T: Scalar + TTScalar + Default, F: Fn(&MultiIndex) -> T,
Expand description

Approximate a function as a tensor train using the TCI1 algorithm (legacy).

Prefer crossinterpolate2 for new code.

§Arguments

  • f – Function to interpolate. Takes &MultiIndex (0-indexed) and returns a scalar.
  • local_dims – Number of values each index can take.
  • first_pivot – Starting multi-index. Must have the same length as local_dims and the function value must be non-zero.
  • options – Algorithm configuration; see TCI1Options.

§Returns

A tuple (tci, ranks, errors):

  • tci: TensorCI1<T> – The interpolation state.
  • ranks: Vec<usize> – Bond dimension after each sweep.
  • errors: Vec<f64> – Error estimate after each sweep.

§Examples

use tensor4all_tensorci::{crossinterpolate1, TCI1Options};
use tensor4all_simplett::AbstractTensorTrain;

let f = |idx: &Vec<usize>| (idx[0] + idx[1] + 1) as f64;
let local_dims = vec![4, 4];
let first_pivot = vec![3, 3];

let (tci, _ranks, _errors) = crossinterpolate1::<f64, _>(
    f,
    local_dims,
    first_pivot,
    TCI1Options { tolerance: 1e-10, ..TCI1Options::default() },
).unwrap();

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