1#ifndef AMREX_FFT_R2C_H_
2#define AMREX_FFT_R2C_H_
3#include <AMReX_Config.H>
14template <
typename T>
class OpenBCSolver;
15template <
typename T>
class Poisson;
16template <
typename T>
class PoissonHybrid;
39template <
typename T = Real, FFT::Direction D = FFT::Direction::both,
bool C = false>
44 using MF = std::conditional_t
45 <
C,
cMF, std::conditional_t<std::is_same_v<T,Real>,
49 template <
typename U>
friend class Poisson;
69 explicit R2C (std::array<int,AMREX_SPACEDIM>
const& domain_size,
103 std::array<int,AMREX_SPACEDIM>
const& local_size);
116 std::pair<std::array<int,AMREX_SPACEDIM>,std::array<int,AMREX_SPACEDIM>>
148 std::array<int,AMREX_SPACEDIM>
const& local_size);
167 std::pair<std::array<int,AMREX_SPACEDIM>,std::array<int,AMREX_SPACEDIM>>
192 std::enable_if_t<DIR == Direction::both, int> = 0>
194 int incomp = 0,
int outcomp = 0)
229 void forward (
MF const& inmf,
cMF& outmf,
int incomp = 0,
int outcomp = 0);
244 template <
typename RT,
typename CT,
Direction DIR=D,
bool CP=
C,
247 && ((
sizeof(RT)*2 ==
sizeof(CT) && !CP) ||
248 (
sizeof(RT) ==
sizeof(CT) && CP)),
int> = 0>
260 template <Direction DIR=D, std::enable_if_t<DIR == Direction::both,
int> = 0>
276 void backward (
cMF const& inmf,
MF& outmf,
int incomp = 0,
int outcomp = 0);
291 template <
typename CT,
typename RT,
Direction DIR=D,
bool CP=
C,
294 && ((
sizeof(RT)*2 ==
sizeof(CT) && !CP) ||
295 (
sizeof(RT) ==
sizeof(CT) && CP)),
int> = 0>
326 template <
typename F>
329 template <
typename F>
334 void prepare_openbc ();
340 void backward_doit (
cMF const& inmf,
MF& outmf,
343 int incomp = 0,
int outcomp = 0);
345 std::pair<Plan<T>,
Plan<T>> make_c2c_plans (
cMF& inout,
int ndims)
const;
347 static Box make_domain_x (
Box const& domain)
362 static Box make_domain_y (
Box const& domain)
367 domain.length(2)-1)),
372 domain.length(2)-1)),
377 static Box make_domain_z (
Box const& domain)
382 domain.length(1)-1)),
387 domain.length(1)-1)),
392 static std::pair<BoxArray,DistributionMapping>
393 make_layout_from_local_domain (std::array<int,AMREX_SPACEDIM>
const& local_start,
394 std::array<int,AMREX_SPACEDIM>
const& local_size);
396 template <
typename FA,
typename RT>
397 std::pair<std::unique_ptr<char,DataDeleter>,std::size_t>
398 install_raw_ptr (FA& fa, RT
const* p);
400 Plan<T> m_fft_fwd_x{};
401 Plan<T> m_fft_bwd_x{};
402 Plan<T> m_fft_fwd_y{};
403 Plan<T> m_fft_bwd_y{};
404 Plan<T> m_fft_fwd_z{};
405 Plan<T> m_fft_bwd_z{};
406 Plan<T> m_fft_fwd_x_half{};
407 Plan<T> m_fft_bwd_x_half{};
412 std::unique_ptr<MultiBlockCommMetaData> m_cmd_x2y;
413 std::unique_ptr<MultiBlockCommMetaData> m_cmd_y2x;
414 std::unique_ptr<MultiBlockCommMetaData> m_cmd_y2z;
415 std::unique_ptr<MultiBlockCommMetaData> m_cmd_z2y;
416 std::unique_ptr<MultiBlockCommMetaData> m_cmd_x2z;
417 std::unique_ptr<MultiBlockCommMetaData> m_cmd_z2x;
418 std::unique_ptr<MultiBlockCommMetaData> m_cmd_x2z_half;
419 std::unique_ptr<MultiBlockCommMetaData> m_cmd_z2x_half;
424 RotateFwd m_dtos_x2z{};
425 RotateBwd m_dtos_z2x{};
433 mutable cMF m_raw_cmf;
435 std::unique_ptr<char,DataDeleter> m_data_1;
436 std::unique_ptr<char,DataDeleter> m_data_2;
439 Box m_spectral_domain_x;
440 Box m_spectral_domain_y;
441 Box m_spectral_domain_z;
443 std::unique_ptr<R2C<T,D,C>> m_r2c_sub;
444 detail::SubHelper m_sub_helper;
448 bool m_do_alld_fft =
false;
449 bool m_slab_decomp =
false;
450 bool m_openbc_half =
false;
453template <
typename T, Direction D,
bool C>
455 : m_real_domain(domain),
456 m_spectral_domain_x(make_domain_x(domain)),
457#if (AMREX_SPACEDIM >= 2)
458 m_spectral_domain_y(make_domain_y(domain)),
459#if (AMREX_SPACEDIM == 3)
460 m_spectral_domain_z(make_domain_z(domain)),
463 m_sub_helper(domain),
468 static_assert(std::is_same_v<float,T> || std::is_same_v<double,T>);
471#if (AMREX_SPACEDIM == 2)
476 int(domain.
length(1) > 1) +
477 int(domain.
length(2) > 1)) >= 2);
482 Box subbox = m_sub_helper.make_box(m_real_domain);
483 if (subbox.
size() != m_real_domain.
size()) {
484 m_r2c_sub = std::make_unique<R2C<T,D,C>>(subbox, m_info);
492#if (AMREX_SPACEDIM == 3)
497 int shortside = m_real_domain.
shortside();
508 m_slab_decomp =
true;
510 m_slab_decomp =
true;
522 m_rx.define(bax, dmx, ncomp, 0,
MFInfo().SetAlloc(
false));
526 for (
auto & b : bl) {
528 b.setBig(0, m_spectral_domain_x.
bigEnd(0));
531 m_cx.
define(cbax, dmx, ncomp, 0,
MFInfo().SetAlloc(
false));
542#if (AMREX_SPACEDIM >= 2)
544 if ((m_real_domain.
length(1) > 1) && !m_slab_decomp && !m_info.
oned_mode)
548 if (cbay.size() == dmx.size()) {
551 cdmy = detail::make_iota_distromap(cbay.size());
553 m_cy.
define(cbay, cdmy, ncomp, 0,
MFInfo().SetAlloc(
false));
557#if (AMREX_SPACEDIM == 3)
558 if (m_real_domain.
length(1) > 1 &&
562 {
false,
true,
true},
true);
564 if (cbaz.size() == dmx.size()) {
566 }
else if (cbaz.size() == cdmy.
size()) {
569 cdmz = detail::make_iota_distromap(cbaz.size());
571 m_cz.
define(cbaz, cdmz, ncomp, 0,
MFInfo().SetAlloc(
false));
577 m_data_1 = detail::make_mfs_share(m_rx, m_cx);
578 m_data_2 = detail::make_mfs_share(m_cz, m_cz);
580 m_data_1 = detail::make_mfs_share(m_rx, m_cz);
581 m_data_2 = detail::make_mfs_share(m_cy, m_cy);
583 if (myproc < m_cx.
size()) {
586 m_cx.
setFab(myproc, FAB(box, ncomp, m_rx[myproc].dataPtr()));
591 m_data_1 = detail::make_mfs_share(m_rx, m_cz);
592 m_data_2 = detail::make_mfs_share(m_cx, m_cx);
594 m_data_1 = detail::make_mfs_share(m_rx, m_cy);
595 m_data_2 = detail::make_mfs_share(m_cx, m_cz);
603#if (AMREX_SPACEDIM >= 2)
604 if (! m_cy.
empty()) {
606 m_cmd_x2y = std::make_unique<MultiBlockCommMetaData>
607 (m_cy, m_spectral_domain_y, m_cx,
IntVect(0), m_dtos_x2y);
608 m_cmd_y2x = std::make_unique<MultiBlockCommMetaData>
609 (m_cx, m_spectral_domain_x, m_cy,
IntVect(0), m_dtos_y2x);
612#if (AMREX_SPACEDIM == 3)
613 if (! m_cz.
empty() ) {
616 m_cmd_x2z = std::make_unique<MultiBlockCommMetaData>
617 (m_cz, m_spectral_domain_z, m_cx,
IntVect(0), m_dtos_x2z);
618 m_cmd_z2x = std::make_unique<MultiBlockCommMetaData>
619 (m_cx, m_spectral_domain_x, m_cz,
IntVect(0), m_dtos_z2x);
622 m_cmd_y2z = std::make_unique<MultiBlockCommMetaData>
623 (m_cz, m_spectral_domain_z, m_cy,
IntVect(0), m_dtos_y2z);
624 m_cmd_z2y = std::make_unique<MultiBlockCommMetaData>
625 (m_cy, m_spectral_domain_y, m_cz,
IntVect(0), m_dtos_z2y);
634 if (myproc < m_rx.size())
637 int ndims = m_slab_decomp ? 2 : 1;
638 std::tie(m_fft_fwd_x, m_fft_bwd_x) = make_c2c_plans(m_cx, ndims);
640 Box const& box = m_rx.box(myproc);
641 auto* pr = m_rx[myproc].dataPtr();
644 m_fft_fwd_x.template init_r2c<Direction::forward>(box, pr, pc, m_slab_decomp, ncomp);
645 m_fft_bwd_x = m_fft_fwd_x;
648 m_fft_fwd_x.template init_r2c<Direction::forward>(box, pr, pc, m_slab_decomp, ncomp);
651 m_fft_bwd_x.template init_r2c<Direction::backward>(box, pr, pc, m_slab_decomp, ncomp);
657#if (AMREX_SPACEDIM >= 2)
658 if (! m_cy.
empty()) {
659 std::tie(m_fft_fwd_y, m_fft_bwd_y) = make_c2c_plans(m_cy,1);
662#if (AMREX_SPACEDIM == 3)
663 if (! m_cz.
empty()) {
664 std::tie(m_fft_fwd_z, m_fft_bwd_z) = make_c2c_plans(m_cz,1);
671 m_data_1 = detail::make_mfs_share(m_rx, m_cx);
672 std::tie(m_fft_fwd_x, m_fft_bwd_x) = make_c2c_plans(m_cx,AMREX_SPACEDIM);
674 m_data_1 = detail::make_mfs_share(m_rx, m_rx);
675 m_data_2 = detail::make_mfs_share(m_cx, m_cx);
677 auto const& len = m_real_domain.
length();
678 auto* pr = (
void*)m_rx[0].dataPtr();
679 auto* pc = (
void*)m_cx[0].dataPtr();
681 m_fft_fwd_x.template init_r2c<Direction::forward>(len, pr, pc,
false, ncomp);
682 m_fft_bwd_x = m_fft_fwd_x;
685 m_fft_fwd_x.template init_r2c<Direction::forward>(len, pr, pc,
false, ncomp);
688 m_fft_bwd_x.template init_r2c<Direction::backward>(len, pr, pc,
false, ncomp);
695template <
typename T, Direction D,
bool C>
700template <
typename T, Direction D,
bool C>
703 if (m_fft_bwd_x.plan != m_fft_fwd_x.plan) {
704 m_fft_bwd_x.destroy();
706 if (m_fft_bwd_y.plan != m_fft_fwd_y.plan) {
707 m_fft_bwd_y.destroy();
709 if (m_fft_bwd_z.plan != m_fft_fwd_z.plan) {
710 m_fft_bwd_z.destroy();
712 m_fft_fwd_x.destroy();
713 m_fft_fwd_y.destroy();
714 m_fft_fwd_z.destroy();
715 if (m_fft_bwd_x_half.plan != m_fft_fwd_x_half.plan) {
716 m_fft_bwd_x_half.destroy();
718 m_fft_fwd_x_half.destroy();
721template <
typename T, Direction D,
bool C>
722std::pair<BoxArray,DistributionMapping>
724 std::array<int,AMREX_SPACEDIM>
const& local_size)
728 Box bx(lo, lo+len-1);
735 pmap.reserve(allboxes.size());
736 for (
int i = 0; i < allboxes.size(); ++i) {
737 if (allboxes[i].ok()) {
741 allboxes.erase(std::remove_if(allboxes.begin(), allboxes.end(),
742 [=] (
Box const& b) { return b.isEmpty(); }),
744 BoxList bl(std::move(allboxes));
745 return std::make_pair(BoxArray(std::move(bl)), DistributionMapping(std::move(pmap)));
747 return std::make_pair(BoxArray(bx), DistributionMapping(Vector<int>({0})));
751template <
typename T, Direction D,
bool C>
753 std::array<int,AMREX_SPACEDIM>
const& local_size)
755 auto const& [ba, dm] = make_layout_from_local_domain(local_start, local_size);
756 m_raw_mf =
MF(ba, dm, m_rx.nComp(), 0,
MFInfo().SetAlloc(
false));
759template <
typename T, Direction D,
bool C>
760std::pair<std::array<int,AMREX_SPACEDIM>,std::array<int,AMREX_SPACEDIM>>
763 m_raw_mf =
MF(m_rx.boxArray(), m_rx.DistributionMap(), m_rx.nComp(), 0,
767 if (myproc < m_rx.size()) {
768 Box const& box = m_rx.box(myproc);
769 return std::make_pair(box.
smallEnd().toArray(),
772 return std::make_pair(std::array<int,AMREX_SPACEDIM>{
AMREX_D_DECL(0,0,0)},
777template <
typename T, Direction D,
bool C>
779 std::array<int,AMREX_SPACEDIM>
const& local_size)
781 auto const& [ba, dm] = make_layout_from_local_domain(local_start, local_size);
782 m_raw_cmf =
cMF(ba, dm, m_rx.nComp(), 0,
MFInfo().SetAlloc(
false));
785template <
typename T, Direction D,
bool C>
786std::pair<std::array<int,AMREX_SPACEDIM>,std::array<int,AMREX_SPACEDIM>>
789 auto const ncomp = m_info.batch_size;
790 auto const& [ba, dm] = getSpectralDataLayout();
795 if (myproc < m_raw_cmf.size()) {
796 Box const& box = m_raw_cmf.box(myproc);
797 return std::make_pair(box.
smallEnd().toArray(), box.
length().toArray());
799 return std::make_pair(std::array<int,AMREX_SPACEDIM>{
AMREX_D_DECL(0,0,0)},
804template <
typename T, Direction D,
bool C>
807 if (
C || m_r2c_sub) {
amrex::Abort(
"R2C: OpenBC not supported with reduced dimensions or complex inputs"); }
809#if (AMREX_SPACEDIM == 3)
810 if (m_do_alld_fft) {
return; }
812 auto const ncomp = m_info.batch_size;
814 if (m_slab_decomp && ! m_fft_fwd_x_half.defined) {
815 auto* fab = detail::get_fab(m_rx);
817 Box bottom_half = m_real_domain;
818 bottom_half.
growHi(2,-m_real_domain.length(2)/2);
819 Box box = fab->box() & bottom_half;
821 auto* pr = fab->dataPtr();
823 detail::get_fab(m_cx)->dataPtr();
825 m_fft_fwd_x_half.template init_r2c<Direction::forward>
826 (box, pr, pc, m_slab_decomp, ncomp);
827 m_fft_bwd_x_half = m_fft_fwd_x_half;
830 m_fft_fwd_x_half.template init_r2c<Direction::forward>
831 (box, pr, pc, m_slab_decomp, ncomp);
834 m_fft_bwd_x_half.template init_r2c<Direction::backward>
835 (box, pr, pc, m_slab_decomp, ncomp);
842 if (m_cmd_x2z && ! m_cmd_x2z_half) {
843 Box bottom_half = m_spectral_domain_z;
846 bottom_half.
growHi(0,-m_spectral_domain_z.length(0)/2);
847 m_cmd_x2z_half = std::make_unique<MultiBlockCommMetaData>
848 (m_cz, bottom_half, m_cx,
IntVect(0), m_dtos_x2z);
851 if (m_cmd_z2x && ! m_cmd_z2x_half) {
852 Box bottom_half = m_spectral_domain_x;
853 bottom_half.
growHi(2,-m_spectral_domain_x.length(2)/2);
854 m_cmd_z2x_half = std::make_unique<MultiBlockCommMetaData>
855 (m_cx, bottom_half, m_cz,
IntVect(0), m_dtos_z2x);
860template <
typename T, Direction D,
bool C>
867 auto const ncomp = m_info.batch_size;
870 if (m_sub_helper.ghost_safe(inmf.nGrowVect())) {
871 m_r2c_sub->forward(m_sub_helper.make_alias_mf(inmf), incomp);
873 MF tmp(inmf.boxArray(), inmf.DistributionMap(), ncomp, 0);
874 tmp.LocalCopy(inmf, incomp, 0, ncomp,
IntVect(0));
875 m_r2c_sub->forward(m_sub_helper.make_alias_mf(tmp),0);
880 if (&m_rx != &inmf) {
881 m_rx.ParallelCopy(inmf, incomp, 0, ncomp);
886 m_fft_fwd_x.template compute_c2c<Direction::forward>();
888 m_fft_fwd_x.template compute_r2c<Direction::forward>();
893 auto& fft_x = m_openbc_half ? m_fft_fwd_x_half : m_fft_fwd_x;
895 fft_x.template compute_c2c<Direction::forward>();
897 fft_x.template compute_r2c<Direction::forward>();
901 ParallelCopy(m_cy, m_cx, *m_cmd_x2y, 0, 0, ncomp, m_dtos_x2y);
903 m_fft_fwd_y.template compute_c2c<Direction::forward>();
906 ParallelCopy(m_cz, m_cy, *m_cmd_y2z, 0, 0, ncomp, m_dtos_y2z);
908#if (AMREX_SPACEDIM == 3)
909 else if ( m_cmd_x2z) {
914 {components, m_dtos_x2z};
915 auto handler = ParallelCopy_nowait(m_cz, m_cx, *m_cmd_x2z_half, packing);
917 Box upper_half = m_spectral_domain_z;
920 upper_half.
growLo (0,-m_spectral_domain_z.length(0)/2);
921 m_cz.setVal(0, upper_half, 0, ncomp);
923 ParallelCopy_finish(m_cz, std::move(handler), *m_cmd_x2z_half, packing);
925 ParallelCopy(m_cz, m_cx, *m_cmd_x2z, 0, 0, ncomp, m_dtos_x2z);
929 m_fft_fwd_z.template compute_c2c<Direction::forward>();
932template <
typename T, Direction D,
bool C>
933template <
typename FA,
typename RT>
934std::pair<std::unique_ptr<char,DataDeleter>,std::size_t>
939 using FAB =
typename FA::FABType::value_type;
940 using T_FAB =
typename FAB::value_type;
941 static_assert(
sizeof(T_FAB) ==
sizeof(RT));
943 auto const ncomp = m_info.batch_size;
944 auto const& ia = fa.IndexArray();
949 if ( ! ia.empty() ) {
951 Box const& box = fa.fabbox(K);
955 sz =
sizeof(T_FAB) * box.
numPts() * ncomp;
958 fa.setFab(K, FAB(box,ncomp,
pp));
962 return std::make_pair(std::unique_ptr<char,DataDeleter>{},std::size_t(0));
964 return std::make_pair(std::unique_ptr<char,DataDeleter>
970template <
typename T, Direction D,
bool C>
971template <
typename RT,
typename CT,
Direction DIR,
bool CP,
974 && ((
sizeof(RT)*2 ==
sizeof(CT) && !CP) ||
975 (
sizeof(RT) ==
sizeof(CT) && CP)),
int> >
978 auto [rdata, rsz] = install_raw_ptr(m_raw_mf, in);
979 auto [cdata, csz] = install_raw_ptr(m_raw_cmf, out);
994template <
typename T, Direction D,
bool C>
995template <Direction DIR, std::enable_if_t<DIR == Direction::both,
int> >
1001template <
typename T, Direction D,
bool C>
1007 auto const ncomp = m_info.batch_size;
1010 if (m_sub_helper.ghost_safe(outmf.nGrowVect())) {
1011 MF submf = m_sub_helper.make_alias_mf(outmf);
1012 IntVect const& subngout = m_sub_helper.make_iv(ngout);
1013 Periodicity const& subperiod = m_sub_helper.make_periodicity(period);
1014 m_r2c_sub->backward_doit(submf, subngout, subperiod, outcomp);
1016 MF tmp(outmf.boxArray(), outmf.DistributionMap(), ncomp,
1017 m_sub_helper.make_safe_ghost(outmf.nGrowVect()));
1018 this->backward_doit(tmp, ngout, period, 0);
1019 outmf.LocalCopy(tmp, 0, outcomp, ncomp, tmp.nGrowVect());
1024 if (m_do_alld_fft) {
1026 m_fft_bwd_x.template compute_c2c<Direction::backward>();
1028 m_fft_bwd_x.template compute_r2c<Direction::backward>();
1030 outmf.ParallelCopy(m_rx, 0, outcomp, ncomp,
IntVect(0),
1035 m_fft_bwd_z.template compute_c2c<Direction::backward>();
1037 ParallelCopy(m_cy, m_cz, *m_cmd_z2y, 0, 0, ncomp, m_dtos_z2y);
1039#if (AMREX_SPACEDIM == 3)
1040 else if ( m_cmd_z2x) {
1041 auto const& cmd = m_openbc_half ? m_cmd_z2x_half : m_cmd_z2x;
1042 ParallelCopy(m_cx, m_cz, *cmd, 0, 0, ncomp, m_dtos_z2x);
1046 m_fft_bwd_y.template compute_c2c<Direction::backward>();
1048 ParallelCopy(m_cx, m_cy, *m_cmd_y2x, 0, 0, ncomp, m_dtos_y2x);
1051 auto& fft_x = m_openbc_half ? m_fft_bwd_x_half : m_fft_bwd_x;
1053 fft_x.template compute_c2c<Direction::backward>();
1055 fft_x.template compute_r2c<Direction::backward>();
1057 outmf.ParallelCopy(m_rx, 0, outcomp, ncomp,
IntVect(0),
1061template <
typename T, Direction D,
bool C>
1062template <
typename CT,
typename RT,
Direction DIR,
bool CP,
1065 && ((
sizeof(RT)*2 ==
sizeof(CT) && !CP) ||
1066 (
sizeof(RT) ==
sizeof(CT) && CP)),
int> >
1069 auto [rdata, rsz] = install_raw_ptr(m_raw_mf, out);
1070 auto [cdata, csz] = install_raw_ptr(m_raw_cmf, in);
1085template <
typename T, Direction D,
bool C>
1092 auto* fab = detail::get_fab(inout);
1093 if (!fab) {
return {fwd, bwd};}
1095 Box const& box = fab->box();
1098 auto const ncomp = m_info.batch_size;
1100#ifdef AMREX_USE_SYCL
1101 fwd.template init_c2c<Direction::forward>(box, pio, ncomp, ndims);
1105 fwd.template init_c2c<Direction::forward>(box, pio, ncomp, ndims);
1108 bwd.template init_c2c<Direction::backward>(box, pio, ncomp, ndims);
1115template <
typename T, Direction D,
bool C>
1116template <
typename F>
1119 if (m_info.twod_mode || m_info.batch_size > 1) {
1121#if (AMREX_SPACEDIM > 1)
1122 }
else if (m_r2c_sub) {
1124#if (AMREX_SPACEDIM == 2)
1126 m_r2c_sub->post_forward_doit_1
1129 post_forward(0, i, 0, sp);
1132 if (m_real_domain.length(0) == 1 && m_real_domain.length(1) == 1) {
1134 m_r2c_sub->post_forward_doit_1
1137 post_forward(0, 0, i, sp);
1139 }
else if (m_real_domain.length(0) == 1 && m_real_domain.length(2) == 1) {
1141 m_r2c_sub->post_forward_doit_1
1144 post_forward(0, i, 0, sp);
1146 }
else if (m_real_domain.length(0) == 1) {
1148 m_r2c_sub->post_forward_doit_1
1151 post_forward(0, i, j, sp);
1153 }
else if (m_real_domain.length(1) == 1) {
1155 m_r2c_sub->post_forward_doit_1
1158 post_forward(i, 0, j, sp);
1161 amrex::Abort(
"R2c::post_forward_doit_0: how did this happen?");
1166 this->post_forward_doit_1(post_forward);
1170template <
typename T, Direction D,
bool C>
1171template <
typename F>
1174 if (m_info.twod_mode || m_info.batch_size > 1) {
1176 }
else if (m_r2c_sub) {
1177 amrex::Abort(
"R2C::post_forward_doit_1: How did this happen?");
1179 if ( ! m_cz.empty()) {
1180 auto* spectral_fab = detail::get_fab(m_cz);
1182 auto const& a = spectral_fab->array();
1186 post_forward(jx,ky,iz,a(iz,jx,ky));
1189 }
else if ( ! m_cy.empty()) {
1190 auto* spectral_fab = detail::get_fab(m_cy);
1192 auto const& a = spectral_fab->array();
1196 post_forward(jx,iy,k,a(iy,jx,k));
1200 auto* spectral_fab = detail::get_fab(m_cx);
1202 auto const& a = spectral_fab->array();
1206 post_forward(i,j,k,a(i,j,k));
1213template <
typename T, Direction D,
bool C>
1216#if (AMREX_SPACEDIM == 3)
1217 if (m_info.twod_mode) {
1218 return T(1)/T(
Long(m_real_domain.length(0)) *
1219 Long(m_real_domain.length(1)));
1223 return T(1)/T(m_real_domain.numPts());
1227template <
typename T, Direction D,
bool C>
1230std::pair<typename R2C<T,D,C>::cMF *,
IntVect>
1233#if (AMREX_SPACEDIM > 1)
1235 auto [cmf, order] = m_r2c_sub->getSpectralData();
1236 return std::make_pair(cmf, m_sub_helper.inverse_order(order));
1239 if (!m_cz.empty()) {
1241 }
else if (!m_cy.empty()) {
1248template <
typename T, Direction D,
bool C>
1255 auto const ncomp = m_info.batch_size;
1259 bool inmf_safe = m_sub_helper.ghost_safe(inmf.nGrowVect());
1260 MF inmf_sub, inmf_tmp;
1263 inmf_sub = m_sub_helper.make_alias_mf(inmf);
1264 incomp_sub = incomp;
1266 inmf_tmp.define(inmf.boxArray(), inmf.DistributionMap(), ncomp, 0);
1267 inmf_tmp.LocalCopy(inmf, incomp, 0, ncomp,
IntVect(0));
1268 inmf_sub = m_sub_helper.make_alias_mf(inmf_tmp);
1272 bool outmf_safe = m_sub_helper.ghost_safe(outmf.
nGrowVect());
1273 cMF outmf_sub, outmf_tmp;
1276 outmf_sub = m_sub_helper.make_alias_mf(outmf);
1277 outcomp_sub = outcomp;
1280 outmf_sub = m_sub_helper.make_alias_mf(outmf_tmp);
1284 m_r2c_sub->forward(inmf_sub, outmf_sub, incomp_sub, outcomp_sub);
1293 if (!m_cz.empty()) {
1296 (outmf, m_spectral_domain_x, m_cz,
IntVect(0), dtos);
1297 ParallelCopy(outmf, m_cz, cmd, 0, outcomp, ncomp, dtos);
1298 }
else if (!m_cy.empty()) {
1300 (outmf, m_spectral_domain_x, m_cy,
IntVect(0), m_dtos_y2x);
1301 ParallelCopy(outmf, m_cy, cmd, 0, outcomp, ncomp, m_dtos_y2x);
1308template <
typename T, Direction D,
bool C>
1316template <
typename T, Direction D,
bool C>
1318 Periodicity const& period,
int incomp,
int outcomp)
1322 auto const ncomp = m_info.batch_size;
1326 bool inmf_safe = m_sub_helper.ghost_safe(inmf.nGrowVect());
1327 cMF inmf_sub, inmf_tmp;
1330 inmf_sub = m_sub_helper.make_alias_mf(inmf);
1331 incomp_sub = incomp;
1333 inmf_tmp.define(inmf.boxArray(), inmf.DistributionMap(), ncomp, 0);
1334 inmf_tmp.LocalCopy(inmf, incomp, 0, ncomp,
IntVect(0));
1335 inmf_sub = m_sub_helper.make_alias_mf(inmf_tmp);
1339 bool outmf_safe = m_sub_helper.ghost_safe(outmf.nGrowVect());
1340 MF outmf_sub, outmf_tmp;
1343 outmf_sub = m_sub_helper.make_alias_mf(outmf);
1344 outcomp_sub = outcomp;
1346 IntVect const& ngtmp = m_sub_helper.make_safe_ghost(outmf.nGrowVect());
1347 outmf_tmp.define(outmf.boxArray(), outmf.DistributionMap(), ncomp, ngtmp);
1348 outmf_sub = m_sub_helper.make_alias_mf(outmf_tmp);
1352 IntVect const& subngout = m_sub_helper.make_iv(ngout);
1353 Periodicity
const& subperiod = m_sub_helper.make_periodicity(period);
1354 m_r2c_sub->backward_doit(inmf_sub, outmf_sub, subngout, subperiod, incomp_sub, outcomp_sub);
1357 outmf.LocalCopy(outmf_tmp, 0, outcomp, ncomp, outmf_tmp.nGrowVect());
1362 if (!m_cz.empty()) {
1364 MultiBlockCommMetaData cmd
1365 (m_cz, m_spectral_domain_z, inmf,
IntVect(0), dtos);
1367 }
else if (!m_cy.empty()) {
1368 MultiBlockCommMetaData cmd
1369 (m_cy, m_spectral_domain_y, inmf,
IntVect(0), m_dtos_x2y);
1370 ParallelCopy(m_cy, inmf, cmd, incomp, 0, ncomp, m_dtos_x2y);
1372 m_cx.ParallelCopy(inmf, incomp, 0, ncomp);
1374 backward_doit(outmf, ngout, period, outcomp);
1378template <
typename T, Direction D,
bool C>
1379std::pair<BoxArray,DistributionMapping>
1382#if (AMREX_SPACEDIM > 1)
1384 auto const& [ba, dm] = m_r2c_sub->getSpectralDataLayout();
1385 return std::make_pair(m_sub_helper.inverse_boxarray(ba), dm);
1389#if (AMREX_SPACEDIM == 3)
1390 if (!m_cz.empty()) {
1391 BoxList bl = m_cz.boxArray().boxList();
1392 for (
auto& b : bl) {
1393 auto lo = b.smallEnd();
1394 auto hi = b.bigEnd();
1395 std::swap(lo[0], lo[1]);
1396 std::swap(lo[1], lo[2]);
1397 std::swap(hi[0], hi[1]);
1398 std::swap(hi[1], hi[2]);
1402 return std::make_pair(
BoxArray(std::move(bl)), m_cz.DistributionMap());
1405#if (AMREX_SPACEDIM >= 2)
1406 if (!m_cy.empty()) {
1407 BoxList bl = m_cy.boxArray().boxList();
1408 for (
auto& b : bl) {
1409 auto lo = b.smallEnd();
1410 auto hi = b.bigEnd();
1411 std::swap(lo[0], lo[1]);
1412 std::swap(hi[0], hi[1]);
1416 return std::make_pair(
BoxArray(std::move(bl)), m_cy.DistributionMap());
1420 return std::make_pair(m_cx.boxArray(), m_cx.DistributionMap());
1425template <
typename T = Real, FFT::Direction D = FFT::Direction::both>
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define AMREX_ALWAYS_ASSERT(EX)
Definition AMReX_BLassert.H:50
#define AMREX_GPU_DEVICE
Definition AMReX_GpuQualifiers.H:18
amrex::ParmParse pp
Input file parser instance for the given namespace.
Definition AMReX_HypreIJIface.cpp:15
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
virtual void * alloc(std::size_t sz)=0
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:567
A class for managing a List of Boxes that share a common IndexType. This class implements operations ...
Definition AMReX_BoxList.H:52
BoxList & shift(int dir, int nzones)
Applies Box::shift(int,int) to each Box in the BoxList.
Definition AMReX_BoxList.cpp:565
__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__ int shortside(int &dir) const noexcept
Return length of shortest side. dir is modified to give direction with shortest side: 0....
Definition AMReX_Box.H:437
__host__ __device__ IntVectND< dim > size() const noexcept
Return the length of the BoxND.
Definition AMReX_Box.H:147
__host__ __device__ BoxND & growLo(int idir, int n_cell=1) noexcept
Grow the BoxND on the low end by n_cell cells in direction idir. NOTE: n_cell negative shrinks the Bo...
Definition AMReX_Box.H:662
__host__ __device__ IndexTypeND< dim > ixType() const noexcept
Return the indexing type.
Definition AMReX_Box.H:135
__host__ __device__ BoxND & growHi(int idir, int n_cell=1) noexcept
Grow the BoxND on the high end by n_cell cells in direction idir. NOTE: n_cell negative shrinks the B...
Definition AMReX_Box.H:673
__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
Long size() const noexcept
Length of the underlying processor map.
Definition AMReX_DistributionMapping.H:129
Definition AMReX_FFT_OpenBCSolver.H:11
3D Poisson solver for periodic, Dirichlet & Neumann boundaries in the first two dimensions,...
Definition AMReX_FFT_Poisson.H:151
Poisson solver for periodic, Dirichlet & Neumann boundaries using FFT.
Definition AMReX_FFT_Poisson.H:61
Parallel Discrete Fourier Transform.
Definition AMReX_FFT_R2C.H:41
std::conditional_t< C, cMF, std::conditional_t< std::is_same_v< T, Real >, MultiFab, FabArray< BaseFab< T > > > > MF
Definition AMReX_FFT_R2C.H:46
R2C & operator=(R2C const &)=delete
void forward(MF const &inmf, cMF &outmf, int incomp=0, int outcomp=0)
Forward transform.
Definition AMReX_FFT_R2C.H:1251
~R2C()
Definition AMReX_FFT_R2C.H:701
void setLocalDomain(std::array< int, 3 > const &local_start, std::array< int, 3 > const &local_size)
Set local domain.
Definition AMReX_FFT_R2C.H:752
R2C(Box const &domain, Info const &info=Info{})
Constructor.
Definition AMReX_FFT_R2C.H:454
void backward(cMF const &inmf, MF &outmf, int incomp=0, int outcomp=0)
Backward transform.
Definition AMReX_FFT_R2C.H:1311
void backward(MF &outmf, int outcomp=0)
Backward transform.
Definition AMReX_FFT_R2C.H:996
void post_forward_doit_1(F const &post_forward)
Definition AMReX_FFT_R2C.H:1172
T scalingFactor() const
Definition AMReX_FFT_R2C.H:1214
void post_forward_doit_0(F const &post_forward)
Definition AMReX_FFT_R2C.H:1117
std::pair< std::array< int, 3 >, std::array< int, 3 > > getLocalDomain() const
Get local domain.
Definition AMReX_FFT_R2C.H:761
void forward(RT const *in, CT *out)
Forward transform.
Definition AMReX_FFT_R2C.H:976
std::pair< BoxArray, DistributionMapping > getSpectralDataLayout() const
Get BoxArray and DistributionMapping for spectral data.
Definition AMReX_FFT_R2C.H:1380
std::pair< cMF *, IntVect > getSpectralData() const
Get the internal spectral data.
R2C(std::array< int, 3 > const &domain_size, Info const &info=Info{})
Constructor.
Definition AMReX_FFT_R2C.H:696
FabArray< BaseFab< GpuComplex< T > > > cMF
Definition AMReX_FFT_R2C.H:43
std::pair< std::array< int, 3 >, std::array< int, 3 > > getLocalSpectralDomain() const
Get local spectral domain.
Definition AMReX_FFT_R2C.H:787
void forwardThenBackward(MF const &inmf, MF &outmf, F const &post_forward, int incomp=0, int outcomp=0)
Forward and then backward transform.
Definition AMReX_FFT_R2C.H:193
void setLocalSpectralDomain(std::array< int, 3 > const &local_start, std::array< int, 3 > const &local_size)
Set local spectral domain.
Definition AMReX_FFT_R2C.H:778
void backward(CT const *in, RT *out)
Backward transform.
Definition AMReX_FFT_R2C.H:1067
void forward(MF const &inmf, int incomp=0)
Forward transform.
Definition AMReX_FFT_R2C.H:863
IntVect nGrowVect() const noexcept
Definition AMReX_FabArrayBase.H:80
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 fabbox(int K) const noexcept
Return the Kth FABs Box in the FabArray. That is, the region the Kth fab is actually defined on.
Definition AMReX_FabArrayBase.cpp:220
const BoxArray & boxArray() const noexcept
Return a constant reference to the BoxArray that defines the valid region associated with this FabArr...
Definition AMReX_FabArrayBase.H:95
void setFab(int boxno, std::unique_ptr< FAB > elem)
Explicitly set the Kth FAB in the FabArray to point to elem.
Definition AMReX_FabArray.H:2355
void ParallelCopy(const FabArray< FAB > &src, const Periodicity &period=Periodicity::NonPeriodic(), CpOp op=FabArrayBase::COPY)
Definition AMReX_FabArray.H:847
typename std::conditional_t< IsBaseFab< BaseFab< GpuComplex< T > > >::value, BaseFab< GpuComplex< T > >, FABType >::value_type value_type
Definition AMReX_FabArray.H:358
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:2171
void LocalCopy(FabArray< SFAB > const &src, int scomp, int dcomp, int ncomp, IntVect const &nghost)
Perform local copy of FabArray data.
Definition AMReX_FabArray.H:1952
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
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
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:243
int MyProc() noexcept
Definition AMReX_ParallelDescriptor.H:128
int NProcs() noexcept
Definition AMReX_ParallelDescriptor.H:255
Definition AMReX_FFT_Helper.H:46
Direction
Definition AMReX_FFT_Helper.H:48
void dtod_memcpy_async(void *p_d_dst, const void *p_d_src, const std::size_t sz) noexcept
Definition AMReX_GpuDevice.H:329
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:263
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
MPI_Comm Communicator() noexcept
Definition AMReX_ParallelDescriptor.H:223
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:27
bool is_aligned(const void *p, std::size_t alignment) noexcept
Definition AMReX_Arena.H:35
BoxArray decompose(Box const &domain, int nboxes, Array< bool, 3 > const &decomp, bool no_overlap)
Decompose domain box into BoxArray.
Definition AMReX_BoxArray.cpp:1940
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:230
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
Arena * The_Arena()
Definition AMReX_Arena.cpp:783
__host__ __device__ constexpr T elemwiseMin(T const &a, T const &b) noexcept
Definition AMReX_Algorithm.H:49
Definition AMReX_FFT_Helper.H:58
bool twod_mode
Definition AMReX_FFT_Helper.H:69
bool oned_mode
We might have a special twod_mode: nx or ny == 1 && nz > 1.
Definition AMReX_FFT_Helper.H:72
int batch_size
Batched FFT size. Only support in R2C, not R2X.
Definition AMReX_FFT_Helper.H:75
DomainStrategy domain_strategy
Domain composition strategy.
Definition AMReX_FFT_Helper.H:60
int nprocs
Max number of processes to use.
Definition AMReX_FFT_Helper.H:78
int pencil_threshold
Definition AMReX_FFT_Helper.H:64
Definition AMReX_FFT_Helper.H:134
std::conditional_t< std::is_same_v< float, T >, cuComplex, cuDoubleComplex > VendorComplex
Definition AMReX_FFT_Helper.H:138
FabArray memory allocation information.
Definition AMReX_FabArray.H:66
MFInfo & SetAlloc(bool a) noexcept
Definition AMReX_FabArray.H:73
This class specializes behaviour on local copies and unpacking receive buffers.
Definition AMReX_NonLocalBC.H:626
Contains information about which components take part of the data transaction.
Definition AMReX_NonLocalBC.H:539
int n_components
Definition AMReX_NonLocalBC.H:542
Communication datatype (note: this structure also works without MPI)
Definition AMReX_ccse-mpi.H:78