1#ifndef AMREX_FFT_POISSON_H_
2#define AMREX_FFT_POISSON_H_
13void fill_physbc (MF& mf, Geometry
const& geom,
14 Array<std::pair<Boundary,Boundary>,AMREX_SPACEDIM>
const& bc);
16inline Geometry shift_geom (Geometry
const& geom)
18 auto const& dom = geom.Domain();
19 auto const& dlo = dom.smallEnd();
23 return Geometry(
amrex::shift(dom,-dlo), geom.ProbDomain(), geom.CoordInt(),
29void shift_mfs (IntVect
const& domain_lo, MF
const& mf1, MF
const& mf2,
30 MF& new_mf1, MF& new_mf2)
33 mf1.DistributionMap() == mf2.DistributionMap());
34 BoxArray ba = mf1.boxArray();
36 new_mf1.define(ba, mf1.DistributionMap(), 1, mf1.nGrowVect(),
37 MFInfo().SetAlloc(
false));
38 new_mf2.define(ba, mf2.DistributionMap(), 1, mf2.nGrowVect(),
39 MFInfo().SetAlloc(
false));
40 for (MFIter mfi(mf1); mfi.isValid(); ++mfi) {
42 fab1.shift(-domain_lo);
43 new_mf1.setFab(mfi, std::move(fab1));
46 fab2.shift(-domain_lo);
47 new_mf2.setFab(mfi, std::move(fab2));
59template <
typename MF = MultiFab>
64 template <
typename FA=MF, std::enable_if_t<IsFabArray_v<FA>,
int> = 0>
66 Array<std::pair<Boundary,Boundary>,AMREX_SPACEDIM>
const& bc)
67 : m_domain_lo(geom.Domain().smallEnd()),
68 m_geom(detail::shift_geom(geom)),
71 bool all_periodic =
true;
72 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
73 all_periodic = all_periodic
78 m_r2c = std::make_unique<R2C<typename MF::value_type>>(m_geom.
Domain());
80 m_r2x = std::make_unique<R2X<typename MF::value_type>> (m_geom.
Domain(), m_bc);
84 template <
typename FA=MF, std::enable_if_t<IsFabArray_v<FA>,
int> = 0>
86 : m_domain_lo(geom.Domain().smallEnd()),
87 m_geom(detail::shift_geom(geom)),
93 m_r2c = std::make_unique<R2C<typename MF::value_type>>(m_geom.
Domain());
106 void solve (MF& a_soln, MF
const& a_rhs);
112 std::unique_ptr<R2X<typename MF::value_type>> m_r2x;
113 std::unique_ptr<R2C<typename MF::value_type>> m_r2c;
116#if (AMREX_SPACEDIM == 3)
121template <
typename MF = MultiFab>
126 template <
typename FA=MF, std::enable_if_t<IsFabArray_v<FA>,
int> = 0>
131 void solve (MF& soln, MF
const& rhs);
149template <
typename MF = MultiFab>
153 using T =
typename MF::value_type;
155 template <
typename FA=MF, std::enable_if_t<IsFabArray_v<FA>,
int> = 0>
157 Array<std::pair<Boundary,Boundary>,AMREX_SPACEDIM>
const& bc)
158 : m_domain_lo(geom.Domain().smallEnd()),
159 m_geom(detail::shift_geom(geom)),
162#if (AMREX_SPACEDIM < 3)
166 bool periodic_xy =
true;
167 for (
int idim = 0; idim < 2; ++idim) {
184 m_r2c = std::make_unique<R2C<typename MF::value_type>>(m_geom.
Domain(),
187 m_r2x = std::make_unique<R2X<typename MF::value_type>> (m_geom.
Domain(),
200 void solve (MF& soln, MF
const& rhs);
204 void solve_2d (MF& a_soln, MF
const& a_rhs);
206 template <
typename TRIA,
typename TRIC>
207 void solve (MF& a_soln, MF
const& a_rhs, TRIA
const& tria, TRIC
const& tric);
210 template <
typename FA,
typename TRIA,
typename TRIC>
211 void solve_z (FA& spmf, TRIA
const& tria, TRIC
const& tric);
222 std::unique_ptr<R2X<typename MF::value_type>> m_r2x;
223 std::unique_ptr<R2C<typename MF::value_type>> m_r2c;
229template <
typename MF>
235 MF
const* rhs = &a_rhs;
237 if (m_domain_lo != 0) {
238 detail::shift_mfs(m_domain_lo, a_soln, a_rhs, solntmp, rhstmp);
243 AMREX_ASSERT(soln->is_cell_centered() && rhs->is_cell_centered());
245 using T =
typename MF::value_type;
248 {
AMREX_D_DECL(Math::pi<T>()/T(m_geom.Domain().length(0)),
249 Math::pi<T>()/T(m_geom.Domain().length(1)),
250 Math::pi<T>()/T(m_geom.Domain().length(2)))};
251 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
257 {
AMREX_D_DECL(T(2)/T(m_geom.CellSize(0)*m_geom.CellSize(0)),
258 T(2)/T(m_geom.CellSize(1)*m_geom.CellSize(1)),
259 T(2)/T(m_geom.CellSize(2)*m_geom.CellSize(2)))};
260 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
261 if (m_geom.Domain().length(idim) == 1) {
265 auto scale = (m_r2x) ? m_r2x->scalingFactor() : m_r2c->scalingFactor();
268 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
287 T b = fac[1]*(j+
offset[1]);,
288 T c = fac[2]*(k+
offset[2]));
290 +dxfac[1]*(std::cos(b)-T(1)),
291 +dxfac[2]*(std::cos(c)-T(1)));
295 spectral_data *=
scale;
301 m_r2x->forwardThenBackward_doit_0(*rhs, *soln, f, ng, m_geom.periodicity());
302 detail::fill_physbc(*soln, m_geom, m_bc);
304 m_r2c->forward(*rhs);
305 m_r2c->post_forward_doit_0(f);
306 m_r2c->backward_doit(*soln, ng, m_geom.periodicity());
310#if (AMREX_SPACEDIM == 3)
312template <
typename MF>
313template <
typename FA, std::enable_if_t<IsFabArray_v<FA>,
int> FOO>
319 m_solver(m_grown_domain)
324template <
typename MF>
327 using T =
typename MF::value_type;
328 auto const& lo = m_grown_domain.smallEnd();
329 auto const dx = T(m_geom.CellSize(0));
330 auto const dy = T(m_geom.CellSize(1));
331 auto const dz = T(m_geom.CellSize(2));
332 auto const gfac = T(1)/T(std::sqrt(T(12)));
334 auto const fac = T(-0.125) * (dx*dy*dz) / (T(4)*Math::pi<T>());
337 auto x = (T(i-lo[0]) - gfac) * dx;
338 auto y = (T(j-lo[1]) - gfac) * dy;
339 auto z = (T(k-lo[2]) - gfac) * dz;
341 for (
int gx = 0; gx < 2; ++gx) {
342 for (
int gy = 0; gy < 2; ++gy) {
343 for (
int gz = 0; gz < 2; ++gz) {
344 auto xg =
x + 2*gx*gfac*dx;
345 auto yg =
y + 2*gy*gfac*dy;
346 auto zg =
z + 2*gz*gfac*dz;
347 r += T(1)/std::sqrt(xg*xg+yg*yg+zg*zg);
353template <
typename MF>
356 AMREX_ASSERT(m_grown_domain.ixType() == soln.ixType() && m_grown_domain.ixType() == rhs.ixType());
357 m_solver.solve(soln, rhs);
363namespace fft_poisson_detail {
364 template <
typename T>
366 [[nodiscard]]
constexpr T operator() (
int,
int,
int)
const
372 template <
typename T>
375 T operator() (
int,
int,
int)
const
382 template <
typename T>
385 T operator() (
int,
int,
int k)
const
387 return (k > 0) ? T(2.0) / (m_dz[k]*(m_dz[k]+m_dz[k-1]))
388 : T(1.0) / (m_dz[k]* m_dz[k]);
393 template <
typename T>
396 T operator() (
int,
int,
int k)
const
398 return (k < m_size-1) ? T(2.0) / (m_dz[k]*(m_dz[k]+m_dz[k+1]))
399 : T(1.0) / (m_dz[k]* m_dz[k]);
407template <
typename MF>
408std::pair<BoxArray,DistributionMapping>
411 if (!m_spmf_r.empty()) {
412 return std::make_pair(m_spmf_r.boxArray(), m_spmf_r.DistributionMap());
414 return std::make_pair(m_spmf_c.boxArray(), m_spmf_c.DistributionMap());
418template <
typename MF>
421#if (AMREX_SPACEDIM == 3)
423 (m_geom.Domain().length(0) > 1 ||
424 m_geom.Domain().length(1) > 1));
427 Box cdomain = m_geom.Domain();
428 if (cdomain.
length(0) > 1) {
435 DistributionMapping dm = detail::make_iota_distromap(cba.size());
436 m_spmf_c.define(cba, dm, 1, 0);
437 }
else if (m_geom.Domain().length(0) > 1 &&
438 m_geom.Domain().length(1) > 1) {
439 if (m_r2x->m_cy.empty()) {
441 {AMREX_D_DECL(true,true,false)});
442 DistributionMapping dm = detail::make_iota_distromap(sba.size());
443 m_spmf_r.define(sba, dm, 1, 0);
445 Box cdomain = m_geom.Domain();
447 cdomain.
setBig(0,cdomain.length(0)/2);
449 cdomain.
setBig(1,cdomain.length(1)/2);
453 DistributionMapping dm = detail::make_iota_distromap(cba.size());
454 m_spmf_c.define(cba, dm, 1, 0);
459 {AMREX_D_DECL(true,true,false)});
460 DistributionMapping dm = detail::make_iota_distromap(sba.size());
461 m_spmf_r.define(sba, dm, 1, 0);
468template <
typename MF>
471 auto delz =
T(m_geom.CellSize(AMREX_SPACEDIM-1));
473 fft_poisson_detail::Tri_Uniform<T>{
T(1)/(delz*delz)},
474 fft_poisson_detail::Tri_Uniform<T>{
T(1)/(delz*delz)});
477template <
typename MF>
480 auto const* pdz = dz.
dataPtr();
482 fft_poisson_detail::TriA<T>{pdz},
483 fft_poisson_detail::TriC<T>{pdz,int(dz.
size())});
486template <
typename MF>
489 AMREX_ASSERT(soln.is_cell_centered() && rhs.is_cell_centered());
494 auto const* pdz = d_dz.
data();
496 auto const* pdz = dz.data();
499 fft_poisson_detail::TriA<T>{pdz},
500 fft_poisson_detail::TriC<T>{pdz,int(dz.
size())});
503template <
typename MF>
506 solve(soln, rhs, fft_poisson_detail::Tri_Zero<T>{}, fft_poisson_detail::Tri_Zero<T>{});
509template <
typename MF>
510template <
typename TRIA,
typename TRIC>
516 AMREX_ASSERT(a_soln.is_cell_centered() && a_rhs.is_cell_centered());
518#if (AMREX_SPACEDIM < 3)
523 MF
const* rhs = &a_rhs;
525 if (m_domain_lo != 0) {
526 detail::shift_mfs(m_domain_lo, a_soln, a_rhs, solntmp, rhstmp);
535 m_r2c->forward(*rhs, m_spmf_c);
536 solve_z(m_spmf_c, tria, tric);
537 m_r2c->backward_doit(m_spmf_c, *soln, ng, m_geom.periodicity());
541 if (m_r2x->m_cy.empty()) {
542 m_r2x->forward(*rhs, m_spmf_r);
543 solve_z(m_spmf_r, tria, tric);
544 m_r2x->backward(m_spmf_r, *soln, ng, m_geom.periodicity());
546 m_r2x->forward(*rhs, m_spmf_c);
547 solve_z(m_spmf_c, tria, tric);
548 m_r2x->backward(m_spmf_c, *soln, ng, m_geom.periodicity());
552 detail::fill_physbc(*soln, m_geom, m_bc);
556template <
typename MF>
557template <
typename FA,
typename TRIA,
typename TRIC>
562#if (AMREX_SPACEDIM < 3)
565 auto facx = Math::pi<T>()/
T(m_geom.Domain().length(0));
566 auto facy = Math::pi<T>()/
T(m_geom.Domain().length(1));
569 auto dxfac =
T(2)/
T(m_geom.CellSize(0)*m_geom.CellSize(0));
570 auto dyfac =
T(2)/
T(m_geom.CellSize(1)*m_geom.CellSize(1));
571 auto scale = (m_r2x) ? m_r2x->scalingFactor() : m_r2c->scalingFactor();
573 if (m_geom.Domain().length(0) == 1) { dxfac = 0; }
574 if (m_geom.Domain().length(1) == 1) { dyfac = 0; }
577 for (
int idim = 0; idim < AMREX_SPACEDIM-1; ++idim) {
578 if (m_geom.Domain().length(idim) > 1) {
598 (std::is_same_v<TRIA,fft_poisson_detail::Tri_Zero<T>> &&
599 std::is_same_v<TRIC,fft_poisson_detail::Tri_Zero<T>>) {
601#if defined(AMREX_USE_OMP) && !defined(AMREX_USE_GPU)
606 auto const& spectral = spmf.array(mfi);
607 auto const& box = mfi.validbox();
612 T k2 = dxfac * (std::cos(a)-
T(1))
613 + dyfac * (std::cos(b)-
T(1));
615 spectral(i,j,k) /= k2;
617 spectral(i,j,k) *=
scale;
624 && zlo_neumann && zhi_neumann;
626 auto nz = m_geom.Domain().length(2);
628#if defined(AMREX_USE_OMP) && !defined(AMREX_USE_GPU)
633 auto const& spectral = spmf.array(mfi);
634 auto const& box = mfi.validbox();
643 auto const& ald = tridiag_workspace.
array(0);
644 auto const& bd = tridiag_workspace.
array(1);
645 auto const& cud = tridiag_workspace.
array(2);
646 auto const& scratch = tridiag_workspace.
array(3);
652 T k2 = dxfac * (std::cos(a)-
T(1))
653 + dyfac * (std::cos(b)-
T(1));
656 for(
int k=0; k < nz; k++) {
659 cud(i,j,k) = tric(i,j,k);
661 bd(i,j,k) = k2 - cud(i,j,k);
663 bd(i,j,k) = k2 - cud(i,j,k) -
T(2.0)*tria(i,j,k);
665 }
else if (k == nz-1) {
666 ald(i,j,k) = tria(i,j,k);
669 bd(i,j,k) = k2 - ald(i,j,k);
670 if (i == 0 && j == 0 && is_singular) {
674 bd(i,j,k) = k2 - ald(i,j,k) -
T(2.0)*tric(i,j,k);
677 ald(i,j,k) = tria(i,j,k);
678 cud(i,j,k) = tric(i,j,k);
679 bd(i,j,k) = k2 -ald(i,j,k)-cud(i,j,k);
683 scratch(i,j,0) = cud(i,j,0)/bd(i,j,0);
684 spectral(i,j,0) = spectral(i,j,0)/bd(i,j,0);
686 for (
int k = 1; k < nz; k++) {
688 scratch(i,j,k) = cud(i,j,k) / (bd(i,j,k) - ald(i,j,k) * scratch(i,j,k-1));
690 spectral(i,j,k) = (spectral(i,j,k) - ald(i,j,k) * spectral(i,j,k - 1))
691 / (bd(i,j,k) - ald(i,j,k) * scratch(i,j,k-1));
694 for (
int k = nz - 2; k >= 0; k--) {
695 spectral(i,j,k) -= scratch(i,j,k) * spectral(i,j,k + 1);
698 for (
int k = 0; k < nz; ++k) {
699 spectral(i,j,k) *=
scale;
715 T k2 = dxfac * (std::cos(a)-
T(1))
716 + dyfac * (std::cos(b)-
T(1));
719 for(
int k=0; k < nz; k++) {
722 cud[k] = tric(i,j,k);
726 bd[k] = k2 - cud[k] -
T(2.0)*tria(i,j,k);
728 }
else if (k == nz-1) {
729 ald[k] = tria(i,j,k);
733 if (i == 0 && j == 0 && is_singular) {
737 bd[k] = k2 - ald[k] -
T(2.0)*tric(i,j,k);
740 ald[k] = tria(i,j,k);
741 cud[k] = tric(i,j,k);
742 bd[k] = k2 -ald[k]-cud[k];
746 scratch[0] = cud[0]/bd[0];
747 spectral(i,j,0) = spectral(i,j,0)/bd[0];
749 for (
int k = 1; k < nz; k++) {
751 scratch[k] = cud[k] / (bd[k] - ald[k] * scratch[k-1]);
753 spectral(i,j,k) = (spectral(i,j,k) - ald[k] * spectral(i,j,k - 1))
754 / (bd[k] - ald[k] * scratch[k-1]);
757 for (
int k = nz - 2; k >= 0; k--) {
758 spectral(i,j,k) -= scratch[k] * spectral(i,j,k + 1);
761 for (
int k = 0; k < nz; ++k) {
762 spectral(i,j,k) *=
scale;
782 Box const& box () const noexcept {
return dbox; }
785template <
typename MF>
786void fill_physbc (MF& mf, Geometry
const& geom,
787 Array<std::pair<Boundary,Boundary>,AMREX_SPACEDIM>
const& bc)
789 using T =
typename MF::value_type;
790 using Tag = FFTPhysBCTag<T>;
793 for (MFIter mfi(mf, MFItInfo{}.DisableDeviceSync()); mfi.isValid(); ++mfi)
795 auto const& box = mfi.fabbox();
796 auto const& arr = mf.array(mfi);
797 for (OrientationIter oit; oit; ++oit) {
798 Orientation face = oit();
799 int idim = face.coordDir();
800 Box b = geom.Domain();
803 b.setRange(idim,geom.Domain().smallEnd(idim)-1);
804 fbc = bc[idim].first;
806 b.setRange(idim,geom.Domain().bigEnd(idim)+1);
807 fbc = bc[idim].second;
810 if (b.ok() && fbc != Boundary::periodic) {
811 tags.push_back({arr, b, fbc, face});
816#if defined(AMREX_USE_GPU)
818 Tag
const& tag)
noexcept
820 auto ntags =
int(tags.size());
822#pragma omp parallel for
824 for (
int itag = 0; itag < ntags; ++itag) {
825 Tag
const& tag = tags[itag];
829 int sgn = tag.face.isLow() ? 1 : -1;
830 IntVect siv = IntVect(AMREX_D_DECL(i,j,k))
831 + sgn * IntVect::TheDimensionVector(tag.face.coordDir());
832 if (tag.bc == Boundary::odd) {
833 tag.dfab(i,j,k) = -tag.dfab(siv);
835 tag.dfab(i,j,k) = tag.dfab(siv);
838#if !defined(AMREX_USE_GPU)
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#define AMREX_ALWAYS_ASSERT(EX)
Definition AMReX_BLassert.H:50
#define AMREX_FORCE_INLINE
Definition AMReX_Extension.H:119
#define AMREX_GPU_DEVICE
Definition AMReX_GpuQualifiers.H:18
#define AMREX_GPU_HOST_DEVICE
Definition AMReX_GpuQualifiers.H:20
Array4< int const > offset
Definition AMReX_HypreMLABecLap.cpp:1089
#define AMREX_D_TERM(a, b, c)
Definition AMReX_SPACE.H:172
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
Array4< T const > array() const noexcept
Definition AMReX_BaseFab.H:381
__host__ __device__ BoxND & setBig(const IntVectND< dim > &bg) noexcept
Redefine the big end of the BoxND.
Definition AMReX_Box.H:487
__host__ __device__ IntVectND< dim > length() const noexcept
Return the length of the BoxND.
Definition AMReX_Box.H:154
A Fortran Array of REALs.
Definition AMReX_FArrayBox.H:231
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
void solve(MF &soln, MF const &rhs)
Definition AMReX_FFT_Poisson.H:469
typename MF::value_type T
Definition AMReX_FFT_Poisson.H:153
void solve_z(FA &spmf, TRIA const &tria, TRIC const &tric)
Definition AMReX_FFT_Poisson.H:558
PoissonHybrid(Geometry const &geom, Array< std::pair< Boundary, Boundary >, 3 > const &bc)
Definition AMReX_FFT_Poisson.H:156
std::pair< BoxArray, DistributionMapping > getSpectralDataLayout() const
Definition AMReX_FFT_Poisson.H:409
void solve_2d(MF &a_soln, MF const &a_rhs)
Definition AMReX_FFT_Poisson.H:504
Poisson solve for Open BC using FFT.
Definition AMReX_FFT_Poisson.H:123
void solve(MF &soln, MF const &rhs)
Definition AMReX_FFT_Poisson.H:354
void define_doit()
Definition AMReX_FFT_Poisson.H:325
PoissonOpenBC(Geometry const &geom, IndexType ixtype=IndexType::TheCellType(), IntVect const &ngrow=IntVect(0))
Definition AMReX_FFT_Poisson.H:314
Poisson solver for periodic, Dirichlet & Neumann boundaries using FFT.
Definition AMReX_FFT_Poisson.H:61
Poisson(Geometry const &geom, Array< std::pair< Boundary, Boundary >, 3 > const &bc)
Definition AMReX_FFT_Poisson.H:65
void solve(MF &a_soln, MF const &a_rhs)
Definition AMReX_FFT_Poisson.H:230
Poisson(Geometry const &geom)
Definition AMReX_FFT_Poisson.H:85
An Array of FortranArrayBox(FAB)-like Objects.
Definition AMReX_FabArray.H:347
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:74
const Box & Domain() const noexcept
Returns our rectangular domain.
Definition AMReX_Geometry.H:211
bool isAllPeriodic() const noexcept
Is domain periodic in all directions?
Definition AMReX_Geometry.H:339
__host__ static __device__ constexpr IndexTypeND< dim > TheCellType() noexcept
This static member function returns an IndexTypeND object of value IndexTypeND::CELL....
Definition AMReX_IndexType.H:152
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:63
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:147
Encapsulation of the Orientation of the Faces of a Box.
Definition AMReX_Orientation.H:29
Dynamically allocated vector for trivially copyable data.
Definition AMReX_PODVector.H:308
size_type size() const noexcept
Definition AMReX_PODVector.H:648
T * dataPtr() noexcept
Definition AMReX_PODVector.H:670
T * data() noexcept
Definition AMReX_PODVector.H:666
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
__host__ __device__ BoxND< dim > grow(const BoxND< dim > &b, int i) noexcept
Grow BoxND in all directions by given amount.
Definition AMReX_Box.H:1280
std::array< T, N > Array
Definition AMReX_Array.H:25
Definition AMReX_FFT_Helper.H:46
Boundary
Definition AMReX_FFT_Helper.H:52
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:263
void htod_memcpy_async(void *p_d, const void *p_h, const std::size_t sz) noexcept
Definition AMReX_GpuDevice.H:301
int NProcsSub() noexcept
number of ranks in current frame
Definition AMReX_ParallelContext.H:74
Definition AMReX_Amr.cpp:49
@ make_alias
Definition AMReX_MakeType.H:7
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:138
__host__ __device__ BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
Return a BoxND with different type.
Definition AMReX_Box.H:1558
__host__ __device__ BoxND< dim > makeSlab(BoxND< dim > const &b, int direction, int slab_index) noexcept
Definition AMReX_Box.H:2123
std::enable_if_t< std::is_integral_v< T > > ParallelFor(TypeList< CTOs... > ctos, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
Definition AMReX_CTOParallelForImpl.H:193
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:27
double second() noexcept
Definition AMReX_Utility.cpp:940
__host__ __device__ BoxND< dim > shift(const BoxND< dim > &b, int dir, int nzones) noexcept
Return a BoxND with indices shifted by nzones in dir direction.
Definition AMReX_Box.H:1495
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
bool TilingIfNotGPU() noexcept
Definition AMReX_MFIter.H:12
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:230
void LoopOnCpu(Dim3 lo, Dim3 hi, F const &f) noexcept
Definition AMReX_Loop.H:365
__host__ __device__ constexpr T elemwiseMin(T const &a, T const &b) noexcept
Definition AMReX_Algorithm.H:49
__host__ __device__ IntVectND< dim > scale(const IntVectND< dim > &p, int s) noexcept
Returns a IntVectND obtained by multiplying each of the components of this IntVectND by s.
Definition AMReX_IntVect.H:1016
Definition AMReX_Array4.H:61
Definition AMReX_FFT_Helper.H:58
Info & setOneDMode(bool x)
Definition AMReX_FFT_Helper.H:83
Info & setTwoDMode(bool x)
Definition AMReX_FFT_Helper.H:82
Fixed-size array that can be used on GPU.
Definition AMReX_Array.H:40