In this library, we provide methods for inverting the slicing operator
Illustration of a forward pass through the slicing operator
Given
This project requires the following Python (version 3.9.21) packages:
/dev requires:
/test additionally requires:
- PyKeOps version 2.3
- simple_torch_NFFT
- tqdm version 4.66.5
- matplotlib version 3.9.2
To install package run:
pip install git+https://github.com/Nicolaj-Rux/slicing-inversionCreate test file test.py with content:
from slicing_inversion.slicing import slicing
S = slicing(F= lambda s: s**2, d=10)
print(S.evaluate())When running python test.py should return
(2.5949782411771594e-06, 6.033632278442383, 14.8980712890625)Import the class slicing as from slicing_inversion.slicing import slicing.
Prepare a torch function F from def F(s): return torch.exp(-s**2) and d=100.
Then S = slicing(F=F, d=d) creates an instance of that class. To compute S.get_matrix().
This step computes the inversion matrix which is independent of S.get_range_coef() to expand method, see section Hyperparamters) and then recover the S.get_domain_coef().
Afterwards the coeficients are stored under S.a.
The complete flow would look like
from slicing_inversion.slicing import slicing
import torch
def F(s): return torch.exp(-s**2)
d = 100
S = slicing(F= lambda s: s**2, d=10)
S.get_matrix()
S.get_range_coef()
S.get_domain_coef()
print(S.a) # Fourier coefficients of fThe class comes with several default hyperparamters
method: int = 0,
T: float = 1.,
L: int = 2**10,
K: int = 2**8,
tau: float = 1e-6,
bs: int = 2**10The first hyperparameter, method, can be chosen as method = 2, a larger regularization around bs can be reduced when using large
The density
Given a function
If
The paper Fast Kernel Summation in High Dimensions via Slicing and Fourier Transforms introduces a method to reduce this quadratic time complexity to linear by slicing the function
Choosing
More information about good choices of slicing directions can be found in FAST SUMMATION OF RADIAL KERNELS VIA QMC SLICING.
The
At first, this seems more complicated. The advantage is that the one-dimensional kernel sums
can be computed in linear time using Fourier methods. Expand
Then, for fixed
Now,
The multiplication of
Since this must be done for each slicing direction, the total time complexity is
Because
As a prerequisite to applying this slicing method, the slicing transformation must be inverted, i.e., finding
Interestingly, (
For some kernels such as the Riesz, Gaussian, or Laplace kernel, the sliced kernel
In this project, we focus on numerical methods to compute
This library was written by Nicolaj Rux in the context of fast kernel summations via slicing. If you find it usefull, please consider to cite
@misc{RHN2025,
title={Numerical Methods for Kernel Slicing},
author={Nicolaj Rux and Johannes Hertrich and Sebastian Neumayer},
year={2025},
eprint={2510.11478},
archivePrefix={arXiv},
primaryClass={math.NA},
url={https://arxiv.org/abs/2510.11478},
}
or
@inproceedings{HJQ2025,
title={Fast Summation of Radial Kernels via {QMC} Slicing},
author={Johannes Hertrich and Tim Jahn and Michael Quellmalz},
booktitle={The Thirteenth International Conference on Learning Representations},
year={2025},
url={https://openreview.net/forum?id=iNmVX9lx9l}
}
or
@article{H2024,
title={Fast Kernel Summation in High Dimensions via Slicing and {F}ourier transforms},
author={Hertrich, Johannes},
journal={SIAM Journal on Mathematics of Data Science},
volume={6},
number={4},
pages={1109--1137},
year={2024}
}

