aboutsummaryrefslogtreecommitdiff
path: root/Source/FieldSolver/SpectralSolver
diff options
context:
space:
mode:
authorGravatar Remi Lehe <remi.lehe@normalesup.org> 2019-04-19 20:03:54 -0700
committerGravatar Remi Lehe <remi.lehe@normalesup.org> 2019-04-23 12:43:53 -0700
commit8bffe04d8eb683230a6b9919438adf9239c87fec (patch)
tree63857f7bb8346b6e04bc8bc96f99b7d72157506c /Source/FieldSolver/SpectralSolver
parent30fb9ade3cfe768b5855e3efbe27a5eb1609fee7 (diff)
downloadWarpX-8bffe04d8eb683230a6b9919438adf9239c87fec.tar.gz
WarpX-8bffe04d8eb683230a6b9919438adf9239c87fec.tar.zst
WarpX-8bffe04d8eb683230a6b9919438adf9239c87fec.zip
Transform rho_old and rho_new
Diffstat (limited to 'Source/FieldSolver/SpectralSolver')
-rw-r--r--Source/FieldSolver/SpectralSolver/SpectralFieldData.H6
-rw-r--r--Source/FieldSolver/SpectralSolver/SpectralFieldData.cpp10
-rw-r--r--Source/FieldSolver/SpectralSolver/SpectralSolver.H10
3 files changed, 16 insertions, 10 deletions
diff --git a/Source/FieldSolver/SpectralSolver/SpectralFieldData.H b/Source/FieldSolver/SpectralSolver/SpectralFieldData.H
index a74bfe22b..2abe81889 100644
--- a/Source/FieldSolver/SpectralSolver/SpectralFieldData.H
+++ b/Source/FieldSolver/SpectralSolver/SpectralFieldData.H
@@ -33,8 +33,10 @@ class SpectralFieldData
SpectralFieldData() = default; // Default constructor
SpectralFieldData& operator=(SpectralFieldData&& field_data) = default; // Default move assignment
~SpectralFieldData();
- void ForwardTransform( const amrex::MultiFab& mf, const int field_index );
- void BackwardTransform( amrex::MultiFab& mf, const int field_index );
+ void ForwardTransform( const amrex::MultiFab& mf,
+ const int field_index, const int i_comp );
+ void BackwardTransform( amrex::MultiFab& mf,
+ const int field_index, const int i_comp );
private:
SpectralField Ex, Ey, Ez, Bx, By, Bz, Jx, Jy, Jz, rho_old, rho_new;
diff --git a/Source/FieldSolver/SpectralSolver/SpectralFieldData.cpp b/Source/FieldSolver/SpectralSolver/SpectralFieldData.cpp
index 66e1e3470..ee6dbff6a 100644
--- a/Source/FieldSolver/SpectralSolver/SpectralFieldData.cpp
+++ b/Source/FieldSolver/SpectralSolver/SpectralFieldData.cpp
@@ -73,7 +73,8 @@ SpectralFieldData::~SpectralFieldData()
* Example: ForwardTransform( Efield_cp[0], SpectralFieldIndex::Ex )
*/
void
-SpectralFieldData::ForwardTransform( const MultiFab& mf, const int field_index )
+SpectralFieldData::ForwardTransform( const MultiFab& mf,
+ const int field_index, const int i_comp )
{
// Loop over boxes
for ( MFIter mfi(mf); mfi.isValid(); ++mfi ){
@@ -91,7 +92,7 @@ SpectralFieldData::ForwardTransform( const MultiFab& mf, const int field_index )
Array4<Complex> tmp_arr = tmpRealField[mfi].array();
ParallelFor( realspace_bx,
[=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
- tmp_arr(i,j,k) = mf_arr(i,j,k);
+ tmp_arr(i,j,k) = mf_arr(i,j,k,i_comp);
});
}
@@ -122,7 +123,8 @@ SpectralFieldData::ForwardTransform( const MultiFab& mf, const int field_index )
/* TODO: Documentation
*/
void
-SpectralFieldData::BackwardTransform( MultiFab& mf, const int field_index )
+SpectralFieldData::BackwardTransform( MultiFab& mf,
+ const int field_index, const int i_comp )
{
// Loop over boxes
for ( MFIter mfi(mf); mfi.isValid(); ++mfi ){
@@ -160,7 +162,7 @@ SpectralFieldData::BackwardTransform( MultiFab& mf, const int field_index )
Array4<const Complex> tmp_arr = tmpRealField[mfi].array();
ParallelFor( realspace_bx,
[=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
- mf_arr(i,j,k) = tmp_arr(i,j,k).real();
+ mf_arr(i,j,k,i_comp) = tmp_arr(i,j,k).real();
});
}
}
diff --git a/Source/FieldSolver/SpectralSolver/SpectralSolver.H b/Source/FieldSolver/SpectralSolver/SpectralSolver.H
index 521e558ba..a471666b9 100644
--- a/Source/FieldSolver/SpectralSolver/SpectralSolver.H
+++ b/Source/FieldSolver/SpectralSolver/SpectralSolver.H
@@ -30,11 +30,13 @@ class SpectralSolver
// - Initialize arrays for fields in Fourier space + FFT plans
field_data = SpectralFieldData( realspace_ba, k_space, dm );
};
- void ForwardTransform( const amrex::MultiFab& mf, const int field_index ){
- field_data.ForwardTransform( mf, field_index );
+ void ForwardTransform( const amrex::MultiFab& mf,
+ const int field_index, const int i_comp=0 ){
+ field_data.ForwardTransform( mf, field_index, i_comp );
};
- void BackwardTransform( amrex::MultiFab& mf, const int field_index ){
- field_data.BackwardTransform( mf, field_index );
+ void BackwardTransform( amrex::MultiFab& mf,
+ const int field_index, const int i_comp=0 ){
+ field_data.BackwardTransform( mf, field_index, i_comp );
};
void pushSpectralFields(){
algorithm.pushSpectralFields( field_data );