aboutsummaryrefslogtreecommitdiff
path: root/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/SpectralBaseAlgorithm.cpp
diff options
context:
space:
mode:
authorGravatar Edoardo Zoni <59625522+EZoni@users.noreply.github.com> 2021-07-15 15:39:34 -0700
committerGravatar GitHub <noreply@github.com> 2021-07-15 15:39:34 -0700
commit300c1659c4bcdae104f828c01de8873743f73d94 (patch)
tree5a220db0c8c74c8d5689955a0d095d6c1cd9d057 /Source/FieldSolver/SpectralSolver/SpectralAlgorithms/SpectralBaseAlgorithm.cpp
parent730e9f416363a6f31a19a06c8c0654aa79b928ce (diff)
downloadWarpX-300c1659c4bcdae104f828c01de8873743f73d94.tar.gz
WarpX-300c1659c4bcdae104f828c01de8873743f73d94.tar.zst
WarpX-300c1659c4bcdae104f828c01de8873743f73d94.zip
Spectral Index: Replace `struct`s of `enum` with Class (#2062)
* Add New Spectral Index Class * Cleaning * Use New Spectral Index Class in PML * Cleaning * Reuse Available Data for divE * Allocate Rho Data Only when Necessary * Cleaning * Fix Bug in RZ Geometry * Revert Commits for Allocation of Rho Data * Cleaning * Update Forward Declaration Header * Do Not Include Unnecessary Header Files * Doxygen * Do Not Use Separate div() Cleaning Flags * SpectralFieldIndex: Add Missing param to Doxygen * Remove Unused getRequiredNumberOfFields
Diffstat (limited to 'Source/FieldSolver/SpectralSolver/SpectralAlgorithms/SpectralBaseAlgorithm.cpp')
-rw-r--r--Source/FieldSolver/SpectralSolver/SpectralAlgorithms/SpectralBaseAlgorithm.cpp22
1 files changed, 12 insertions, 10 deletions
diff --git a/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/SpectralBaseAlgorithm.cpp b/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/SpectralBaseAlgorithm.cpp
index 4445705cb..92e3fcf94 100644
--- a/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/SpectralBaseAlgorithm.cpp
+++ b/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/SpectralBaseAlgorithm.cpp
@@ -29,11 +29,13 @@ using namespace amrex;
*/
SpectralBaseAlgorithm::SpectralBaseAlgorithm(const SpectralKSpace& spectral_kspace,
const amrex::DistributionMapping& dm,
+ const SpectralFieldIndex& spectral_index,
const int norder_x, const int norder_y,
const int norder_z, const bool nodal,
const amrex::IntVect& fill_guards):
- // Compute and assign the modified k vectors
m_fill_guards(fill_guards),
+ m_spectral_index(spectral_index),
+ // Compute and assign the modified k vectors
modified_kx_vec(spectral_kspace.getModifiedKComponent(dm,0,norder_x,nodal)),
#if (AMREX_SPACEDIM==3)
modified_ky_vec(spectral_kspace.getModifiedKComponent(dm,1,norder_y,nodal)),
@@ -57,12 +59,12 @@ SpectralBaseAlgorithm::ComputeSpectralDivE (
const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
amrex::MultiFab& divE )
{
- using Idx = SpectralFieldIndex;
+ const SpectralFieldIndex& Idx = m_spectral_index;
// Forward Fourier transform of E
- field_data.ForwardTransform(lev, *Efield[0], Idx::Ex, 0 );
- field_data.ForwardTransform(lev, *Efield[1], Idx::Ey, 0 );
- field_data.ForwardTransform(lev, *Efield[2], Idx::Ez, 0 );
+ field_data.ForwardTransform(lev, *Efield[0], Idx.Ex, 0 );
+ field_data.ForwardTransform(lev, *Efield[1], Idx.Ey, 0 );
+ field_data.ForwardTransform(lev, *Efield[2], Idx.Ez, 0 );
const amrex::IntVect& fill_guards = m_fill_guards;
@@ -85,9 +87,9 @@ SpectralBaseAlgorithm::ComputeSpectralDivE (
[=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
{
// Shortcuts for the components of E
- const Complex Ex = fields(i,j,k,Idx::Ex);
- const Complex Ey = fields(i,j,k,Idx::Ey);
- const Complex Ez = fields(i,j,k,Idx::Ez);
+ const Complex Ex = fields(i,j,k,Idx.Ex);
+ const Complex Ey = fields(i,j,k,Idx.Ey);
+ const Complex Ez = fields(i,j,k,Idx.Ez);
// k vector values
const Real kx = modified_kx_arr[i];
#if (AMREX_SPACEDIM==3)
@@ -100,10 +102,10 @@ SpectralBaseAlgorithm::ComputeSpectralDivE (
const Complex I = Complex{0,1};
// div(E) in Fourier space
- fields(i,j,k,Idx::divE) = I*(kx*Ex+ky*Ey+kz*Ez);
+ fields(i,j,k,Idx.divE) = I*(kx*Ex+ky*Ey+kz*Ez);
});
}
// Backward Fourier transform
- field_data.BackwardTransform(lev, divE, Idx::divE, 0, fill_guards);
+ field_data.BackwardTransform(lev, divE, Idx.divE, 0, fill_guards);
}