aboutsummaryrefslogtreecommitdiff
path: root/Source/FieldSolver/SpectralSolver/SpectralFieldDataRZ.cpp
diff options
context:
space:
mode:
authorGravatar Axel Huebl <axel.huebl@plasma.ninja> 2023-05-31 17:41:59 -0700
committerGravatar GitHub <noreply@github.com> 2023-06-01 02:41:59 +0200
commit5f053ec122bcf751bcc2fb41660e43ac120c26c0 (patch)
treec955ca293739a7fc89247a67a6ebf6504445b16a /Source/FieldSolver/SpectralSolver/SpectralFieldDataRZ.cpp
parent3c341f88be150c427511e1d05c9c528497087757 (diff)
downloadWarpX-5f053ec122bcf751bcc2fb41660e43ac120c26c0.tar.gz
WarpX-5f053ec122bcf751bcc2fb41660e43ac120c26c0.tar.zst
WarpX-5f053ec122bcf751bcc2fb41660e43ac120c26c0.zip
Fix: PSATD+RZ in SP (#3956)
* Fix: PSATD+RZ in SP Fix the FFTW logic for RZ simulations using the PSATD field solver in single precision. Previously, the API was using the double precision interfaces for FFTW unconditionally. * Simplify Plans Reuse AnyFFT defines.
Diffstat (limited to 'Source/FieldSolver/SpectralSolver/SpectralFieldDataRZ.cpp')
-rw-r--r--Source/FieldSolver/SpectralSolver/SpectralFieldDataRZ.cpp39
1 files changed, 33 insertions, 6 deletions
diff --git a/Source/FieldSolver/SpectralSolver/SpectralFieldDataRZ.cpp b/Source/FieldSolver/SpectralSolver/SpectralFieldDataRZ.cpp
index b13c3121a..eecae3037 100644
--- a/Source/FieldSolver/SpectralSolver/SpectralFieldDataRZ.cpp
+++ b/Source/FieldSolver/SpectralSolver/SpectralFieldDataRZ.cpp
@@ -7,10 +7,14 @@
#include "SpectralFieldDataRZ.H"
#include "Utils/WarpXUtil.H"
+#include "FieldSolver/SpectralSolver/AnyFFT.H"
#include "WarpX.H"
#include <ablastr/warn_manager/WarnManager.H>
+#include <AMReX_Config.H>
+
+
using amrex::operator""_rt;
/* \brief Initialize fields in spectral space, and FFT plans
@@ -162,21 +166,31 @@ SpectralFieldDataRZ::SpectralFieldDataRZ (const int lev,
howmany_dims[1].os = 1;
forward_plan[mfi] =
// Note that AMReX FAB are Fortran-order.
- fftw_plan_guru_dft(1, // int rank
+#ifdef AMREX_USE_FLOAT
+ fftwf_plan_guru_dft
+#else
+ fftw_plan_guru_dft
+#endif
+ (1, // int rank
dims,
2, // int howmany_rank,
howmany_dims,
- reinterpret_cast<fftw_complex*>(tempHTransformed[mfi].dataPtr()), // fftw_complex *in
- reinterpret_cast<fftw_complex*>(tmpSpectralField[mfi].dataPtr()), // fftw_complex *out
+ reinterpret_cast<AnyFFT::Complex*>(tempHTransformed[mfi].dataPtr()), // complex *in
+ reinterpret_cast<AnyFFT::Complex*>(tmpSpectralField[mfi].dataPtr()), // complex *out
FFTW_FORWARD, // int sign
FFTW_ESTIMATE); // unsigned flags
backward_plan[mfi] =
- fftw_plan_guru_dft(1, // int rank
+#ifdef AMREX_USE_FLOAT
+ fftwf_plan_guru_dft
+#else
+ fftw_plan_guru_dft
+#endif
+ (1, // int rank
dims,
2, // int howmany_rank,
howmany_dims,
- reinterpret_cast<fftw_complex*>(tmpSpectralField[mfi].dataPtr()), // fftw_complex *in
- reinterpret_cast<fftw_complex*>(tempHTransformed[mfi].dataPtr()), // fftw_complex *out
+ reinterpret_cast<AnyFFT::Complex*>(tmpSpectralField[mfi].dataPtr()), // complex *in
+ reinterpret_cast<AnyFFT::Complex*>(tempHTransformed[mfi].dataPtr()), // complex *out
FFTW_BACKWARD, // int sign
FFTW_ESTIMATE); // unsigned flags
#endif
@@ -201,8 +215,13 @@ SpectralFieldDataRZ::~SpectralFieldDataRZ()
rocfft_plan_destroy(backward_plan[mfi]);
#else
// Destroy FFTW plans.
+# ifdef AMREX_USE_FLOAT
+ fftwf_destroy_plan(forward_plan[mfi]);
+ fftwf_destroy_plan(backward_plan[mfi]);
+# else
fftw_destroy_plan(forward_plan[mfi]);
fftw_destroy_plan(backward_plan[mfi]);
+# endif
#endif
}
}
@@ -280,7 +299,11 @@ SpectralFieldDataRZ::FABZForwardTransform (amrex::MFIter const & mfi, amrex::Box
amrex::The_Arena()->free(buffer);
result = rocfft_execution_info_destroy(execinfo);
#else
+# ifdef AMREX_USE_FLOAT
+ fftwf_execute(forward_plan[mfi]);
+# else
fftw_execute(forward_plan[mfi]);
+# endif
#endif
// Copy the spectral-space field `tmpSpectralField` to the appropriate
@@ -393,7 +416,11 @@ SpectralFieldDataRZ::FABZBackwardTransform (amrex::MFIter const & mfi, amrex::Bo
amrex::The_Arena()->free(buffer);
result = rocfft_execution_info_destroy(execinfo);
#else
+# ifdef AMREX_USE_FLOAT
+ fftwf_execute(backward_plan[mfi]);
+# else
fftw_execute(backward_plan[mfi]);
+# endif
#endif
// Copy the interleaved complex to the split complex.