aboutsummaryrefslogtreecommitdiff
path: root/Source/FieldSolver/SpectralSolver
diff options
context:
space:
mode:
authorGravatar Remi Lehe <remi.lehe@normalesup.org> 2019-04-24 10:35:05 -0700
committerGravatar Remi Lehe <remi.lehe@normalesup.org> 2019-04-24 10:35:05 -0700
commit8e5df6f8622a90c29ecfd4dad27971be76971cd9 (patch)
treeb22f580ee6acabcf0c2f9cfed8fc7c957e5cab04 /Source/FieldSolver/SpectralSolver
parent4a3265daf691f0a8fa7461322c7299016fd706bf (diff)
downloadWarpX-8e5df6f8622a90c29ecfd4dad27971be76971cd9.tar.gz
WarpX-8e5df6f8622a90c29ecfd4dad27971be76971cd9.tar.zst
WarpX-8e5df6f8622a90c29ecfd4dad27971be76971cd9.zip
Correct out-of-bounds with real-space box
Diffstat (limited to 'Source/FieldSolver/SpectralSolver')
-rw-r--r--Source/FieldSolver/SpectralSolver/SpectralFieldData.cpp6
1 files changed, 2 insertions, 4 deletions
diff --git a/Source/FieldSolver/SpectralSolver/SpectralFieldData.cpp b/Source/FieldSolver/SpectralSolver/SpectralFieldData.cpp
index 15652addc..291fe945e 100644
--- a/Source/FieldSolver/SpectralSolver/SpectralFieldData.cpp
+++ b/Source/FieldSolver/SpectralSolver/SpectralFieldData.cpp
@@ -237,13 +237,11 @@ SpectralFieldData::BackwardTransform( MultiFab& mf,
#endif
// Copy the temporary field `tmpRealField` to the real-space field `mf`
- // in the *valid* part of the domain (guard cells are filled later,
- // through guard cell exchange)
{
- const Box realspace_valid_bx = mfi.validbox();
+ const Box realspace_bx = tmpRealField[mfi].box();
Array4<Real> mf_arr = mf[mfi].array();
Array4<const Complex> tmp_arr = tmpRealField[mfi].array();
- ParallelFor( realspace_valid_bx,
+ ParallelFor( realspace_bx,
[=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
mf_arr(i,j,k,i_comp) = tmp_arr(i,j,k).real();
});