13#include <KokkosFFT.hpp>
14#include <Kokkos_Core.hpp>
19
20
21
22
23template <
typename Dim>
27
28
29
30
37
38
39
40
47
48
49
50
55namespace ddc::detail::fft {
64struct RealType<Kokkos::complex<T>>
70using real_type_t =
typename RealType<T>::type;
74struct is_complex : std::false_type
79struct is_complex<Kokkos::complex<T>> : std::true_type
84constexpr bool is_complex_v = is_complex<T>::value;
87
88
89
90
98template <
typename... DDimX>
99KokkosFFT::axis_type<
sizeof...(DDimX)> axes()
101 return KokkosFFT::axis_type<
sizeof...(DDimX)> {
102 static_cast<
int>(
ddc::type_seq_rank_v<DDimX,
ddc::detail::TypeSeq<DDimX...>>)...};
105KokkosFFT::Normalization ddc_fft_normalization_to_kokkos_fft(
114 explicit ScaleFn(T coef)
noexcept : m_coef(std::move(coef)) {}
117 [[nodiscard]] KOKKOS_FUNCTION U operator()(U
const& value)
const noexcept
119 return m_coef * value;
124Real forward_full_norm_coef(
DiscreteDomain<DDim>
const& ddom)
noexcept
126 return rlength(ddom) / Kokkos::sqrt(2 * Kokkos::numbers::pi_v<Real>)
127 / (ddom.extents() - 1).value();
131Real backward_full_norm_coef(
DiscreteDomain<DDim>
const& ddom)
noexcept
133 return 1 / (forward_full_norm_coef(ddom) * ddom.extents().value());
141 typename MemorySpace,
147 ExecSpace
const& exec_space,
150 KwArgsImpl
const& kwargs)
153 std::is_same_v<real_type_t<Tin>,
float> || std::is_same_v<real_type_t<Tin>,
double>,
154 "Base type of Tin (and Tout) must be float or double.");
156 std::is_same_v<real_type_t<Tin>, real_type_t<Tout>>,
157 "Types Tin and Tout must be based on same type (float or double)");
159 Kokkos::SpaceAccessibility<ExecSpace, MemorySpace>::accessible,
160 "MemorySpace has to be accessible for ExecutionSpace.");
163 ddc::detail::mdspan_to_kokkos_element_t<Tin,
sizeof...(DDimIn)>,
164 ddc::detail::mdspan_to_kokkos_layout_t<LayoutIn>,
165 MemorySpace>
const in_view
166 = in.allocation_kokkos_view();
168 ddc::detail::mdspan_to_kokkos_element_t<Tout,
sizeof...(DDimIn)>,
169 ddc::detail::mdspan_to_kokkos_layout_t<LayoutOut>,
170 MemorySpace>
const out_view
171 = out.allocation_kokkos_view();
172 KokkosFFT::Normalization
const kokkos_fft_normalization
173 = ddc_fft_normalization_to_kokkos_fft(kwargs.normalization);
176 if constexpr (std::is_same_v<Tin, Tout>) {
183 kokkos_fft_normalization);
190 kokkos_fft_normalization);
194 if constexpr (is_complex_v<Tout>) {
201 kokkos_fft_normalization);
209 kokkos_fft_normalization);
218 norm_coef = (forward_full_norm_coef(DiscreteDomain<DDimIn>(ddom_in)) * ...);
221 norm_coef = (backward_full_norm_coef(DiscreteDomain<DDimOut>(ddom_out)) * ...);
224 ddc::parallel_transform(exec_space, out, ScaleFn<real_type_t<Tout>>(norm_coef));
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251template <
typename DDimFx,
typename DDimX>
256 is_uniform_point_sampling_v<DDimX>,
257 "DDimX dimension must derive from UniformPointSampling");
259 is_periodic_sampling_v<DDimFx>,
260 "DDimFx dimension must derive from PeriodicSampling");
261 using CDimFx =
typename DDimFx::continuous_dimension_type;
262 using CDimX =
typename DDimX::continuous_dimension_type;
264 std::is_same_v<CDimFx,
ddc::Fourier<CDimX>>,
265 "DDimX and DDimFx dimensions must be defined over the same continuous dimension");
267 DiscreteVectorElement
const nx = get<DDimX>(x_mesh.extents());
268 double const lx =
ddc::rlength(x_mesh);
269 auto [impl, ddom] = DDimFx::
template init<DDimFx>(
270 ddc::Coordinate<CDimFx>(0),
271 ddc::Coordinate<CDimFx>(2 * (nx - 1) * (nx - 1) / (nx * lx) * Kokkos::numbers::pi),
274 return std::move(impl);
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292template <
typename... DDimFx,
typename... DDimX>
296 (is_uniform_point_sampling_v<DDimX> && ...),
297 "DDimX dimensions should derive from UniformPointSampling");
299 (is_periodic_sampling_v<DDimFx> && ...),
300 "DDimFx dimensions should derive from PeriodicPointSampling");
303 detail::array(extents).back() = detail::array(extents).back() / 2 + 1;
306 ddc::DiscreteElement<DDimFx>(0),
311
312
313
314
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
348 typename MemorySpace,
352 ExecSpace
const& exec_space,
358 std::is_same_v<LayoutIn, Kokkos::layout_right>
359 && std::is_same_v<LayoutOut, Kokkos::layout_right>,
360 "Layouts must be right-handed");
362 (is_uniform_point_sampling_v<DDimX> && ...),
363 "DDimX dimensions should derive from UniformPointSampling");
365 (is_periodic_sampling_v<DDimFx> && ...),
366 "DDimFx dimensions should derive from PeriodicPointSampling");
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
401 typename MemorySpace,
405 ExecSpace
const& exec_space,
411 std::is_same_v<LayoutIn, Kokkos::layout_right>
412 && std::is_same_v<LayoutOut, Kokkos::layout_right>,
413 "Layouts must be right-handed");
415 (is_uniform_point_sampling_v<DDimX> && ...),
416 "DDimX dimensions should derive from UniformPointSampling");
418 (is_periodic_sampling_v<DDimFx> && ...),
419 "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.