16#ifndef CABANA_GRID_FASTFOURIERTRANSFORM_HPP
17#define CABANA_GRID_FASTFOURIERTRANSFORM_HPP
22#include <Kokkos_Core.hpp>
23#include <Kokkos_Profiling_ScopedRegion.hpp>
25#include <heffte_fft3d.h>
64struct FFTBackendDefault
71template <
class ArrayEntity,
class ArrayMesh,
class ArrayMemorySpace,
72 class ArrayScalar,
class Entity,
class Mesh,
class MemorySpace,
73 class Scalar,
typename SFINAE =
void>
76 static_assert( std::is_same<ArrayEntity, Entity>::value,
77 "Array entity type must match FFT entity type." );
78 static_assert( std::is_same<ArrayMesh, Mesh>::value,
79 "Array mesh type must match FFT mesh type." );
80 static_assert( std::is_same<ArrayMemorySpace, MemorySpace>::value,
81 "Array memory space must match FFT memory space." );
85template <
class ArrayEntity,
class ArrayMesh,
class ArrayMemorySpace,
86 class ArrayScalar,
class Entity,
class Mesh,
class MemorySpace,
89 ArrayEntity, ArrayMesh, ArrayMemorySpace, ArrayScalar, Entity, Mesh,
91 typename std::enable_if<
92 std::is_same<ArrayEntity, Entity>::value &&
93 std::is_same<ArrayMesh, Mesh>::value &&
94 std::is_same<ArrayMemorySpace, MemorySpace>::value>::type>
95 :
public std::true_type
133 FFTcomm = FFTCommPattern::alltoallv;
135 FFTcomm = FFTCommPattern::p2p;
170template <
class EntityType,
class MeshType,
class Scalar,
class MemorySpace,
184 static_assert( Kokkos::is_memory_space<MemorySpace>() );
209 const auto& global_grid = layout.
localGrid()->globalGrid();
218 int local_num_entity =
231 throw std::logic_error(
232 "Cabana::Grid::Experimental::FastFourierTransform::"
233 "checkArrayDofs: Only 1 complex value per entity allowed in "
242 template <
class Array_t,
class ScaleType>
244 const Array_t& x,
const ScaleType scaling,
245 typename std::enable_if<
248 typename Array_t::entity_type,
typename Array_t::mesh_type,
249 typename Array_t::memory_space,
typename Array_t::value_type,
253 Kokkos::Profiling::ScopedRegion region(
"Cabana::FFT::forward" );
256 static_cast<Derived*
>( this )->forwardImpl( x, scaling );
264 template <
class Array_t,
class ScaleType>
266 const Array_t& x,
const ScaleType scaling,
267 typename std::enable_if<
270 typename Array_t::entity_type,
typename Array_t::mesh_type,
271 typename Array_t::memory_space,
typename Array_t::value_type,
275 Kokkos::Profiling::ScopedRegion region(
"Cabana::FFT::reverse" );
278 static_cast<Derived*
>( this )->reverseImpl( x, scaling );
284 template <
class ExecutionSpace,
class IndexSpaceType,
class LViewType,
286 std::enable_if_t<3 == NSD, void>
287 copyToLocal( ExecutionSpace exec_space,
const IndexSpaceType own_space,
288 LViewType& l_view,
const LGViewType lg_view )
290 Kokkos::parallel_for(
291 "Cabana::Grid::FastFourierTransform::copyTo",
293 KOKKOS_LAMBDA(
const int i,
const int j,
const int k ) {
294 auto iw = i - own_space.min( Dim::I );
295 auto jw = j - own_space.min( Dim::J );
296 auto kw = k - own_space.min( Dim::K );
297 l_view( iw, jw, kw, 0 ) = lg_view( i, j, k, 0 );
298 l_view( iw, jw, kw, 1 ) = lg_view( i, j, k, 1 );
305 template <
class ExecutionSpace,
class IndexSpaceType,
class LViewType,
307 std::enable_if_t<2 == NSD, void>
308 copyToLocal( ExecutionSpace space,
const IndexSpaceType own_space,
309 LViewType& l_view,
const LGViewType lg_view )
311 Kokkos::parallel_for(
312 "Cabana::Grid::FastFourierTransform::copyTo",
314 KOKKOS_LAMBDA(
const int i,
const int j ) {
315 auto iw = i - own_space.min( Dim::I );
316 auto jw = j - own_space.min( Dim::J );
317 l_view( iw, jw, 0 ) = lg_view( i, j, 0 );
318 l_view( iw, jw, 1 ) = lg_view( i, j, 1 );
325 template <
class ExecutionSpace,
class IndexSpaceType,
class LViewType,
327 std::enable_if_t<3 == NSD, void>
329 const LViewType l_view, LGViewType& lg_view )
331 Kokkos::parallel_for(
332 "Cabana::Grid::FastFourierTransform::copyFrom",
334 KOKKOS_LAMBDA(
const int i,
const int j,
const int k ) {
335 auto iw = i - own_space.min( Dim::I );
336 auto jw = j - own_space.min( Dim::J );
337 auto kw = k - own_space.min( Dim::K );
338 lg_view( i, j, k, 0 ) = l_view( iw, jw, kw, 0 );
339 lg_view( i, j, k, 1 ) = l_view( iw, jw, kw, 1 );
346 template <
class ExecutionSpace,
class IndexSpaceType,
class LViewType,
348 std::enable_if_t<2 == NSD, void>
350 const LViewType l_view, LGViewType& lg_view )
352 Kokkos::parallel_for(
353 "Cabana::Grid::FastFourierTransform::copyFrom",
355 KOKKOS_LAMBDA(
const int i,
const int j ) {
356 auto iw = i - own_space.min( Dim::I );
357 auto jw = j - own_space.min( Dim::J );
358 lg_view( i, j, 0 ) = l_view( iw, jw, 0 );
359 lg_view( i, j, 1 ) = l_view( iw, jw, 1 );
370template <
class ExecutionSpace,
class BackendType>
371struct HeffteBackendTraits
374#ifdef Heffte_ENABLE_MKL
375template <
class ExecutionSpace>
376struct HeffteBackendTraits<ExecutionSpace, FFTBackendMKL>
378 using backend_type = heffte::backend::mkl;
381#ifdef Heffte_ENABLE_FFTW
382template <
class ExecutionSpace>
383struct HeffteBackendTraits<ExecutionSpace, FFTBackendFFTW>
385 using backend_type = heffte::backend::fftw;
388#ifdef Heffte_ENABLE_FFTW
389template <
class ExecutionSpace>
390struct HeffteBackendTraits<ExecutionSpace, Impl::FFTBackendDefault>
392 using backend_type = heffte::backend::fftw;
395#ifdef Heffte_ENABLE_MKL
396template <
class ExecutionSpace>
397struct HeffteBackendTraits<ExecutionSpace, Impl::FFTBackendDefault>
399 using backend_type = heffte::backend::mkl;
402throw std::runtime_error(
"Cabana::Grid::FastFourierTransform: Must enable at "
403 "least one heFFTe host backend." );
406#ifdef Heffte_ENABLE_CUDA
407#ifdef KOKKOS_ENABLE_CUDA
409struct HeffteBackendTraits<Kokkos::Cuda, Impl::FFTBackendDefault>
411 using backend_type = heffte::backend::cufft;
415#ifdef Heffte_ENABLE_ROCM
416#ifdef KOKKOS_ENABLE_HIP
418struct HeffteBackendTraits<Kokkos::Experimental::HIP, Impl::FFTBackendDefault>
420 using backend_type = heffte::backend::rocfft;
424#ifdef Heffte_ENABLE_ONEAPI
425#ifdef KOKKOS_ENABLE_SYCL
427struct HeffteBackendTraits<Kokkos::Experimental::SYCL, Impl::FFTBackendDefault>
429 using backend_type = heffte::backend::onemkl;
434template <
class ScaleType>
435struct HeffteScalingTraits
439struct HeffteScalingTraits<FFTScaleNone>
441 static const auto scaling_type = heffte::scale::none;
444struct HeffteScalingTraits<FFTScaleFull>
446 static const auto scaling_type = heffte::scale::full;
449struct HeffteScalingTraits<FFTScaleSymmetric>
451 static const auto scaling_type = heffte::scale::symmetric;
454#ifdef KOKKOS_ENABLE_SYCL
456template <
class ExecSpace,
class HeffteBackendType>
457auto createHeffteFft3d(
458 ExecSpace exec_space, HeffteBackendType, heffte::box3d<> inbox,
459 heffte::box3d<> outbox, MPI_Comm comm, heffte::plan_options params,
460 typename std::enable_if<
461 std::is_same<HeffteBackendType, heffte::backend::onemkl>::value,
466 sycl::queue& q = exec_space.sycl_queue();
467 return std::make_shared<heffte::fft3d<HeffteBackendType>>( q, inbox, outbox,
472template <
class ExecSpace,
class HeffteBackendType>
473auto createHeffteFft3d(
474 ExecSpace, HeffteBackendType, heffte::box3d<> inbox, heffte::box3d<> outbox,
475 MPI_Comm comm, heffte::plan_options params,
476 typename std::enable_if<
477 std::is_same<HeffteBackendType, heffte::backend::fftw>::value ||
478 std::is_same<HeffteBackendType, heffte::backend::mkl>::value ||
479 std::is_same<HeffteBackendType, heffte::backend::cufft>::value ||
480 std::is_same<HeffteBackendType, heffte::backend::rocfft>::value,
485 return std::make_shared<heffte::fft3d<HeffteBackendType>>( inbox, outbox,
496template <
class EntityType,
class MeshType,
class Scalar,
class MemorySpace,
497 class ExecSpace,
class BackendType>
500 EntityType, MeshType, Scalar, MemorySpace,
501 HeffteFastFourierTransform<EntityType, MeshType, Scalar, MemorySpace,
502 ExecSpace, BackendType>>
510 static_assert( Kokkos::is_memory_space<MemorySpace>() );
540 EntityType, MeshType, Scalar, MemorySpace,
542 MemorySpace, ExecSpace, BackendType>>(
547 heffte::box3d<> inbox = { this->global_low, this->global_high };
548 heffte::box3d<> outbox = { this->global_low, this->global_high };
550 heffte::plan_options heffte_params =
551 heffte::default_options<heffte_backend_type>();
555 case Cabana::Grid::Experimental::FFTCommPattern::p2p:
556 heffte_params.algorithm = heffte::reshape_algorithm::p2p;
558 case Cabana::Grid::Experimental::FFTCommPattern::alltoallv:
559 heffte_params.algorithm = heffte::reshape_algorithm::alltoallv;
561 case Cabana::Grid::Experimental::FFTCommPattern::alltoall:
562 heffte_params.algorithm = heffte::reshape_algorithm::alltoall;
564 case Cabana::Grid::Experimental::FFTCommPattern::p2p_plined:
565 heffte_params.algorithm = heffte::reshape_algorithm::p2p_plined;
568 heffte_params.algorithm = heffte::reshape_algorithm::alltoallv;
572 heffte_params.use_pencils = params.
getPencils();
573 heffte_params.use_reorder = params.
getReorder();
577 _fft = Impl::createHeffteFft3d(
579 layout.
localGrid()->globalGrid().comm(), heffte_params );
580 long long fftsize = std::max( _fft->size_outbox(), _fft->size_inbox() );
585 if ( fftsize < (
int)entity_space.size() )
586 throw std::logic_error(
587 "Cabana::Grid::Experimental::HeffteFastFourierTransform: "
588 "Expected FFT allocation size smaller "
589 "than local grid size" );
591 _fft_work = Kokkos::View<Scalar*, memory_space>(
592 Kokkos::ViewAllocateWithoutInitializing(
"fft_work" ),
594 _workspace = Kokkos::View<Scalar* [2], memory_space>(
595 Kokkos::ViewAllocateWithoutInitializing(
"workspace" ),
596 _fft->size_workspace() );
604 template <
class Array_t,
class ScaleType>
607 compute( x, 1, Impl::HeffteScalingTraits<ScaleType>().scaling_type );
615 template <
class Array_t,
class ScaleType>
618 compute( x, -1, Impl::HeffteScalingTraits<ScaleType>().scaling_type );
627 template <
class Array_t>
628 void compute(
const Array_t& x,
const int flag,
const heffte::scale scale )
632 x.layout()->localGrid()->indexSpace(
Own(), EntityType(),
Local() );
635 local_view_space, _fft_work.data() );
639 auto localghost_view = x.view();
641 this->
copyToLocal( heffte_execution_space, own_space, local_view,
647 reinterpret_cast<std::complex<Scalar>*
>( _fft_work.data() ),
648 reinterpret_cast<std::complex<Scalar>*
>( _fft_work.data() ),
649 reinterpret_cast<std::complex<Scalar>*
>( _workspace.data() ),
652 else if ( flag == -1 )
655 reinterpret_cast<std::complex<Scalar>*
>( _fft_work.data() ),
656 reinterpret_cast<std::complex<Scalar>*
>( _fft_work.data() ),
657 reinterpret_cast<std::complex<Scalar>*
>( _workspace.data() ),
662 throw std::logic_error(
663 "Cabana::Grid::Experimental::HeffteFastFourierTransform::"
664 "compute: Only 1:forward and -1:backward are allowed as "
669 this->
copyFromLocal( heffte_execution_space, own_space, local_view,
675 std::shared_ptr<heffte::fft3d<heffte_backend_type>> _fft;
676 Kokkos::View<Scalar*, memory_space> _fft_work;
677 Kokkos::View<Scalar* [2], memory_space> _workspace;
687template <
class Scalar,
class MemorySpace,
class BackendType,
class EntityType,
688 class MeshType,
class ExecSpace>
694 EntityType, MeshType, Scalar, MemorySpace, ExecSpace, BackendType>>(
695 exec_space, layout, params );
702template <
class Scalar,
class MemorySpace,
class EntityType,
class MeshType,
709 Impl::FFTBackendDefault>(
710 exec_space, layout, params );
717template <
class Scalar,
class MemorySpace,
class BackendType,
class EntityType,
718 class MeshType,
class ExecSpace>
722 using heffte_backend_type =
723 typename Impl::HeffteBackendTraits<ExecSpace,
724 BackendType>::backend_type;
727 const heffte::plan_options heffte_params =
728 heffte::default_options<heffte_backend_type>();
731 params.
setPencils( heffte_params.use_pencils );
732 params.
setReorder( heffte_params.use_reorder );
735 EntityType, MeshType, Scalar, MemorySpace, ExecSpace, BackendType>>(
736 exec_space, layout, params );
743template <
class Scalar,
class MemorySpace,
class EntityType,
class MeshType,
749 Scalar, MemorySpace, Impl::FFTBackendDefault, EntityType, MeshType>(
750 exec_space, layout );
756template <
class Scalar,
class MemorySpace,
class BackendType,
class EntityType,
762 using exec_space =
typename MemorySpace::execution_space;
764 EntityType, MeshType>(
765 exec_space{}, layout, params );
771template <
class Scalar,
class MemorySpace,
class EntityType,
class MeshType>
776 using exec_space =
typename MemorySpace::execution_space;
778 Scalar, MemorySpace, Impl::FFTBackendDefault, EntityType, MeshType>(
779 exec_space{}, layout, params );
785template <
class Scalar,
class MemorySpace,
class BackendType,
class EntityType,
790 using exec_space =
typename MemorySpace::execution_space;
792 EntityType, MeshType>( exec_space{},
799template <
class Scalar,
class MemorySpace,
class EntityType,
class MeshType>
803 using exec_space =
typename MemorySpace::execution_space;
805 Scalar, MemorySpace, Impl::FFTBackendDefault, EntityType, MeshType>(
806 exec_space{}, layout );
Kokkos::View< Scalar *, Params... > createView(const std::string &label, const IndexSpace< 1 > &index_space)
Given an index space create a view over the extent of that index space.
Definition Cabana_Grid_IndexSpace.hpp:237
Kokkos::RangePolicy< ExecutionSpace > createExecutionPolicy(const IndexSpace< 1 > &index_space, const ExecutionSpace &exec_space)
Create a multi-dimensional execution policy over an index space.
Definition Cabana_Grid_IndexSpace.hpp:175
IndexSpace< N+1 > appendDimension(const IndexSpace< N > &index_space, const long size)
Given an N-dimensional index space append an additional dimension with the given size.
Definition Cabana_Grid_IndexSpace.hpp:444
Entity layout for array data on the local mesh.
Definition Cabana_Grid_Array.hpp:46
int dofsPerEntity() const
Get the number of degrees-of-freedom on each grid entity.
Definition Cabana_Grid_Array.hpp:77
const std::shared_ptr< LocalGrid< MeshType > > localGrid() const
Get the local grid over which this layout is defined.
Definition Cabana_Grid_Array.hpp:71
Core: particle data structures and algorithms.
Definition Cabana_AoSoA.hpp:36
Tag specifying FFTW backend for FFT (host default).
Definition Cabana_Grid_FastFourierTransform.hpp:54
Tag specifying MKL backend for FFT.
Definition Cabana_Grid_FastFourierTransform.hpp:58
Tag for full scaling of FFT.
Definition Cabana_Grid_FastFourierTransform.hpp:41
Tag for no scaling of FFT.
Definition Cabana_Grid_FastFourierTransform.hpp:45
Tag for symmetric scaling of FFT.
Definition Cabana_Grid_FastFourierTransform.hpp:49
Matching Array static type checker.
Definition Cabana_Grid_FastFourierTransform.hpp:75
Local index tag.
Definition Cabana_Grid_Types.hpp:208
Owned decomposition tag.
Definition Cabana_Grid_Types.hpp:190
Definition Cabana_Grid_Array.hpp:324