diff options
author | 2021-07-15 15:39:34 -0700 | |
---|---|---|
committer | 2021-07-15 15:39:34 -0700 | |
commit | 300c1659c4bcdae104f828c01de8873743f73d94 (patch) | |
tree | 5a220db0c8c74c8d5689955a0d095d6c1cd9d057 /Source/FieldSolver/SpectralSolver/SpectralAlgorithms/PsatdAlgorithmRZ.cpp | |
parent | 730e9f416363a6f31a19a06c8c0654aa79b928ce (diff) | |
download | WarpX-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/PsatdAlgorithmRZ.cpp')
-rw-r--r-- | Source/FieldSolver/SpectralSolver/SpectralAlgorithms/PsatdAlgorithmRZ.cpp | 58 |
1 files changed, 30 insertions, 28 deletions
diff --git a/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/PsatdAlgorithmRZ.cpp b/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/PsatdAlgorithmRZ.cpp index cb055a27b..6a3c86bf6 100644 --- a/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/PsatdAlgorithmRZ.cpp +++ b/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/PsatdAlgorithmRZ.cpp @@ -17,12 +17,13 @@ using amrex::operator""_rt; /* \brief Initialize coefficients for the update equation */ PsatdAlgorithmRZ::PsatdAlgorithmRZ (SpectralKSpaceRZ const & spectral_kspace, amrex::DistributionMapping const & dm, + const SpectralFieldIndex& spectral_index, int const n_rz_azimuthal_modes, int const norder_z, bool const nodal, amrex::Real const dt, bool const update_with_rho) // Initialize members of base class - : SpectralBaseAlgorithmRZ(spectral_kspace, dm, - norder_z, nodal), + : SpectralBaseAlgorithmRZ(spectral_kspace, dm, spectral_index, norder_z, nodal), + m_spectral_index(spectral_index), m_dt(dt), m_update_with_rho(update_with_rho) { @@ -53,6 +54,8 @@ PsatdAlgorithmRZ::pushSpectralFields(SpectralFieldDataRZ & f) coefficients_initialized = true; } + const SpectralFieldIndex& Idx = m_spectral_index; + // Loop over boxes for (amrex::MFIter mfi(f.fields); mfi.isValid(); ++mfi){ @@ -82,18 +85,17 @@ PsatdAlgorithmRZ::pushSpectralFields(SpectralFieldDataRZ & f) { // All of the fields of each mode are grouped together - using Idx = SpectralFieldIndex; - auto const Ep_m = Idx::Ex + Idx::n_fields*mode; - auto const Em_m = Idx::Ey + Idx::n_fields*mode; - auto const Ez_m = Idx::Ez + Idx::n_fields*mode; - auto const Bp_m = Idx::Bx + Idx::n_fields*mode; - auto const Bm_m = Idx::By + Idx::n_fields*mode; - auto const Bz_m = Idx::Bz + Idx::n_fields*mode; - auto const Jp_m = Idx::Jx + Idx::n_fields*mode; - auto const Jm_m = Idx::Jy + Idx::n_fields*mode; - auto const Jz_m = Idx::Jz + Idx::n_fields*mode; - auto const rho_old_m = Idx::rho_old + Idx::n_fields*mode; - auto const rho_new_m = Idx::rho_new + Idx::n_fields*mode; + auto const Ep_m = Idx.Ex + Idx.n_fields*mode; + auto const Em_m = Idx.Ey + Idx.n_fields*mode; + auto const Ez_m = Idx.Ez + Idx.n_fields*mode; + auto const Bp_m = Idx.Bx + Idx.n_fields*mode; + auto const Bm_m = Idx.By + Idx.n_fields*mode; + auto const Bz_m = Idx.Bz + Idx.n_fields*mode; + auto const Jp_m = Idx.Jx + Idx.n_fields*mode; + auto const Jm_m = Idx.Jy + Idx.n_fields*mode; + auto const Jz_m = Idx.Jz + Idx.n_fields*mode; + auto const rho_old_m = Idx.rho_old + Idx.n_fields*mode; + auto const rho_new_m = Idx.rho_new + Idx.n_fields*mode; // Record old values of the fields to be updated Complex const Ep_old = fields(i,j,k,Ep_m); @@ -223,15 +225,15 @@ PsatdAlgorithmRZ::CurrentCorrection (const int lev, // Profiling WARPX_PROFILE( "PsatdAlgorithmRZ::CurrentCorrection" ); - using Idx = SpectralFieldIndex; + const SpectralFieldIndex& Idx = m_spectral_index; // Forward Fourier transform of J and rho field_data.ForwardTransform( lev, - *current[0], Idx::Jx, - *current[1], Idx::Jy); - field_data.ForwardTransform( lev, *current[2], Idx::Jz, 0); - field_data.ForwardTransform( lev, *rho, Idx::rho_old, 0 ); - field_data.ForwardTransform( lev, *rho, Idx::rho_new, 1 ); + *current[0], Idx.Jx, + *current[1], Idx.Jy); + field_data.ForwardTransform( lev, *current[2], Idx.Jz, 0); + field_data.ForwardTransform( lev, *rho, Idx.rho_old, 0 ); + field_data.ForwardTransform( lev, *rho, Idx.rho_new, 1 ); // Loop over boxes for (amrex::MFIter mfi(field_data.fields); mfi.isValid(); ++mfi){ @@ -256,11 +258,11 @@ PsatdAlgorithmRZ::CurrentCorrection (const int lev, [=] AMREX_GPU_DEVICE(int i, int j, int k, int mode) noexcept { // All of the fields of each mode are grouped together - auto const Jp_m = Idx::Jx + Idx::n_fields*mode; - auto const Jm_m = Idx::Jy + Idx::n_fields*mode; - auto const Jz_m = Idx::Jz + Idx::n_fields*mode; - auto const rho_old_m = Idx::rho_old + Idx::n_fields*mode; - auto const rho_new_m = Idx::rho_new + Idx::n_fields*mode; + auto const Jp_m = Idx.Jx + Idx.n_fields*mode; + auto const Jm_m = Idx.Jy + Idx.n_fields*mode; + auto const Jz_m = Idx.Jz + Idx.n_fields*mode; + auto const rho_old_m = Idx.rho_old + Idx.n_fields*mode; + auto const rho_new_m = Idx.rho_new + Idx.n_fields*mode; // Shortcuts for the values of J and rho Complex const Jp = fields(i,j,k,Jp_m); @@ -292,10 +294,10 @@ PsatdAlgorithmRZ::CurrentCorrection (const int lev, // Backward Fourier transform of J field_data.BackwardTransform( lev, - *current[0], Idx::Jx, - *current[1], Idx::Jy); + *current[0], Idx.Jx, + *current[1], Idx.Jy); field_data.BackwardTransform( lev, - *current[2], Idx::Jz, 0 ); + *current[2], Idx.Jz, 0 ); } void |