Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_FFT_R2X.H
Go to the documentation of this file.
1#ifndef AMREX_FFT_R2X_H_
2#define AMREX_FFT_R2X_H_
3#include <AMReX_Config.H>
4
5#include <AMReX_MultiFab.H>
6#include <AMReX_FFT_Helper.H>
7#include <algorithm>
8#include <numeric>
9#include <tuple>
10
11namespace amrex::FFT
12{
13
20template <typename T> class Poisson;
21template <typename T> class PoissonHybrid;
22
30template <typename T = Real>
31class R2X
32{
33public:
34 using MF = std::conditional_t<std::is_same_v<T,Real>,
37
38 template <typename U> friend class Poisson;
39 template <typename U> friend class PoissonHybrid;
40
48 R2X (Box const& domain,
49 Array<std::pair<Boundary,Boundary>,AMREX_SPACEDIM> const& bc,
50 Info const& info = Info{});
51
55 ~R2X ();
56
57 R2X (R2X const&) = delete;
58 R2X (R2X &&) = delete;
59 R2X& operator= (R2X const&) = delete;
60 R2X& operator= (R2X &&) = delete;
61
67 [[nodiscard]] T scalingFactor () const;
68
69 template <typename F>
77 void forwardThenBackward (MF const& inmf, MF& outmf, F const& post_forward);
78
79 // public for cuda
80 template <int dim, typename FAB, typename F>
87 void post_forward_doit (FAB* fab, F const& f);
88
89 // private function made public for cuda
90 template <typename F>
100 void forwardThenBackward_doit_0 (MF const& inmf, MF& outmf, F const& post_forward,
101 IntVect const& ngout = IntVect(0),
102 Periodicity const& period = Periodicity::NonPeriodic());
103 template <typename F>
113 void forwardThenBackward_doit_1 (MF const& inmf, MF& outmf, F const& post_forward,
114 IntVect const& ngout = IntVect(0),
115 Periodicity const& period = Periodicity::NonPeriodic());
116
117private:
118
125 void forward (MF const& inmf, MF& outmf);
132 void forward (MF const& inmf, cMF& outmf);
138 void forward (MF const& inmf);
147 void backward (MF const& inmf, MF& outmf, IntVect const& ngout,
148 Periodicity const& period);
157 void backward (cMF const& inmf, MF& outmf, IntVect const& ngout,
158 Periodicity const& period);
162 void backward ();
163
164 Box m_dom_0;
165 Array<std::pair<Boundary,Boundary>,AMREX_SPACEDIM> m_bc;
166
167 Plan<T> m_fft_fwd_x{};
168 Plan<T> m_fft_bwd_x{};
169 Plan<T> m_fft_fwd_y{};
170 Plan<T> m_fft_bwd_y{};
171 Plan<T> m_fft_fwd_z{};
172 Plan<T> m_fft_bwd_z{};
173
174 std::unique_ptr<MultiBlockCommMetaData> m_cmd_cx2cy;
175 std::unique_ptr<MultiBlockCommMetaData> m_cmd_rx2ry;
176 std::unique_ptr<MultiBlockCommMetaData> m_cmd_cy2cz;
177 std::unique_ptr<MultiBlockCommMetaData> m_cmd_ry2rz;
178
179 std::unique_ptr<MultiBlockCommMetaData> m_cmd_cy2cx;
180 std::unique_ptr<MultiBlockCommMetaData> m_cmd_ry2rx;
181 std::unique_ptr<MultiBlockCommMetaData> m_cmd_cz2cy;
182 std::unique_ptr<MultiBlockCommMetaData> m_cmd_rz2ry;
183
184 Swap01 m_dtos_x2y{};
185 Swap01 m_dtos_y2x{};
186 Swap02 m_dtos_y2z{};
187 Swap02 m_dtos_z2y{};
188
189 MF m_rx;
190 MF m_ry;
191 MF m_rz;
192 cMF m_cx;
193 cMF m_cy;
194 cMF m_cz;
195
196 std::unique_ptr<char,DataDeleter> m_data_1;
197 std::unique_ptr<char,DataDeleter> m_data_2;
198
199 Box m_dom_rx;
200 Box m_dom_ry;
201 Box m_dom_rz;
202 Box m_dom_cx;
203 Box m_dom_cy;
204 Box m_dom_cz;
205
206 std::unique_ptr<R2X<T>> m_r2x_sub;
207 detail::SubHelper m_sub_helper;
208
209 Info m_info;
210};
211
212template <typename T>
213R2X<T>::R2X (Box const& domain,
214 Array<std::pair<Boundary,Boundary>,AMREX_SPACEDIM> const& bc,
215 Info const& info)
216 : m_dom_0(domain),
217 m_bc(bc),
218 m_sub_helper(domain),
219 m_info(info)
220{
221 BL_PROFILE("FFT::R2X");
222
223 static_assert(std::is_same_v<float,T> || std::is_same_v<double,T>);
224
225 AMREX_ALWAYS_ASSERT((m_dom_0.numPts() > 1) && (m_info.batch_size == 1));
227 m_dom_0.smallEnd() == IntVect(AMREX_D_DECL(0,0,0)),
228 "FFT::R2X currently requires domain.smallEnd() == 0");
229#if (AMREX_SPACEDIM == 2)
231#else
232 if (m_info.twod_mode) {
233 AMREX_ALWAYS_ASSERT((int(domain.length(0) > 1) +
234 int(domain.length(1) > 1) +
235 int(domain.length(2) > 1)) >= 2);
236 }
237#endif
238
239 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
240 if (bc[idim].first == Boundary::periodic ||
241 bc[idim].second == Boundary::periodic) {
242 AMREX_ALWAYS_ASSERT(bc[idim].first == bc[idim].second);
243 }
244 }
245
246 {
247 Box subbox = m_sub_helper.make_box(m_dom_0);
248 if (subbox.size() != m_dom_0.size()) {
249 m_r2x_sub = std::make_unique<R2X<T>>
250 (subbox, m_sub_helper.make_array(bc), info);
251 return;
252 }
253 }
254
255 int myproc = ParallelContext::MyProcSub();
256 int nprocs = std::min(ParallelContext::NProcsSub(), m_info.nprocs);
257
258 //
259 // make data containers
260 //
261
262 m_dom_rx = m_dom_0;
263 auto bax = amrex::decompose(m_dom_rx, nprocs, {AMREX_D_DECL(false,true,true)});
264 DistributionMapping dmx = detail::make_iota_distromap(bax.size());
265 m_rx.define(bax, dmx, 1, 0, MFInfo().SetAlloc(false));
266
267 // x-direction
268 if (bc[0].first == Boundary::periodic) {
269 // x-fft: r2c(m_rx->m_cx)
270 m_dom_cx = Box(IntVect(0), IntVect(AMREX_D_DECL(domain.length(0)/2,
271 domain.bigEnd(1),
272 domain.bigEnd(2))));
273 BoxList bl = bax.boxList();
274 for (auto & b : bl) {
275 b.setBig(0, m_dom_cx.bigEnd(0));
276 }
277 BoxArray cbax(std::move(bl));
278 m_cx.define(cbax, dmx, 1, 0, MFInfo().SetAlloc(false));
279 } // else: x-fft: r2r(m_rx)
280
281#if (AMREX_SPACEDIM >= 2)
282 if (domain.length(1) > 1 && !m_info.oned_mode) {
283 if (! m_cx.empty()) {
284 // copy(m_cx->m_cy)
285 m_dom_cy = Box(IntVect(0), IntVect(AMREX_D_DECL(m_dom_cx.bigEnd(1),
286 m_dom_cx.bigEnd(0),
287 m_dom_cx.bigEnd(2))));
288 auto ba = amrex::decompose(m_dom_cy, nprocs, {AMREX_D_DECL(false,true,true)});
290 if (ba.size() == m_cx.size()) {
291 dm = m_cx.DistributionMap();
292 } else {
293 dm = detail::make_iota_distromap(ba.size());
294 }
295 m_cy.define(ba, dm, 1, 0, MFInfo().SetAlloc(false));
296 // if bc[1] is periodic:
297 // c2c(m_cy->m_cy)
298 // else:
299 // r2r(m_cy.re) & r2r(m_cy.im)
300 } else {
301 // copy(m_rx->m_ry)
302 m_dom_ry = Box(IntVect(0), IntVect(AMREX_D_DECL(m_dom_rx.bigEnd(1),
303 m_dom_rx.bigEnd(0),
304 m_dom_rx.bigEnd(2))));
305 auto ba = amrex::decompose(m_dom_ry, nprocs, {AMREX_D_DECL(false,true,true)});
307 if (ba.size() == m_rx.size()) {
308 dm = m_rx.DistributionMap();
309 } else {
310 dm = detail::make_iota_distromap(ba.size());
311 }
312 m_ry.define(ba, dm, 1, 0, MFInfo().SetAlloc(false));
313 // if bc[1] is periodic:
314 // r2c(m_ry->m_cy)
315 // else:
316 // r2r(m_ry)
317 if (bc[1].first == Boundary::periodic) {
318 m_dom_cy = Box(IntVect(0), IntVect(AMREX_D_DECL(m_dom_ry.length(0)/2,
319 m_dom_ry.bigEnd(1),
320 m_dom_ry.bigEnd(2))));
321 BoxList bl = ba.boxList();
322 for (auto & b : bl) {
323 b.setBig(0, m_dom_cy.bigEnd(0));
324 }
325 BoxArray cba(std::move(bl));
326 m_cy.define(cba, dm, 1, 0, MFInfo().SetAlloc(false));
327 }
328 }
329 }
330#endif
331
332#if (AMREX_SPACEDIM == 3)
333 if (domain.length(2) > 1 && !m_info.twod_mode) {
334 if (! m_cy.empty()) {
335 // copy(m_cy, m_cz)
336 m_dom_cz = Box(IntVect(0), IntVect(AMREX_D_DECL(m_dom_cy.bigEnd(2),
337 m_dom_cy.bigEnd(1),
338 m_dom_cy.bigEnd(0))));
339 auto ba = amrex::decompose(m_dom_cz, nprocs, {AMREX_D_DECL(false,true,true)});
341 if (ba.size() == m_cy.size()) {
342 dm = m_cy.DistributionMap();
343 } else {
344 dm = detail::make_iota_distromap(ba.size());
345 }
346 m_cz.define(ba, dm, 1, 0, MFInfo().SetAlloc(false));
347 // if bc[2] is periodic:
348 // c2c(m_cz->m_cz)
349 // else:
350 // r2r(m_cz.re) & r2r(m_cz.im)
351 } else {
352 // copy(m_ry, m_rz)
353 m_dom_rz = Box(IntVect(0), IntVect(AMREX_D_DECL(m_dom_ry.bigEnd(2),
354 m_dom_ry.bigEnd(1),
355 m_dom_ry.bigEnd(0))));
356 auto ba = amrex::decompose(m_dom_rz, nprocs, {AMREX_D_DECL(false,true,true)});
358 if (ba.size() == m_ry.size()) {
359 dm = m_ry.DistributionMap();
360 } else {
361 dm = detail::make_iota_distromap(ba.size());
362 }
363 m_rz.define(ba, dm, 1, 0, MFInfo().SetAlloc(false));
364 // if bc[2] is periodic:
365 // r2c(m_rz->m_cz)
366 // else:
367 // r2r(m_rz)
368 if (bc[2].first == Boundary::periodic) {
369 m_dom_cz = Box(IntVect(0), IntVect(AMREX_D_DECL(m_dom_rz.length(0)/2,
370 m_dom_rz.bigEnd(1),
371 m_dom_rz.bigEnd(2))));
372 BoxList bl = ba.boxList();
373 for (auto & b : bl) {
374 b.setBig(0, m_dom_cz.bigEnd(0));
375 }
376 BoxArray cba(std::move(bl));
377 m_cz.define(cba, dm, 1, 0, MFInfo().SetAlloc(false));
378 }
379 }
380 }
381#endif
382
383 // There are several different execution paths.
384 //
385 // (1) x-r2c(m_rx->m_cx), copy(m_cx->m_cy), y-fft(m_cy),
386 // copy(m_cy->m_cz), z-fft(m_cz)
387 // In this case, we have m_rx, m_cx, m_cy, & m_cz.
388 // we can alias(m_rx,m_cy) and alias(m_cx,m_cz).
389 //
390 // (2) x_r2r(m_rx), copy(m_rx->m_ry), y-r2c(m_ry->m_cy),
391 // copy(m_cy->m_cz), z-fft(m_cz)
392 // In this case, we have m_rx, m_ry, m_cy, & m_cz.
393 // We can alias(m_rx,m_cy) and alias(m_ry,m_cz).
394 //
395 // (3) x_r2r(m_rx), copy(m_rx->m_ry), y-r2r(m_ry),
396 // copy(m_ry->m_rz), z-r2c(m_rz->m_rz)
397 // In this case, we have m_rx, m_ry, m_rz, & m_cz
398 // We can alias(m_rx,m_rz) and alias(m_ry,m_cz)
399 //
400 // (4) x_r2r(m_rx), copy(m_rx->m_ry), y-r2r(m_ry),
401 // copy(m_ry->m_rz), z-r2r(m_rz)
402 // In this case, we have m_rx, m_ry, & m_rz.
403 // We can alias(m_rx,m_rz).
404
405 if (! m_cx.empty()) {
406 m_data_1 = detail::make_mfs_share(m_rx, m_cy);
407 m_data_2 = detail::make_mfs_share(m_cx, m_cz);
408 } else if (! m_cy.empty()) {
409 m_data_1 = detail::make_mfs_share(m_rx, m_cy);
410 m_data_2 = detail::make_mfs_share(m_ry, m_cz);
411 } else if (! m_cz.empty()) {
412 m_data_1 = detail::make_mfs_share(m_rx, m_rz);
413 m_data_2 = detail::make_mfs_share(m_ry, m_cz);
414 } else {
415 m_data_1 = detail::make_mfs_share(m_rx, m_rz);
416 m_data_2 = detail::make_mfs_share(m_ry, m_cz); // It's okay m_cz is empty.
417 }
418
419 //
420 // make copiers
421 //
422
423#if (AMREX_SPACEDIM >= 2)
424 if (!m_cy.empty() || !m_ry.empty()) {
425 if (! m_cx.empty()) {
426 // copy(m_cx->m_cy)
427 m_cmd_cx2cy = std::make_unique<MultiBlockCommMetaData>
428 (m_cy, m_dom_cy, m_cx, IntVect(0), m_dtos_x2y);
429 m_cmd_cy2cx = std::make_unique<MultiBlockCommMetaData>
430 (m_cx, m_dom_cx, m_cy, IntVect(0), m_dtos_y2x);
431 } else {
432 // copy(m_rx->m_ry)
433 m_cmd_rx2ry = std::make_unique<MultiBlockCommMetaData>
434 (m_ry, m_dom_ry, m_rx, IntVect(0), m_dtos_x2y);
435 m_cmd_ry2rx = std::make_unique<MultiBlockCommMetaData>
436 (m_rx, m_dom_rx, m_ry, IntVect(0), m_dtos_y2x);
437 }
438 }
439#endif
440
441#if (AMREX_SPACEDIM == 3)
442 if (!m_cz.empty() || !m_rz.empty()) {
443 if (! m_cy.empty()) {
444 // copy(m_cy, m_cz)
445 m_cmd_cy2cz = std::make_unique<MultiBlockCommMetaData>
446 (m_cz, m_dom_cz, m_cy, IntVect(0), m_dtos_y2z);
447 m_cmd_cz2cy = std::make_unique<MultiBlockCommMetaData>
448 (m_cy, m_dom_cy, m_cz, IntVect(0), m_dtos_z2y);
449 } else {
450 // copy(m_ry, m_rz)
451 m_cmd_ry2rz = std::make_unique<MultiBlockCommMetaData>
452 (m_rz, m_dom_rz, m_ry, IntVect(0), m_dtos_y2z);
453 m_cmd_rz2ry = std::make_unique<MultiBlockCommMetaData>
454 (m_ry, m_dom_ry, m_rz, IntVect(0), m_dtos_z2y);
455 }
456 }
457#endif
458
459 //
460 // make plans
461 //
462
463 using VendorComplex = typename Plan<T>::VendorComplex;
464
465 if (myproc < m_rx.size())
466 {
467 Box const& box = m_rx.box(myproc);
468 auto* pf = m_rx[myproc].dataPtr();
469 if (bc[0].first == Boundary::periodic) {
470 auto* pb = (VendorComplex*) m_cx[myproc].dataPtr();
471 m_fft_fwd_x.template init_r2c<Direction::forward>(box, pf, pb);
472#if defined(AMREX_USE_SYCL)
473 m_fft_bwd_x = m_fft_fwd_x;
474#else
475 m_fft_bwd_x.template init_r2c<Direction::backward>(box, pf, pb);
476#endif
477 } else {
478 m_fft_fwd_x.template init_r2r<Direction::forward>(box, pf, bc[0]);
479#if defined(AMREX_USE_GPU)
480 if ((bc[0].first == Boundary::even && bc[0].second == Boundary::odd) ||
481 (bc[0].first == Boundary::odd && bc[0].second == Boundary::even)) {
482 m_fft_bwd_x = m_fft_fwd_x;
483 } else
484#endif
485 {
486 m_fft_bwd_x.template init_r2r<Direction::backward>(box, pf, bc[0]);
487 }
488 }
489 }
490
491#if (AMREX_SPACEDIM >= 2)
492 if (m_ry.empty() && m_bc[1].first == Boundary::periodic) {
493 if (myproc < m_cy.size()) {
494 Box const& box = m_cy.box(myproc);
495 auto* p = (VendorComplex *)m_cy[myproc].dataPtr();
496 m_fft_fwd_y.template init_c2c<Direction::forward>(box, p);
497#if defined(AMREX_USE_SYCL)
498 m_fft_bwd_y = m_fft_fwd_y;
499#else
500 m_fft_bwd_y.template init_c2c<Direction::backward>(box, p);
501#endif
502 }
503 } else if (!m_ry.empty() && m_bc[1].first == Boundary::periodic) {
504 if (myproc < m_ry.size()) {
505 Box const& box = m_ry.box(myproc);
506 auto* pr = m_ry[myproc].dataPtr();
507 auto* pc = (VendorComplex*)m_cy[myproc].dataPtr();
508 m_fft_fwd_y.template init_r2c<Direction::forward>(box, pr, pc);
509#if defined(AMREX_USE_SYCL)
510 m_fft_bwd_y = m_fft_fwd_y;
511#else
512 m_fft_bwd_y.template init_r2c<Direction::backward>(box, pr, pc);
513#endif
514 }
515 } else if (!m_cy.empty()) {
516 if (myproc < m_cy.size()) {
517 Box const& box = m_cy.box(myproc);
518 auto* p = (VendorComplex*) m_cy[myproc].dataPtr();
519 m_fft_fwd_y.template init_r2r<Direction::forward>(box, p, bc[1]);
520#if defined(AMREX_USE_GPU)
521 if ((bc[1].first == Boundary::even && bc[1].second == Boundary::odd) ||
522 (bc[1].first == Boundary::odd && bc[1].second == Boundary::even)) {
523 m_fft_bwd_y = m_fft_fwd_y;
524 } else
525#endif
526 {
527 m_fft_bwd_y.template init_r2r<Direction::backward>(box, p, bc[1]);
528 }
529 }
530 } else {
531 if (myproc < m_ry.size()) {
532 Box const& box = m_ry.box(myproc);
533 auto* p = m_ry[myproc].dataPtr();
534 m_fft_fwd_y.template init_r2r<Direction::forward>(box, p, bc[1]);
535#if defined(AMREX_USE_GPU)
536 if ((bc[1].first == Boundary::even && bc[1].second == Boundary::odd) ||
537 (bc[1].first == Boundary::odd && bc[1].second == Boundary::even)) {
538 m_fft_bwd_y = m_fft_fwd_y;
539 } else
540#endif
541 {
542 m_fft_bwd_y.template init_r2r<Direction::backward>(box, p, bc[1]);
543 }
544 }
545 }
546#endif
547
548#if (AMREX_SPACEDIM == 3)
549 if (m_rz.empty() && m_bc[2].first == Boundary::periodic) {
550 if (myproc < m_cz.size()) {
551 Box const& box = m_cz.box(myproc);
552 auto* p = (VendorComplex*)m_cz[myproc].dataPtr();
553 m_fft_fwd_z.template init_c2c<Direction::forward>(box, p);
554#if defined(AMREX_USE_SYCL)
555 m_fft_bwd_z = m_fft_fwd_z;
556#else
557 m_fft_bwd_z.template init_c2c<Direction::backward>(box, p);
558#endif
559 }
560 } else if (!m_rz.empty() && m_bc[2].first == Boundary::periodic) {
561 if (myproc < m_rz.size()) {
562 Box const& box = m_rz.box(myproc);
563 auto* pr = m_rz[myproc].dataPtr();
564 auto* pc = (VendorComplex*)m_cz[myproc].dataPtr();
565 m_fft_fwd_z.template init_r2c<Direction::forward>(box, pr, pc);
566#if defined(AMREX_USE_SYCL)
567 m_fft_bwd_z = m_fft_fwd_z;
568#else
569 m_fft_bwd_z.template init_r2c<Direction::backward>(box, pr, pc);
570#endif
571 }
572 } else if (!m_cz.empty()) {
573 if (myproc < m_cz.size()) {
574 Box const& box = m_cz.box(myproc);
575 auto* p = (VendorComplex*) m_cz[myproc].dataPtr();
576 m_fft_fwd_z.template init_r2r<Direction::forward>(box, p, bc[2]);
577#if defined(AMREX_USE_GPU)
578 if ((bc[2].first == Boundary::even && bc[2].second == Boundary::odd) ||
579 (bc[2].first == Boundary::odd && bc[2].second == Boundary::even)) {
580 m_fft_bwd_z = m_fft_fwd_z;
581 } else
582#endif
583 {
584 m_fft_bwd_z.template init_r2r<Direction::backward>(box, p, bc[2]);
585 }
586 }
587 } else {
588 if (myproc < m_rz.size()) {
589 Box const& box = m_rz.box(myproc);
590 auto* p = m_rz[myproc].dataPtr();
591 m_fft_fwd_z.template init_r2r<Direction::forward>(box, p, bc[2]);
592#if defined(AMREX_USE_GPU)
593 if ((bc[2].first == Boundary::even && bc[2].second == Boundary::odd) ||
594 (bc[2].first == Boundary::odd && bc[2].second == Boundary::even)) {
595 m_fft_bwd_z = m_fft_fwd_z;
596 } else
597#endif
598 {
599 m_fft_bwd_z.template init_r2r<Direction::backward>(box, p, bc[2]);
600 }
601 }
602 }
603#endif
604}
605
606template <typename T>
608{
609 if (m_fft_bwd_x.plan != m_fft_fwd_x.plan) {
610 m_fft_bwd_x.destroy();
611 }
612 if (m_fft_bwd_y.plan != m_fft_fwd_y.plan) {
613 m_fft_bwd_y.destroy();
614 }
615 if (m_fft_bwd_z.plan != m_fft_fwd_z.plan) {
616 m_fft_bwd_z.destroy();
617 }
618 m_fft_fwd_x.destroy();
619 m_fft_fwd_y.destroy();
620 m_fft_fwd_z.destroy();
621}
622
623template <typename T>
625{
626 Long r = 1;
627 int ndims = m_info.twod_mode ? AMREX_SPACEDIM-1 : AMREX_SPACEDIM;
628#if (AMREX_SPACEDIM == 3)
629 if (m_info.twod_mode && m_dom_0.length(2) == 1) { ndims = 1; };
630#endif
631 for (int idim = 0; idim < ndims; ++idim) {
632 r *= m_dom_0.length(idim);
633 if (m_bc[idim].first != Boundary::periodic && (m_dom_0.length(idim) > 1)) {
634 r *= 2;
635 }
636 }
637 return T(1)/T(r);
638}
639
640template <typename T>
641template <typename F>
642void R2X<T>::forwardThenBackward (MF const& inmf, MF& outmf, F const& post_forward)
643{
645 !m_info.twod_mode,
646 "FFT::R2X::forwardThenBackward(post_forward) does not support twod_mode yet");
647 forwardThenBackward_doit_0(inmf, outmf, post_forward);
648}
649
650template <typename T>
651template <typename F>
652void R2X<T>::forwardThenBackward_doit_0 (MF const& inmf, MF& outmf,
653 F const& post_forward,
654 IntVect const& ngout,
655 Periodicity const& period)
656{
657 BL_PROFILE("FFT::R2X::forwardbackward_0");
658
659 if (m_r2x_sub) {
660 bool inmf_safe = m_sub_helper.ghost_safe(inmf.nGrowVect());
661 MF inmf_sub, inmf_tmp;
662 if (inmf_safe) {
663 inmf_sub = m_sub_helper.make_alias_mf(inmf);
664 } else {
665 inmf_tmp.define(inmf.boxArray(), inmf.DistributionMap(), 1, 0);
666 inmf_tmp.LocalCopy(inmf, 0, 0, 1, IntVect(0));
667 inmf_sub = m_sub_helper.make_alias_mf(inmf_tmp);
668 }
669
670 bool outmf_safe = m_sub_helper.ghost_safe(outmf.nGrowVect());
671 MF outmf_sub, outmf_tmp;
672 if (outmf_safe) {
673 outmf_sub = m_sub_helper.make_alias_mf(outmf);
674 } else {
675 IntVect const& ngtmp = m_sub_helper.make_safe_ghost(outmf.nGrowVect());
676 outmf_tmp.define(outmf.boxArray(), outmf.DistributionMap(), 1, ngtmp);
677 outmf_sub = m_sub_helper.make_alias_mf(outmf_tmp);
678 }
679
680 IntVect const& subngout = m_sub_helper.make_iv(ngout);
681 Periodicity const& subperiod = m_sub_helper.make_periodicity(period);
682 GpuArray<int,3> const& order = m_sub_helper.xyz_order();
683 m_r2x_sub->forwardThenBackward_doit_1
684 (inmf_sub, outmf_sub,
685 [=] AMREX_GPU_DEVICE (int i, int j, int k, auto& sp)
686 {
687 GpuArray<int,3> idx{i,j,k};
688 post_forward(idx[order[0]], idx[order[1]], idx[order[2]], sp);
689 },
690 subngout, subperiod);
691
692 if (!outmf_safe) {
693 outmf.LocalCopy(outmf_tmp, 0, 0, 1, outmf_tmp.nGrowVect());
694 }
695 }
696 else
697 {
698 this->forwardThenBackward_doit_1(inmf, outmf, post_forward, ngout, period);
699 }
700}
701
702template <typename T>
703template <typename F>
704void R2X<T>::forwardThenBackward_doit_1 (MF const& inmf, MF& outmf,
705 F const& post_forward,
706 IntVect const& ngout,
707 Periodicity const& period)
708{
709 BL_PROFILE("FFT::R2X::forwardbackward_1");
710
711 if (m_r2x_sub) {
712 amrex::Abort("R2X::forwardThenBackward_doit_1: How did this happen?");
713 }
714 else
715 {
716 this->forward(inmf);
717
718 // post-forward
719
720 int actual_dim = AMREX_SPACEDIM;
721#if (AMREX_SPACEDIM >= 2)
722 if (m_dom_0.length(1) == 1) { actual_dim = 1; }
723#endif
724#if (AMREX_SPACEDIM == 3)
725 if ((m_dom_0.length(2) == 1) && (m_dom_0.length(1) > 1)) { actual_dim = 2; }
726#endif
727
728 if (actual_dim == 1) {
729 if (m_cx.empty()) {
730 post_forward_doit<0>(detail::get_fab(m_rx), post_forward);
731 } else {
732 post_forward_doit<0>(detail::get_fab(m_cx), post_forward);
733 }
734 }
735#if (AMREX_SPACEDIM >= 2)
736 else if (actual_dim == 2) {
737 if (m_cy.empty()) {
738 post_forward_doit<1>(detail::get_fab(m_ry), post_forward);
739 } else {
740 post_forward_doit<1>(detail::get_fab(m_cy), post_forward);
741 }
742 }
743#endif
744#if (AMREX_SPACEDIM == 3)
745 else if (actual_dim == 3) {
746 if (m_cz.empty()) {
747 post_forward_doit<2>(detail::get_fab(m_rz), post_forward);
748 } else {
749 post_forward_doit<2>(detail::get_fab(m_cz), post_forward);
750 }
751 }
752#endif
753
754 this->backward();
755
756 outmf.ParallelCopy(m_rx, 0, 0, 1, IntVect(0),
757 amrex::elemwiseMin(ngout,outmf.nGrowVect()), period);
758 }
759}
760
761template <typename T>
762void R2X<T>::forward (MF const& inmf)
763{
764 BL_PROFILE("FFT::R2X::forward");
765
766 if (m_r2x_sub) {
767 if (m_sub_helper.ghost_safe(inmf.nGrowVect())) {
768 m_r2x_sub->forward(m_sub_helper.make_alias_mf(inmf));
769 } else {
770 MF tmp(inmf.boxArray(), inmf.DistributionMap(), 1, 0);
771 tmp.LocalCopy(inmf, 0, 0, 1, IntVect(0));
772 m_r2x_sub->forward(m_sub_helper.make_alias_mf(tmp));
773 }
774 return;
775 }
776
777 m_rx.ParallelCopy(inmf, 0, 0, 1);
778 if (m_bc[0].first == Boundary::periodic) {
779 m_fft_fwd_x.template compute_r2c<Direction::forward>();
780 } else {
781 m_fft_fwd_x.template compute_r2r<Direction::forward>();
782 }
783
784#if (AMREX_SPACEDIM >= 2)
785 if ( m_cmd_cx2cy) {
786 ParallelCopy(m_cy, m_cx, *m_cmd_cx2cy, 0, 0, 1, m_dtos_x2y);
787 } else if ( m_cmd_rx2ry) {
788 ParallelCopy(m_ry, m_rx, *m_cmd_rx2ry, 0, 0, 1, m_dtos_x2y);
789 }
790 if (m_bc[1].first != Boundary::periodic)
791 {
792 m_fft_fwd_y.template compute_r2r<Direction::forward>();
793 }
794 else if (m_bc[0].first == Boundary::periodic)
795 {
796 m_fft_fwd_y.template compute_c2c<Direction::forward>();
797 }
798 else
799 {
800 m_fft_fwd_y.template compute_r2c<Direction::forward>();
801 }
802#endif
803
804#if (AMREX_SPACEDIM == 3)
805 if ( m_cmd_cy2cz) {
806 ParallelCopy(m_cz, m_cy, *m_cmd_cy2cz, 0, 0, 1, m_dtos_y2z);
807 } else if ( m_cmd_ry2rz) {
808 ParallelCopy(m_rz, m_ry, *m_cmd_ry2rz, 0, 0, 1, m_dtos_y2z);
809 }
810 if (m_bc[2].first != Boundary::periodic)
811 {
812 m_fft_fwd_z.template compute_r2r<Direction::forward>();
813 }
814 else if (m_bc[0].first == Boundary::periodic ||
815 m_bc[1].first == Boundary::periodic)
816 {
817 m_fft_fwd_z.template compute_c2c<Direction::forward>();
818 }
819 else
820 {
821 m_fft_fwd_z.template compute_r2c<Direction::forward>();
822 }
823#endif
824}
825
826template <typename T>
827void R2X<T>::forward (MF const& inmf, MF& outmf)
828{
829 if (m_r2x_sub)
830 {
831 bool inmf_safe = m_sub_helper.ghost_safe(inmf.nGrowVect());
832 MF inmf_sub, inmf_tmp;
833 if (inmf_safe) {
834 inmf_sub = m_sub_helper.make_alias_mf(inmf);
835 } else {
836 inmf_tmp.define(inmf.boxArray(), inmf.DistributionMap(), 1, 0);
837 inmf_tmp.LocalCopy(inmf, 0, 0, 1, IntVect(0));
838 inmf_sub = m_sub_helper.make_alias_mf(inmf_tmp);
839 }
840
841 bool outmf_safe = m_sub_helper.ghost_safe(outmf.nGrowVect());
842 MF outmf_sub, outmf_tmp;
843 if (outmf_safe) {
844 outmf_sub = m_sub_helper.make_alias_mf(outmf);
845 } else {
846 outmf_tmp.define(outmf.boxArray(), outmf.DistributionMap(), 1, 0);
847 outmf_sub = m_sub_helper.make_alias_mf(outmf_tmp);
848 }
849
850 m_r2x_sub->forward(inmf_sub, outmf_sub);
851
852 if (!outmf_safe) {
853 outmf.LocalCopy(outmf_tmp, 0, 0, 1, IntVect(0));
854 }
855 }
856 else
857 {
858 this->forward(inmf);
859
860#if (AMREX_SPACEDIM == 3)
861 if (m_info.twod_mode) {
862 if (m_cy.empty() && !m_ry.empty()) {
863 ParallelCopy(outmf, m_dom_rx, m_ry, 0, 0, 1, IntVect(0), Swap01{});
864 } else if (m_ry.empty() && m_cy.empty() && m_cx.empty()) {
865 outmf.ParallelCopy(m_rx, 0, 0, 1);
866 } else {
867 amrex::Abort("R2X::forward(MF,MF): How did this happen?");
868 }
869 } else
870#endif
871 {
873 amrex::Abort("R2X::forward(MF,MF): TODO");
874 }
875 }
876}
877
878template <typename T>
879void R2X<T>::forward (MF const& inmf, cMF& outmf)
880{
881 if (m_r2x_sub)
882 {
883 bool inmf_safe = m_sub_helper.ghost_safe(inmf.nGrowVect());
884 MF inmf_sub, inmf_tmp;
885 if (inmf_safe) {
886 inmf_sub = m_sub_helper.make_alias_mf(inmf);
887 } else {
888 inmf_tmp.define(inmf.boxArray(), inmf.DistributionMap(), 1, 0);
889 inmf_tmp.LocalCopy(inmf, 0, 0, 1, IntVect(0));
890 inmf_sub = m_sub_helper.make_alias_mf(inmf_tmp);
891 }
892
893 bool outmf_safe = m_sub_helper.ghost_safe(outmf.nGrowVect());
894 cMF outmf_sub, outmf_tmp;
895 if (outmf_safe) {
896 outmf_sub = m_sub_helper.make_alias_mf(outmf);
897 } else {
898 outmf_tmp.define(outmf.boxArray(), outmf.DistributionMap(), 1, 0);
899 outmf_sub = m_sub_helper.make_alias_mf(outmf_tmp);
900 }
901
902 m_r2x_sub->forward(inmf_sub, outmf_sub);
903
904 if (!outmf_safe) {
905 outmf.LocalCopy(outmf_tmp, 0, 0, 1, IntVect(0));
906 }
907 }
908 else
909 {
910 this->forward(inmf);
911
912#if (AMREX_SPACEDIM == 3)
913 if (m_info.twod_mode) {
914 if (!m_cy.empty()) {
915 auto lo = m_dom_cy.smallEnd();
916 auto hi = m_dom_cy.bigEnd();
917 std::swap(lo[0],lo[1]);
918 std::swap(hi[0],hi[1]);
919 Box dom(lo,hi);
920 ParallelCopy(outmf, dom, m_cy, 0, 0, 1, IntVect(0), Swap01{});
921 } else if (m_ry.empty() && m_cy.empty() && !m_cx.empty()) {
922 outmf.ParallelCopy(m_cx, 0, 0, 1);
923 } else {
924 amrex::Abort("R2X::forward(MF,cMF): How did this happen?");
925 }
926 } else
927#endif
928 {
930 amrex::Abort("R2X::forward(MF,cMF): TODO");
931 }
932 }
933}
934
935template <typename T>
936void R2X<T>::backward ()
937{
938 BL_PROFILE("FFT::R2X::backward");
939
940 AMREX_ALWAYS_ASSERT(m_r2x_sub == nullptr);
941
942#if (AMREX_SPACEDIM == 3)
943 if (m_bc[2].first != Boundary::periodic)
944 {
945 m_fft_bwd_z.template compute_r2r<Direction::backward>();
946 }
947 else if (m_bc[0].first == Boundary::periodic ||
948 m_bc[1].first == Boundary::periodic)
949 {
950 m_fft_bwd_z.template compute_c2c<Direction::backward>();
951 }
952 else
953 {
954 m_fft_bwd_z.template compute_r2c<Direction::backward>();
955 }
956 if ( m_cmd_cz2cy) {
957 ParallelCopy(m_cy, m_cz, *m_cmd_cz2cy, 0, 0, 1, m_dtos_z2y);
958 } else if ( m_cmd_rz2ry) {
959 ParallelCopy(m_ry, m_rz, *m_cmd_rz2ry, 0, 0, 1, m_dtos_z2y);
960 }
961#endif
962
963#if (AMREX_SPACEDIM >= 2)
964 if (m_bc[1].first != Boundary::periodic)
965 {
966 m_fft_bwd_y.template compute_r2r<Direction::backward>();
967 }
968 else if (m_bc[0].first == Boundary::periodic)
969 {
970 m_fft_bwd_y.template compute_c2c<Direction::backward>();
971 }
972 else
973 {
974 m_fft_bwd_y.template compute_r2c<Direction::backward>();
975 }
976 if ( m_cmd_cy2cx) {
977 ParallelCopy(m_cx, m_cy, *m_cmd_cy2cx, 0, 0, 1, m_dtos_y2x);
978 } else if ( m_cmd_ry2rx) {
979 ParallelCopy(m_rx, m_ry, *m_cmd_ry2rx, 0, 0, 1, m_dtos_y2x);
980 }
981#endif
982
983 if (m_bc[0].first == Boundary::periodic) {
984 m_fft_bwd_x.template compute_r2c<Direction::backward>();
985 } else {
986 m_fft_bwd_x.template compute_r2r<Direction::backward>();
987 }
988}
989
990template <typename T>
991void R2X<T>::backward (MF const& inmf, MF& outmf, IntVect const& ngout,
992 Periodicity const& period)
993{
994 if (m_r2x_sub)
995 {
996 bool inmf_safe = m_sub_helper.ghost_safe(inmf.nGrowVect());
997 MF inmf_sub, inmf_tmp;
998 if (inmf_safe) {
999 inmf_sub = m_sub_helper.make_alias_mf(inmf);
1000 } else {
1001 inmf_tmp.define(inmf.boxArray(), inmf.DistributionMap(), 1, 0);
1002 inmf_tmp.LocalCopy(inmf, 0, 0, 1, IntVect(0));
1003 inmf_sub = m_sub_helper.make_alias_mf(inmf_tmp);
1004 }
1005
1006 bool outmf_safe = m_sub_helper.ghost_safe(outmf.nGrowVect());
1007 MF outmf_sub, outmf_tmp;
1008 if (outmf_safe) {
1009 outmf_sub = m_sub_helper.make_alias_mf(outmf);
1010 } else {
1011 IntVect const& ngtmp = m_sub_helper.make_safe_ghost(outmf.nGrowVect());
1012 outmf_tmp.define(outmf.boxArray(), outmf.DistributionMap(), 1, ngtmp);
1013 outmf_sub = m_sub_helper.make_alias_mf(outmf_tmp);
1014 }
1015
1016 IntVect const& subngout = m_sub_helper.make_iv(ngout);
1017 Periodicity const& subperiod = m_sub_helper.make_periodicity(period);
1018 m_r2x_sub->backward(inmf_sub, outmf_sub, subngout, subperiod);
1019
1020 if (!outmf_safe) {
1021 outmf.LocalCopy(outmf_tmp, 0, 0, 1, outmf_tmp.nGrowVect());
1022 }
1023 }
1024 else
1025 {
1026#if (AMREX_SPACEDIM == 3)
1027 if (m_info.twod_mode) {
1028 if (m_cy.empty() && !m_ry.empty()) {
1029 ParallelCopy(m_ry, m_dom_ry, inmf, 0, 0, 1, IntVect(0), Swap01{});
1030 } else if (m_ry.empty() && m_cy.empty() && m_cx.empty()) {
1031 m_rx.ParallelCopy(inmf, 0, 0, 1);
1032 } else {
1033 amrex::Abort("R2X::backward(MF,MF): How did this happen?");
1034 }
1035 } else
1036#endif
1037 {
1038 amrex::ignore_unused(inmf,outmf,ngout,period);
1039 amrex::Abort("R2X::backward(MF,MF): TODO");
1040 }
1041
1042 this->backward();
1043
1044 outmf.ParallelCopy(m_rx, 0, 0, 1, IntVect(0),
1045 amrex::elemwiseMin(ngout,outmf.nGrowVect()), period);
1046 }
1047}
1048
1049template <typename T>
1050void R2X<T>::backward (cMF const& inmf, MF& outmf, IntVect const& ngout,
1051 Periodicity const& period)
1052{
1053 if (m_r2x_sub)
1054 {
1055 bool inmf_safe = m_sub_helper.ghost_safe(inmf.nGrowVect());
1056 cMF inmf_sub, inmf_tmp;
1057 if (inmf_safe) {
1058 inmf_sub = m_sub_helper.make_alias_mf(inmf);
1059 } else {
1060 inmf_tmp.define(inmf.boxArray(), inmf.DistributionMap(), 1, 0);
1061 inmf_tmp.LocalCopy(inmf, 0, 0, 1, IntVect(0));
1062 inmf_sub = m_sub_helper.make_alias_mf(inmf_tmp);
1063 }
1064
1065 bool outmf_safe = m_sub_helper.ghost_safe(outmf.nGrowVect());
1066 MF outmf_sub, outmf_tmp;
1067 if (outmf_safe) {
1068 outmf_sub = m_sub_helper.make_alias_mf(outmf);
1069 } else {
1070 IntVect const& ngtmp = m_sub_helper.make_safe_ghost(outmf.nGrowVect());
1071 outmf_tmp.define(outmf.boxArray(), outmf.DistributionMap(), 1, ngtmp);
1072 outmf_sub = m_sub_helper.make_alias_mf(outmf_tmp);
1073 }
1074
1075 IntVect const& subngout = m_sub_helper.make_iv(ngout);
1076 Periodicity const& subperiod = m_sub_helper.make_periodicity(period);
1077 m_r2x_sub->backward(inmf_sub, outmf_sub, subngout, subperiod);
1078
1079 if (!outmf_safe) {
1080 outmf.LocalCopy(outmf_tmp, 0, 0, 1, outmf_tmp.nGrowVect());
1081 }
1082 }
1083 else
1084 {
1085#if (AMREX_SPACEDIM == 3)
1086 if (m_info.twod_mode) {
1087 if (!m_cy.empty()) {
1088 ParallelCopy(m_cy, m_dom_cy, inmf, 0, 0, 1, IntVect(0), Swap01{});
1089 } else if (m_ry.empty() && m_cy.empty() && !m_cx.empty()) {
1090 m_cx.ParallelCopy(inmf, 0, 0, 1);
1091 } else {
1092 amrex::Abort("R2X::backward(cMF,MF): How did this happen?");
1093 }
1094 } else
1095#endif
1096 {
1097 amrex::ignore_unused(inmf,outmf,ngout,period);
1098 amrex::Abort("R2X::backward(cMF,MF): TODO");
1099 }
1100
1101 this->backward();
1102
1103 outmf.ParallelCopy(m_rx, 0, 0, 1, IntVect(0),
1104 amrex::elemwiseMin(ngout,outmf.nGrowVect()), period);
1105 }
1106}
1107
1108template <typename T>
1109template <int dim, typename FAB, typename F>
1110void R2X<T>::post_forward_doit (FAB* fab, F const& f)
1111{
1112 if (m_info.twod_mode) {
1113 amrex::Abort("xxxxx post_forward_doit: todo");
1114 }
1115 if (fab) {
1116 auto const& a = fab->array();
1117 ParallelForOMP(fab->box(),
1118 [f=f,a=a] AMREX_GPU_DEVICE (int i, int j, int k)
1119 {
1120 if constexpr (dim == 0) {
1121 f(i,j,k,a(i,j,k));
1122 } else if constexpr (dim == 1) {
1123 f(j,i,k,a(i,j,k));
1124 } else {
1125 f(j,k,i,a(i,j,k));
1126 }
1127 });
1128 }
1129}
1130
1131}
1132
1133#endif
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define AMREX_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition AMReX_BLassert.H:49
#define AMREX_ALWAYS_ASSERT(EX)
Definition AMReX_BLassert.H:50
#define AMREX_GPU_DEVICE
Definition AMReX_GpuQualifiers.H:18
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:568
A class for managing a List of Boxes that share a common IndexType. This class implements operations ...
Definition AMReX_BoxList.H:52
__host__ __device__ const IntVectND< dim > & bigEnd() const &noexcept
Return the inclusive upper bound of the box.
Definition AMReX_Box.H:123
__host__ __device__ Long numPts() const noexcept
Return the number of points contained in the BoxND.
Definition AMReX_Box.H:356
__host__ __device__ IntVectND< dim > length() const noexcept
Return the length of the BoxND.
Definition AMReX_Box.H:154
__host__ __device__ IntVectND< dim > size() const noexcept
Return the length of the BoxND.
Definition AMReX_Box.H:147
__host__ __device__ const IntVectND< dim > & smallEnd() const &noexcept
Return the inclusive lower bound of the box.
Definition AMReX_Box.H:111
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
3D Poisson solver for periodic, Dirichlet & Neumann boundaries in the first two dimensions,...
Definition AMReX_FFT_Poisson.H:187
Poisson solver for periodic, Dirichlet & Neumann boundaries using FFT.
Definition AMReX_FFT_Poisson.H:67
Discrete Fourier Transform.
Definition AMReX_FFT_R2X.H:32
void post_forward_doit(FAB *fab, F const &f)
Apply a user functor to a FAB after the forward transform along dimension dim.
Definition AMReX_FFT_R2X.H:1110
std::conditional_t< std::is_same_v< T, Real >, MultiFab, FabArray< BaseFab< T > > > MF
Definition AMReX_FFT_R2X.H:35
~R2X()
Destroy any FFT resources held by this object.
Definition AMReX_FFT_R2X.H:607
void forwardThenBackward(MF const &inmf, MF &outmf, F const &post_forward)
Execute forward transform, apply post_forward, then backward transform.
Definition AMReX_FFT_R2X.H:642
R2X(Box const &domain, Array< std::pair< Boundary, Boundary >, 3 > const &bc, Info const &info=Info{})
Build an FFT plan for the given domain and boundary types.
Definition AMReX_FFT_R2X.H:213
FabArray< BaseFab< GpuComplex< T > > > cMF
Definition AMReX_FFT_R2X.H:36
R2X(R2X &&)=delete
T scalingFactor() const
Normalization applied after a forward/backward sequence.
Definition AMReX_FFT_R2X.H:624
R2X(R2X const &)=delete
R2X & operator=(R2X const &)=delete
void forwardThenBackward_doit_1(MF const &inmf, MF &outmf, F const &post_forward, IntVect const &ngout=IntVect(0), Periodicity const &period=Periodicity::NonPeriodic())
CUDA-visible helper that operates on the complex FAB path.
Definition AMReX_FFT_R2X.H:704
void forwardThenBackward_doit_0(MF const &inmf, MF &outmf, F const &post_forward, IntVect const &ngout=IntVect(0), Periodicity const &period=Periodicity::NonPeriodic())
CUDA-visible helper that performs the forward/modify/backward cycle.
Definition AMReX_FFT_R2X.H:652
int size() const noexcept
Return the number of FABs in the FabArray.
Definition AMReX_FabArrayBase.H:110
const DistributionMapping & DistributionMap() const noexcept
Return constant reference to associated DistributionMapping.
Definition AMReX_FabArrayBase.H:131
bool empty() const noexcept
Definition AMReX_FabArrayBase.H:89
Box box(int K) const noexcept
Return the Kth Box in the BoxArray. That is, the valid region of the Kth grid.
Definition AMReX_FabArrayBase.H:101
An Array of FortranArrayBox(FAB)-like Objects.
Definition AMReX_FabArray.H:349
void define(const BoxArray &bxs, const DistributionMapping &dm, int nvar, int ngrow, const MFInfo &info=MFInfo(), const FabFactory< FAB > &factory=DefaultFabFactory< FAB >())
Define this FabArray identically to that performed by the constructor having an analogous function si...
Definition AMReX_FabArray.H:2173
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:40
This provides length of period for periodic domains. 0 means it is not periodic in that direction....
Definition AMReX_Periodicity.H:17
static const Periodicity & NonPeriodic() noexcept
Definition AMReX_Periodicity.cpp:52
amrex_long Long
Definition AMReX_INT.H:30
void ParallelForOMP(T n, L const &f) noexcept
Performance-portable kernel launch function with optional OpenMP threading.
Definition AMReX_GpuLaunch.H:319
std::array< T, N > Array
Definition AMReX_Array.H:26
Definition AMReX_FFT_Helper.H:52
int MyProcSub() noexcept
my sub-rank in current frame
Definition AMReX_ParallelContext.H:76
int NProcsSub() noexcept
number of ranks in current frame
Definition AMReX_ParallelContext.H:74
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:139
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30
double second() noexcept
Definition AMReX_Utility.cpp:940
BoxArray decompose(Box const &domain, int nboxes, Array< bool, 3 > const &decomp, bool no_overlap)
Decompose domain box into BoxArray.
Definition AMReX_BoxArray.cpp:1947
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:240
void ParallelCopy(MF &dst, MF const &src, int scomp, int dcomp, int ncomp, IntVect const &ng_src=IntVect(0), IntVect const &ng_dst=IntVect(0), Periodicity const &period=Periodicity::NonPeriodic())
dst = src w/ MPI communication
Definition AMReX_FabArrayUtility.H:2019
__host__ __device__ constexpr T elemwiseMin(T const &a, T const &b) noexcept
Return the element-wise minimum of the given values for types like XDim3.
Definition AMReX_Algorithm.H:62
Definition AMReX_FFT_Helper.H:64
bool twod_mode
Definition AMReX_FFT_Helper.H:75
bool oned_mode
We might have a special twod_mode: nx or ny == 1 && nz > 1.
Definition AMReX_FFT_Helper.H:78
int batch_size
Batched FFT size. Only support in R2C, not R2X.
Definition AMReX_FFT_Helper.H:81
int nprocs
Max number of processes to use.
Definition AMReX_FFT_Helper.H:84
Definition AMReX_FFT_Helper.H:188
std::conditional_t< std::is_same_v< float, T >, cuComplex, cuDoubleComplex > VendorComplex
Definition AMReX_FFT_Helper.H:192
Fixed-size array that can be used on GPU.
Definition AMReX_Array.H:43
FabArray memory allocation information.
Definition AMReX_FabArray.H:66