11#include "./InternalHeaderCheck.h"
28inline T *fftw_cast(
const T *p) {
29 return const_cast<T *
>(p);
32inline fftw_complex *fftw_cast(
const std::complex<double> *p) {
33 return const_cast<fftw_complex *
>(
reinterpret_cast<const fftw_complex *
>(p));
36inline fftwf_complex *fftw_cast(
const std::complex<float> *p) {
37 return const_cast<fftwf_complex *
>(
reinterpret_cast<const fftwf_complex *
>(p));
40inline fftwl_complex *fftw_cast(
const std::complex<long double> *p) {
41 return const_cast<fftwl_complex *
>(
reinterpret_cast<const fftwl_complex *
>(p));
48struct fftw_plan<float> {
49 typedef float scalar_type;
50 typedef fftwf_complex complex_type;
51 std::shared_ptr<fftwf_plan_s> m_plan;
52 fftw_plan() =
default;
54 void set_plan(fftwf_plan p) { m_plan.reset(p, fftwf_destroy_plan); }
55 inline void fwd(complex_type *dst, complex_type *src,
int nfft) {
56 if (m_plan == NULL) set_plan(fftwf_plan_dft_1d(nfft, src, dst, FFTW_FORWARD, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT));
57 fftwf_execute_dft(m_plan.get(), src, dst);
59 inline void inv(complex_type *dst, complex_type *src,
int nfft) {
60 if (m_plan == NULL) set_plan(fftwf_plan_dft_1d(nfft, src, dst, FFTW_BACKWARD, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT));
61 fftwf_execute_dft(m_plan.get(), src, dst);
63 inline void fwd(complex_type *dst, scalar_type *src,
int nfft) {
64 if (m_plan == NULL) set_plan(fftwf_plan_dft_r2c_1d(nfft, src, dst, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT));
65 fftwf_execute_dft_r2c(m_plan.get(), src, dst);
67 inline void inv(scalar_type *dst, complex_type *src,
int nfft) {
68 if (m_plan == NULL) set_plan(fftwf_plan_dft_c2r_1d(nfft, src, dst, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT));
69 fftwf_execute_dft_c2r(m_plan.get(), src, dst);
72 inline void fwd2(complex_type *dst, complex_type *src,
int n0,
int n1) {
74 set_plan(fftwf_plan_dft_2d(n0, n1, src, dst, FFTW_FORWARD, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT));
75 fftwf_execute_dft(m_plan.get(), src, dst);
77 inline void inv2(complex_type *dst, complex_type *src,
int n0,
int n1) {
79 set_plan(fftwf_plan_dft_2d(n0, n1, src, dst, FFTW_BACKWARD, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT));
80 fftwf_execute_dft(m_plan.get(), src, dst);
84struct fftw_plan<double> {
85 typedef double scalar_type;
86 typedef fftw_complex complex_type;
87 std::shared_ptr<fftw_plan_s> m_plan;
88 fftw_plan() =
default;
90 void set_plan(::fftw_plan p) { m_plan.reset(p, fftw_destroy_plan); }
91 inline void fwd(complex_type *dst, complex_type *src,
int nfft) {
92 if (m_plan == NULL) set_plan(fftw_plan_dft_1d(nfft, src, dst, FFTW_FORWARD, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT));
93 fftw_execute_dft(m_plan.get(), src, dst);
95 inline void inv(complex_type *dst, complex_type *src,
int nfft) {
96 if (m_plan == NULL) set_plan(fftw_plan_dft_1d(nfft, src, dst, FFTW_BACKWARD, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT));
97 fftw_execute_dft(m_plan.get(), src, dst);
99 inline void fwd(complex_type *dst, scalar_type *src,
int nfft) {
100 if (m_plan == NULL) set_plan(fftw_plan_dft_r2c_1d(nfft, src, dst, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT));
101 fftw_execute_dft_r2c(m_plan.get(), src, dst);
103 inline void inv(scalar_type *dst, complex_type *src,
int nfft) {
104 if (m_plan == NULL) set_plan(fftw_plan_dft_c2r_1d(nfft, src, dst, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT));
105 fftw_execute_dft_c2r(m_plan.get(), src, dst);
107 inline void fwd2(complex_type *dst, complex_type *src,
int n0,
int n1) {
108 if (m_plan == NULL) set_plan(fftw_plan_dft_2d(n0, n1, src, dst, FFTW_FORWARD, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT));
109 fftw_execute_dft(m_plan.get(), src, dst);
111 inline void inv2(complex_type *dst, complex_type *src,
int n0,
int n1) {
113 set_plan(fftw_plan_dft_2d(n0, n1, src, dst, FFTW_BACKWARD, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT));
114 fftw_execute_dft(m_plan.get(), src, dst);
118struct fftw_plan<long double> {
119 typedef long double scalar_type;
120 typedef fftwl_complex complex_type;
121 std::shared_ptr<fftwl_plan_s> m_plan;
122 fftw_plan() =
default;
124 void set_plan(fftwl_plan p) { m_plan.reset(p, fftwl_destroy_plan); }
125 inline void fwd(complex_type *dst, complex_type *src,
int nfft) {
126 if (m_plan == NULL) set_plan(fftwl_plan_dft_1d(nfft, src, dst, FFTW_FORWARD, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT));
127 fftwl_execute_dft(m_plan.get(), src, dst);
129 inline void inv(complex_type *dst, complex_type *src,
int nfft) {
130 if (m_plan == NULL) set_plan(fftwl_plan_dft_1d(nfft, src, dst, FFTW_BACKWARD, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT));
131 fftwl_execute_dft(m_plan.get(), src, dst);
133 inline void fwd(complex_type *dst, scalar_type *src,
int nfft) {
134 if (m_plan == NULL) set_plan(fftwl_plan_dft_r2c_1d(nfft, src, dst, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT));
135 fftwl_execute_dft_r2c(m_plan.get(), src, dst);
137 inline void inv(scalar_type *dst, complex_type *src,
int nfft) {
138 if (m_plan == NULL) set_plan(fftwl_plan_dft_c2r_1d(nfft, src, dst, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT));
139 fftwl_execute_dft_c2r(m_plan.get(), src, dst);
141 inline void fwd2(complex_type *dst, complex_type *src,
int n0,
int n1) {
143 set_plan(fftwl_plan_dft_2d(n0, n1, src, dst, FFTW_FORWARD, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT));
144 fftwl_execute_dft(m_plan.get(), src, dst);
146 inline void inv2(complex_type *dst, complex_type *src,
int n0,
int n1) {
148 set_plan(fftwl_plan_dft_2d(n0, n1, src, dst, FFTW_BACKWARD, FFTW_ESTIMATE | FFTW_PRESERVE_INPUT));
149 fftwl_execute_dft(m_plan.get(), src, dst);
153template <
typename Scalar_>
155 typedef Scalar_ Scalar;
156 typedef std::complex<Scalar> Complex;
158 inline void clear() { m_plans.clear(); }
161 inline void fwd(Complex *dst,
const Complex *src,
int nfft) {
162 get_plan(nfft,
false, dst, src).fwd(fftw_cast(dst), fftw_cast(src), nfft);
166 inline void fwd(Complex *dst,
const Scalar *src,
int nfft) {
167 get_plan(nfft,
false, dst, src).fwd(fftw_cast(dst), fftw_cast(src), nfft);
171 inline void fwd2(Complex *dst,
const Complex *src,
int n0,
int n1) {
172 get_plan(n0, n1,
false, dst, src).fwd2(fftw_cast(dst), fftw_cast(src), n0, n1);
176 inline void inv(Complex *dst,
const Complex *src,
int nfft) {
177 get_plan(nfft,
true, dst, src).inv(fftw_cast(dst), fftw_cast(src), nfft);
181 inline void inv(Scalar *dst,
const Complex *src,
int nfft) {
182 get_plan(nfft,
true, dst, src).inv(fftw_cast(dst), fftw_cast(src), nfft);
186 inline void inv2(Complex *dst,
const Complex *src,
int n0,
int n1) {
187 get_plan(n0, n1,
true, dst, src).inv2(fftw_cast(dst), fftw_cast(src), n0, n1);
191 typedef fftw_plan<Scalar> PlanData;
193 typedef Eigen::numext::int64_t int64_t;
195 typedef std::map<int64_t, PlanData> PlanMap;
199 inline PlanData &get_plan(
int nfft,
bool inverse,
void *dst,
const void *src) {
200 bool inplace = (dst == src);
201 bool aligned = ((
reinterpret_cast<size_t>(src) & 15) | (
reinterpret_cast<size_t>(dst) & 15)) == 0;
202 int64_t key = ((nfft << 3) | (
inverse << 2) | (inplace << 1) | aligned) << 1;
206 inline PlanData &get_plan(
int n0,
int n1,
bool inverse,
void *dst,
const void *src) {
207 bool inplace = (dst == src);
208 bool aligned = ((
reinterpret_cast<size_t>(src) & 15) | (
reinterpret_cast<size_t>(dst) & 15)) == 0;
209 int64_t key = (((((int64_t)n0) << 30) | (n1 << 3) | (
inverse << 2) | (inplace << 1) | aligned) << 1) + 1;
Namespace containing all symbols from the Eigen library.
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_inverse_op< typename Derived::Scalar >, const Derived > inverse(const Eigen::ArrayBase< Derived > &x)