1#ifndef AMREX_FFT_POISSON_H_
2#define AMREX_FFT_POISSON_H_
19void fill_physbc (MF& mf, Geometry
const& geom,
20 Array<std::pair<Boundary,Boundary>,AMREX_SPACEDIM>
const& bc);
22inline Geometry shift_geom (Geometry
const& geom)
24 auto const& dom = geom.Domain();
25 auto const& dlo = dom.smallEnd();
29 return Geometry(
amrex::shift(dom,-dlo), geom.ProbDomain(), geom.CoordInt(),
35void shift_mfs (IntVect
const& domain_lo, MF
const& mf1, MF
const& mf2,
36 MF& new_mf1, MF& new_mf2)
39 mf1.DistributionMap() == mf2.DistributionMap());
40 BoxArray ba = mf1.boxArray();
42 new_mf1.define(ba, mf1.DistributionMap(), 1, mf1.nGrowVect(),
43 MFInfo().SetAlloc(
false));
44 new_mf2.define(ba, mf2.DistributionMap(), 1, mf2.nGrowVect(),
45 MFInfo().SetAlloc(
false));
46 for (MFIter mfi(mf1); mfi.isValid(); ++mfi) {
48 fab1.shift(-domain_lo);
49 new_mf1.setFab(mfi, std::move(fab1));
52 fab2.shift(-domain_lo);
53 new_mf2.setFab(mfi, std::move(fab2));
65template <
typename MF = MultiFab>
76 template <
typename FA=MF, std::enable_if_t<IsFabArray_v<FA>,
int> = 0>
78 Array<std::pair<Boundary,Boundary>,AMREX_SPACEDIM>
const& bc)
79 : m_domain_lo(geom.Domain().smallEnd()),
80 m_geom(detail::shift_geom(geom)),
83 bool all_periodic =
true;
84 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
85 all_periodic = all_periodic
90 m_r2c = std::make_unique<R2C<typename MF::value_type>>(m_geom.
Domain());
92 m_r2x = std::make_unique<R2X<typename MF::value_type>> (m_geom.
Domain(), m_bc);
101 template <
typename FA=MF, std::enable_if_t<IsFabArray_v<FA>,
int> = 0>
103 : m_domain_lo(geom.Domain().smallEnd()),
104 m_geom(detail::shift_geom(geom)),
110 m_r2c = std::make_unique<R2C<typename MF::value_type>>(m_geom.
Domain());
126 void solve (MF& a_soln, MF
const& a_rhs);
132 std::unique_ptr<R2X<typename MF::value_type>> m_r2x;
133 std::unique_ptr<R2C<typename MF::value_type>> m_r2c;
136#if (AMREX_SPACEDIM == 3)
141template <
typename MF = MultiFab>
153 template <
typename FA=MF, std::enable_if_t<IsFabArray_v<FA>,
int> = 0>
164 void solve (MF& soln, MF
const& rhs);
185template <
typename MF = MultiFab>
189 using T =
typename MF::value_type;
197 template <
typename FA=MF, std::enable_if_t<IsFabArray_v<FA>,
int> = 0>
199 Array<std::pair<Boundary,Boundary>,AMREX_SPACEDIM>
const& bc)
200 : m_domain_lo(geom.Domain().smallEnd()),
201 m_geom(detail::shift_geom(geom)),
204#if (AMREX_SPACEDIM < 3)
208 bool periodic_xy =
true;
209 for (
int idim = 0; idim < 2; ++idim) {
226 m_r2c = std::make_unique<R2C<typename MF::value_type>>(m_geom.
Domain(),
229 m_r2x = std::make_unique<R2X<typename MF::value_type>> (m_geom.
Domain(),
244 void solve (MF& soln, MF
const& rhs);
268 void solve_2d (MF& a_soln, MF
const& a_rhs);
278 template <
typename TRIA,
typename TRIC>
279 void solve (MF& a_soln, MF
const& a_rhs, TRIA
const& tria, TRIC
const& tric);
289 template <
typename FA,
typename TRIA,
typename TRIC>
290 void solve_z (FA& spmf, TRIA
const& tria, TRIC
const& tric);
306 std::unique_ptr<R2X<typename MF::value_type>> m_r2x;
307 std::unique_ptr<R2C<typename MF::value_type>> m_r2c;
313template <
typename MF>
319 MF
const* rhs = &a_rhs;
321 if (m_domain_lo != 0) {
322 detail::shift_mfs(m_domain_lo, a_soln, a_rhs, solntmp, rhstmp);
327 AMREX_ASSERT(soln->is_cell_centered() && rhs->is_cell_centered());
329 using T =
typename MF::value_type;
332 {
AMREX_D_DECL(Math::pi<T>()/T(m_geom.Domain().length(0)),
333 Math::pi<T>()/T(m_geom.Domain().length(1)),
334 Math::pi<T>()/T(m_geom.Domain().length(2)))};
335 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
341 {
AMREX_D_DECL(T(2)/T(m_geom.CellSize(0)*m_geom.CellSize(0)),
342 T(2)/T(m_geom.CellSize(1)*m_geom.CellSize(1)),
343 T(2)/T(m_geom.CellSize(2)*m_geom.CellSize(2)))};
344 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
345 if (m_geom.Domain().length(idim) == 1) {
349 auto scale = (m_r2x) ? m_r2x->scalingFactor() : m_r2c->scalingFactor();
352 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
371 T b = fac[1]*(j+
offset[1]);,
372 T c = fac[2]*(k+
offset[2]));
374 +dxfac[1]*(std::cos(b)-T(1)),
375 +dxfac[2]*(std::cos(c)-T(1)));
379 spectral_data *=
scale;
385 m_r2x->forwardThenBackward_doit_0(*rhs, *soln, f, ng, m_geom.periodicity());
386 detail::fill_physbc(*soln, m_geom, m_bc);
388 m_r2c->forward(*rhs);
389 m_r2c->post_forward_doit_0(f);
390 m_r2c->backward_doit(*soln, ng, m_geom.periodicity());
394#if (AMREX_SPACEDIM == 3)
396template <
typename MF>
397template <
typename FA, std::enable_if_t<IsFabArray_v<FA>,
int> FOO>
403 m_solver(m_grown_domain)
408template <
typename MF>
411 using T =
typename MF::value_type;
412 auto const& lo = m_grown_domain.smallEnd();
413 auto const dx = T(m_geom.CellSize(0));
414 auto const dy = T(m_geom.CellSize(1));
415 auto const dz = T(m_geom.CellSize(2));
416 auto const gfac = T(1)/T(std::sqrt(T(12)));
418 auto const fac = T(-0.125) * (dx*dy*dz) / (T(4)*Math::pi<T>());
421 auto x = (T(i-lo[0]) - gfac) * dx;
422 auto y = (T(j-lo[1]) - gfac) * dy;
423 auto z = (T(k-lo[2]) - gfac) * dz;
425 for (
int gx = 0; gx < 2; ++gx) {
426 for (
int gy = 0; gy < 2; ++gy) {
427 for (
int gz = 0; gz < 2; ++gz) {
428 auto xg =
x + 2*gx*gfac*dx;
429 auto yg =
y + 2*gy*gfac*dy;
430 auto zg =
z + 2*gz*gfac*dz;
431 r += T(1)/std::sqrt(xg*xg+yg*yg+zg*zg);
437template <
typename MF>
440 AMREX_ASSERT(m_grown_domain.ixType() == soln.ixType() && m_grown_domain.ixType() == rhs.ixType());
441 m_solver.solve(soln, rhs);
447namespace fft_poisson_detail {
448 template <
typename T>
450 [[nodiscard]]
constexpr T operator() (
int,
int,
int)
const
456 template <
typename T>
459 T operator() (
int,
int,
int)
const
466 template <
typename T>
469 T operator() (
int,
int,
int k)
const
471 return (k > 0) ? T(2.0) / (m_dz[k]*(m_dz[k]+m_dz[k-1]))
472 : T(1.0) / (m_dz[k]* m_dz[k]);
477 template <
typename T>
480 T operator() (
int,
int,
int k)
const
482 return (k < m_size-1) ? T(2.0) / (m_dz[k]*(m_dz[k]+m_dz[k+1]))
483 : T(1.0) / (m_dz[k]* m_dz[k]);
491template <
typename MF>
492std::pair<BoxArray,DistributionMapping>
495 if (!m_spmf_r.empty()) {
496 return std::make_pair(m_spmf_r.boxArray(), m_spmf_r.DistributionMap());
498 return std::make_pair(m_spmf_c.boxArray(), m_spmf_c.DistributionMap());
502template <
typename MF>
505#if (AMREX_SPACEDIM == 3)
507 (m_geom.Domain().length(0) > 1 ||
508 m_geom.Domain().length(1) > 1));
511 Box cdomain = m_geom.Domain();
512 if (cdomain.
length(0) > 1) {
519 DistributionMapping dm = detail::make_iota_distromap(cba.size());
520 m_spmf_c.define(cba, dm, 1, 0);
521 }
else if (m_geom.Domain().length(0) > 1 &&
522 m_geom.Domain().length(1) > 1) {
523 if (m_r2x->m_cy.empty()) {
525 {AMREX_D_DECL(true,true,false)});
526 DistributionMapping dm = detail::make_iota_distromap(sba.size());
527 m_spmf_r.define(sba, dm, 1, 0);
529 Box cdomain = m_geom.Domain();
531 cdomain.
setBig(0,cdomain.length(0)/2);
533 cdomain.
setBig(1,cdomain.length(1)/2);
537 DistributionMapping dm = detail::make_iota_distromap(cba.size());
538 m_spmf_c.define(cba, dm, 1, 0);
543 {AMREX_D_DECL(true,true,false)});
544 DistributionMapping dm = detail::make_iota_distromap(sba.size());
545 m_spmf_r.define(sba, dm, 1, 0);
552template <
typename MF>
555 auto delz =
T(m_geom.CellSize(AMREX_SPACEDIM-1));
557 fft_poisson_detail::Tri_Uniform<T>{
T(1)/(delz*delz)},
558 fft_poisson_detail::Tri_Uniform<T>{
T(1)/(delz*delz)});
561template <
typename MF>
564 auto const* pdz = dz.
dataPtr();
566 fft_poisson_detail::TriA<T>{pdz},
567 fft_poisson_detail::TriC<T>{pdz,
int(dz.
size())});
570template <
typename MF>
573 AMREX_ASSERT(soln.is_cell_centered() && rhs.is_cell_centered());
578 auto const* pdz = d_dz.
data();
580 auto const* pdz = dz.data();
583 fft_poisson_detail::TriA<T>{pdz},
584 fft_poisson_detail::TriC<T>{pdz,
int(dz.
size())});
587template <
typename MF>
590 solve(soln, rhs, fft_poisson_detail::Tri_Zero<T>{}, fft_poisson_detail::Tri_Zero<T>{});
593template <
typename MF>
594template <
typename TRIA,
typename TRIC>
600 AMREX_ASSERT(a_soln.is_cell_centered() && a_rhs.is_cell_centered());
602#if (AMREX_SPACEDIM < 3)
607 MF
const* rhs = &a_rhs;
609 if (m_domain_lo != 0) {
610 detail::shift_mfs(m_domain_lo, a_soln, a_rhs, solntmp, rhstmp);
619 m_r2c->forward(*rhs, m_spmf_c);
620 solve_z(m_spmf_c, tria, tric);
621 m_r2c->backward_doit(m_spmf_c, *soln, ng, m_geom.periodicity());
625 if (m_r2x->m_cy.empty()) {
626 m_r2x->forward(*rhs, m_spmf_r);
627 solve_z(m_spmf_r, tria, tric);
628 m_r2x->backward(m_spmf_r, *soln, ng, m_geom.periodicity());
630 m_r2x->forward(*rhs, m_spmf_c);
631 solve_z(m_spmf_c, tria, tric);
632 m_r2x->backward(m_spmf_c, *soln, ng, m_geom.periodicity());
636 detail::fill_physbc(*soln, m_geom, m_bc);
640template <
typename MF>
641template <
typename FA,
typename TRIA,
typename TRIC>
646#if (AMREX_SPACEDIM < 3)
649 auto facx = Math::pi<T>()/
T(m_geom.Domain().length(0));
650 auto facy = Math::pi<T>()/
T(m_geom.Domain().length(1));
653 auto dxfac =
T(2)/
T(m_geom.CellSize(0)*m_geom.CellSize(0));
654 auto dyfac =
T(2)/
T(m_geom.CellSize(1)*m_geom.CellSize(1));
655 auto scale = (m_r2x) ? m_r2x->scalingFactor() : m_r2c->scalingFactor();
657 if (m_geom.Domain().length(0) == 1) { dxfac = 0; }
658 if (m_geom.Domain().length(1) == 1) { dyfac = 0; }
661 for (
int idim = 0; idim < AMREX_SPACEDIM-1; ++idim) {
662 if (m_geom.Domain().length(idim) > 1) {
682 (std::is_same_v<TRIA,fft_poisson_detail::Tri_Zero<T>> &&
683 std::is_same_v<TRIC,fft_poisson_detail::Tri_Zero<T>>) {
685#if defined(AMREX_USE_OMP) && !defined(AMREX_USE_GPU)
690 auto const& spectral = spmf.array(mfi);
691 auto const& box = mfi.validbox();
696 T k2 = dxfac * (std::cos(a)-
T(1))
697 + dyfac * (std::cos(b)-
T(1));
699 spectral(i,j,k) /= k2;
701 spectral(i,j,k) *=
scale;
708 && zlo_neumann && zhi_neumann;
710 auto nz = m_geom.Domain().length(2);
712#if defined(AMREX_USE_OMP) && !defined(AMREX_USE_GPU)
717 auto const& spectral = spmf.array(mfi);
718 auto const& box = mfi.validbox();
727 auto const& ald = tridiag_workspace.
array(0);
728 auto const& bd = tridiag_workspace.
array(1);
729 auto const& cud = tridiag_workspace.
array(2);
730 auto const& scratch = tridiag_workspace.
array(3);
736 T k2 = dxfac * (std::cos(a)-
T(1))
737 + dyfac * (std::cos(b)-
T(1));
740 for(
int k=0; k < nz; k++) {
743 cud(i,j,k) = tric(i,j,k);
745 bd(i,j,k) = k2 - cud(i,j,k);
747 bd(i,j,k) = k2 - cud(i,j,k) -
T(2.0)*tria(i,j,k);
749 }
else if (k == nz-1) {
750 ald(i,j,k) = tria(i,j,k);
753 bd(i,j,k) = k2 - ald(i,j,k);
754 if (i == 0 && j == 0 && is_singular) {
758 bd(i,j,k) = k2 - ald(i,j,k) -
T(2.0)*tric(i,j,k);
761 ald(i,j,k) = tria(i,j,k);
762 cud(i,j,k) = tric(i,j,k);
763 bd(i,j,k) = k2 -ald(i,j,k)-cud(i,j,k);
767 scratch(i,j,0) = cud(i,j,0)/bd(i,j,0);
768 spectral(i,j,0) = spectral(i,j,0)/bd(i,j,0);
770 for (
int k = 1; k < nz; k++) {
772 scratch(i,j,k) = cud(i,j,k) / (bd(i,j,k) - ald(i,j,k) * scratch(i,j,k-1));
774 spectral(i,j,k) = (spectral(i,j,k) - ald(i,j,k) * spectral(i,j,k - 1))
775 / (bd(i,j,k) - ald(i,j,k) * scratch(i,j,k-1));
778 for (
int k = nz - 2; k >= 0; k--) {
779 spectral(i,j,k) -= scratch(i,j,k) * spectral(i,j,k + 1);
782 for (
int k = 0; k < nz; ++k) {
783 spectral(i,j,k) *=
scale;
799 T k2 = dxfac * (std::cos(a)-
T(1))
800 + dyfac * (std::cos(b)-
T(1));
803 for(
int k=0; k < nz; k++) {
806 cud[k] = tric(i,j,k);
810 bd[k] = k2 - cud[k] -
T(2.0)*tria(i,j,k);
812 }
else if (k == nz-1) {
813 ald[k] = tria(i,j,k);
817 if (i == 0 && j == 0 && is_singular) {
821 bd[k] = k2 - ald[k] -
T(2.0)*tric(i,j,k);
824 ald[k] = tria(i,j,k);
825 cud[k] = tric(i,j,k);
826 bd[k] = k2 -ald[k]-cud[k];
830 scratch[0] = cud[0]/bd[0];
831 spectral(i,j,0) = spectral(i,j,0)/bd[0];
833 for (
int k = 1; k < nz; k++) {
835 scratch[k] = cud[k] / (bd[k] - ald[k] * scratch[k-1]);
837 spectral(i,j,k) = (spectral(i,j,k) - ald[k] * spectral(i,j,k - 1))
838 / (bd[k] - ald[k] * scratch[k-1]);
841 for (
int k = nz - 2; k >= 0; k--) {
842 spectral(i,j,k) -= scratch[k] * spectral(i,j,k + 1);
845 for (
int k = 0; k < nz; ++k) {
846 spectral(i,j,k) *=
scale;
866 Box const& box () const noexcept {
return dbox; }
869template <
typename MF>
870void fill_physbc (MF& mf, Geometry
const& geom,
871 Array<std::pair<Boundary,Boundary>,AMREX_SPACEDIM>
const& bc)
873 using T =
typename MF::value_type;
874 using Tag = FFTPhysBCTag<T>;
877 for (MFIter mfi(mf, MFItInfo{}.DisableDeviceSync()); mfi.isValid(); ++mfi)
879 auto const& box = mfi.fabbox();
880 auto const& arr = mf.array(mfi);
881 for (OrientationIter oit; oit; ++oit) {
882 Orientation face = oit();
883 int idim = face.coordDir();
884 Box b = geom.Domain();
887 b.setRange(idim,geom.Domain().smallEnd(idim)-1);
888 fbc = bc[idim].first;
890 b.setRange(idim,geom.Domain().bigEnd(idim)+1);
891 fbc = bc[idim].second;
894 if (b.ok() && fbc != Boundary::periodic) {
895 tags.push_back({arr, b, fbc, face});
900#if defined(AMREX_USE_GPU)
902 Tag
const& tag)
noexcept
904 auto ntags =
int(tags.size());
906#pragma omp parallel for
908 for (
int itag = 0; itag < ntags; ++itag) {
909 Tag
const& tag = tags[itag];
913 int sgn = tag.face.isLow() ? 1 : -1;
914 IntVect siv = IntVect(AMREX_D_DECL(i,j,k))
915 + sgn * IntVect::TheDimensionVector(tag.face.coordDir());
916 if (tag.bc == Boundary::odd) {
917 tag.dfab(i,j,k) = -tag.dfab(siv);
919 tag.dfab(i,j,k) = tag.dfab(siv);
922#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:382
__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
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
void solve(MF &soln, MF const &rhs)
Solve del dot grad soln = rhs for uniform spacing in all directions.
Definition AMReX_FFT_Poisson.H:553
typename MF::value_type T
Definition AMReX_FFT_Poisson.H:189
void solve_z(FA &spmf, TRIA const &tria, TRIC const &tric)
CUDA helper that applies the supplied tridiagonal operator along z.
Definition AMReX_FFT_Poisson.H:642
PoissonHybrid(Geometry const &geom, Array< std::pair< Boundary, Boundary >, 3 > const &bc)
Construct the hybrid Poisson solver (mixed BCs in z with optional nonuniform spacing).
Definition AMReX_FFT_Poisson.H:198
std::pair< BoxArray, DistributionMapping > getSpectralDataLayout() const
Layout information for spectral storage used by the hybrid solver.
Definition AMReX_FFT_Poisson.H:493
void solve_2d(MF &a_soln, MF const &a_rhs)
Solve an independent 2-D Poisson problem on every z-plane.
Definition AMReX_FFT_Poisson.H:588
Poisson solve for Open BC using FFT.
Definition AMReX_FFT_Poisson.H:143
void solve(MF &soln, MF const &rhs)
Solve the open-boundary Poisson problem.
Definition AMReX_FFT_Poisson.H:438
void define_doit()
Initialize the discretized Green's function cache (public for CUDA kernels).
Definition AMReX_FFT_Poisson.H:409
PoissonOpenBC(Geometry const &geom, IndexType ixtype=IndexType::TheCellType(), IntVect const &ngrow=IntVect(0))
Build an open-boundary FFT solver over the grown domain.
Definition AMReX_FFT_Poisson.H:398
Poisson solver for periodic, Dirichlet & Neumann boundaries using FFT.
Definition AMReX_FFT_Poisson.H:67
Poisson(Geometry const &geom, Array< std::pair< Boundary, Boundary >, 3 > const &bc)
Construct a Poisson solver with explicit boundary types.
Definition AMReX_FFT_Poisson.H:77
void solve(MF &a_soln, MF const &a_rhs)
Solve del dot grad soln = rhs.
Definition AMReX_FFT_Poisson.H:314
Poisson(Geometry const &geom)
Construct a purely periodic Poisson solver.
Definition AMReX_FFT_Poisson.H:102
An Array of FortranArrayBox(FAB)-like Objects.
Definition AMReX_FabArray.H:349
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:85
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:169
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:26
Definition AMReX_FFT_Helper.H:52
Boundary
Definition AMReX_FFT_Helper.H:58
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:139
__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:30
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:1947
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
__host__ __device__ constexpr 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:1013
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
const int[]
Definition AMReX_BLProfiler.cpp:1664
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
Return the element-wise minimum of the given values for types like XDim3.
Definition AMReX_Algorithm.H:62
A multidimensional array accessor.
Definition AMReX_Array4.H:283
Definition AMReX_FFT_Helper.H:64
Info & setOneDMode(bool x)
Flag the degenerate 2-D mode (nx==1 or ny==1) that still batches along z.
Definition AMReX_FFT_Helper.H:117
Info & setTwoDMode(bool x)
Restrict transforms to the first two dimensions (3-D problems only).
Definition AMReX_FFT_Helper.H:106
Fixed-size array that can be used on GPU.
Definition AMReX_Array.H:43