1#ifndef AMREX_FFT_R2C_H_
2#define AMREX_FFT_R2C_H_
3#include <AMReX_Config.H>
20template <
typename T>
class OpenBCSolver;
21template <
typename T>
class Poisson;
22template <
typename T>
class PoissonHybrid;
45template <
typename T = Real, FFT::Direction D = FFT::Direction::both,
bool C = false>
50 using MF = std::conditional_t
51 <
C,
cMF, std::conditional_t<std::is_same_v<T,Real>,
55 template <
typename U>
friend class Poisson;
75 explicit R2C (std::array<int,AMREX_SPACEDIM>
const& domain_size,
109 std::array<int,AMREX_SPACEDIM>
const& local_size);
124 std::pair<std::array<int,AMREX_SPACEDIM>,std::array<int,AMREX_SPACEDIM>>
156 std::array<int,AMREX_SPACEDIM>
const& local_size);
176 std::pair<std::array<int,AMREX_SPACEDIM>,std::array<int,AMREX_SPACEDIM>>
201 std::enable_if_t<DIR == Direction::both, int> = 0>
203 int incomp = 0,
int outcomp = 0)
207 "FFT::R2C::forwardThenBackward(post_forward) currently supports only !twod_mode and batch_size==1");
241 void forward (
MF const& inmf,
cMF& outmf,
int incomp = 0,
int outcomp = 0);
259 template <
typename RT,
typename CT,
Direction DIR=D,
bool CP=
C,
262 && ((
sizeof(RT)*2 ==
sizeof(CT) && !CP) ||
263 (
sizeof(RT) ==
sizeof(CT) && CP)),
int> = 0>
275 template <Direction DIR=D, std::enable_if_t<DIR == Direction::both,
int> = 0>
291 void backward (
cMF const& inmf,
MF& outmf,
int incomp = 0,
int outcomp = 0);
309 template <
typename CT,
typename RT,
Direction DIR=D,
bool CP=
C,
312 && ((
sizeof(RT)*2 ==
sizeof(CT) && !CP) ||
313 (
sizeof(RT) ==
sizeof(CT) && CP)),
int> = 0>
356 template <
typename F>
359 template <
typename F>
369 void prepare_openbc ();
375 void backward_doit (
cMF const& inmf,
MF& outmf,
378 int incomp = 0,
int outcomp = 0);
380 std::pair<Plan<T>,
Plan<T>> make_c2c_plans (
cMF& inout,
int ndims)
const;
382 static Box make_domain_x (
Box const& domain)
397 static Box make_domain_y (
Box const& domain)
402 domain.length(2)-1)),
407 domain.length(2)-1)),
412 static Box make_domain_z (
Box const& domain)
417 domain.length(1)-1)),
422 domain.length(1)-1)),
427 static std::pair<BoxArray,DistributionMapping>
428 make_layout_from_local_domain (std::array<int,AMREX_SPACEDIM>
const& local_start,
429 std::array<int,AMREX_SPACEDIM>
const& local_size);
431 template <
typename FA,
typename RT>
432 std::pair<std::unique_ptr<char,DataDeleter>,std::size_t>
433 install_raw_ptr (FA& fa, RT
const* p);
435 Plan<T> m_fft_fwd_x{};
436 Plan<T> m_fft_bwd_x{};
437 Plan<T> m_fft_fwd_y{};
438 Plan<T> m_fft_bwd_y{};
439 Plan<T> m_fft_fwd_z{};
440 Plan<T> m_fft_bwd_z{};
441 Plan<T> m_fft_fwd_x_half{};
442 Plan<T> m_fft_bwd_x_half{};
447 std::unique_ptr<MultiBlockCommMetaData> m_cmd_x2y;
448 std::unique_ptr<MultiBlockCommMetaData> m_cmd_y2x;
449 std::unique_ptr<MultiBlockCommMetaData> m_cmd_y2z;
450 std::unique_ptr<MultiBlockCommMetaData> m_cmd_z2y;
451 std::unique_ptr<MultiBlockCommMetaData> m_cmd_x2z;
452 std::unique_ptr<MultiBlockCommMetaData> m_cmd_z2x;
453 std::unique_ptr<MultiBlockCommMetaData> m_cmd_x2z_half;
454 std::unique_ptr<MultiBlockCommMetaData> m_cmd_z2x_half;
459 RotateFwd m_dtos_x2z{};
460 RotateBwd m_dtos_z2x{};
468 mutable cMF m_raw_cmf;
470 std::unique_ptr<char,DataDeleter> m_data_1;
471 std::unique_ptr<char,DataDeleter> m_data_2;
474 Box m_spectral_domain_x;
475 Box m_spectral_domain_y;
476 Box m_spectral_domain_z;
478 std::unique_ptr<R2C<T,D,C>> m_r2c_sub;
479 detail::SubHelper m_sub_helper;
483 bool m_do_alld_fft =
false;
484 bool m_slab_decomp =
false;
485 bool m_openbc_half =
false;
488template <
typename T, Direction D,
bool C>
490 : m_real_domain(domain),
491 m_spectral_domain_x(make_domain_x(domain)),
492#if (AMREX_SPACEDIM >= 2)
493 m_spectral_domain_y(make_domain_y(domain)),
494#if (AMREX_SPACEDIM == 3)
495 m_spectral_domain_z(make_domain_z(domain)),
498 m_sub_helper(domain),
503 static_assert(std::is_same_v<float,T> || std::is_same_v<double,T>);
506#if (AMREX_SPACEDIM == 2)
511 int(domain.
length(1) > 1) +
512 int(domain.
length(2) > 1)) >= 2);
517 Box subbox = m_sub_helper.make_box(m_real_domain);
518 if (subbox.
size() != m_real_domain.
size()) {
519 m_r2c_sub = std::make_unique<R2C<T,D,C>>(subbox, m_info);
527#if (AMREX_SPACEDIM == 3)
532 int shortside = m_real_domain.
shortside();
543 m_slab_decomp =
true;
545 m_slab_decomp =
true;
556 m_rx.define(bax, dmx, ncomp, 0,
MFInfo().SetAlloc(
false));
560 for (
auto & b : bl) {
562 b.setBig(0, m_spectral_domain_x.
bigEnd(0));
565 m_cx.
define(cbax, dmx, ncomp, 0,
MFInfo().SetAlloc(
false));
577#if (AMREX_SPACEDIM >= 2)
579 if ((m_real_domain.
length(1) > 1) && !m_slab_decomp && !m_info.
oned_mode)
583 if (cbay.size() == dmx.size()) {
586 cdmy = detail::make_iota_distromap(cbay.size());
588 m_cy.
define(cbay, cdmy, ncomp, 0,
MFInfo().SetAlloc(
false));
592#if (AMREX_SPACEDIM == 3)
594 m_real_domain.
length(1) > 1 &&
595 m_real_domain.
length(2) > 1)
598 {
false,
true,
true},
true);
600 if (cbaz.size() == dmx.size()) {
602 }
else if (cbaz.size() == cdmy.
size()) {
605 cdmz = detail::make_iota_distromap(cbaz.size());
607 m_cz.
define(cbaz, cdmz, ncomp, 0,
MFInfo().SetAlloc(
false));
613 m_data_1 = detail::make_mfs_share(m_rx, m_cx);
614 m_data_2 = detail::make_mfs_share(m_cz, m_cz);
616 m_data_1 = detail::make_mfs_share(m_rx, m_cz);
617 m_data_2 = detail::make_mfs_share(m_cy, m_cy);
619 if (myproc < m_cx.
size()) {
622 m_cx.
setFab(myproc, FAB(box, ncomp, m_rx[myproc].dataPtr()));
627 m_data_1 = detail::make_mfs_share(m_rx, m_cz);
628 m_data_2 = detail::make_mfs_share(m_cx, m_cx);
630 m_data_1 = detail::make_mfs_share(m_rx, m_cy);
631 m_data_2 = detail::make_mfs_share(m_cx, m_cz);
639#if (AMREX_SPACEDIM >= 2)
640 if (! m_cy.
empty()) {
642 m_cmd_x2y = std::make_unique<MultiBlockCommMetaData>
643 (m_cy, m_spectral_domain_y, m_cx,
IntVect(0), m_dtos_x2y);
644 m_cmd_y2x = std::make_unique<MultiBlockCommMetaData>
645 (m_cx, m_spectral_domain_x, m_cy,
IntVect(0), m_dtos_y2x);
648#if (AMREX_SPACEDIM == 3)
649 if (! m_cz.
empty() ) {
652 m_cmd_x2z = std::make_unique<MultiBlockCommMetaData>
653 (m_cz, m_spectral_domain_z, m_cx,
IntVect(0), m_dtos_x2z);
654 m_cmd_z2x = std::make_unique<MultiBlockCommMetaData>
655 (m_cx, m_spectral_domain_x, m_cz,
IntVect(0), m_dtos_z2x);
658 m_cmd_y2z = std::make_unique<MultiBlockCommMetaData>
659 (m_cz, m_spectral_domain_z, m_cy,
IntVect(0), m_dtos_y2z);
660 m_cmd_z2y = std::make_unique<MultiBlockCommMetaData>
661 (m_cy, m_spectral_domain_y, m_cz,
IntVect(0), m_dtos_z2y);
670 if (myproc < m_rx.size())
673 int ndims = m_slab_decomp ? 2 : 1;
674 std::tie(m_fft_fwd_x, m_fft_bwd_x) = make_c2c_plans(m_cx, ndims);
676 Box const& box = m_rx.box(myproc);
677 auto* pr = m_rx[myproc].dataPtr();
680 m_fft_fwd_x.template init_r2c<Direction::forward>(box, pr, pc, m_slab_decomp, ncomp);
681 m_fft_bwd_x = m_fft_fwd_x;
684 m_fft_fwd_x.template init_r2c<Direction::forward>(box, pr, pc, m_slab_decomp, ncomp);
687 m_fft_bwd_x.template init_r2c<Direction::backward>(box, pr, pc, m_slab_decomp, ncomp);
693#if (AMREX_SPACEDIM >= 2)
694 if (! m_cy.
empty()) {
695 std::tie(m_fft_fwd_y, m_fft_bwd_y) = make_c2c_plans(m_cy,1);
698#if (AMREX_SPACEDIM == 3)
699 if (! m_cz.
empty()) {
700 std::tie(m_fft_fwd_z, m_fft_bwd_z) = make_c2c_plans(m_cz,1);
707 m_data_1 = detail::make_mfs_share(m_rx, m_cx);
708 std::tie(m_fft_fwd_x, m_fft_bwd_x) = make_c2c_plans(m_cx,AMREX_SPACEDIM);
710 m_data_1 = detail::make_mfs_share(m_rx, m_rx);
711 m_data_2 = detail::make_mfs_share(m_cx, m_cx);
713 auto const& len = m_real_domain.
length();
714 auto* pr = (
void*)m_rx[0].dataPtr();
715 auto* pc = (
void*)m_cx[0].dataPtr();
717 m_fft_fwd_x.template init_r2c<Direction::forward>(len, pr, pc,
false, ncomp);
718 m_fft_bwd_x = m_fft_fwd_x;
721 m_fft_fwd_x.template init_r2c<Direction::forward>(len, pr, pc,
false, ncomp);
724 m_fft_bwd_x.template init_r2c<Direction::backward>(len, pr, pc,
false, ncomp);
731template <
typename T, Direction D,
bool C>
736template <
typename T, Direction D,
bool C>
739 if (m_fft_bwd_x.plan != m_fft_fwd_x.plan) {
740 m_fft_bwd_x.destroy();
742 if (m_fft_bwd_y.plan != m_fft_fwd_y.plan) {
743 m_fft_bwd_y.destroy();
745 if (m_fft_bwd_z.plan != m_fft_fwd_z.plan) {
746 m_fft_bwd_z.destroy();
748 m_fft_fwd_x.destroy();
749 m_fft_fwd_y.destroy();
750 m_fft_fwd_z.destroy();
751 if (m_fft_bwd_x_half.plan != m_fft_fwd_x_half.plan) {
752 m_fft_bwd_x_half.destroy();
754 m_fft_fwd_x_half.destroy();
757template <
typename T, Direction D,
bool C>
758std::pair<BoxArray,DistributionMapping>
760 std::array<int,AMREX_SPACEDIM>
const& local_size)
764 Box bx(lo, lo+len-1);
771 pmap.reserve(allboxes.size());
772 for (
int i = 0; i < allboxes.size(); ++i) {
773 if (allboxes[i].ok()) {
777 allboxes.erase(std::remove_if(allboxes.begin(), allboxes.end(),
778 [=] (
Box const& b) { return b.isEmpty(); }),
780 BoxList bl(std::move(allboxes));
781 return std::make_pair(BoxArray(std::move(bl)), DistributionMapping(std::move(pmap)));
783 return std::make_pair(BoxArray(bx), DistributionMapping(Vector<int>({0})));
787template <
typename T, Direction D,
bool C>
789 std::array<int,AMREX_SPACEDIM>
const& local_size)
791 auto const& [ba, dm] = make_layout_from_local_domain(local_start, local_size);
792 m_raw_mf =
MF(ba, dm, m_rx.nComp(), 0,
MFInfo().SetAlloc(
false));
795template <
typename T, Direction D,
bool C>
796std::pair<std::array<int,AMREX_SPACEDIM>,std::array<int,AMREX_SPACEDIM>>
799 m_raw_mf =
MF(m_rx.boxArray(), m_rx.DistributionMap(), m_rx.nComp(), 0,
803 if (myproc < m_rx.size()) {
804 Box const& box = m_rx.box(myproc);
805 return std::make_pair(box.
smallEnd().toArray(),
808 return std::make_pair(std::array<int,AMREX_SPACEDIM>{
AMREX_D_DECL(0,0,0)},
813template <
typename T, Direction D,
bool C>
815 std::array<int,AMREX_SPACEDIM>
const& local_size)
817 auto const& [ba, dm] = make_layout_from_local_domain(local_start, local_size);
818 m_raw_cmf =
cMF(ba, dm, m_rx.nComp(), 0,
MFInfo().SetAlloc(
false));
821template <
typename T, Direction D,
bool C>
822std::pair<std::array<int,AMREX_SPACEDIM>,std::array<int,AMREX_SPACEDIM>>
825 auto const ncomp = m_info.batch_size;
826 auto const& [ba, dm] = getSpectralDataLayout();
831 if (myproc < m_raw_cmf.size()) {
832 Box const& box = m_raw_cmf.box(myproc);
833 return std::make_pair(box.
smallEnd().toArray(), box.
length().toArray());
835 return std::make_pair(std::array<int,AMREX_SPACEDIM>{
AMREX_D_DECL(0,0,0)},
840template <
typename T, Direction D,
bool C>
843 if (
C || m_r2c_sub) {
amrex::Abort(
"R2C: OpenBC not supported with reduced dimensions or complex inputs"); }
845#if (AMREX_SPACEDIM == 3)
846 if (m_do_alld_fft) {
return; }
848 auto const ncomp = m_info.batch_size;
850 if (m_slab_decomp && ! m_fft_fwd_x_half.defined) {
851 auto* fab = detail::get_fab(m_rx);
853 Box bottom_half = m_real_domain;
854 bottom_half.
growHi(2,-m_real_domain.length(2)/2);
855 Box box = fab->box() & bottom_half;
857 auto* pr = fab->dataPtr();
859 detail::get_fab(m_cx)->dataPtr();
861 m_fft_fwd_x_half.template init_r2c<Direction::forward>
862 (box, pr, pc, m_slab_decomp, ncomp);
863 m_fft_bwd_x_half = m_fft_fwd_x_half;
866 m_fft_fwd_x_half.template init_r2c<Direction::forward>
867 (box, pr, pc, m_slab_decomp, ncomp);
870 m_fft_bwd_x_half.template init_r2c<Direction::backward>
871 (box, pr, pc, m_slab_decomp, ncomp);
878 if (m_cmd_x2z && ! m_cmd_x2z_half) {
879 Box bottom_half = m_spectral_domain_z;
882 bottom_half.
growHi(0,-m_spectral_domain_z.length(0)/2);
883 m_cmd_x2z_half = std::make_unique<MultiBlockCommMetaData>
884 (m_cz, bottom_half, m_cx,
IntVect(0), m_dtos_x2z);
887 if (m_cmd_z2x && ! m_cmd_z2x_half) {
888 Box bottom_half = m_spectral_domain_x;
889 bottom_half.
growHi(2,-m_spectral_domain_x.length(2)/2);
890 m_cmd_z2x_half = std::make_unique<MultiBlockCommMetaData>
891 (m_cx, bottom_half, m_cz,
IntVect(0), m_dtos_z2x);
896template <
typename T, Direction D,
bool C>
903 auto const ncomp = m_info.batch_size;
906 if (m_sub_helper.ghost_safe(inmf.nGrowVect())) {
907 m_r2c_sub->forward(m_sub_helper.make_alias_mf(inmf), incomp);
909 MF tmp(inmf.boxArray(), inmf.DistributionMap(), ncomp, 0);
910 tmp.LocalCopy(inmf, incomp, 0, ncomp,
IntVect(0));
911 m_r2c_sub->forward(m_sub_helper.make_alias_mf(tmp),0);
916 if (&m_rx != &inmf) {
917 m_rx.ParallelCopy(inmf, incomp, 0, ncomp);
922 m_fft_fwd_x.template compute_c2c<Direction::forward>();
924 m_fft_fwd_x.template compute_r2c<Direction::forward>();
929 auto& fft_x = m_openbc_half ? m_fft_fwd_x_half : m_fft_fwd_x;
931 fft_x.template compute_c2c<Direction::forward>();
933 fft_x.template compute_r2c<Direction::forward>();
937 ParallelCopy(m_cy, m_cx, *m_cmd_x2y, 0, 0, ncomp, m_dtos_x2y);
939 m_fft_fwd_y.template compute_c2c<Direction::forward>();
942 ParallelCopy(m_cz, m_cy, *m_cmd_y2z, 0, 0, ncomp, m_dtos_y2z);
944#if (AMREX_SPACEDIM == 3)
945 else if ( m_cmd_x2z) {
950 {components, m_dtos_x2z};
951 auto handler = ParallelCopy_nowait(m_cz, m_cx, *m_cmd_x2z_half, packing);
953 Box upper_half = m_spectral_domain_z;
956 upper_half.
growLo (0,-m_spectral_domain_z.length(0)/2);
957 m_cz.setVal(0, upper_half, 0, ncomp);
959 ParallelCopy_finish(m_cz, std::move(handler), *m_cmd_x2z_half, packing);
961 ParallelCopy(m_cz, m_cx, *m_cmd_x2z, 0, 0, ncomp, m_dtos_x2z);
965 m_fft_fwd_z.template compute_c2c<Direction::forward>();
968template <
typename T, Direction D,
bool C>
969template <
typename FA,
typename RT>
970std::pair<std::unique_ptr<char,DataDeleter>,std::size_t>
975 using FAB =
typename FA::FABType::value_type;
976 using T_FAB =
typename FAB::value_type;
977 static_assert(
sizeof(T_FAB) ==
sizeof(RT));
979 auto const ncomp = m_info.batch_size;
980 auto const& ia = fa.IndexArray();
985 if ( ! ia.empty() ) {
987 Box const& box = fa.fabbox(K);
991 sz =
sizeof(T_FAB) * box.
numPts() * ncomp;
994 fa.setFab(K, FAB(box,ncomp,
pp));
998 return std::make_pair(std::unique_ptr<char,DataDeleter>{},std::size_t(0));
1000 return std::make_pair(std::unique_ptr<char,DataDeleter>
1006template <
typename T, Direction D,
bool C>
1007template <
typename RT,
typename CT,
Direction DIR,
bool CP,
1010 && ((
sizeof(RT)*2 ==
sizeof(CT) && !CP) ||
1011 (
sizeof(RT) ==
sizeof(CT) && CP)),
int> >
1014 auto [rdata, rsz] = install_raw_ptr(m_raw_mf, in);
1015 auto [cdata, csz] = install_raw_ptr(m_raw_cmf, out);
1030template <
typename T, Direction D,
bool C>
1031template <Direction DIR, std::enable_if_t<DIR == Direction::both,
int> >
1037template <
typename T, Direction D,
bool C>
1043 auto const ncomp = m_info.batch_size;
1046 if (m_sub_helper.ghost_safe(outmf.nGrowVect())) {
1047 MF submf = m_sub_helper.make_alias_mf(outmf);
1048 IntVect const& subngout = m_sub_helper.make_iv(ngout);
1049 Periodicity const& subperiod = m_sub_helper.make_periodicity(period);
1050 m_r2c_sub->backward_doit(submf, subngout, subperiod, outcomp);
1052 MF tmp(outmf.boxArray(), outmf.DistributionMap(), ncomp,
1053 m_sub_helper.make_safe_ghost(outmf.nGrowVect()));
1054 this->backward_doit(tmp, ngout, period, 0);
1055 outmf.LocalCopy(tmp, 0, outcomp, ncomp, tmp.nGrowVect());
1060 if (m_do_alld_fft) {
1062 m_fft_bwd_x.template compute_c2c<Direction::backward>();
1064 m_fft_bwd_x.template compute_r2c<Direction::backward>();
1066 outmf.ParallelCopy(m_rx, 0, outcomp, ncomp,
IntVect(0),
1071 m_fft_bwd_z.template compute_c2c<Direction::backward>();
1073 ParallelCopy(m_cy, m_cz, *m_cmd_z2y, 0, 0, ncomp, m_dtos_z2y);
1075#if (AMREX_SPACEDIM == 3)
1076 else if ( m_cmd_z2x) {
1077 auto const& cmd = m_openbc_half ? m_cmd_z2x_half : m_cmd_z2x;
1078 ParallelCopy(m_cx, m_cz, *cmd, 0, 0, ncomp, m_dtos_z2x);
1082 m_fft_bwd_y.template compute_c2c<Direction::backward>();
1084 ParallelCopy(m_cx, m_cy, *m_cmd_y2x, 0, 0, ncomp, m_dtos_y2x);
1087 auto& fft_x = m_openbc_half ? m_fft_bwd_x_half : m_fft_bwd_x;
1089 fft_x.template compute_c2c<Direction::backward>();
1091 fft_x.template compute_r2c<Direction::backward>();
1093 outmf.ParallelCopy(m_rx, 0, outcomp, ncomp,
IntVect(0),
1097template <
typename T, Direction D,
bool C>
1098template <
typename CT,
typename RT,
Direction DIR,
bool CP,
1101 && ((
sizeof(RT)*2 ==
sizeof(CT) && !CP) ||
1102 (
sizeof(RT) ==
sizeof(CT) && CP)),
int> >
1105 auto [rdata, rsz] = install_raw_ptr(m_raw_mf, out);
1106 auto [cdata, csz] = install_raw_ptr(m_raw_cmf, in);
1121template <
typename T, Direction D,
bool C>
1128 auto* fab = detail::get_fab(inout);
1129 if (!fab) {
return {fwd, bwd};}
1131 Box const& box = fab->box();
1134 auto const ncomp = m_info.batch_size;
1136#ifdef AMREX_USE_SYCL
1137 fwd.template init_c2c<Direction::forward>(box, pio, ncomp, ndims);
1141 fwd.template init_c2c<Direction::forward>(box, pio, ncomp, ndims);
1144 bwd.template init_c2c<Direction::backward>(box, pio, ncomp, ndims);
1151template <
typename T, Direction D,
bool C>
1152template <
typename F>
1155 if (m_info.twod_mode || m_info.batch_size > 1) {
1157#if (AMREX_SPACEDIM > 1)
1158 }
else if (m_r2c_sub) {
1160#if (AMREX_SPACEDIM == 2)
1162 m_r2c_sub->post_forward_doit_1
1165 post_forward(0, i, 0, sp);
1168 if (m_real_domain.length(0) == 1 && m_real_domain.length(1) == 1) {
1170 m_r2c_sub->post_forward_doit_1
1173 post_forward(0, 0, i, sp);
1175 }
else if (m_real_domain.length(0) == 1 && m_real_domain.length(2) == 1) {
1177 m_r2c_sub->post_forward_doit_1
1180 post_forward(0, i, 0, sp);
1182 }
else if (m_real_domain.length(0) == 1) {
1184 m_r2c_sub->post_forward_doit_1
1187 post_forward(0, i, j, sp);
1189 }
else if (m_real_domain.length(1) == 1) {
1191 m_r2c_sub->post_forward_doit_1
1194 post_forward(i, 0, j, sp);
1197 amrex::Abort(
"R2c::post_forward_doit_0: how did this happen?");
1202 this->post_forward_doit_1(post_forward);
1206template <
typename T, Direction D,
bool C>
1207template <
typename F>
1210 if (m_info.twod_mode || m_info.batch_size > 1) {
1212 }
else if (m_r2c_sub) {
1213 amrex::Abort(
"R2C::post_forward_doit_1: How did this happen?");
1215 if ( ! m_cz.empty()) {
1216 auto* spectral_fab = detail::get_fab(m_cz);
1218 auto const& a = spectral_fab->array();
1222 post_forward(jx,ky,iz,a(iz,jx,ky));
1225 }
else if ( ! m_cy.empty()) {
1226 auto* spectral_fab = detail::get_fab(m_cy);
1228 auto const& a = spectral_fab->array();
1232 post_forward(jx,iy,k,a(iy,jx,k));
1236 auto* spectral_fab = detail::get_fab(m_cx);
1238 auto const& a = spectral_fab->array();
1242 post_forward(i,j,k,a(i,j,k));
1249template <
typename T, Direction D,
bool C>
1252#if (AMREX_SPACEDIM == 3)
1253 if (m_info.oned_mode && !m_info.twod_mode) {
1254 return T(1)/T(
Long(m_real_domain.length(0)));
1255 }
else if (m_info.twod_mode) {
1256 return T(1)/T(
Long(m_real_domain.length(0)) *
1257 Long(m_real_domain.length(1)));
1259#elif (AMREX_SPACEDIM == 2)
1260 if (m_info.oned_mode) {
1261 return T(1)/T(
Long(m_real_domain.length(0)));
1265 return T(1)/T(m_real_domain.numPts());
1269template <
typename T, Direction D,
bool C>
1272std::pair<typename R2C<T,D,C>::cMF *,
IntVect>
1275#if (AMREX_SPACEDIM > 1)
1277 auto [cmf, order] = m_r2c_sub->getSpectralData();
1278 return std::make_pair(cmf, m_sub_helper.inverse_order(order));
1281 if (!m_cz.empty()) {
1283 }
else if (!m_cy.empty()) {
1290template <
typename T, Direction D,
bool C>
1297 auto const ncomp = m_info.batch_size;
1301 bool inmf_safe = m_sub_helper.ghost_safe(inmf.nGrowVect());
1302 MF inmf_sub, inmf_tmp;
1305 inmf_sub = m_sub_helper.make_alias_mf(inmf);
1306 incomp_sub = incomp;
1308 inmf_tmp.define(inmf.boxArray(), inmf.DistributionMap(), ncomp, 0);
1309 inmf_tmp.LocalCopy(inmf, incomp, 0, ncomp,
IntVect(0));
1310 inmf_sub = m_sub_helper.make_alias_mf(inmf_tmp);
1314 bool outmf_safe = m_sub_helper.ghost_safe(outmf.
nGrowVect());
1315 cMF outmf_sub, outmf_tmp;
1318 outmf_sub = m_sub_helper.make_alias_mf(outmf);
1319 outcomp_sub = outcomp;
1322 outmf_sub = m_sub_helper.make_alias_mf(outmf_tmp);
1326 m_r2c_sub->forward(inmf_sub, outmf_sub, incomp_sub, outcomp_sub);
1335 if (!m_cz.empty()) {
1338 (outmf, m_spectral_domain_x, m_cz,
IntVect(0), dtos);
1339 ParallelCopy(outmf, m_cz, cmd, 0, outcomp, ncomp, dtos);
1340 }
else if (!m_cy.empty()) {
1342 (outmf, m_spectral_domain_x, m_cy,
IntVect(0), m_dtos_y2x);
1343 ParallelCopy(outmf, m_cy, cmd, 0, outcomp, ncomp, m_dtos_y2x);
1350template <
typename T, Direction D,
bool C>
1358template <
typename T, Direction D,
bool C>
1360 Periodicity const& period,
int incomp,
int outcomp)
1364 auto const ncomp = m_info.batch_size;
1368 bool inmf_safe = m_sub_helper.ghost_safe(inmf.nGrowVect());
1369 cMF inmf_sub, inmf_tmp;
1372 inmf_sub = m_sub_helper.make_alias_mf(inmf);
1373 incomp_sub = incomp;
1375 inmf_tmp.define(inmf.boxArray(), inmf.DistributionMap(), ncomp, 0);
1376 inmf_tmp.LocalCopy(inmf, incomp, 0, ncomp,
IntVect(0));
1377 inmf_sub = m_sub_helper.make_alias_mf(inmf_tmp);
1381 bool outmf_safe = m_sub_helper.ghost_safe(outmf.nGrowVect());
1382 MF outmf_sub, outmf_tmp;
1385 outmf_sub = m_sub_helper.make_alias_mf(outmf);
1386 outcomp_sub = outcomp;
1388 IntVect const& ngtmp = m_sub_helper.make_safe_ghost(outmf.nGrowVect());
1389 outmf_tmp.define(outmf.boxArray(), outmf.DistributionMap(), ncomp, ngtmp);
1390 outmf_sub = m_sub_helper.make_alias_mf(outmf_tmp);
1394 IntVect const& subngout = m_sub_helper.make_iv(ngout);
1395 Periodicity
const& subperiod = m_sub_helper.make_periodicity(period);
1396 m_r2c_sub->backward_doit(inmf_sub, outmf_sub, subngout, subperiod, incomp_sub, outcomp_sub);
1399 outmf.LocalCopy(outmf_tmp, 0, outcomp, ncomp, outmf_tmp.nGrowVect());
1404 if (!m_cz.empty()) {
1406 MultiBlockCommMetaData cmd
1407 (m_cz, m_spectral_domain_z, inmf,
IntVect(0), dtos);
1409 }
else if (!m_cy.empty()) {
1410 MultiBlockCommMetaData cmd
1411 (m_cy, m_spectral_domain_y, inmf,
IntVect(0), m_dtos_x2y);
1412 ParallelCopy(m_cy, inmf, cmd, incomp, 0, ncomp, m_dtos_x2y);
1414 m_cx.ParallelCopy(inmf, incomp, 0, ncomp);
1416 backward_doit(outmf, ngout, period, outcomp);
1420template <
typename T, Direction D,
bool C>
1421std::pair<BoxArray,DistributionMapping>
1424#if (AMREX_SPACEDIM > 1)
1426 auto const& [ba, dm] = m_r2c_sub->getSpectralDataLayout();
1427 return std::make_pair(m_sub_helper.inverse_boxarray(ba), dm);
1431#if (AMREX_SPACEDIM == 3)
1432 if (!m_cz.empty()) {
1433 BoxList bl = m_cz.boxArray().boxList();
1434 for (
auto& b : bl) {
1435 auto lo = b.smallEnd();
1436 auto hi = b.bigEnd();
1437 std::swap(lo[0], lo[1]);
1438 std::swap(lo[1], lo[2]);
1439 std::swap(hi[0], hi[1]);
1440 std::swap(hi[1], hi[2]);
1444 return std::make_pair(
BoxArray(std::move(bl)), m_cz.DistributionMap());
1447#if (AMREX_SPACEDIM >= 2)
1448 if (!m_cy.empty()) {
1449 BoxList bl = m_cy.boxArray().boxList();
1450 for (
auto& b : bl) {
1451 auto lo = b.smallEnd();
1452 auto hi = b.bigEnd();
1453 std::swap(lo[0], lo[1]);
1454 std::swap(hi[0], hi[1]);
1458 return std::make_pair(
BoxArray(std::move(bl)), m_cy.DistributionMap());
1462 return std::make_pair(m_cx.boxArray(), m_cx.DistributionMap());
1467template <
typename T = Real, FFT::Direction D = FFT::Direction::both>
#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
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:568
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
Convolution-based solver for open boundary conditions using Green's functions.
Definition AMReX_FFT_OpenBCSolver.H:24
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
Parallel Discrete Fourier Transform.
Definition AMReX_FFT_R2C.H:47
std::conditional_t< C, cMF, std::conditional_t< std::is_same_v< T, Real >, MultiFab, FabArray< BaseFab< T > > > > MF
Definition AMReX_FFT_R2C.H:52
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:1293
~R2C()
Definition AMReX_FFT_R2C.H:737
void setLocalDomain(std::array< int, 3 > const &local_start, std::array< int, 3 > const &local_size)
Set local domain.
Definition AMReX_FFT_R2C.H:788
R2C(Box const &domain, Info const &info=Info{})
Constructor.
Definition AMReX_FFT_R2C.H:489
void backward(cMF const &inmf, MF &outmf, int incomp=0, int outcomp=0)
Backward transform.
Definition AMReX_FFT_R2C.H:1353
void backward(MF &outmf, int outcomp=0)
Backward transform.
Definition AMReX_FFT_R2C.H:1032
void post_forward_doit_1(F const &post_forward)
CUDA-visible helper that redistributes and applies post_forward for the batched layout.
Definition AMReX_FFT_R2C.H:1208
T scalingFactor() const
Scaling factor. If the data goes through forward and then backward, the result multiplied by the scal...
Definition AMReX_FFT_R2C.H:1250
void post_forward_doit_0(F const &post_forward)
CUDA-visible hook that walks internal spectral data and applies post_forward.
Definition AMReX_FFT_R2C.H:1153
std::pair< std::array< int, 3 >, std::array< int, 3 > > getLocalDomain() const
Get local domain.
Definition AMReX_FFT_R2C.H:797
void forward(RT const *in, CT *out)
Forward transform.
Definition AMReX_FFT_R2C.H:1012
std::pair< BoxArray, DistributionMapping > getSpectralDataLayout() const
Get BoxArray and DistributionMapping for spectral data.
Definition AMReX_FFT_R2C.H:1422
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:732
FabArray< BaseFab< GpuComplex< T > > > cMF
Definition AMReX_FFT_R2C.H:49
std::pair< std::array< int, 3 >, std::array< int, 3 > > getLocalSpectralDomain() const
Get local spectral domain.
Definition AMReX_FFT_R2C.H:823
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:202
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:814
void backward(CT const *in, RT *out)
Backward transform.
Definition AMReX_FFT_R2C.H:1103
void forward(MF const &inmf, int incomp=0)
Forward transform.
Definition AMReX_FFT_R2C.H:899
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:217
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:2358
void ParallelCopy(const FabArray< FAB > &src, const Periodicity &period=Periodicity::NonPeriodic(), CpOp op=FabArrayBase::COPY)
Definition AMReX_FabArray.H:850
typename std::conditional_t< IsBaseFab< BaseFab< GpuComplex< T > > >::value, BaseFab< GpuComplex< T > >, FABType >::value_type value_type
Definition AMReX_FabArray.H:361
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:2174
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:1955
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:319
Arena * The_Arena()
Definition AMReX_Arena.cpp:805
int NProcs() noexcept
Definition AMReX_ParallelDescriptor.H:255
Definition AMReX_FFT_Helper.H:52
Direction
Definition AMReX_FFT_Helper.H:54
void dtod_memcpy_async(void *p_d_dst, const void *p_d_src, const std::size_t sz) noexcept
Definition AMReX_GpuDevice.H:449
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:310
MPI_Comm CommunicatorSub() noexcept
sub-communicator for current frame
Definition AMReX_ParallelContext.H:70
int MyProcSub() noexcept
my sub-rank in current frame
Definition AMReX_ParallelContext.H:76
int local_to_global_rank(int rank) noexcept
translate between local rank and global rank
Definition AMReX_ParallelContext.H:95
int NProcsSub() noexcept
number of ranks in current frame
Definition AMReX_ParallelContext.H:74
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30
bool is_aligned(const void *p, std::size_t alignment) noexcept
Return whether the address p is aligned to alignment bytes.
Definition AMReX_Arena.H:39
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
Definition AMReX_FFT_Helper.H:84
int batch_size
Batched FFT size. Only support in R2C, not R2X.
Definition AMReX_FFT_Helper.H:87
DomainStrategy domain_strategy
Domain composition strategy.
Definition AMReX_FFT_Helper.H:66
int nprocs
Max number of processes to use.
Definition AMReX_FFT_Helper.H:90
int pencil_threshold
Definition AMReX_FFT_Helper.H:70
Definition AMReX_FFT_Helper.H:194
std::conditional_t< std::is_same_v< float, T >, cuComplex, cuDoubleComplex > VendorComplex
Definition AMReX_FFT_Helper.H:198
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