IPPL (Independent Parallel Particle Layer)
IPPL
Loading...
Searching...
No Matches
FFT.hpp
Go to the documentation of this file.
1//
2// Class FFT
3// The FFT class performs complex-to-complex,
4// real-to-complex on IPPL Fields.
5// FFT is templated on the type of transform to be performed,
6// the dimensionality of the Field to transform, and the
7// floating-point precision type of the Field (float or double).
8// Currently, we use heffte for taking the transforms and the class FFT
9// serves as an interface between IPPL and heffte. In making this interface,
10// we have referred Cabana library.
11// https://github.com/ECP-copa/Cabana.
12//
13// Copyright (c) 2021, Sriramkrishnan Muralikrishnan,
14// Paul Scherrer Institut, Villigen PSI, Switzerland
15// All rights reserved
16//
17// This file is part of IPPL.
18//
22
23#include "Utility/IpplTimings.h"
24
25#include "Field/BareField.h"
26
28
29namespace ippl {
30
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;
35
36 const NDIndex<Dim> lDom = layout.getLocalNDIndex();
37 domainToBounds(lDom, low, high);
38
39 heffte::box3d<long long> inbox = {low, high};
40 heffte::box3d<long long> outbox = {low, high};
41
42 setup(inbox, outbox, params);
43 }
44
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) {
49 low.fill(0);
50 high.fill(0);
51
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);
59 }
60 }
61
65 template <typename Field, template <typename...> class FFT, typename Backend, typename T>
66 void FFTBase<Field, FFT, Backend, T>::setup(const heffte::box3d<long long>& inbox,
67 const heffte::box3d<long long>& outbox,
68 const ParameterList& params) {
69 heffte::plan_options heffteOptions = heffte::default_options<heffteBackend>();
70
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");
76#endif
77
78 switch (params.get<int>("comm")) {
79 case a2a:
80 heffteOptions.algorithm = heffte::reshape_algorithm::alltoall;
81 break;
82 case a2av:
83 heffteOptions.algorithm = heffte::reshape_algorithm::alltoallv;
84 break;
85 case p2p:
86 heffteOptions.algorithm = heffte::reshape_algorithm::p2p;
87 break;
88 case p2p_pl:
89 heffteOptions.algorithm = heffte::reshape_algorithm::p2p_plined;
90 break;
91 default:
92 throw IpplException("FFT::setup", "Unrecognized heffte communication type");
93 }
94 }
95
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);
99 } else {
100 heffte_m = std::make_shared<FFT<heffteBackend, long long>>(
101 inbox, outbox, params.get<int>("r2c_direction"), Comm->getCommunicator(),
102 heffteOptions);
103 }
104
105 // heffte::gpu::device_set(Comm->rank() % heffte::gpu::device_count());
106 if (workspace_m.size() < heffte_m->size_workspace()) {
107 workspace_m = workspace_t(heffte_m->size_workspace());
108 }
109 }
110
111 template <typename ComplexField>
113 this->transform(FORWARD, f);
114 this->transform(BACKWARD, f);
115 }
116
117 template <typename ComplexField>
119 static_assert(Dim == 2 || Dim == 3, "heFFTe only supports 2D and 3D");
120
121 auto fview = f.getView();
122 const int nghost = f.getNghost();
123
131 auto& tempField = this->tempField;
132 if (tempField.size() != f.getOwned().size()) {
133 tempField = detail::shrinkView("tempField", fview, nghost);
134 }
135
136 using index_array_type = typename RangePolicy<Dim>::index_array_type;
138 "copy from Kokkos FFT", getRangePolicy(fview, nghost),
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());
142 });
143
144 if (direction == FORWARD) {
145 this->heffte_m->forward(tempField.data(), tempField.data(), this->workspace_m.data(),
146 heffte::scale::full);
147 } else if (direction == BACKWARD) {
148 this->heffte_m->backward(tempField.data(), tempField.data(), this->workspace_m.data(),
149 heffte::scale::none);
150 } else {
151 throw std::logic_error("Only 1:forward and -1:backward are allowed as directions");
152 }
153
155 "copy to Kokkos FFT", getRangePolicy(fview, nghost),
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();
159 });
160 }
161
162 //========================================================================
163 // FFT RCTransform Constructors
164 //========================================================================
165
170
171 template <typename RealField>
172 FFT<RCTransform, RealField>::FFT(const Layout_t& layoutInput, const Layout_t& layoutOutput,
173 const ParameterList& params) {
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;
183
184 const NDIndex<Dim>& lDomInput = layoutInput.getLocalNDIndex();
185 const NDIndex<Dim>& lDomOutput = layoutOutput.getLocalNDIndex();
186
187 this->domainToBounds(lDomInput, lowInput, highInput);
188 this->domainToBounds(lDomOutput, lowOutput, highOutput);
189
190 heffte::box3d<long long> inbox = {lowInput, highInput};
191 heffte::box3d<long long> outbox = {lowOutput, highOutput};
192
193 this->setup(inbox, outbox, params);
194 }
195
196 template <typename RealField>
198 this->transform(FORWARD, f, g);
199 this->transform(BACKWARD, f, g);
200 }
201
202 template <typename RealField>
204 ComplexField& g) {
205 static_assert(Dim == 2 || Dim == 3, "heFFTe only supports 2D and 3D");
206
207 auto fview = f.getView();
208 auto gview = g.getView();
209 const int nghostf = f.getNghost();
210 const int nghostg = g.getNghost();
211
219 auto& tempFieldf = this->tempField;
220 auto& tempFieldg = this->tempFieldComplex;
221 if (tempFieldf.size() != f.getOwned().size()) {
222 tempFieldf = detail::shrinkView("tempFieldf", fview, nghostf);
223 }
224 if (tempFieldg.size() != g.getOwned().size()) {
225 tempFieldg = detail::shrinkView("tempFieldg", gview, nghostg);
226 }
227
228 using index_array_type = typename RangePolicy<Dim>::index_array_type;
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);
233 });
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());
239 });
240
241 if (direction == FORWARD) {
242 this->heffte_m->forward(tempFieldf.data(), tempFieldg.data(), this->workspace_m.data(),
243 heffte::scale::full);
244 } else if (direction == BACKWARD) {
245 this->heffte_m->backward(tempFieldg.data(), tempFieldf.data(), this->workspace_m.data(),
246 heffte::scale::none);
247 } else {
248 throw std::logic_error("Only 1:forward and -1:backward are allowed as directions");
249 }
250
252 "copy to Kokkos f field FFT", getRangePolicy(fview, nghostf),
253 KOKKOS_LAMBDA(const index_array_type& args) {
254 apply(fview, args) = apply(tempFieldf, args - nghostf);
255 });
256
258 "copy to Kokkos g field FFT", getRangePolicy(gview, nghostg),
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();
262 });
263 }
264
265 template <typename Field>
267 this->transform(FORWARD, f);
268 this->transform(BACKWARD, f);
269 }
270
271 template <typename Field>
273 static_assert(Dim == 2 || Dim == 3, "heFFTe only supports 2D and 3D");
274#ifdef Heffte_ENABLE_FFTW
275 if (direction == FORWARD) {
276 f = f / 8.0;
277 }
278#endif
279
280 auto fview = f.getView();
281 const int nghost = f.getNghost();
282
290 auto& tempField = this->tempField;
291 if (tempField.size() != f.getOwned().size()) {
292 tempField = detail::shrinkView("tempField", fview, nghost);
293 }
294
295 using index_array_type = typename RangePolicy<Dim>::index_array_type;
297 "copy from Kokkos FFT", getRangePolicy(fview, nghost),
298 KOKKOS_LAMBDA(const index_array_type& args) {
299 apply(tempField, args - nghost) = apply(fview, args);
300 });
301
302 if (direction == FORWARD) {
303 this->heffte_m->forward(tempField.data(), tempField.data(), this->workspace_m.data(),
304 heffte::scale::full);
305 } else if (direction == BACKWARD) {
306 this->heffte_m->backward(tempField.data(), tempField.data(), this->workspace_m.data(),
307 heffte::scale::none);
308 } else {
309 throw std::logic_error("Only 1:forward and -1:backward are allowed as directions");
310 }
311
313 "copy to Kokkos FFT", getRangePolicy(fview, nghost),
314 KOKKOS_LAMBDA(const index_array_type& args) {
315 apply(fview, args) = apply(tempField, args - nghost);
316 });
317#ifdef Heffte_ENABLE_FFTW
318 if (direction == BACKWARD) {
319 f = f * 8.0;
320 }
321#endif
322 }
323
324 template <typename Field>
326 this->transform(FORWARD, f);
327 this->transform(BACKWARD, f);
328 }
329
330 template <typename Field>
332 static_assert(Dim == 2 || Dim == 3, "heFFTe only supports 2D and 3D");
333#ifdef Heffte_ENABLE_FFTW
334 if (direction == FORWARD) {
335 f = f / 8.0;
336 }
337#endif
338
339 auto fview = f.getView();
340 const int nghost = f.getNghost();
341
349 auto& tempField = this->tempField;
350 if (tempField.size() != f.getOwned().size()) {
351 tempField = detail::shrinkView("tempField", fview, nghost);
352 }
353
354 using index_array_type = typename RangePolicy<Dim>::index_array_type;
356 "copy from Kokkos FFT", getRangePolicy(fview, nghost),
357 KOKKOS_LAMBDA(const index_array_type& args) {
358 apply(tempField, args - nghost) = apply(fview, args);
359 });
360
361 if (direction == FORWARD) {
362 this->heffte_m->forward(tempField.data(), tempField.data(), this->workspace_m.data(),
363 heffte::scale::full);
364 } else if (direction == BACKWARD) {
365 this->heffte_m->backward(tempField.data(), tempField.data(), this->workspace_m.data(),
366 heffte::scale::none);
367 } else {
368 throw std::logic_error("Only 1:forward and -1:backward are allowed as directions");
369 }
370
372 "copy to Kokkos FFT", getRangePolicy(fview, nghost),
373 KOKKOS_LAMBDA(const index_array_type& args) {
374 apply(fview, args) = apply(tempField, args - nghost);
375 });
376#ifdef Heffte_ENABLE_FFTW
377 if (direction == BACKWARD) {
378 f = f * 8.0;
379 }
380#endif
381 }
382
383 template <typename Field>
385 this->transform(FORWARD, f);
386 this->transform(BACKWARD, f);
387 }
388
389 template <typename Field>
391 static_assert(Dim == 2 || Dim == 3, "heFFTe only supports 2D and 3D");
392
397#ifdef Heffte_ENABLE_FFTW
398 if (direction == FORWARD) {
399 f = f / 8.0;
400 }
401#endif
402
403 auto fview = f.getView();
404 const int nghost = f.getNghost();
405
413 auto& tempField = this->tempField;
414 if (tempField.size() != f.getOwned().size()) {
415 tempField = detail::shrinkView("tempField", fview, nghost);
416 }
417
418 using index_array_type = typename RangePolicy<Dim>::index_array_type;
420 "copy from Kokkos FFT", getRangePolicy(fview, nghost),
421 KOKKOS_LAMBDA(const index_array_type& args) {
422 apply(tempField, args - nghost) = apply(fview, args);
423 });
424
425 if (direction == FORWARD) {
426 this->heffte_m->forward(tempField.data(), tempField.data(), this->workspace_m.data(),
427 heffte::scale::full);
428 } else if (direction == BACKWARD) {
429 this->heffte_m->backward(tempField.data(), tempField.data(), this->workspace_m.data(),
430 heffte::scale::none);
431 } else {
432 throw std::logic_error("Only 1:forward and -1:backward are allowed as directions");
433 }
434
436 "copy to Kokkos FFT", getRangePolicy(fview, nghost),
437 KOKKOS_LAMBDA(const index_array_type& args) {
438 apply(fview, args) = apply(tempField, args - nghost);
439 });
440
445#ifdef Heffte_ENABLE_FFTW
446 if (direction == BACKWARD) {
447 f = f * 8.0;
448 }
449#endif
450 }
451
452} // namespace ippl
453
454// vi: set et ts=4 sw=4 sts=4:
455// Local Variables:
456// mode:c
457// c-basic-offset: 4
458// indent-tabs-mode: nil
459// require-final-newline: nil
460// End:
constexpr KOKKOS_INLINE_FUNCTION auto first()
Definition AbsorbingBC.h:10
Definition Archive.h:20
TransformDirection
Definition FFT.h:62
@ FORWARD
Definition FFT.h:63
@ BACKWARD
Definition FFT.h:64
@ p2p_pl
Definition FFT.h:59
@ a2av
Definition FFT.h:56
@ p2p
Definition FFT.h:58
@ a2a
Definition FFT.h:57
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
Definition Ippl.h:22
decltype(auto) shrinkView(std::string label, const View &view, int nghost)
Definition ViewUtils.h:122
typename FFT< heffteBackend >::template buffer_container< BufferType > workspace_t
Definition FFT.h:149
FieldLayout< Dim > Layout_t
Definition FFT.h:150
void domainToBounds(const NDIndex< Dim > &domain, std::array< long long, 3 > &low, std::array< long long, 3 > &high)
Definition FFT.hpp:46
workspace_t workspace_m
Definition FFT.h:164
FFTBase()=default
static constexpr unsigned Dim
Definition FFT.h:145
void setup(const heffte::box3d< long long > &inbox, const heffte::box3d< long long > &outbox, const ParameterList &params)
Definition FFT.hpp:66
std::shared_ptr< FFT< heffteBackend, long long > > heffte_m
Definition FFT.h:163
void transform(TransformDirection direction, ComplexField &f)
Definition FFT.hpp:118
static constexpr unsigned Dim
Definition FFT.h:192
static constexpr unsigned Dim
Definition FFT.h:222
typename Field< Complex_t, Dim, typename RealField::Mesh_t, typename RealField::Centering_t, typename RealField::execution_space >::uniform_type ComplexField
Definition FFT.h:229
void transform(TransformDirection direction, RealField &f, ComplexField &g)
Definition FFT.hpp:203
Base::template temp_view_type< ComplexField > tempFieldComplex
Definition FFT.h:256
FFT(const Layout_t &layoutInput, const Layout_t &layoutOutput, const ParameterList &params)
Definition FFT.hpp:172
static constexpr unsigned Dim
Definition FFT.h:264
void transform(TransformDirection direction, Field &f)
Definition FFT.hpp:272
void transform(TransformDirection direction, Field &f)
Definition FFT.hpp:331
static constexpr unsigned Dim
Definition FFT.h:289
void transform(TransformDirection direction, Field &f)
Definition FFT.hpp:390
static constexpr unsigned Dim
Definition FFT.h:314
int getNghost() const
Definition BareField.h:125
view_type & getView()
Definition BareField.h:169
const Domain_t & getOwned() const
Definition BareField.h:117
const NDIndex_t & getLocalNDIndex() const
KOKKOS_INLINE_FUNCTION unsigned size() const noexcept
Definition NDIndex.hpp:43
KOKKOS_INLINE_FUNCTION Vector< int, Dim > first() const
Definition NDIndex.hpp:170
KOKKOS_INLINE_FUNCTION Vector< size_t, Dim > length() const
Definition NDIndex.hpp:162
::ippl::Vector< index_type, Dim > index_array_type
T get(const std::string &key) const