diff options
author | 2020-02-04 17:06:35 -0800 | |
---|---|---|
committer | 2020-02-04 17:06:46 -0800 | |
commit | 13fb97ac903245f7232f8153d467e233fc71a2eb (patch) | |
tree | 99fc97094507c78192b78d196b90511c89213c84 /Source/WarpX.cpp | |
parent | aa8e26d4295b02864daf2279adbfb13e5b0f1f26 (diff) | |
download | WarpX-13fb97ac903245f7232f8153d467e233fc71a2eb.tar.gz WarpX-13fb97ac903245f7232f8153d467e233fc71a2eb.tar.zst WarpX-13fb97ac903245f7232f8153d467e233fc71a2eb.zip |
Implement 'warpx_build_buffer_masks' in C++.
Diffstat (limited to 'Source/WarpX.cpp')
-rw-r--r-- | Source/WarpX.cpp | 46 |
1 files changed, 42 insertions, 4 deletions
diff --git a/Source/WarpX.cpp b/Source/WarpX.cpp index 3fe61274a..907176908 100644 --- a/Source/WarpX.cpp +++ b/Source/WarpX.cpp @@ -1281,16 +1281,54 @@ WarpX::BuildBufferMasks () for (MFIter mfi(*bmasks, true); mfi.isValid(); ++mfi) { const Box& tbx = mfi.growntilebox(); - warpx_build_buffer_masks (BL_TO_FORTRAN_BOX(tbx), - BL_TO_FORTRAN_ANYD((*bmasks)[mfi]), - BL_TO_FORTRAN_ANYD(tmp[mfi]), - &ngbuffer); + //warpx_build_buffer_masks(BL_TO_FORTRAN_BOX(tbx), + // BL_TO_FORTRAN_ANYD((*bmasks)[mfi]), + // BL_TO_FORTRAN_ANYD(tmp[mfi]), + // &ngbuffer); + Array4<int> msk = (*bmasks)[mfi].array(); + Array4<int> gmsk = tmp[mfi].array(); + BuildBufferMasksInBox( tbx, msk, gmsk, ngbuffer ); } } } } } +void +WarpX::BuildBufferMasksInBox ( const Box &tbx, Array4<int> &msk, const Array4<int> gmsk, const int ng ) +{ + const auto lo = amrex::lbound( tbx ); + const auto hi = amrex::ubound( tbx ); +#if (AMREX_SPACEDIM == 2) + int k = lo.z; + for (int j = lo.y; j <= hi.y; ++j) { + for (int i = lo.x; i <= hi.x; ++i) { + for (int jj = j-ng; jj <= j+ng; ++jj) { + for (int ii = i-ng; ii <= i+ng; ++ii) { + if ( gmsk(ii,jj,k) == 0 ) msk(i,j,k) = 0; + else msk(i,j,k) = 1; + } + } + } + } +#elif (AMREX_SPACEDIM == 3) + for (int k = lo.z; k <= hi.z; ++k) { + for (int j = lo.y; j <= hi.y; ++j) { + for (int i = lo.x; i <= hi.x; ++i) { + for (int kk = k-ng; kk <= k+ng; ++kk) { + for (int jj = j-ng; jj <= j+ng; ++jj) { + for (int ii = i-ng; ii <= i+ng; ++ii) { + if ( gmsk(ii,jj,k) == 0 ) msk(i,j,k) = 0; + else msk(i,j,k) = 1; + } + } + } + } + } + } +#endif +} + const iMultiFab* WarpX::CurrentBufferMasks (int lev) { |