Page MenuHomec4science

fftw_interface.hh
No OneTemporary

File Metadata

Created
Sat, Apr 27, 15:36

fftw_interface.hh

/**
* @file
* LICENSE
*
* Copyright (©) 2016-2021 EPFL (École Polytechnique Fédérale de Lausanne),
* Laboratory (LSMS - Laboratoire de Simulation en Mécanique des Solides)
*
* This program is free software: you can redistribute it and/or modify
* it under the terms of the GNU Affero General Public License as published
* by the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* This program is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU Affero General Public License for more details.
*
* You should have received a copy of the GNU Affero General Public License
* along with this program. If not, see <https://www.gnu.org/licenses/>.
*
*/
/* -------------------------------------------------------------------------- */
#ifndef FFTW_INTERFACE
#define FFTW_INTERFACE
/* -------------------------------------------------------------------------- */
#include <cstddef>
#include <functional>
#include <numeric>
#include <utility>
#include <fftw3.h>
#ifdef TAMAAS_USE_MPI
#include <fftw3-mpi.h>
#endif
namespace fftw {
template <typename T>
struct helper;
template <>
struct helper<double> {
using complex = fftw_complex;
using plan = fftw_plan;
static auto alloc_real(std::size_t size) { return fftw_alloc_real(size); }
static auto alloc_complex(std::size_t size) {
return fftw_alloc_complex(size);
}
};
template <>
struct helper<long double> {
using complex = fftwl_complex;
using plan = fftwl_plan;
static auto alloc_real(std::size_t size) { return fftwl_alloc_real(size); }
static auto alloc_complex(std::size_t size) {
return fftwl_alloc_complex(size);
}
};
template <typename T>
inline auto free(T* ptr) {
fftw_free(ptr);
}
inline auto destroy(fftw_plan plan) { fftw_destroy_plan(plan); }
inline auto destroy(fftwl_plan plan) { fftwl_destroy_plan(plan); }
inline auto init_threads() { return fftw_init_threads(); }
inline auto plan_with_nthreads(int nthreads) {
return fftw_plan_with_nthreads(nthreads);
}
inline auto cleanup_threads() { return fftw_cleanup_threads(); }
/// Holder type for fftw plans
template <typename T>
struct plan {
typename helper<T>::plan _plan;
/// Create from plan
explicit plan(typename helper<T>::plan _plan = nullptr) : _plan(_plan) {}
/// Move constructor to avoid accidental plan destruction
plan(plan&& o) : _plan(std::exchange(o._plan, nullptr)) {}
/// Move operator
plan& operator=(plan&& o) {
_plan = std::exchange(o._plan, nullptr);
return *this;
}
/// Destroy plan
~plan() noexcept {
if (_plan)
fftw::destroy(_plan);
}
/// For seamless use with fftw api
operator typename helper<T>::plan() const { return _plan; }
};
/// RAII helper for fftw_free
template <typename T>
struct ptr {
T* _ptr;
~ptr() noexcept {
if (_ptr)
fftw::free(_ptr);
}
operator T*() { return _ptr; }
};
/* -------------------------------------------------------------------------- */
inline auto plan_many_forward(int rank, const int* n, int howmany, double* in,
const int* inembed, int istride, int idist,
fftw_complex* out, const int* onembed,
int ostride, int odist, unsigned flags) {
return fftw_plan_many_dft_r2c(rank, n, howmany, in, inembed, istride, idist,
out, onembed, ostride, odist, flags);
}
inline auto plan_many_backward(int rank, const int* n, int howmany,
fftw_complex* in, const int* inembed,
int istride, int idist, double* out,
const int* onembed, int ostride, int odist,
unsigned flags) {
return fftw_plan_many_dft_c2r(rank, n, howmany, in, inembed, istride, idist,
out, onembed, ostride, odist, flags);
}
inline auto plan_1d_forward(int n, double* in, fftw_complex* out,
unsigned flags) {
return fftw_plan_dft_r2c_1d(n, in, out, flags);
}
inline auto plan_1d_backward(int n, fftw_complex* in, double* out,
unsigned flags) {
return fftw_plan_dft_c2r_1d(n, in, out, flags);
}
inline auto plan_2d_forward(int n0, int n1, double* in, fftw_complex* out,
unsigned flags) {
return fftw_plan_dft_r2c_2d(n0, n1, in, out, flags);
}
inline auto plan_2d_backward(int n0, int n1, fftw_complex* out, double* in,
unsigned flags) {
return fftw_plan_dft_c2r_2d(n0, n1, out, in, flags);
}
inline auto execute(fftw_plan plan) { fftw_execute(plan); }
inline auto execute(fftw_plan plan, double* in, fftw_complex* out) {
fftw_execute_dft_r2c(plan, in, out);
}
inline auto execute(fftw_plan plan, fftw_complex* in, double* out) {
fftw_execute_dft_c2r(plan, in, out);
}
/* -------------------------------------------------------------------------- */
inline auto plan_many_forward(int rank, const int* n, int howmany,
long double* in, const int* inembed, int istride,
int idist, fftwl_complex* out, const int* onembed,
int ostride, int odist, unsigned flags) {
return fftwl_plan_many_dft_r2c(rank, n, howmany, in, inembed, istride, idist,
out, onembed, ostride, odist, flags);
}
inline auto plan_many_backward(int rank, const int* n, int howmany,
fftwl_complex* in, const int* inembed,
int istride, int idist, long double* out,
const int* onembed, int ostride, int odist,
unsigned flags) {
return fftwl_plan_many_dft_c2r(rank, n, howmany, in, inembed, istride, idist,
out, onembed, ostride, odist, flags);
}
inline auto plan_1d_forward(int n, long double* in, fftwl_complex* out,
unsigned flags) {
return fftwl_plan_dft_r2c_1d(n, in, out, flags);
}
inline auto plan_1d_backward(int n, fftwl_complex* in, long double* out,
unsigned flags) {
return fftwl_plan_dft_c2r_1d(n, in, out, flags);
}
inline auto plan_2d_forward(int n0, int n1, long double* in, fftwl_complex* out,
unsigned flags) {
return fftwl_plan_dft_r2c_2d(n0, n1, in, out, flags);
}
inline auto plan_2d_backward(int n0, int n1, fftwl_complex* out,
long double* in, unsigned flags) {
return fftwl_plan_dft_c2r_2d(n0, n1, out, in, flags);
}
inline auto execute(fftwl_plan plan) { fftwl_execute(plan); }
inline auto execute(fftwl_plan plan, long double* in, fftwl_complex* out) {
fftwl_execute_dft_r2c(plan, in, out);
}
inline auto execute(fftwl_plan plan, fftwl_complex* in, long double* out) {
fftwl_execute_dft_c2r(plan, in, out);
}
} // namespace fftw
#endif // FFTW_INTERFACE

Event Timeline