1#ifndef AMREX_FFT_POISSON_H_
2#define AMREX_FFT_POISSON_H_
13 Array<std::pair<Boundary,Boundary>,AMREX_SPACEDIM>
const& bc);
17 auto const& dom = geom.
Domain();
29 MF& new_mf1, MF& new_mf2)
32 mf1.DistributionMap() == mf2.DistributionMap());
35 new_mf1.
define(ba, mf1.DistributionMap(), 1, mf1.nGrowVect(),
37 new_mf2.define(ba, mf2.DistributionMap(), 1, mf2.nGrowVect(),
41 fab1.shift(-domain_lo);
42 new_mf1.setFab(mfi, std::move(fab1));
45 fab2.shift(-domain_lo);
46 new_mf2.setFab(mfi, std::move(fab2));
56template <
typename MF = MultiFab>
61 template <
typename FA=MF, std::enable_if_t<IsFabArray_v<FA>,
int> = 0>
63 Array<std::pair<Boundary,Boundary>,AMREX_SPACEDIM>
const& bc)
68 bool all_periodic =
true;
69 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
70 all_periodic = all_periodic
81 template <
typename FA=MF, std::enable_if_t<IsFabArray_v<FA>,
int> = 0>
103 void solve (MF& a_soln, MF
const& a_rhs);
109 std::unique_ptr<R2X<typename MF::value_type>>
m_r2x;
110 std::unique_ptr<R2C<typename MF::value_type>>
m_r2c;
113#if (AMREX_SPACEDIM == 3)
117template <
typename MF = MultiFab>
122 template <
typename FA=MF, std::enable_if_t<IsFabArray_v<FA>,
int> = 0>
127 void solve (MF& soln, MF
const& rhs);
144template <
typename MF = MultiFab>
148 using T =
typename MF::value_type;
150 template <
typename FA=MF, std::enable_if_t<IsFabArray_v<FA>,
int> = 0>
152 Array<std::pair<Boundary,Boundary>,AMREX_SPACEDIM>
const& bc)
157#if (AMREX_SPACEDIM < 3)
161 bool periodic_xy =
true;
162 for (
int idim = 0; idim < 2; ++idim) {
195 void solve (MF& soln, MF
const& rhs);
199 void solve_2d (MF& a_soln, MF
const& a_rhs);
201 template <
typename TRIA,
typename TRIC>
202 void solve (MF& a_soln, MF
const& a_rhs, TRIA
const& tria, TRIC
const& tric);
205 template <
typename FA,
typename TRIA,
typename TRIC>
206 void solve_z (FA& spmf, TRIA
const& tria, TRIC
const& tric);
217 std::unique_ptr<R2X<typename MF::value_type>>
m_r2x;
218 std::unique_ptr<R2C<typename MF::value_type>>
m_r2c;
224template <
typename MF>
230 MF
const* rhs = &a_rhs;
232 if (m_domain_lo != 0) {
238 AMREX_ASSERT(soln->is_cell_centered() && rhs->is_cell_centered());
240 using T =
typename MF::value_type;
243 {
AMREX_D_DECL(Math::pi<T>()/T(m_geom.Domain().length(0)),
244 Math::pi<T>()/T(m_geom.Domain().length(1)),
245 Math::pi<T>()/T(m_geom.Domain().length(2)))};
246 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
252 {
AMREX_D_DECL(T(2)/T(m_geom.CellSize(0)*m_geom.CellSize(0)),
253 T(2)/T(m_geom.CellSize(1)*m_geom.CellSize(1)),
254 T(2)/T(m_geom.CellSize(2)*m_geom.CellSize(2)))};
255 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
256 if (m_geom.Domain().length(idim) == 1) {
260 auto scale = (m_r2x) ? m_r2x->scalingFactor() : m_r2c->scalingFactor();
263 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
283 T c = fac[2]*(k+
offset[2]));
285 +dxfac[1]*(std::cos(
b)-T(1)),
286 +dxfac[2]*(std::cos(c)-T(1)));
290 spectral_data *=
scale;
296 m_r2x->forwardThenBackward_doit_0(*rhs, *soln, f, ng, m_geom.periodicity());
299 m_r2c->forward(*rhs);
300 m_r2c->post_forward_doit_0(f);
301 m_r2c->backward_doit(*soln, ng, m_geom.periodicity());
305#if (AMREX_SPACEDIM == 3)
307template <
typename MF>
308template <
typename FA, std::enable_if_t<IsFabArray_v<FA>,
int> FOO>
314 m_solver(m_grown_domain)
319template <
typename MF>
322 using T =
typename MF::value_type;
323 auto const& lo = m_grown_domain.smallEnd();
324 auto const dx = T(m_geom.CellSize(0));
325 auto const dy = T(m_geom.CellSize(1));
326 auto const dz = T(m_geom.CellSize(2));
327 auto const gfac = T(1)/T(std::sqrt(T(12)));
329 auto const fac = T(-0.125) * (dx*dy*dz) / (T(4)*Math::pi<T>());
332 auto x = (T(i-lo[0]) - gfac) * dx;
333 auto y = (T(j-lo[1]) - gfac) * dy;
334 auto z = (T(k-lo[2]) - gfac) * dz;
336 for (
int gx = 0; gx < 2; ++gx) {
337 for (
int gy = 0; gy < 2; ++gy) {
338 for (
int gz = 0; gz < 2; ++gz) {
339 auto xg =
x + 2*gx*gfac*dx;
340 auto yg =
y + 2*gy*gfac*dy;
341 auto zg =
z + 2*gz*gfac*dz;
342 r += T(1)/std::sqrt(xg*xg+yg*yg+zg*zg);
348template <
typename MF>
351 AMREX_ASSERT(m_grown_domain.ixType() == soln.ixType() && m_grown_domain.ixType() == rhs.ixType());
352 m_solver.solve(soln, rhs);
357namespace fft_poisson_detail {
358 template <
typename T>
366 template <
typename T>
376 template <
typename T>
387 template <
typename T>
400template <
typename MF>
401std::pair<BoxArray,DistributionMapping>
404 if (!m_spmf_r.empty()) {
405 return std::make_pair(m_spmf_r.boxArray(), m_spmf_r.DistributionMap());
407 return std::make_pair(m_spmf_c.boxArray(), m_spmf_c.DistributionMap());
411template <
typename MF>
414#if (AMREX_SPACEDIM == 3)
416 (m_geom.Domain().length(0) > 1 ||
417 m_geom.Domain().length(1) > 1));
420 Box cdomain = m_geom.Domain();
421 if (cdomain.
length(0) > 1) {
429 m_spmf_c.define(cba, dm, 1, 0);
430 }
else if (m_geom.Domain().length(0) > 1 &&
431 m_geom.Domain().length(1) > 1) {
432 if (m_r2x->m_cy.empty()) {
434 {AMREX_D_DECL(true,true,false)});
436 m_spmf_r.define(sba, dm, 1, 0);
438 Box cdomain = m_geom.Domain();
447 m_spmf_c.define(cba, dm, 1, 0);
452 {AMREX_D_DECL(true,true,false)});
454 m_spmf_r.define(sba, dm, 1, 0);
461template <
typename MF>
464 auto delz =
T(m_geom.CellSize(AMREX_SPACEDIM-1));
470template <
typename MF>
473 auto const* pdz = dz.
dataPtr();
479template <
typename MF>
482 AMREX_ASSERT(soln.is_cell_centered() && rhs.is_cell_centered());
487 auto const* pdz = d_dz.
data();
489 auto const* pdz = dz.data();
496template <
typename MF>
502template <
typename MF>
503template <
typename TRIA,
typename TRIC>
509 AMREX_ASSERT(a_soln.is_cell_centered() && a_rhs.is_cell_centered());
511#if (AMREX_SPACEDIM < 3)
516 MF
const* rhs = &a_rhs;
518 if (m_domain_lo != 0) {
528 m_r2c->forward(*rhs, m_spmf_c);
529 solve_z(m_spmf_c, tria, tric);
530 m_r2c->backward_doit(m_spmf_c, *soln, ng, m_geom.periodicity());
534 if (m_r2x->m_cy.empty()) {
535 m_r2x->forward(*rhs, m_spmf_r);
536 solve_z(m_spmf_r, tria, tric);
537 m_r2x->backward(m_spmf_r, *soln, ng, m_geom.periodicity());
539 m_r2x->forward(*rhs, m_spmf_c);
540 solve_z(m_spmf_c, tria, tric);
541 m_r2x->backward(m_spmf_c, *soln, ng, m_geom.periodicity());
549template <
typename MF>
550template <
typename FA,
typename TRIA,
typename TRIC>
555#if (AMREX_SPACEDIM < 3)
558 auto facx = Math::pi<T>()/
T(m_geom.Domain().length(0));
559 auto facy = Math::pi<T>()/
T(m_geom.Domain().length(1));
562 auto dxfac =
T(2)/
T(m_geom.CellSize(0)*m_geom.CellSize(0));
563 auto dyfac =
T(2)/
T(m_geom.CellSize(1)*m_geom.CellSize(1));
564 auto scale = (m_r2x) ? m_r2x->scalingFactor() : m_r2c->scalingFactor();
566 if (m_geom.Domain().length(0) == 1) { dxfac = 0; }
567 if (m_geom.Domain().length(1) == 1) { dyfac = 0; }
570 for (
int idim = 0; idim < AMREX_SPACEDIM-1; ++idim) {
571 if (m_geom.Domain().length(idim) > 1) {
591 (std::is_same_v<TRIA,fft_poisson_detail::Tri_Zero<T>> &&
592 std::is_same_v<TRIC,fft_poisson_detail::Tri_Zero<T>>) {
594#if defined(AMREX_USE_OMP) && !defined(AMREX_USE_GPU)
599 auto const& spectral = spmf.array(mfi);
600 auto const& box = mfi.validbox();
605 T k2 = dxfac * (std::cos(a)-
T(1))
606 + dyfac * (std::cos(
b)-
T(1));
608 spectral(i,j,k) /= k2;
610 spectral(i,j,k) *=
scale;
617 && zlo_neumann && zhi_neumann;
619 auto nz = m_geom.Domain().length(2);
621#if defined(AMREX_USE_OMP) && !defined(AMREX_USE_GPU)
626 auto const& spectral = spmf.array(mfi);
627 auto const& box = mfi.validbox();
636 auto const& ald = tridiag_workspace.
array(0);
637 auto const& bd = tridiag_workspace.
array(1);
638 auto const& cud = tridiag_workspace.
array(2);
639 auto const& scratch = tridiag_workspace.
array(3);
645 T k2 = dxfac * (std::cos(a)-
T(1))
646 + dyfac * (std::cos(
b)-
T(1));
649 for(
int k=0; k < nz; k++) {
652 cud(i,j,k) = tric(i,j,k);
654 bd(i,j,k) = k2 - cud(i,j,k);
656 bd(i,j,k) = k2 - cud(i,j,k) -
T(2.0)*tria(i,j,k);
658 }
else if (k == nz-1) {
659 ald(i,j,k) = tria(i,j,k);
662 bd(i,j,k) = k2 - ald(i,j,k);
663 if (i == 0 && j == 0 && is_singular) {
667 bd(i,j,k) = k2 - ald(i,j,k) -
T(2.0)*tric(i,j,k);
670 ald(i,j,k) = tria(i,j,k);
671 cud(i,j,k) = tric(i,j,k);
672 bd(i,j,k) = k2 -ald(i,j,k)-cud(i,j,k);
676 scratch(i,j,0) = cud(i,j,0)/bd(i,j,0);
677 spectral(i,j,0) = spectral(i,j,0)/bd(i,j,0);
679 for (
int k = 1; k < nz; k++) {
681 scratch(i,j,k) = cud(i,j,k) / (bd(i,j,k) - ald(i,j,k) * scratch(i,j,k-1));
683 spectral(i,j,k) = (spectral(i,j,k) - ald(i,j,k) * spectral(i,j,k - 1))
684 / (bd(i,j,k) - ald(i,j,k) * scratch(i,j,k-1));
687 for (
int k = nz - 2; k >= 0; k--) {
688 spectral(i,j,k) -= scratch(i,j,k) * spectral(i,j,k + 1);
691 for (
int k = 0; k < nz; ++k) {
692 spectral(i,j,k) *=
scale;
708 T k2 = dxfac * (std::cos(a)-
T(1))
709 + dyfac * (std::cos(
b)-
T(1));
712 for(
int k=0; k < nz; k++) {
715 cud[k] = tric(i,j,k);
719 bd[k] = k2 - cud[k] -
T(2.0)*tria(i,j,k);
721 }
else if (k == nz-1) {
722 ald[k] = tria(i,j,k);
726 if (i == 0 && j == 0 && is_singular) {
730 bd[k] = k2 - ald[k] -
T(2.0)*tric(i,j,k);
733 ald[k] = tria(i,j,k);
734 cud[k] = tric(i,j,k);
735 bd[k] = k2 -ald[k]-cud[k];
739 scratch[0] = cud[0]/bd[0];
740 spectral(i,j,0) = spectral(i,j,0)/bd[0];
742 for (
int k = 1; k < nz; k++) {
744 scratch[k] = cud[k] / (bd[k] - ald[k] * scratch[k-1]);
746 spectral(i,j,k) = (spectral(i,j,k) - ald[k] * spectral(i,j,k - 1))
747 / (bd[k] - ald[k] * scratch[k-1]);
750 for (
int k = nz - 2; k >= 0; k--) {
751 spectral(i,j,k) -= scratch[k] * spectral(i,j,k + 1);
754 for (
int k = 0; k < nz; ++k) {
755 spectral(i,j,k) *=
scale;
777template <
typename MF>
779 Array<std::pair<Boundary,Boundary>,AMREX_SPACEDIM>
const& bc)
781 using T =
typename MF::value_type;
787 auto const& box = mfi.fabbox();
788 auto const& arr = mf.array(mfi);
796 fbc = bc[idim].first;
799 fbc = bc[idim].second;
803 tags.push_back({arr,
b, fbc, face});
808#if defined(AMREX_USE_GPU)
810 Tag
const& tag)
noexcept
812 auto ntags =
int(tags.
size());
814#pragma omp parallel for
816 for (
int itag = 0; itag < ntags; ++itag) {
817 Tag
const& tag = tags[itag];
821 int sgn = tag.face.isLow() ? 1 : -1;
822 IntVect siv = IntVect(AMREX_D_DECL(i,j,k))
823 + sgn * IntVect::TheDimensionVector(tag.face.coordDir());
824 if (tag.bc == Boundary::odd) {
825 tag.dfab(i,j,k) = -tag.dfab(siv);
827 tag.dfab(i,j,k) = tag.dfab(siv);
830#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:375
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:551
void define(const Box &bx)
Initialize the BoxArray from a single box. It is an error if the BoxArray has already been initialize...
Definition AMReX_BoxArray.cpp:345
BoxArray & shift(int dir, int nzones)
Apply Box::shift(int,int) to each Box in the BoxArray.
Definition AMReX_BoxArray.cpp:840
__host__ __device__ const IntVectND< dim > & bigEnd() const &noexcept
Get the bigend.
Definition AMReX_Box.H:119
__host__ __device__ BoxND & setBig(const IntVectND< dim > &bg) noexcept
Redefine the big end of the BoxND.
Definition AMReX_Box.H:477
__host__ __device__ IntVectND< dim > length() const noexcept
Return the length of the BoxND.
Definition AMReX_Box.H:149
__host__ __device__ BoxND & setRange(int dir, int sm_index, int n_cells=1) noexcept
Set the entire range in a given direction, starting at sm_index with length n_cells....
Definition AMReX_Box.H:1064
__host__ __device__ const IntVectND< dim > & smallEnd() const &noexcept
Get the smallend of the BoxND.
Definition AMReX_Box.H:108
int CoordInt() const noexcept
Returns the CoordType as an int.
Definition AMReX_CoordSys.H:44
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:41
A Fortran Array of REALs.
Definition AMReX_FArrayBox.H:229
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
void solve(MF &soln, MF const &rhs)
Definition AMReX_FFT_Poisson.H:462
std::unique_ptr< R2C< typename MF::value_type > > m_r2c
Definition AMReX_FFT_Poisson.H:218
typename MF::value_type T
Definition AMReX_FFT_Poisson.H:148
cMF m_spmf_c
Definition AMReX_FFT_Poisson.H:221
void solve_z(FA &spmf, TRIA const &tria, TRIC const &tric)
Definition AMReX_FFT_Poisson.H:551
MF m_spmf_r
Definition AMReX_FFT_Poisson.H:219
std::unique_ptr< R2X< typename MF::value_type > > m_r2x
Definition AMReX_FFT_Poisson.H:217
PoissonHybrid(Geometry const &geom, Array< std::pair< Boundary, Boundary >, 3 > const &bc)
Definition AMReX_FFT_Poisson.H:151
Geometry m_geom
Definition AMReX_FFT_Poisson.H:215
IntVect m_domain_lo
Definition AMReX_FFT_Poisson.H:214
std::pair< BoxArray, DistributionMapping > getSpectralDataLayout() const
Definition AMReX_FFT_Poisson.H:402
Array< std::pair< Boundary, Boundary >, 3 > m_bc
Definition AMReX_FFT_Poisson.H:216
void build_spmf()
Definition AMReX_FFT_Poisson.H:412
void solve_2d(MF &a_soln, MF const &a_rhs)
Definition AMReX_FFT_Poisson.H:497
Poisson solve for Open BC using FFT.
Definition AMReX_FFT_Poisson.H:119
Box m_grown_domain
Definition AMReX_FFT_Poisson.H:133
IntVect m_ngrow
Definition AMReX_FFT_Poisson.H:134
Geometry m_geom
Definition AMReX_FFT_Poisson.H:132
void solve(MF &soln, MF const &rhs)
Definition AMReX_FFT_Poisson.H:349
OpenBCSolver< typename MF::value_type > m_solver
Definition AMReX_FFT_Poisson.H:135
void define_doit()
Definition AMReX_FFT_Poisson.H:320
PoissonOpenBC(Geometry const &geom, IndexType ixtype=IndexType::TheCellType(), IntVect const &ngrow=IntVect(0))
Definition AMReX_FFT_Poisson.H:309
Poisson solver for periodic, Dirichlet & Neumann boundaries using FFT.
Definition AMReX_FFT_Poisson.H:58
Poisson(Geometry const &geom, Array< std::pair< Boundary, Boundary >, 3 > const &bc)
Definition AMReX_FFT_Poisson.H:62
IntVect m_domain_lo
Definition AMReX_FFT_Poisson.H:106
std::unique_ptr< R2X< typename MF::value_type > > m_r2x
Definition AMReX_FFT_Poisson.H:109
Geometry m_geom
Definition AMReX_FFT_Poisson.H:107
void solve(MF &a_soln, MF const &a_rhs)
Definition AMReX_FFT_Poisson.H:225
Poisson(Geometry const &geom)
Definition AMReX_FFT_Poisson.H:82
std::unique_ptr< R2C< typename MF::value_type > > m_r2c
Definition AMReX_FFT_Poisson.H:110
Array< std::pair< Boundary, Boundary >, 3 > m_bc
Definition AMReX_FFT_Poisson.H:108
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:73
const Box & Domain() const noexcept
Returns our rectangular domain.
Definition AMReX_Geometry.H:210
const RealBox & ProbDomain() const noexcept
Returns the problem domain.
Definition AMReX_Geometry.H:170
bool isAllPeriodic() const noexcept
Is domain periodic in all directions?
Definition AMReX_Geometry.H:338
bool isPeriodic(int dir) const noexcept
Is the domain periodic in the specified direction?
Definition AMReX_Geometry.H:331
__host__ static __device__ constexpr IndexTypeND< dim > TheCellType() noexcept
This static member function returns an IndexTypeND object of value IndexTypeND::CELL....
Definition AMReX_IndexType.H:149
Definition AMReX_MFIter.H:57
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:141
An Iterator over the Orientation of Faces of a Box.
Definition AMReX_Orientation.H:135
Encapsulation of the Orientation of the Faces of a Box.
Definition AMReX_Orientation.H:29
__host__ __device__ bool isLow() const noexcept
Returns true if Orientation is low.
Definition AMReX_Orientation.H:89
__host__ __device__ int coordDir() const noexcept
Returns the coordinate direction.
Definition AMReX_Orientation.H:83
Definition AMReX_PODVector.H:297
size_type size() const noexcept
Definition AMReX_PODVector.H:637
T * dataPtr() noexcept
Definition AMReX_PODVector.H:659
T * data() noexcept
Definition AMReX_PODVector.H:655
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
void fill_physbc(MF &mf, Geometry const &geom, Array< std::pair< Boundary, Boundary >, 3 > const &bc)
Definition AMReX_FFT_Poisson.H:778
Geometry shift_geom(Geometry const &geom)
Definition AMReX_FFT_Poisson.H:15
DistributionMapping make_iota_distromap(Long n)
Definition AMReX_FFT.cpp:88
void shift_mfs(IntVect const &domain_lo, MF const &mf1, MF const &mf2, MF &new_mf1, MF &new_mf2)
Definition AMReX_FFT_Poisson.H:28
Definition AMReX_FFT.cpp:7
Boundary
Definition AMReX_FFT_Helper.H:52
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:260
void htod_memcpy_async(void *p_d, const void *p_h, const std::size_t sz) noexcept
Definition AMReX_GpuDevice.H:289
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
Returns a BoxND with different type.
Definition AMReX_Box.H:1453
__host__ __device__ BoxND< dim > makeSlab(BoxND< dim > const &b, int direction, int slab_index) noexcept
Definition AMReX_Box.H:2016
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:191
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:1390
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
bool TilingIfNotGPU() noexcept
Definition AMReX_MFIter.H:12
__host__ __device__ BoxND< dim > grow(const BoxND< dim > &b, int i) noexcept
Grow BoxND in all directions by given amount.
Definition AMReX_Box.H:1229
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:355
__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:1011
std::array< T, N > Array
Definition AMReX_Array.H:24
Definition AMReX_FabArrayCommI.H:1000
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
Definition AMReX_FFT_Poisson.H:767
Orientation face
Definition AMReX_FFT_Poisson.H:771
__host__ __device__ Box const & box() const noexcept
Definition AMReX_FFT_Poisson.H:774
Boundary bc
Definition AMReX_FFT_Poisson.H:770
Array4< T > dfab
Definition AMReX_FFT_Poisson.H:768
Box dbox
Definition AMReX_FFT_Poisson.H:769
Definition AMReX_FFT_Poisson.H:377
T const * m_dz
Definition AMReX_FFT_Poisson.H:384
__device__ T operator()(int, int, int k) const
Definition AMReX_FFT_Poisson.H:379
Definition AMReX_FFT_Poisson.H:388
int m_size
Definition AMReX_FFT_Poisson.H:396
T const * m_dz
Definition AMReX_FFT_Poisson.H:395
__device__ T operator()(int, int, int k) const
Definition AMReX_FFT_Poisson.H:390
Definition AMReX_FFT_Poisson.H:359
constexpr T operator()(int, int, int) const
Definition AMReX_FFT_Poisson.H:360
Definition AMReX_Array.H:34
FabArray memory allocation information.
Definition AMReX_FabArray.H:66
Definition AMReX_MFIter.H:20
MFItInfo & DisableDeviceSync() noexcept
Definition AMReX_MFIter.H:38