Skip to content

This issue was moved to a discussion.

You can continue the conversation there. Go to discussion →

New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Exploit structure of upsampled input data when doing FFTs #304

Closed
blackwer opened this issue Jun 28, 2023 · 6 comments
Closed

Exploit structure of upsampled input data when doing FFTs #304

blackwer opened this issue Jun 28, 2023 · 6 comments

Comments

@blackwer
Copy link
Member

blackwer commented Jun 28, 2023

There may be other ways than pruning to speed up the FFT. In a unifom-to-nonuniform 2D NUFFT with, say, an oversampling factor of 2, you only need to transform over the first axis for half of the array, since the other half of the array contains vectors that are identically zero. This saves 25% of FFT cost without having to use complicated hand-crafted algorithms. In 3D, you have to transform a quarter of the array along the first axis and half of the array along the second axis, saving an even larger fraction of FFT time

Originally posted by @mreineck in #293 (comment)

I've implemented a POC for the suggestion by @mreineck that shows typical speedups ~1.5 on 2D data vs the naive way with MKL. I'm not doing upsampling quite the way it's done in finufft, but the answer should correct up to a phase shift as far as I can tell. I don't have a lot of cycles to work on this right now, but this seems like a very easy way to get a huge performance bump in higher dimensions.

#include <cstring>
#include <cstdio>
#include <random>
#include <tuple>
#include <utility>

#include <fftw3.h>
#include <omp.h>

std::pair<fftw_complex *, fftw_complex *> initialize_arrays(int n_inner[2], int n_tot[2]) {
    std::random_device rand_dev;
    std::mt19937 generator(1);
    std::uniform_real_distribution<double> distr(-1.0, 1.0);

    int N = n_tot[0] * n_tot[1];
    fftw_complex *in = (fftw_complex *)fftw_malloc(sizeof(fftw_complex) * N);
    fftw_complex *out = (fftw_complex *)fftw_malloc(sizeof(fftw_complex) * N);

    std::memset(in, 0, n_tot[0] * n_tot[1] * sizeof(fftw_complex));
    for (int i = 0; i < n_inner[0]; ++i)
        for (int j = 0; j < n_inner[1]; ++j)
            for (int k = 0; k < 2; ++k)
                in[i * n_tot[0] + j][k] = distr(generator);

    std::memset(out, 0, n_tot[0] * n_tot[1] * sizeof(fftw_complex));
    return std::make_pair(in, out);
}

double eval_naive(int n[2], fftw_complex *in, fftw_complex *out) {
    int n_threads;
#pragma omp parallel
    n_threads = omp_get_num_threads();
    fftw_plan_with_nthreads(n_threads);
    fftw_plan p = fftw_plan_dft(2, n, in, out, FFTW_FORWARD, FFTW_MEASURE);
    double t_elapsed = omp_get_wtime();

    fftw_execute(p);
    fftw_destroy_plan(p);
    return omp_get_wtime() - t_elapsed;
}

double eval_clever(int n_actual[2], int n_upsampled[2], fftw_complex *in, fftw_complex *out) {
    int n_threads;
#pragma omp parallel
    n_threads = omp_get_num_threads();
    fftw_plan_with_nthreads(n_threads);

    fftw_plan p_cols;
    {
        int rank = 1;               // we're doing 1d transforms of columns first
        int n[] = {n_upsampled[0]}; // each column has n_upsampled[0] rows
        int howmany = n_actual[1];  // we only need to tranform the original columns
        int idist = 1, odist = 1;
        int istride = n_upsampled[1], ostride = n_upsampled[1];
        int *oembed = NULL, *iembed = NULL;
        p_cols = fftw_plan_many_dft(rank, n, howmany, in, iembed, istride, idist, out, oembed, ostride, odist,
                                    FFTW_FORWARD, FFTW_MEASURE);
    }
    fftw_plan p_rows;
    {
        int rank = 1;               // we're doing 1d transforms of rows now
        int n[] = {n_upsampled[1]}; // each row has n_upsampled[1] columns to transform
        int howmany = n_upsampled[0];  // we need to tranform ALL rows
        int idist = n_upsampled[1], odist = n_upsampled[1];
        int istride = 1, ostride = 1;
        int *oembed = NULL, *iembed = NULL;
        p_rows = fftw_plan_many_dft(rank, n, howmany, out, iembed, istride, idist, out, oembed, ostride, odist,
                                    FFTW_FORWARD, FFTW_MEASURE);
    }
    double t_elapsed = omp_get_wtime();

    fftw_execute(p_cols);
    fftw_execute(p_rows);


    fftw_destroy_plan(p_cols);
    fftw_destroy_plan(p_rows);
    return omp_get_wtime() - t_elapsed;
}

