aboutsummaryrefslogtreecommitdiff
path: root/Source/FieldSolver/SpectralSolver
diff options
context:
space:
mode:
authorGravatar Revathi Jambunathan <revanathan@login2.summit.olcf.ornl.gov> 2019-06-06 19:19:57 -0400
committerGravatar Revathi Jambunathan <revanathan@login2.summit.olcf.ornl.gov> 2019-06-06 19:19:57 -0400
commit2569bcd08921227bedc7ccdb1018a5614ab31610 (patch)
tree88dcc1af385ea676782bb96327192679d9b2432a /Source/FieldSolver/SpectralSolver
parent984fb82044026463892db86894a68ed0343acba2 (diff)
downloadWarpX-2569bcd08921227bedc7ccdb1018a5614ab31610.tar.gz
WarpX-2569bcd08921227bedc7ccdb1018a5614ab31610.tar.zst
WarpX-2569bcd08921227bedc7ccdb1018a5614ab31610.zip
Included revisions as suggested in the PR
Diffstat (limited to 'Source/FieldSolver/SpectralSolver')
-rw-r--r--Source/FieldSolver/SpectralSolver/SpectralAlgorithms/PsatdAlgorithm.H4
-rw-r--r--Source/FieldSolver/SpectralSolver/SpectralFieldData.H1
-rw-r--r--Source/FieldSolver/SpectralSolver/SpectralFieldData.cpp16
3 files changed, 9 insertions, 12 deletions
diff --git a/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/PsatdAlgorithm.H b/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/PsatdAlgorithm.H
index 0a8614e54..12718e38b 100644
--- a/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/PsatdAlgorithm.H
+++ b/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/PsatdAlgorithm.H
@@ -1,8 +1,6 @@
#ifndef WARPX_PSATD_ALGORITHM_H_
#define WARPX_PSATD_ALGORITHM_H_
-//#include <SpectralKSpace.H>
-//#include <SpectralFieldData.H>
#include <SpectralBaseAlgorithm.H>
/* \brief Class that updates the field in spectral space
@@ -10,7 +8,6 @@
*/
class PsatdAlgorithm : public SpectralBaseAlgorithm
{
- using SpectralCoefficients = amrex::FabArray< amrex::BaseFab <amrex::Real> >;
public:
PsatdAlgorithm(const SpectralKSpace& spectral_kspace,
@@ -25,7 +22,6 @@ class PsatdAlgorithm : public SpectralBaseAlgorithm
void pushSpectralFields(SpectralFieldData& f) const override final;
private:
- // Modified finite-order vectors
SpectralCoefficients C_coef, S_ck_coef, X1_coef, X2_coef, X3_coef;
};
diff --git a/Source/FieldSolver/SpectralSolver/SpectralFieldData.H b/Source/FieldSolver/SpectralSolver/SpectralFieldData.H
index 1d64817ef..7954414b8 100644
--- a/Source/FieldSolver/SpectralSolver/SpectralFieldData.H
+++ b/Source/FieldSolver/SpectralSolver/SpectralFieldData.H
@@ -25,7 +25,6 @@ class SpectralFieldData
// (plans are only initialized for the boxes that are owned by
// the local MPI rank)
#ifdef AMREX_USE_GPU
- // Add cuFFT-specific code
using FFTplans = amrex::LayoutData<cufftHandle>;
#else
using FFTplans = amrex::LayoutData<fftw_plan>;
diff --git a/Source/FieldSolver/SpectralSolver/SpectralFieldData.cpp b/Source/FieldSolver/SpectralSolver/SpectralFieldData.cpp
index 6eeb266e7..a2b695568 100644
--- a/Source/FieldSolver/SpectralSolver/SpectralFieldData.cpp
+++ b/Source/FieldSolver/SpectralSolver/SpectralFieldData.cpp
@@ -53,9 +53,11 @@ SpectralFieldData::SpectralFieldData( const BoxArray& realspace_ba,
// 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
+ // Create cuFFT plans
// Creating 3D plan for real to complex -- double precision
// Assuming CUDA is used for programming GPU
+ // Note that D2Z is inherently forward plan
+ // and Z2D is inherently backward plan
cufftResult result;
#if (AMREX_SPACEDIM == 3)
result = cufftPlan3d( &forward_plan[mfi], fft_size[2],
@@ -70,8 +72,6 @@ SpectralFieldData::SpectralFieldData( const BoxArray& realspace_ba,
amrex::Print() << " cufftplan3d backward failed! \n";
}
#else
- // Add 2D cuFFT-spacific code for D2Z
- // Note that D2Z is inherently forward plan
result = cufftPlan2d( &forward_plan[mfi], fft_size[1],
fft_size[0], CUFFT_D2Z );
if ( result != CUFFT_SUCCESS ) {
@@ -117,7 +117,7 @@ SpectralFieldData::~SpectralFieldData()
if (tmpRealField.size() > 0){
for ( MFIter mfi(tmpRealField); mfi.isValid(); ++mfi ){
#ifdef AMREX_USE_GPU
- // Add cuFFT-specific code
+ // Destroy cuFFT plans
cufftDestroy( forward_plan[mfi] );
cufftDestroy( backward_plan[mfi] );
#else
@@ -168,8 +168,9 @@ SpectralFieldData::ForwardTransform( const MultiFab& mf,
// Perform Fourier transform from `tmpRealField` to `tmpSpectralField`
#ifdef AMREX_USE_GPU
- // Add cuFFT-specific code ; make sure that this is done on the same
- // GPU stream as the above copy
+ // Perform Fast Fourier Transform on GPU using cuFFT
+ // make sure that this is done on the same
+ // GPU stream as the above copy
cufftResult result;
cudaStream_t stream = amrex::Gpu::Device::cudaStream();
cufftSetStream ( forward_plan[mfi], stream);
@@ -270,7 +271,8 @@ SpectralFieldData::BackwardTransform( MultiFab& mf,
// Perform Fourier transform from `tmpSpectralField` to `tmpRealField`
#ifdef AMREX_USE_GPU
- // Add cuFFT-specific code ; make sure that this is done on the same
+ // Perform Fast Fourier Transform on GPU using cuFFT.
+ // make sure that this is done on the same
// GPU stream as the above copy
cufftResult result;
cudaStream_t stream = amrex::Gpu::Device::cudaStream();