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;
38template <
typename T = Real, FFT::Direction D = FFT::Direction::both,
bool C = false>
43 using MF = std::conditional_t
44 <
C,
cMF, std::conditional_t<std::is_same_v<T,Real>,
48 template <
typename U>
friend class Poisson;
68 explicit R2C (std::array<int,AMREX_SPACEDIM>
const& domain_size,
102 std::array<int,AMREX_SPACEDIM>
const& local_size);
115 std::pair<std::array<int,AMREX_SPACEDIM>,std::array<int,AMREX_SPACEDIM>>
147 std::array<int,AMREX_SPACEDIM>
const& local_size);
166 std::pair<std::array<int,AMREX_SPACEDIM>,std::array<int,AMREX_SPACEDIM>>
191 std::enable_if_t<DIR == Direction::both, int> = 0>
193 int incomp = 0,
int outcomp = 0)
228 void forward (
MF const& inmf,
cMF& outmf,
int incomp = 0,
int outcomp = 0);
243 template <
typename RT,
typename CT,
Direction DIR=D,
bool CP=
C,
246 && ((
sizeof(RT)*2 ==
sizeof(CT) && !CP) ||
247 (
sizeof(RT) ==
sizeof(CT) && CP)),
int> = 0>
259 template <Direction DIR=D, std::enable_if_t<DIR == Direction::both,
int> = 0>
275 void backward (
cMF const& inmf,
MF& outmf,
int incomp = 0,
int outcomp = 0);
290 template <
typename CT,
typename RT,
Direction DIR=D,
bool CP=
C,
293 && ((
sizeof(RT)*2 ==
sizeof(CT) && !CP) ||
294 (
sizeof(RT) ==
sizeof(CT) && CP)),
int> = 0>
325 template <
typename F>
328 template <
typename F>
342 int incomp = 0,
int outcomp = 0);
391 static std::pair<BoxArray,DistributionMapping>
393 std::array<int,AMREX_SPACEDIM>
const& local_size);
395 template <
typename FA,
typename RT>
396 std::pair<std::unique_ptr<char,DataDeleter>,std::size_t>
452template <
typename T, Direction D,
bool C>
454 : m_real_domain(domain),
455 m_spectral_domain_x(make_domain_x(domain)),
456#if (AMREX_SPACEDIM >= 2)
457 m_spectral_domain_y(make_domain_y(domain)),
458#if (AMREX_SPACEDIM == 3)
459 m_spectral_domain_z(make_domain_z(domain)),
462 m_sub_helper(domain),
467 static_assert(std::is_same_v<float,T> || std::is_same_v<double,T>);
470#if (AMREX_SPACEDIM == 2)
475 int(domain.
length(1) > 1) +
476 int(domain.
length(2) > 1)) >= 2);
491#if (AMREX_SPACEDIM == 3)
521 m_rx.define(bax, dmx, ncomp, 0,
MFInfo().SetAlloc(
false));
525 for (
auto &
b : bl) {
541#if (AMREX_SPACEDIM >= 2)
547 if (cbay.size() == dmx.size()) {
556#if (AMREX_SPACEDIM == 3)
561 {
false,
true,
true},
true);
563 if (cbaz.size() == dmx.size()) {
565 }
else if (cbaz.size() == cdmy.
size()) {
602#if (AMREX_SPACEDIM >= 2)
605 m_cmd_x2y = std::make_unique<MultiBlockCommMetaData>
607 m_cmd_y2x = std::make_unique<MultiBlockCommMetaData>
611#if (AMREX_SPACEDIM == 3)
615 m_cmd_x2z = std::make_unique<MultiBlockCommMetaData>
617 m_cmd_z2x = std::make_unique<MultiBlockCommMetaData>
621 m_cmd_y2z = std::make_unique<MultiBlockCommMetaData>
623 m_cmd_z2y = std::make_unique<MultiBlockCommMetaData>
633 if (myproc <
m_rx.size())
639 Box const& box =
m_rx.box(myproc);
640 auto* pr =
m_rx[myproc].dataPtr();
656#if (AMREX_SPACEDIM >= 2)
661#if (AMREX_SPACEDIM == 3)
677 auto* pr = (
void*)
m_rx[0].dataPtr();
678 auto* pc = (
void*)
m_cx[0].dataPtr();
680 m_fft_fwd_x.template init_r2c<Direction::forward>(len, pr, pc,
false, ncomp);
684 m_fft_fwd_x.template init_r2c<Direction::forward>(len, pr, pc,
false, ncomp);
687 m_fft_bwd_x.template init_r2c<Direction::backward>(len, pr, pc,
false, ncomp);
694template <
typename T, Direction D,
bool C>
699template <
typename T, Direction D,
bool C>
702 if (m_fft_bwd_x.plan != m_fft_fwd_x.plan) {
703 m_fft_bwd_x.destroy();
705 if (m_fft_bwd_y.plan != m_fft_fwd_y.plan) {
706 m_fft_bwd_y.destroy();
708 if (m_fft_bwd_z.plan != m_fft_fwd_z.plan) {
709 m_fft_bwd_z.destroy();
711 m_fft_fwd_x.destroy();
712 m_fft_fwd_y.destroy();
713 m_fft_fwd_z.destroy();
714 if (m_fft_bwd_x_half.plan != m_fft_fwd_x_half.plan) {
715 m_fft_bwd_x_half.destroy();
717 m_fft_fwd_x_half.destroy();
720template <
typename T, Direction D,
bool C>
721std::pair<BoxArray,DistributionMapping>
723 std::array<int,AMREX_SPACEDIM>
const& local_size)
727 Box bx(lo, lo+len-1);
734 pmap.reserve(allboxes.
size());
735 for (
int i = 0; i < allboxes.
size(); ++i) {
736 if (allboxes[i].ok()) {
740 allboxes.erase(std::remove_if(allboxes.begin(), allboxes.end(),
741 [=] (
Box const&
b) { return b.isEmpty(); }),
743 BoxList bl(std::move(allboxes));
750template <
typename T, Direction D,
bool C>
752 std::array<int,AMREX_SPACEDIM>
const& local_size)
754 auto const& [ba, dm] = make_layout_from_local_domain(local_start, local_size);
755 m_raw_mf =
MF(ba, dm, m_rx.nComp(), 0,
MFInfo().SetAlloc(
false));
758template <
typename T, Direction D,
bool C>
759std::pair<std::array<int,AMREX_SPACEDIM>,std::array<int,AMREX_SPACEDIM>>
762 m_raw_mf =
MF(m_rx.boxArray(), m_rx.DistributionMap(), m_rx.nComp(), 0,
766 if (myproc < m_rx.size()) {
767 Box const& box = m_rx.box(myproc);
768 return std::make_pair(box.
smallEnd().toArray(),
771 return std::make_pair(std::array<int,AMREX_SPACEDIM>{
AMREX_D_DECL(0,0,0)},
776template <
typename T, Direction D,
bool C>
778 std::array<int,AMREX_SPACEDIM>
const& local_size)
780 auto const& [ba, dm] = make_layout_from_local_domain(local_start, local_size);
781 m_raw_cmf =
cMF(ba, dm, m_rx.nComp(), 0,
MFInfo().SetAlloc(
false));
784template <
typename T, Direction D,
bool C>
785std::pair<std::array<int,AMREX_SPACEDIM>,std::array<int,AMREX_SPACEDIM>>
788 auto const ncomp = m_info.batch_size;
789 auto const& [ba, dm] = getSpectralDataLayout();
794 if (myproc < m_raw_cmf.size()) {
795 Box const& box = m_raw_cmf.box(myproc);
796 return std::make_pair(box.
smallEnd().toArray(), box.
length().toArray());
798 return std::make_pair(std::array<int,AMREX_SPACEDIM>{
AMREX_D_DECL(0,0,0)},
803template <
typename T, Direction D,
bool C>
806 if (
C || m_r2c_sub) {
amrex::Abort(
"R2C: OpenBC not supported with reduced dimensions or complex inputs"); }
808#if (AMREX_SPACEDIM == 3)
809 if (m_do_alld_fft) {
return; }
811 auto const ncomp = m_info.batch_size;
813 if (m_slab_decomp && ! m_fft_fwd_x_half.defined) {
816 Box bottom_half = m_real_domain;
817 bottom_half.
growHi(2,-m_real_domain.length(2)/2);
818 Box box = fab->box() & bottom_half;
820 auto* pr = fab->dataPtr();
824 m_fft_fwd_x_half.template init_r2c<Direction::forward>
825 (box, pr, pc, m_slab_decomp, ncomp);
826 m_fft_bwd_x_half = m_fft_fwd_x_half;
829 m_fft_fwd_x_half.template init_r2c<Direction::forward>
830 (box, pr, pc, m_slab_decomp, ncomp);
833 m_fft_bwd_x_half.template init_r2c<Direction::backward>
834 (box, pr, pc, m_slab_decomp, ncomp);
841 if (m_cmd_x2z && ! m_cmd_x2z_half) {
842 Box bottom_half = m_spectral_domain_z;
845 bottom_half.
growHi(0,-m_spectral_domain_z.length(0)/2);
846 m_cmd_x2z_half = std::make_unique<MultiBlockCommMetaData>
847 (m_cz, bottom_half, m_cx,
IntVect(0), m_dtos_x2z);
850 if (m_cmd_z2x && ! m_cmd_z2x_half) {
851 Box bottom_half = m_spectral_domain_x;
852 bottom_half.
growHi(2,-m_spectral_domain_x.length(2)/2);
853 m_cmd_z2x_half = std::make_unique<MultiBlockCommMetaData>
854 (m_cx, bottom_half, m_cz,
IntVect(0), m_dtos_z2x);
859template <
typename T, Direction D,
bool C>
866 auto const ncomp = m_info.batch_size;
869 if (m_sub_helper.ghost_safe(inmf.nGrowVect())) {
870 m_r2c_sub->forward(m_sub_helper.make_alias_mf(inmf), incomp);
872 MF tmp(inmf.boxArray(), inmf.DistributionMap(), ncomp, 0);
873 tmp.LocalCopy(inmf, incomp, 0, ncomp,
IntVect(0));
874 m_r2c_sub->forward(m_sub_helper.make_alias_mf(tmp),0);
879 if (&m_rx != &inmf) {
880 m_rx.ParallelCopy(inmf, incomp, 0, ncomp);
885 m_fft_fwd_x.template compute_c2c<Direction::forward>();
887 m_fft_fwd_x.template compute_r2c<Direction::forward>();
892 auto& fft_x = m_openbc_half ? m_fft_fwd_x_half : m_fft_fwd_x;
894 fft_x.template compute_c2c<Direction::forward>();
896 fft_x.template compute_r2c<Direction::forward>();
900 ParallelCopy(m_cy, m_cx, *m_cmd_x2y, 0, 0, ncomp, m_dtos_x2y);
902 m_fft_fwd_y.template compute_c2c<Direction::forward>();
905 ParallelCopy(m_cz, m_cy, *m_cmd_y2z, 0, 0, ncomp, m_dtos_y2z);
907#if (AMREX_SPACEDIM == 3)
908 else if ( m_cmd_x2z) {
913 {components, m_dtos_x2z};
914 auto handler = ParallelCopy_nowait(m_cz, m_cx, *m_cmd_x2z_half, packing);
916 Box upper_half = m_spectral_domain_z;
919 upper_half.
growLo (0,-m_spectral_domain_z.length(0)/2);
920 m_cz.setVal(0, upper_half, 0, ncomp);
922 ParallelCopy_finish(m_cz, std::move(handler), *m_cmd_x2z_half, packing);
924 ParallelCopy(m_cz, m_cx, *m_cmd_x2z, 0, 0, ncomp, m_dtos_x2z);
928 m_fft_fwd_z.template compute_c2c<Direction::forward>();
931template <
typename T, Direction D,
bool C>
932template <
typename FA,
typename RT>
933std::pair<std::unique_ptr<char,DataDeleter>,std::size_t>
938 using FAB =
typename FA::FABType::value_type;
939 using T_FAB =
typename FAB::value_type;
940 static_assert(
sizeof(T_FAB) ==
sizeof(RT));
942 auto const ncomp = m_info.batch_size;
943 auto const& ia = fa.IndexArray();
948 if ( ! ia.empty() ) {
950 Box const& box = fa.fabbox(K);
954 sz =
sizeof(T_FAB) * box.
numPts() * ncomp;
957 fa.setFab(K, FAB(box,ncomp,
pp));
961 return std::make_pair(std::unique_ptr<char,DataDeleter>{},std::size_t(0));
963 return std::make_pair(std::unique_ptr<char,DataDeleter>
969template <
typename T, Direction D,
bool C>
970template <
typename RT,
typename CT,
Direction DIR,
bool CP,
973 && ((
sizeof(RT)*2 ==
sizeof(CT) && !CP) ||
974 (
sizeof(RT) ==
sizeof(CT) && CP)),
int> >
977 auto [rdata, rsz] = install_raw_ptr(m_raw_mf, in);
978 auto [cdata, csz] = install_raw_ptr(m_raw_cmf, out);
993template <
typename T, Direction D,
bool C>
994template <Direction DIR, std::enable_if_t<DIR == Direction::both,
int> >
1000template <
typename T, Direction D,
bool C>
1006 auto const ncomp = m_info.batch_size;
1009 if (m_sub_helper.ghost_safe(outmf.nGrowVect())) {
1010 MF submf = m_sub_helper.make_alias_mf(outmf);
1011 IntVect const& subngout = m_sub_helper.make_iv(ngout);
1012 Periodicity const& subperiod = m_sub_helper.make_periodicity(period);
1013 m_r2c_sub->backward_doit(submf, subngout, subperiod, outcomp);
1015 MF tmp(outmf.boxArray(), outmf.DistributionMap(), ncomp,
1016 m_sub_helper.make_safe_ghost(outmf.nGrowVect()));
1017 this->backward_doit(tmp, ngout, period, 0);
1018 outmf.LocalCopy(tmp, 0, outcomp, ncomp, tmp.nGrowVect());
1023 if (m_do_alld_fft) {
1025 m_fft_bwd_x.template compute_c2c<Direction::backward>();
1027 m_fft_bwd_x.template compute_r2c<Direction::backward>();
1029 outmf.ParallelCopy(m_rx, 0, outcomp, ncomp,
IntVect(0),
1034 m_fft_bwd_z.template compute_c2c<Direction::backward>();
1036 ParallelCopy(m_cy, m_cz, *m_cmd_z2y, 0, 0, ncomp, m_dtos_z2y);
1038#if (AMREX_SPACEDIM == 3)
1039 else if ( m_cmd_z2x) {
1040 auto const& cmd = m_openbc_half ? m_cmd_z2x_half : m_cmd_z2x;
1041 ParallelCopy(m_cx, m_cz, *cmd, 0, 0, ncomp, m_dtos_z2x);
1045 m_fft_bwd_y.template compute_c2c<Direction::backward>();
1047 ParallelCopy(m_cx, m_cy, *m_cmd_y2x, 0, 0, ncomp, m_dtos_y2x);
1050 auto& fft_x = m_openbc_half ? m_fft_bwd_x_half : m_fft_bwd_x;
1052 fft_x.template compute_c2c<Direction::backward>();
1054 fft_x.template compute_r2c<Direction::backward>();
1056 outmf.ParallelCopy(m_rx, 0, outcomp, ncomp,
IntVect(0),
1060template <
typename T, Direction D,
bool C>
1061template <
typename CT,
typename RT,
Direction DIR,
bool CP,
1064 && ((
sizeof(RT)*2 ==
sizeof(CT) && !CP) ||
1065 (
sizeof(RT) ==
sizeof(CT) && CP)),
int> >
1068 auto [rdata, rsz] = install_raw_ptr(m_raw_mf, out);
1069 auto [cdata, csz] = install_raw_ptr(m_raw_cmf, in);
1084template <
typename T, Direction D,
bool C>
1092 if (!fab) {
return {fwd, bwd};}
1094 Box const& box = fab->box();
1097 auto const ncomp = m_info.batch_size;
1099#ifdef AMREX_USE_SYCL
1100 fwd.template init_c2c<Direction::forward>(box, pio, ncomp, ndims);
1104 fwd.template init_c2c<Direction::forward>(box, pio, ncomp, ndims);
1107 bwd.template init_c2c<Direction::backward>(box, pio, ncomp, ndims);
1114template <
typename T, Direction D,
bool C>
1115template <
typename F>
1118 if (m_info.twod_mode || m_info.batch_size > 1) {
1120#if (AMREX_SPACEDIM > 1)
1121 }
else if (m_r2c_sub) {
1123#if (AMREX_SPACEDIM == 2)
1125 m_r2c_sub->post_forward_doit_1
1128 post_forward(0, i, 0, sp);
1131 if (m_real_domain.length(0) == 1 && m_real_domain.length(1) == 1) {
1133 m_r2c_sub->post_forward_doit_1
1136 post_forward(0, 0, i, sp);
1138 }
else if (m_real_domain.length(0) == 1 && m_real_domain.length(2) == 1) {
1140 m_r2c_sub->post_forward_doit_1
1143 post_forward(0, i, 0, sp);
1145 }
else if (m_real_domain.length(0) == 1) {
1147 m_r2c_sub->post_forward_doit_1
1150 post_forward(0, i, j, sp);
1152 }
else if (m_real_domain.length(1) == 1) {
1154 m_r2c_sub->post_forward_doit_1
1157 post_forward(i, 0, j, sp);
1160 amrex::Abort(
"R2c::post_forward_doit_0: how did this happen?");
1165 this->post_forward_doit_1(post_forward);
1169template <
typename T, Direction D,
bool C>
1170template <
typename F>
1173 if (m_info.twod_mode || m_info.batch_size > 1) {
1175 }
else if (m_r2c_sub) {
1176 amrex::Abort(
"R2C::post_forward_doit_1: How did this happen?");
1178 if ( ! m_cz.empty()) {
1181 auto const& a = spectral_fab->array();
1185 post_forward(jx,ky,iz,a(iz,jx,ky));
1188 }
else if ( ! m_cy.empty()) {
1191 auto const& a = spectral_fab->array();
1195 post_forward(jx,iy,k,a(iy,jx,k));
1201 auto const& a = spectral_fab->array();
1205 post_forward(i,j,k,a(i,j,k));
1212template <
typename T, Direction D,
bool C>
1215#if (AMREX_SPACEDIM == 3)
1216 if (m_info.twod_mode) {
1217 return T(1)/T(Long(m_real_domain.length(0)) *
1218 Long(m_real_domain.length(1)));
1222 return T(1)/T(m_real_domain.numPts());
1226template <
typename T, Direction D,
bool C>
1229std::pair<typename R2C<T,D,C>::cMF *,
IntVect>
1232#if (AMREX_SPACEDIM > 1)
1234 auto [cmf, order] = m_r2c_sub->getSpectralData();
1235 return std::make_pair(cmf, m_sub_helper.inverse_order(order));
1238 if (!m_cz.empty()) {
1240 }
else if (!m_cy.empty()) {
1247template <
typename T, Direction D,
bool C>
1254 auto const ncomp = m_info.batch_size;
1258 bool inmf_safe = m_sub_helper.ghost_safe(inmf.nGrowVect());
1259 MF inmf_sub, inmf_tmp;
1262 inmf_sub = m_sub_helper.make_alias_mf(inmf);
1263 incomp_sub = incomp;
1265 inmf_tmp.define(inmf.boxArray(), inmf.DistributionMap(), ncomp, 0);
1266 inmf_tmp.LocalCopy(inmf, incomp, 0, ncomp,
IntVect(0));
1267 inmf_sub = m_sub_helper.make_alias_mf(inmf_tmp);
1271 bool outmf_safe = m_sub_helper.ghost_safe(outmf.
nGrowVect());
1272 cMF outmf_sub, outmf_tmp;
1275 outmf_sub = m_sub_helper.make_alias_mf(outmf);
1276 outcomp_sub = outcomp;
1279 outmf_sub = m_sub_helper.make_alias_mf(outmf_tmp);
1283 m_r2c_sub->forward(inmf_sub, outmf_sub, incomp_sub, outcomp_sub);
1292 if (!m_cz.empty()) {
1295 (outmf, m_spectral_domain_x, m_cz,
IntVect(0), dtos);
1296 ParallelCopy(outmf, m_cz, cmd, 0, outcomp, ncomp, dtos);
1297 }
else if (!m_cy.empty()) {
1299 (outmf, m_spectral_domain_x, m_cy,
IntVect(0), m_dtos_y2x);
1300 ParallelCopy(outmf, m_cy, cmd, 0, outcomp, ncomp, m_dtos_y2x);
1307template <
typename T, Direction D,
bool C>
1315template <
typename T, Direction D,
bool C>
1317 Periodicity const& period,
int incomp,
int outcomp)
1321 auto const ncomp = m_info.batch_size;
1325 bool inmf_safe = m_sub_helper.ghost_safe(inmf.
nGrowVect());
1326 cMF inmf_sub, inmf_tmp;
1329 inmf_sub = m_sub_helper.make_alias_mf(inmf);
1330 incomp_sub = incomp;
1334 inmf_sub = m_sub_helper.make_alias_mf(inmf_tmp);
1338 bool outmf_safe = m_sub_helper.ghost_safe(outmf.nGrowVect());
1339 MF outmf_sub, outmf_tmp;
1342 outmf_sub = m_sub_helper.make_alias_mf(outmf);
1343 outcomp_sub = outcomp;
1345 IntVect const& ngtmp = m_sub_helper.make_safe_ghost(outmf.nGrowVect());
1346 outmf_tmp.define(outmf.boxArray(), outmf.DistributionMap(), ncomp, ngtmp);
1347 outmf_sub = m_sub_helper.make_alias_mf(outmf_tmp);
1351 IntVect const& subngout = m_sub_helper.make_iv(ngout);
1352 Periodicity const& subperiod = m_sub_helper.make_periodicity(period);
1353 m_r2c_sub->backward_doit(inmf_sub, outmf_sub, subngout, subperiod, incomp_sub, outcomp_sub);
1356 outmf.LocalCopy(outmf_tmp, 0, outcomp, ncomp, outmf_tmp.nGrowVect());
1361 if (!m_cz.empty()) {
1364 (m_cz, m_spectral_domain_z, inmf,
IntVect(0), dtos);
1366 }
else if (!m_cy.empty()) {
1368 (m_cy, m_spectral_domain_y, inmf,
IntVect(0), m_dtos_x2y);
1369 ParallelCopy(m_cy, inmf, cmd, incomp, 0, ncomp, m_dtos_x2y);
1371 m_cx.ParallelCopy(inmf, incomp, 0, ncomp);
1373 backward_doit(outmf, ngout, period, outcomp);
1377template <
typename T, Direction D,
bool C>
1378std::pair<BoxArray,DistributionMapping>
1381#if (AMREX_SPACEDIM > 1)
1383 auto const& [ba, dm] = m_r2c_sub->getSpectralDataLayout();
1384 return std::make_pair(m_sub_helper.inverse_boxarray(ba), dm);
1388#if (AMREX_SPACEDIM == 3)
1389 if (!m_cz.empty()) {
1390 BoxList bl = m_cz.boxArray().boxList();
1391 for (
auto&
b : bl) {
1392 auto lo =
b.smallEnd();
1393 auto hi =
b.bigEnd();
1394 std::swap(lo[0], lo[1]);
1395 std::swap(lo[1], lo[2]);
1396 std::swap(hi[0], hi[1]);
1397 std::swap(hi[1], hi[2]);
1401 return std::make_pair(
BoxArray(std::move(bl)), m_cz.DistributionMap());
1404#if (AMREX_SPACEDIM >= 2)
1405 if (!m_cy.empty()) {
1406 BoxList bl = m_cy.boxArray().boxList();
1407 for (
auto&
b : bl) {
1408 auto lo =
b.smallEnd();
1409 auto hi =
b.bigEnd();
1410 std::swap(lo[0], lo[1]);
1411 std::swap(hi[0], hi[1]);
1415 return std::make_pair(
BoxArray(std::move(bl)), m_cy.DistributionMap());
1419 return std::make_pair(m_cx.boxArray(), m_cx.DistributionMap());
1424template <
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:551
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
Get the bigend.
Definition AMReX_Box.H:119
__host__ __device__ Long numPts() const noexcept
Returns the number of points contained in the BoxND.
Definition AMReX_Box.H:349
__host__ __device__ IntVectND< dim > length() const noexcept
Return the length of the BoxND.
Definition AMReX_Box.H:149
__host__ __device__ int shortside(int &dir) const noexcept
Returns length of shortest side. dir is modified to give direction with shortest side: 0....
Definition AMReX_Box.H:430
__host__ __device__ IntVectND< dim > size() const noexcept
Return the length of the BoxND.
Definition AMReX_Box.H:142
__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:651
__host__ __device__ IndexTypeND< dim > ixType() const noexcept
Returns the indexing type.
Definition AMReX_Box.H:130
__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:662
__host__ __device__ bool ok() const noexcept
Checks if it is a proper BoxND (including a valid type).
Definition AMReX_Box.H:203
__host__ __device__ const IntVectND< dim > & smallEnd() const &noexcept
Get the smallend of the BoxND.
Definition AMReX_Box.H:108
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:41
Long size() const noexcept
Length of the underlying processor map.
Definition AMReX_DistributionMapping.H:127
Definition AMReX_FFT_OpenBCSolver.H:11
3D Poisson solver for periodic, Dirichlet & Neumann boundaries in the first two dimensions,...
Definition AMReX_FFT_Poisson.H:146
Poisson solver for periodic, Dirichlet & Neumann boundaries using FFT.
Definition AMReX_FFT_Poisson.H:58
Parallel Discrete Fourier Transform.
Definition AMReX_FFT_R2C.H:40
std::unique_ptr< MultiBlockCommMetaData > m_cmd_x2z_half
Definition AMReX_FFT_R2C.H:417
Info m_info
Definition AMReX_FFT_R2C.H:445
std::conditional_t< C, cMF, std::conditional_t< std::is_same_v< T, Real >, MultiFab, FabArray< BaseFab< T > > > > MF
Definition AMReX_FFT_R2C.H:45
R2C & operator=(R2C const &)=delete
void backward_doit(cMF const &inmf, MF &outmf, IntVect const &ngout=IntVect(0), Periodicity const &period=Periodicity::NonPeriodic(), int incomp=0, int outcomp=0)
Definition AMReX_FFT_R2C.H:1316
void backward_doit(MF &outmf, IntVect const &ngout=IntVect(0), Periodicity const &period=Periodicity::NonPeriodic(), int outcomp=0)
Definition AMReX_FFT_R2C.H:1001
Box m_spectral_domain_y
Definition AMReX_FFT_R2C.H:439
cMF m_cy
Definition AMReX_FFT_R2C.H:428
Plan< T > m_fft_fwd_y
Definition AMReX_FFT_R2C.H:401
MF m_raw_mf
Definition AMReX_FFT_R2C.H:431
Plan< T > m_fft_fwd_z
Definition AMReX_FFT_R2C.H:403
void forward(MF const &inmf, cMF &outmf, int incomp=0, int outcomp=0)
Forward transform.
Definition AMReX_FFT_R2C.H:1250
~R2C()
Definition AMReX_FFT_R2C.H:700
std::unique_ptr< MultiBlockCommMetaData > m_cmd_z2x_half
Definition AMReX_FFT_R2C.H:418
void setLocalDomain(std::array< int, 3 > const &local_start, std::array< int, 3 > const &local_size)
Set local domain.
Definition AMReX_FFT_R2C.H:751
std::unique_ptr< MultiBlockCommMetaData > m_cmd_x2y
Definition AMReX_FFT_R2C.H:411
static Box make_domain_x(Box const &domain)
Definition AMReX_FFT_R2C.H:346
Box m_spectral_domain_x
Definition AMReX_FFT_R2C.H:438
static std::pair< BoxArray, DistributionMapping > make_layout_from_local_domain(std::array< int, 3 > const &local_start, std::array< int, 3 > const &local_size)
Definition AMReX_FFT_R2C.H:722
R2C(Box const &domain, Info const &info=Info{})
Constructor.
Definition AMReX_FFT_R2C.H:453
cMF m_cx
Definition AMReX_FFT_R2C.H:427
Plan< T > m_fft_bwd_x
Definition AMReX_FFT_R2C.H:400
std::unique_ptr< char, DataDeleter > m_data_2
Definition AMReX_FFT_R2C.H:435
static Box make_domain_y(Box const &domain)
Definition AMReX_FFT_R2C.H:361
void backward(cMF const &inmf, MF &outmf, int incomp=0, int outcomp=0)
Backward transform.
Definition AMReX_FFT_R2C.H:1310
Swap02 m_dtos_y2z
Definition AMReX_FFT_R2C.H:421
void backward(MF &outmf, int outcomp=0)
Backward transform.
Definition AMReX_FFT_R2C.H:995
bool m_do_alld_fft
Definition AMReX_FFT_R2C.H:447
void post_forward_doit_1(F const &post_forward)
Definition AMReX_FFT_R2C.H:1171
Swap01 m_dtos_x2y
Definition AMReX_FFT_R2C.H:419
std::unique_ptr< MultiBlockCommMetaData > m_cmd_x2z
Definition AMReX_FFT_R2C.H:415
RotateFwd m_dtos_x2z
Definition AMReX_FFT_R2C.H:423
std::unique_ptr< MultiBlockCommMetaData > m_cmd_y2z
Definition AMReX_FFT_R2C.H:413
Plan< T > m_fft_bwd_x_half
Definition AMReX_FFT_R2C.H:406
T scalingFactor() const
Definition AMReX_FFT_R2C.H:1213
Box m_spectral_domain_z
Definition AMReX_FFT_R2C.H:440
cMF m_cz
Definition AMReX_FFT_R2C.H:429
std::unique_ptr< MultiBlockCommMetaData > m_cmd_z2y
Definition AMReX_FFT_R2C.H:414
std::unique_ptr< R2C< T, D, C > > m_r2c_sub
Definition AMReX_FFT_R2C.H:442
void post_forward_doit_0(F const &post_forward)
Definition AMReX_FFT_R2C.H:1116
detail::SubHelper m_sub_helper
Definition AMReX_FFT_R2C.H:443
Swap01 m_dtos_y2x
Definition AMReX_FFT_R2C.H:420
std::pair< std::array< int, 3 >, std::array< int, 3 > > getLocalDomain() const
Get local domain.
Definition AMReX_FFT_R2C.H:760
std::pair< Plan< T >, Plan< T > > make_c2c_plans(cMF &inout, int ndims) const
Definition AMReX_FFT_R2C.H:1086
std::pair< std::unique_ptr< char, DataDeleter >, std::size_t > install_raw_ptr(FA &fa, RT const *p)
Definition AMReX_FFT_R2C.H:934
void forward(RT const *in, CT *out)
Forward transform.
Definition AMReX_FFT_R2C.H:975
Plan< T > m_fft_bwd_z
Definition AMReX_FFT_R2C.H:404
Plan< T > m_fft_fwd_x_half
Definition AMReX_FFT_R2C.H:405
std::unique_ptr< MultiBlockCommMetaData > m_cmd_z2x
Definition AMReX_FFT_R2C.H:416
bool m_openbc_half
Definition AMReX_FFT_R2C.H:449
std::pair< BoxArray, DistributionMapping > getSpectralDataLayout() const
Get BoxArray and DistributionMapping for spectral data.
Definition AMReX_FFT_R2C.H:1379
void prepare_openbc()
Definition AMReX_FFT_R2C.H:804
static Box make_domain_z(Box const &domain)
Definition AMReX_FFT_R2C.H:376
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:695
FabArray< BaseFab< GpuComplex< T > > > cMF
Definition AMReX_FFT_R2C.H:42
std::pair< std::array< int, 3 >, std::array< int, 3 > > getLocalSpectralDomain() const
Get local spectral domain.
Definition AMReX_FFT_R2C.H:786
std::unique_ptr< MultiBlockCommMetaData > m_cmd_y2x
Definition AMReX_FFT_R2C.H:412
Box m_real_domain
Definition AMReX_FFT_R2C.H:437
bool m_slab_decomp
Definition AMReX_FFT_R2C.H:448
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:192
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:777
MF m_rx
Definition AMReX_FFT_R2C.H:426
cMF m_raw_cmf
Definition AMReX_FFT_R2C.H:432
Plan< T > m_fft_bwd_y
Definition AMReX_FFT_R2C.H:402
void backward(CT const *in, RT *out)
Backward transform.
Definition AMReX_FFT_R2C.H:1066
std::unique_ptr< char, DataDeleter > m_data_1
Definition AMReX_FFT_R2C.H:434
Plan< T > m_fft_fwd_x
Definition AMReX_FFT_R2C.H:399
RotateBwd m_dtos_z2x
Definition AMReX_FFT_R2C.H:424
void forward(MF const &inmf, int incomp=0)
Forward transform.
Definition AMReX_FFT_R2C.H:862
Swap02 m_dtos_z2y
Definition AMReX_FFT_R2C.H:422
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:2312
void ParallelCopy(const FabArray< FAB > &src, const Periodicity &period=Periodicity::NonPeriodic(), CpOp op=FabArrayBase::COPY)
Definition AMReX_FabArray.H:845
typename std::conditional_t< IsBaseFab< BaseFab< GpuComplex< T > > >::value, BaseFab< GpuComplex< T > >, FABType >::value_type value_type
Definition AMReX_FabArray.H:356
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:2128
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:1909
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:38
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
Long size() const noexcept
Definition AMReX_Vector.H:53
std::unique_ptr< char, DataDeleter > make_mfs_share(FA1 &fa1, FA2 &fa2)
Definition AMReX_FFT_Helper.H:1387
FA::FABType::value_type * get_fab(FA &fa)
Definition AMReX_FFT_Helper.H:1376
DistributionMapping make_iota_distromap(Long n)
Definition AMReX_FFT.cpp:88
Definition AMReX_FFT.cpp:7
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:317
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:260
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:210
int MyProc() noexcept
return the rank number local to the current Parallel Context
Definition AMReX_ParallelDescriptor.H:125
int NProcs() noexcept
return the number of MPI ranks local to the current Parallel Context
Definition AMReX_ParallelDescriptor.H:243
void ParallelForOMP(T n, L const &f) noexcept
Definition AMReX_GpuLaunch.H:249
BoxND< 3 > Box
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:1931
IntVectND< 3 > IntVect
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:1967
Arena * The_Arena()
Definition AMReX_Arena.cpp:705
__host__ __device__ constexpr T elemwiseMin(T const &a, T const &b) noexcept
Definition AMReX_Algorithm.H:49
Definition AMReX_DataAllocator.H:29
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:130
std::conditional_t< std::is_same_v< float, T >, cuComplex, cuDoubleComplex > VendorComplex
Definition AMReX_FFT_Helper.H:134
Definition AMReX_FFT_Helper.H:1499
Definition AMReX_FFT_Helper.H:1474
Definition AMReX_FFT_Helper.H:1428
Definition AMReX_FFT_Helper.H:1451
Definition AMReX_FFT_Helper.H:1526
Box make_box(Box const &box) const
Definition AMReX_FFT.cpp:142
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