Module anova_func: construct TT-tensor of interpolation coefs by TT-ANOVA

Package teneva, module anova_func: Functional ANOVA in the TT-format.

This module contains the function “anova_func” which computes the TT-approximation of Chebyshev interpolation coefficients’ tensor, using given train dataset with random points.




teneva.anova_func.anova_func(X_trn, y_trn, n, a=-1.0, b=1.0, lamb=1e-07, e=1e-08)[source]

Build functional TT-tensor by TT-ANOVA from the given random points.

Parameters:
  • X_trn (np.ndarray) – train points in the form of array of the shape [samples, d].

  • y_trn (np.ndarray) – values of the tensor for train points X_trn in the form of array of the shape [samples].

  • n (int) – mode size of the tensor (i.e, the maximum power of the interpolation coefficient plus 1).

  • a (float, list, np.ndarray) – grid lower bounds for each dimension (list or np.ndarray of length d or float).

  • b (float, list, np.ndarray) – grid upper bounds for each dimension (list or np.ndarray of length d or float).

  • lamb (float) – regularization parameter.

  • e (float) – optional truncation accuracy for resulting TT-tensor.

Returns:

TT-tensor, which represents the interpolation coefficients in the TT-format.

Return type:

list

Note

A class “ANOVA_func” that represents a wider set of methods for working with this decomposition is also available; see “anova_func.py”.

Examples:

d = 5                           # Dimension of the function
a = [-5., -4., -3., -2., -1.]   # Lower bounds for spatial grid
b = [+6., +3., +3., +1., +2.]   # Upper bounds for spatial grid
n = 3                           # Shape of interpolation tensor
m     = 1.E+4  # Number of calls to target function
e     = 1.E-8  # Truncation accuracy

We set the target function (the function takes as input a set of multidimensional points X of the shape [samples, dimension]):

def func(X):
    return np.sum(X, axis=1)

We prepare train data from the random distribution:

X_trn = teneva.sample_rand_poi(a, b, m)
y_trn = func(X_trn)

We prepare test data from random points:

# Number of test points:
m_tst = int(1.E+4)

# Random points:
X_tst = teneva.sample_rand_poi(a, b, m_tst)

# Function values for the test points:
y_tst = func(X_tst)

We build the TT-tensor of interpolation coefficients:

t = tpc()
A = teneva.anova_func(X_trn, y_trn, n, a, b, 1.E-5, e=e)
t = tpc() - t

print(f'Build time     : {t:-10.2f}')

# >>> ----------------------------------------
# >>> Output:

# Build time     :       0.01
#

And now we can check the result:

# Compute approximation in train points:
y_our = teneva.func_get(X_trn, A, a, b)

# Accuracy of the result for train points:
e_trn = np.linalg.norm(y_our - y_trn)
e_trn /= np.linalg.norm(y_trn)

# Compute approximation in test points:
y_our = teneva.func_get(X_tst, A, a, b)

# Accuracy of the result for test points:
e_tst = np.linalg.norm(y_our - y_tst)
e_tst /= np.linalg.norm(y_tst)

print(f'Error on train : {e_trn:-10.2e}')
print(f'Error on test  : {e_tst:-10.2e}')

# >>> ----------------------------------------
# >>> Output:

# Error on train :   2.87e-02
# Error on test  :   2.78e-02
#