int main(int argc, char *argv[]) {
    constexpr int n_runs = 10;
    int n_per_dim = 500;
    if (argc > 1)
        n_per_dim = std::atof(argv[1]);
    constexpr int dim = 2;
    constexpr double upsamp = 2.0;

    int n_inner[dim];
    int n_upsampled[dim];
    for (int i = 0; i < dim; ++i) {
        n_inner[i] = n_per_dim;
        n_upsampled[i] = upsamp * n_inner[i];
    }
    int N_UPSAMP = n_upsampled[0] * n_upsampled[1];

    auto [in_naive, out_naive] = initialize_arrays(n_inner, n_upsampled);
    auto [in_clever, out_clever] = initialize_arrays(n_inner, n_upsampled);

    double t_naive = 0.0;
    for (int i = 0; i < n_runs; ++i) {
        std::memset(out_naive, 0, N_UPSAMP * sizeof(fftw_complex));
        t_naive += eval_naive(n_upsampled, in_naive, out_naive);
    }

    double t_clever =  0.0;
    for (int i = 0; i < n_runs; ++i) {
        std::memset(out_clever, 0, N_UPSAMP * sizeof(fftw_complex));
        t_clever += eval_clever(n_inner, n_upsampled, in_clever, out_clever);
    }

    printf("%.4f\t%.4f\n", t_naive, t_naive / n_runs);
    printf("%.4f\t%.4f\n", t_clever, t_clever / n_runs);
    printf("%.4f\n", t_naive / t_clever);

    fftw_free(in_naive);
    fftw_free(in_clever);
    fftw_free(out_naive);
    fftw_free(out_clever);
    return 0;
}
% g++ upsample_test.cpp -lmkl_rt -fopenmp -std=c++17 -g -O3
% OMP_NUM_THREADS=1 taskset -c 0 ./a.out 2000
1.6956  0.1696
1.1748  0.1175
1.4433
@blackwer blackwer changed the title There may be other ways than pruning to speed up the FFT. In a unifom-to-nonuniform 2D NUFFT with, say, an oversampling factor of 2, you only need to transform over the first axis for half of the array, since the other half of the array contains vectors that are identically zero. This saves 25% of FFT cost without having to use complicated hand-crafted algorithms. In 3D, you have to transform a quarter of the array along the first axis and half of the array along the second axis, saving an even larger fraction of FFT time. Exploit structure of upsampled input data when doing FFTs Jun 28, 2023
@mreineck
Copy link
Collaborator

I pushed a demo for the partial FFTs as part of the #287 PR.

@turbotage
Copy link

Any news here? Seems like an amazing performance boost! Are there any limitations as to why this shouldn't work well both for cufinufft and finufft or certain upsampling factors?

@ahbarnett
Copy link
Collaborator

ahbarnett commented Dec 29, 2023 via email

@DiamonDinoia
Copy link
Collaborator

I think #463 achieves this for DUCC. I do not think supporting FFTW in the same way is worth.

@ahbarnett
Copy link
Collaborator

ahbarnett commented Jul 17, 2024 via email

@DiamonDinoia
Copy link
Collaborator

To clarify. I do not thing adding support for sparsity in fftw is worth the maintainability tradeoff. I am not advocating for dropping fftw. This can be reopened in the future in case there a strong need, or can be moved to discussions for future reference.

@flatironinstitute flatironinstitute locked and limited conversation to collaborators Aug 1, 2024
@ahbarnett ahbarnett converted this issue into discussion #510 Aug 1, 2024

This issue was moved to a discussion.

You can continue the conversation there. Go to discussion →

Projects
None yet
Development

No branches or pull requests

5 participants