31 template <
typename Field,
template <
typename...>
class FFT,
typename Backend,
typename T>
33 std::array<long long, 3> low;
34 std::array<long long, 3> high;
39 heffte::box3d<long long> inbox = {low, high};
40 heffte::box3d<long long> outbox = {low, high};
42 setup(inbox, outbox, params);
45 template <
typename Field,
template <
typename...>
class FFT,
typename Backend,
typename T>
47 std::array<long long, 3>& low,
48 std::array<long long, 3>& high) {
56 for (
size_t d = 0; d <
Dim; ++d) {
57 low[d] =
static_cast<long long>(domain[d].
first());
58 high[d] =
static_cast<long long>(domain[d].
length() + domain[d].
first() - 1);
65 template <
typename Field,
template <
typename...>
class FFT,
typename Backend,
typename T>
67 const heffte::box3d<long long>& outbox,
69 heffte::plan_options heffteOptions = heffte::default_options<heffteBackend>();
71 if (!params.
get<
bool>(
"use_heffte_defaults")) {
72 heffteOptions.use_pencils = params.
get<
bool>(
"use_pencils");
73 heffteOptions.use_reorder = params.
get<
bool>(
"use_reorder");
74#ifdef Heffte_ENABLE_GPU
75 heffteOptions.use_gpu_aware = params.
get<
bool>(
"use_gpu_aware");
78 switch (params.
get<
int>(
"comm")) {
80 heffteOptions.algorithm = heffte::reshape_algorithm::alltoall;
83 heffteOptions.algorithm = heffte::reshape_algorithm::alltoallv;
86 heffteOptions.algorithm = heffte::reshape_algorithm::p2p;
89 heffteOptions.algorithm = heffte::reshape_algorithm::p2p_plined;
92 throw IpplException(
"FFT::setup",
"Unrecognized heffte communication type");
96 if constexpr (std::is_same_v<FFT<heffteBackend>, heffte::fft3d<heffteBackend>>) {
97 heffte_m = std::make_shared<FFT<heffteBackend, long long>>(
98 inbox, outbox,
Comm->getCommunicator(), heffteOptions);
100 heffte_m = std::make_shared<FFT<heffteBackend, long long>>(
101 inbox, outbox, params.
get<
int>(
"r2c_direction"),
Comm->getCommunicator(),
111 template <
typename ComplexField>
117 template <
typename ComplexField>
119 static_assert(
Dim == 2 ||
Dim == 3,
"heFFTe only supports 2D and 3D");
121 auto fview = f.getView();
122 const int nghost = f.getNghost();
131 auto& tempField = this->tempField;
132 if (tempField.size() != f.getOwned().size()) {
139 KOKKOS_LAMBDA(
const index_array_type& args) {
140 apply(tempField, args - nghost).real(
apply(fview, args).real());
141 apply(tempField, args - nghost).imag(
apply(fview, args).imag());
145 this->heffte_m->forward(tempField.data(), tempField.data(), this->workspace_m.data(),
146 heffte::scale::full);
148 this->heffte_m->backward(tempField.data(), tempField.data(), this->workspace_m.data(),
149 heffte::scale::none);
151 throw std::logic_error(
"Only 1:forward and -1:backward are allowed as directions");
156 KOKKOS_LAMBDA(
const index_array_type& args) {
157 apply(fview, args).real() =
apply(tempField, args - nghost).real();
158 apply(fview, args).imag() =
apply(tempField, args - nghost).imag();
171 template <
typename RealField>
179 std::array<long long, 3> lowInput;
180 std::array<long long, 3> highInput;
181 std::array<long long, 3> lowOutput;
182 std::array<long long, 3> highOutput;
184 const NDIndex<Dim>& lDomInput = layoutInput.getLocalNDIndex();
185 const NDIndex<Dim>& lDomOutput = layoutOutput.getLocalNDIndex();
187 this->domainToBounds(lDomInput, lowInput, highInput);
188 this->domainToBounds(lDomOutput, lowOutput, highOutput);
190 heffte::box3d<long long> inbox = {lowInput, highInput};
191 heffte::box3d<long long> outbox = {lowOutput, highOutput};
193 this->setup(inbox, outbox, params);
196 template <
typename RealField>
202 template <
typename RealField>
205 static_assert(
Dim == 2 ||
Dim == 3,
"heFFTe only supports 2D and 3D");
207 auto fview = f.getView();
208 auto gview = g.getView();
209 const int nghostf = f.getNghost();
210 const int nghostg = g.getNghost();
219 auto& tempFieldf = this->tempField;
221 if (tempFieldf.size() != f.getOwned().size()) {
224 if (tempFieldg.size() != g.getOwned().size()) {
230 "copy from Kokkos f field in FFT",
getRangePolicy(fview, nghostf),
231 KOKKOS_LAMBDA(
const index_array_type& args) {
232 apply(tempFieldf, args - nghostf) =
apply(fview, args);
235 "copy from Kokkos g field in FFT",
getRangePolicy(gview, nghostg),
236 KOKKOS_LAMBDA(
const index_array_type& args) {
237 apply(tempFieldg, args - nghostg).real(
apply(gview, args).real());
238 apply(tempFieldg, args - nghostg).imag(
apply(gview, args).imag());
242 this->heffte_m->forward(tempFieldf.data(), tempFieldg.data(), this->workspace_m.data(),
243 heffte::scale::full);
245 this->heffte_m->backward(tempFieldg.data(), tempFieldf.data(), this->workspace_m.data(),
246 heffte::scale::none);
248 throw std::logic_error(
"Only 1:forward and -1:backward are allowed as directions");
253 KOKKOS_LAMBDA(
const index_array_type& args) {
254 apply(fview, args) =
apply(tempFieldf, args - nghostf);
259 KOKKOS_LAMBDA(
const index_array_type& args) {
260 apply(gview, args).real() =
apply(tempFieldg, args - nghostg).real();
261 apply(gview, args).imag() =
apply(tempFieldg, args - nghostg).imag();
265 template <
typename Field>
271 template <
typename Field>
273 static_assert(
Dim == 2 ||
Dim == 3,
"heFFTe only supports 2D and 3D");
274#ifdef Heffte_ENABLE_FFTW
290 auto& tempField = this->tempField;
298 KOKKOS_LAMBDA(
const index_array_type& args) {
299 apply(tempField, args - nghost) =
apply(fview, args);
303 this->heffte_m->forward(tempField.data(), tempField.data(), this->workspace_m.data(),
304 heffte::scale::full);
306 this->heffte_m->backward(tempField.data(), tempField.data(), this->workspace_m.data(),
307 heffte::scale::none);
309 throw std::logic_error(
"Only 1:forward and -1:backward are allowed as directions");
314 KOKKOS_LAMBDA(
const index_array_type& args) {
315 apply(fview, args) =
apply(tempField, args - nghost);
317#ifdef Heffte_ENABLE_FFTW
324 template <
typename Field>
330 template <
typename Field>
332 static_assert(
Dim == 2 ||
Dim == 3,
"heFFTe only supports 2D and 3D");
333#ifdef Heffte_ENABLE_FFTW
349 auto& tempField = this->tempField;
357 KOKKOS_LAMBDA(
const index_array_type& args) {
358 apply(tempField, args - nghost) =
apply(fview, args);
362 this->heffte_m->forward(tempField.data(), tempField.data(), this->workspace_m.data(),
363 heffte::scale::full);
365 this->heffte_m->backward(tempField.data(), tempField.data(), this->workspace_m.data(),
366 heffte::scale::none);
368 throw std::logic_error(
"Only 1:forward and -1:backward are allowed as directions");
373 KOKKOS_LAMBDA(
const index_array_type& args) {
374 apply(fview, args) =
apply(tempField, args - nghost);
376#ifdef Heffte_ENABLE_FFTW
383 template <
typename Field>
389 template <
typename Field>
391 static_assert(
Dim == 2 ||
Dim == 3,
"heFFTe only supports 2D and 3D");
397#ifdef Heffte_ENABLE_FFTW
413 auto& tempField = this->tempField;
421 KOKKOS_LAMBDA(
const index_array_type& args) {
422 apply(tempField, args - nghost) =
apply(fview, args);
426 this->heffte_m->forward(tempField.data(), tempField.data(), this->workspace_m.data(),
427 heffte::scale::full);
429 this->heffte_m->backward(tempField.data(), tempField.data(), this->workspace_m.data(),
430 heffte::scale::none);
432 throw std::logic_error(
"Only 1:forward and -1:backward are allowed as directions");
437 KOKKOS_LAMBDA(
const index_array_type& args) {
438 apply(fview, args) =
apply(tempField, args - nghost);
445#ifdef Heffte_ENABLE_FFTW
constexpr KOKKOS_INLINE_FUNCTION auto first()
KOKKOS_INLINE_FUNCTION constexpr decltype(auto) apply(const View &view, const Coords &coords)
RangePolicy< View::rank, typenameView::execution_space, PolicyArgs... >::policy_type getRangePolicy(const View &view, int shift=0)
void parallel_for(const std::string &name, const ExecPolicy &policy, const FunctorType &functor)
std::unique_ptr< mpi::Communicator > Comm
decltype(auto) shrinkView(std::string label, const View &view, int nghost)
typename FFT< heffteBackend >::template buffer_container< BufferType > workspace_t
FieldLayout< Dim > Layout_t
void domainToBounds(const NDIndex< Dim > &domain, std::array< long long, 3 > &low, std::array< long long, 3 > &high)
static constexpr unsigned Dim
void setup(const heffte::box3d< long long > &inbox, const heffte::box3d< long long > &outbox, const ParameterList ¶ms)
std::shared_ptr< FFT< heffteBackend, long long > > heffte_m
void transform(TransformDirection direction, ComplexField &f)
static constexpr unsigned Dim
static constexpr unsigned Dim
typename Field< Complex_t, Dim, typename RealField::Mesh_t, typename RealField::Centering_t, typename RealField::execution_space >::uniform_type ComplexField
void transform(TransformDirection direction, RealField &f, ComplexField &g)
Base::template temp_view_type< ComplexField > tempFieldComplex
FFT(const Layout_t &layoutInput, const Layout_t &layoutOutput, const ParameterList ¶ms)
static constexpr unsigned Dim
void transform(TransformDirection direction, Field &f)
void transform(TransformDirection direction, Field &f)
static constexpr unsigned Dim
void transform(TransformDirection direction, Field &f)
static constexpr unsigned Dim
const Domain_t & getOwned() const
const NDIndex_t & getLocalNDIndex() const
KOKKOS_INLINE_FUNCTION unsigned size() const noexcept
KOKKOS_INLINE_FUNCTION Vector< int, Dim > first() const
KOKKOS_INLINE_FUNCTION Vector< size_t, Dim > length() const
::ippl::Vector< index_type, Dim > index_array_type
T get(const std::string &key) const