14#include <KokkosFFT.hpp>
15#include <Kokkos_Core.hpp>
20
21
22
23
24template <
typename Dim>
28
29
30
31
38
39
40
41
48
49
50
51
56namespace ddc::detail::fft {
65struct RealType<Kokkos::complex<T>>
71using real_type_t =
typename RealType<T>::type;
75struct is_complex : std::false_type
80struct is_complex<Kokkos::complex<T>> : std::true_type
85constexpr bool is_complex_v = is_complex<T>::value;
88
89
90
91
99template <
typename... DDimX>
100KokkosFFT::axis_type<
sizeof...(DDimX)> axes()
102 return KokkosFFT::axis_type<
sizeof...(DDimX)> {
103 static_cast<
int>(
ddc::type_seq_rank_v<DDimX,
ddc::detail::TypeSeq<DDimX...>>)...};
106inline KokkosFFT::Normalization ddc_fft_normalization_to_kokkos_fft(
111 return KokkosFFT::Normalization::none;
115 return KokkosFFT::Normalization::forward;
119 return KokkosFFT::Normalization::backward;
123 return KokkosFFT::Normalization::ortho;
126 throw std::runtime_error(
"ddc::FFT_Normalization not handled");
135 explicit ScaleFn(T coef)
noexcept : m_coef(std::move(coef)) {}
138 [[nodiscard]] KOKKOS_FUNCTION U operator()(U
const& value)
const noexcept
140 return m_coef * value;
145Real forward_full_norm_coef(
DiscreteDomain<DDim>
const& ddom)
noexcept
147 return rlength(ddom) / Kokkos::sqrt(2 * Kokkos::numbers::pi_v<Real>)
148 / (ddom.extents() - 1).value();
152Real backward_full_norm_coef(
DiscreteDomain<DDim>
const& ddom)
noexcept
154 return 1 / (forward_full_norm_coef(ddom) * ddom.extents().value());
162 typename MemorySpace,
168 ExecSpace
const& exec_space,
171 KwArgsImpl
const& kwargs)
174 std::is_same_v<real_type_t<Tin>,
float> || std::is_same_v<real_type_t<Tin>,
double>,
175 "Base type of Tin (and Tout) must be float or double.");
177 std::is_same_v<real_type_t<Tin>, real_type_t<Tout>>,
178 "Types Tin and Tout must be based on same type (float or double)");
180 Kokkos::SpaceAccessibility<ExecSpace, MemorySpace>::accessible,
181 "MemorySpace has to be accessible for ExecutionSpace.");
184 ddc::detail::mdspan_to_kokkos_element_t<Tin,
sizeof...(DDimIn)>,
185 ddc::detail::mdspan_to_kokkos_layout_t<LayoutIn>,
186 MemorySpace>
const in_view
187 = in.allocation_kokkos_view();
189 ddc::detail::mdspan_to_kokkos_element_t<Tout,
sizeof...(DDimIn)>,
190 ddc::detail::mdspan_to_kokkos_layout_t<LayoutOut>,
191 MemorySpace>
const out_view
192 = out.allocation_kokkos_view();
193 KokkosFFT::Normalization
const kokkos_fft_normalization
194 = ddc_fft_normalization_to_kokkos_fft(kwargs.normalization);
197 if constexpr (std::is_same_v<Tin, Tout>) {
204 kokkos_fft_normalization);
211 kokkos_fft_normalization);
215 if constexpr (is_complex_v<Tout>) {
216 assert(kwargs.direction == ddc::FFT_Direction::FORWARD);
222 kokkos_fft_normalization);
224 assert(kwargs.direction == ddc::FFT_Direction::BACKWARD);
230 kokkos_fft_normalization);
239 norm_coef = (forward_full_norm_coef(DiscreteDomain<DDimIn>(ddom_in)) * ...);
242 norm_coef = (backward_full_norm_coef(DiscreteDomain<DDimOut>(ddom_out)) * ...);
245 ddc::parallel_transform(exec_space, out, ScaleFn<real_type_t<Tout>>(norm_coef));
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272template <
typename DDimFx,
typename DDimX>
277 is_uniform_point_sampling_v<DDimX>,
278 "DDimX dimension must derive from UniformPointSampling");
280 is_periodic_sampling_v<DDimFx>,
281 "DDimFx dimension must derive from PeriodicSampling");
282 using CDimFx =
typename DDimFx::continuous_dimension_type;
283 using CDimX =
typename DDimX::continuous_dimension_type;
285 std::is_same_v<CDimFx,
ddc::Fourier<CDimX>>,
286 "DDimX and DDimFx dimensions must be defined over the same continuous dimension");
288 DiscreteVectorElement
const nx = get<DDimX>(x_mesh.extents());
289 double const lx =
ddc::rlength(x_mesh);
290 auto [impl, ddom] = DDimFx::
template init<DDimFx>(
291 ddc::Coordinate<CDimFx>(0),
292 ddc::Coordinate<CDimFx>(2 * (nx - 1) * (nx - 1) / (nx * lx) * Kokkos::numbers::pi),
295 return std::move(impl);
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313template <
typename... DDimFx,
typename... DDimX>
317 (is_uniform_point_sampling_v<DDimX> && ...),
318 "DDimX dimensions should derive from UniformPointSampling");
320 (is_periodic_sampling_v<DDimFx> && ...),
321 "DDimFx dimensions should derive from PeriodicPointSampling");
324 detail::array(extents).back() = detail::array(extents).back() / 2 + 1;
327 ddc::DiscreteElement<DDimFx>(0),
332
333
334
335
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
369 typename MemorySpace,
373 ExecSpace
const& exec_space,
379 std::is_same_v<LayoutIn, Kokkos::layout_right>
380 && std::is_same_v<LayoutOut, Kokkos::layout_right>,
381 "Layouts must be right-handed");
383 (is_uniform_point_sampling_v<DDimX> && ...),
384 "DDimX dimensions should derive from UniformPointSampling");
386 (is_periodic_sampling_v<DDimFx> && ...),
387 "DDimFx dimensions should derive from PeriodicPointSampling");
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
422 typename MemorySpace,
426 ExecSpace
const& exec_space,
432 std::is_same_v<LayoutIn, Kokkos::layout_right>
433 && std::is_same_v<LayoutOut, Kokkos::layout_right>,
434 "Layouts must be right-handed");
436 (is_uniform_point_sampling_v<DDimX> && ...),
437 "DDimX dimensions should derive from UniformPointSampling");
439 (is_periodic_sampling_v<DDimFx> && ...),
440 "DDimFx dimensions should derive from PeriodicPointSampling");
friend class DiscreteDomain
KOKKOS_FUNCTION constexpr bool operator!=(DiscreteVector< OTags... > const &rhs) const noexcept
The top-level namespace of DDC.
ddc::FFT_Normalization normalization
Enum member to identify the type of normalization performed.
void ifft(ExecSpace const &exec_space, ddc::ChunkSpan< Tout, ddc::DiscreteDomain< DDimX... >, LayoutOut, MemorySpace > out, ddc::ChunkSpan< Tin, ddc::DiscreteDomain< DDimFx... >, LayoutIn, MemorySpace > in, ddc::kwArgs_fft kwargs={ddc::FFT_Normalization::OFF})
Perform an inverse Fast Fourier Transform.
void fft(ExecSpace const &exec_space, ddc::ChunkSpan< Tout, ddc::DiscreteDomain< DDimFx... >, LayoutOut, MemorySpace > out, ddc::ChunkSpan< Tin, ddc::DiscreteDomain< DDimX... >, LayoutIn, MemorySpace > in, ddc::kwArgs_fft kwargs={ddc::FFT_Normalization::OFF})
Perform a direct Fast Fourier Transform.
FFT_Normalization
A named argument to choose the type of normalization of the FFT.
@ BACKWARD
No normalization for forward FFT, multiply by 1/N for backward FFT.
@ OFF
No normalization. Un-normalized FFT is sum_j f(x_j)*e^-ikx_j.
@ ORTHO
Multiply by 1/sqrt(N)
@ FULL
Multiply by dx/sqrt(2*pi) for forward FFT and dk/sqrt(2*pi) for backward FFT.
@ FORWARD
Multiply by 1/N for forward FFT, no normalization for backward FFT.
ddc::DiscreteDomain< DDimFx... > fourier_mesh(ddc::DiscreteDomain< DDimX... > x_mesh, bool C2C)
Get the Fourier mesh.
DDimFx::template Impl< DDimFx, Kokkos::HostSpace > init_fourier_space(ddc::DiscreteDomain< DDimX > x_mesh)
Initialize a Fourier discrete dimension.
FFT_Direction
A named argument to choose the direction of the FFT.
@ BACKWARD
Backward, corresponds to inverse FFT up to normalization.
@ FORWARD
Forward, corresponds to direct FFT up to normalization.
A structure embedding the configuration of the exposed FFT function with the type of normalization.