$darkmode
Eigen-unsupported  5.0.1-dev
fftw_impl.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2009 Mark Borgerding mark a borgerding net
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9 
10 // IWYU pragma: private
11 #include "./InternalHeaderCheck.h"
12 
13 #include <memory>
14 
15 namespace Eigen {
16 
17 namespace internal {
18 
19 // FFTW uses non-const arguments
20 // so we must use ugly const_cast calls for all the args it uses
21 //
22 // This should be safe as long as
23 // 1. we use FFTW_ESTIMATE for all our planning
24 // see the FFTW docs section 4.3.2 "Planner Flags"
25 // 2. fftw_complex is compatible with std::complex
26 // This assumes std::complex<T> layout is array of size 2 with real,imag
27 template <typename T>
28 inline T *fftw_cast(const T *p) {
29  return const_cast<T *>(p);
30 }
31 
32 inline fftw_complex *fftw_cast(const std::complex<double> *p) {
33  return const_cast<fftw_complex *>(reinterpret_cast<const fftw_complex *>(p));
34 }
35 
36 inline fftwf_complex *fftw_cast(const std::complex<float> *p) {
37  return const_cast<fftwf_complex *>(reinterpret_cast<const fftwf_complex *>(p));
38 }
39 
40 inline fftwl_complex *fftw_cast(const std::complex<long double> *p) {
41  return const_cast<fftwl_complex *>(reinterpret_cast<const fftwl_complex *>(p));
42 }
43 
44 template <typename T>
45 struct fftw_plan {};
46 
47 template <>
48 struct 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;
53 
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);
58  }
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);
62  }
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);
66  }
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);
70  }
71 
72  inline void fwd2(complex_type *dst, complex_type *src, int n0, int n1) {
73  if (m_plan == NULL)
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);
76  }
77  inline void inv2(complex_type *dst, complex_type *src, int n0, int n1) {
78  if (m_plan == NULL)
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);
81  }
82 };
83 template <>
84 struct 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;
89 
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);
94  }
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);
98  }
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);
102  }
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);
106  }
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);
110  }
111  inline void inv2(complex_type *dst, complex_type *src, int n0, int n1) {
112  if (m_plan == NULL)
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);
115  }
116 };
117 template <>
118 struct 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;
123 
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);
128  }
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);
132  }
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);
136  }
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);
140  }
141  inline void fwd2(complex_type *dst, complex_type *src, int n0, int n1) {
142  if (m_plan == NULL)
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);
145  }
146  inline void inv2(complex_type *dst, complex_type *src, int n0, int n1) {
147  if (m_plan == NULL)
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);
150  }
151 };
152 
153 template <typename Scalar_>
154 struct fftw_impl {
155  typedef Scalar_ Scalar;
156  typedef std::complex<Scalar> Complex;
157 
158  inline void clear() { m_plans.clear(); }
159 
160  // complex-to-complex forward FFT
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);
163  }
164 
165  // real-to-complex forward FFT
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);
168  }
169 
170  // 2-d complex-to-complex
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);
173  }
174 
175  // inverse complex-to-complex
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);
178  }
179 
180  // half-complex to scalar
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);
183  }
184 
185  // 2-d complex-to-complex
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);
188  }
189 
190  protected:
191  typedef fftw_plan<Scalar> PlanData;
192 
193  typedef Eigen::numext::int64_t int64_t;
194 
195  typedef std::map<int64_t, PlanData> PlanMap;
196 
197  PlanMap m_plans;
198 
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;
203  return m_plans[key];
204  }
205 
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;
210  return m_plans[key];
211  }
212 };
213 
214 } // end namespace internal
215 
216 } // end namespace Eigen
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)