Skip to content

Commit

Permalink
Fix FFT size
Browse files Browse the repository at this point in the history
  • Loading branch information
RemiLehe committed May 1, 2019
1 parent 52cf650 commit 581510f
Showing 1 changed file with 8 additions and 5 deletions.
13 changes: 8 additions & 5 deletions Source/FieldSolver/SpectralSolver/SpectralFieldData.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -48,27 +48,30 @@ SpectralFieldData::SpectralFieldData( const BoxArray& realspace_ba,
// Loop over boxes and allocate the corresponding plan
// for each box owned by the local MPI proc
for ( MFIter mfi(spectralspace_ba, dm); mfi.isValid(); ++mfi ){
Box bx = spectralspace_ba[mfi];
// Note: the size of the real-space box and spectral-space box
// differ when using real-to-complex FFT. When initializing
// the FFT plan, the valid dimensions are those of the real-space box.
IntVect fft_size = realspace_ba[mfi].length();
#ifdef AMREX_USE_GPU
// Add cuFFT-specific code
#else
// Create FFTW plans
forward_plan[mfi] =
// Swap dimensions: AMReX FAB are Fortran-order but FFTW is C-order
#if (AMREX_SPACEDIM == 3)
fftw_plan_dft_r2c_3d( bx.length(2), bx.length(1), bx.length(0),
fftw_plan_dft_r2c_3d( fft_size[2], fft_size[1], fft_size[0],
#else
fftw_plan_dft_r2c_2d( bx.length(1), bx.length(0),
fftw_plan_dft_r2c_2d( fft_size[1], fft_size[0],
#endif
tmpRealField[mfi].dataPtr(),
reinterpret_cast<fftw_complex*>( tmpSpectralField[mfi].dataPtr() ),
FFTW_ESTIMATE );
backward_plan[mfi] =
// Swap dimensions: AMReX FAB are Fortran-order but FFTW is C-order
#if (AMREX_SPACEDIM == 3)
fftw_plan_dft_c2r_3d( bx.length(2), bx.length(1), bx.length(0),
fftw_plan_dft_c2r_3d( fft_size[2], fft_size[1], fft_size[0],
#else
fftw_plan_dft_c2r_2d( bx.length(1), bx.length(0),
fftw_plan_dft_c2r_2d( fft_size[1], fft_size[0],
#endif
reinterpret_cast<fftw_complex*>( tmpSpectralField[mfi].dataPtr() ),
tmpRealField[mfi].dataPtr(),
Expand Down

0 comments on commit 581510f

Please sign in to comment.