aboutsummaryrefslogtreecommitdiff
path: root/Source/FieldSolver/SpectralSolver/SpectralFieldData.cpp
diff options
context:
space:
mode:
authorGravatar Remi Lehe <remi.lehe@normalesup.org> 2019-04-19 20:23:10 -0700
committerGravatar Remi Lehe <remi.lehe@normalesup.org> 2019-04-23 12:43:53 -0700
commit0b77b789043461c77780ab56ab46a1e564b98d6e (patch)
tree852698700497c3a367b4752cefd2b75c82aff206 /Source/FieldSolver/SpectralSolver/SpectralFieldData.cpp
parent8bffe04d8eb683230a6b9919438adf9239c87fec (diff)
downloadWarpX-0b77b789043461c77780ab56ab46a1e564b98d6e.tar.gz
WarpX-0b77b789043461c77780ab56ab46a1e564b98d6e.tar.zst
WarpX-0b77b789043461c77780ab56ab46a1e564b98d6e.zip
Normalize inverse FFT
Diffstat (limited to 'Source/FieldSolver/SpectralSolver/SpectralFieldData.cpp')
-rw-r--r--Source/FieldSolver/SpectralSolver/SpectralFieldData.cpp5
1 files changed, 4 insertions, 1 deletions
diff --git a/Source/FieldSolver/SpectralSolver/SpectralFieldData.cpp b/Source/FieldSolver/SpectralSolver/SpectralFieldData.cpp
index ee6dbff6a..595cbac16 100644
--- a/Source/FieldSolver/SpectralSolver/SpectralFieldData.cpp
+++ b/Source/FieldSolver/SpectralSolver/SpectralFieldData.cpp
@@ -154,15 +154,18 @@ SpectralFieldData::BackwardTransform( MultiFab& mf,
// The copy does *not* fill the *last* point of `mf`
// in any direction that has *nodal* index type (but this point is
// in the guard cells and will be filled by guard cell exchange)
+ // Normalize (divide by 1/N) since the FFT result in a factor N
{
Box bx = mf[mfi].box();
const Box realspace_bx = bx.enclosedCells(); // discards last point in each nodal direction
AMREX_ALWAYS_ASSERT( realspace_bx == tmpRealField[mfi].box() );
Array4<Real> mf_arr = mf[mfi].array();
Array4<const Complex> tmp_arr = tmpRealField[mfi].array();
+ // For normalization: divide by the number of points in the box
+ const Real inv_N = 1./(bx.length(0)*bx.length(1)*bx.length(2));
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();
+ mf_arr(i,j,k,i_comp) = inv_N*tmp_arr(i,j,k).real();
});
}
}