aboutsummaryrefslogtreecommitdiff
path: root/Source/FieldSolver/SpectralSolver/SpectralSolverRZ.H
diff options
context:
space:
mode:
authorGravatar Remi Lehe <remi.lehe@normalesup.org> 2020-09-10 13:36:50 -0700
committerGravatar GitHub <noreply@github.com> 2020-09-10 13:36:50 -0700
commitebde54faa8bcae2f1b37b81270a8d0a64bf58a98 (patch)
tree73e95be4a0d9b3bc2ad0014cbcf28487c84f829d /Source/FieldSolver/SpectralSolver/SpectralSolverRZ.H
parent614dc2962f9b5b576b4f734532c969b89f1316c0 (diff)
downloadWarpX-ebde54faa8bcae2f1b37b81270a8d0a64bf58a98.tar.gz
WarpX-ebde54faa8bcae2f1b37b81270a8d0a64bf58a98.tar.zst
WarpX-ebde54faa8bcae2f1b37b81270a8d0a64bf58a98.zip
SpectralSolverRZ.H: Move Implementations to Source File (#1290)
* Moved routines from SpectralSolverRZ.H to .cpp * Correct missing Vay deposition routine * Apply suggestions from code review Co-authored-by: Axel Huebl <axel.huebl@plasma.ninja> * Apply suggestions from code review * Apply suggestions from code review Co-authored-by: Edoardo Zoni <59625522+EZoni@users.noreply.github.com> Co-authored-by: Dave Grote <dpgrote@lbl.gov> Co-authored-by: Axel Huebl <axel.huebl@plasma.ninja> Co-authored-by: Edoardo Zoni <59625522+EZoni@users.noreply.github.com>
Diffstat (limited to 'Source/FieldSolver/SpectralSolver/SpectralSolverRZ.H')
-rw-r--r--Source/FieldSolver/SpectralSolver/SpectralSolverRZ.H50
1 files changed, 12 insertions, 38 deletions
diff --git a/Source/FieldSolver/SpectralSolver/SpectralSolverRZ.H b/Source/FieldSolver/SpectralSolver/SpectralSolverRZ.H
index 122e980ac..5bb6d7f90 100644
--- a/Source/FieldSolver/SpectralSolver/SpectralSolverRZ.H
+++ b/Source/FieldSolver/SpectralSolver/SpectralSolverRZ.H
@@ -36,49 +36,30 @@ class SpectralSolverRZ
/* \brief Transform the component `i_comp` of MultiFab `field_mf`
* to spectral space, and store the corresponding result internally
* (in the spectral field specified by `field_index`) */
+
void ForwardTransform (amrex::MultiFab const & field_mf,
int const field_index,
- int const i_comp=0 ) {
- BL_PROFILE("SpectralSolverRZ::ForwardTransform");
- field_data.ForwardTransform(field_mf, field_index, i_comp);
- }
+ int const i_comp=0);
/* \brief Transform the two MultiFabs `field_mf1` and `field_mf2`
* to spectral space, and store the corresponding results internally
* (in the spectral field specified by `field_index1` and `field_index2`) */
void ForwardTransform (amrex::MultiFab const & field_mf1, int const field_index1,
- amrex::MultiFab const & field_mf2, int const field_index2) {
- BL_PROFILE("SpectralSolverRZ::ForwardTransform");
- field_data.ForwardTransform(field_mf1, field_index1,
- field_mf2, field_index2);
- }
+ amrex::MultiFab const & field_mf2, int const field_index2);
/* \brief Transform spectral field specified by `field_index` back to
* real space, and store it in the component `i_comp` of `field_mf` */
void BackwardTransform (amrex::MultiFab& field_mf,
int const field_index,
- int const i_comp=0) {
- BL_PROFILE("SpectralSolverRZ::BackwardTransform");
- field_data.BackwardTransform(field_mf, field_index, i_comp);
- }
+ int const i_comp=0);
/* \brief Transform spectral fields specified by `field_index1` and `field_index2`
* back to real space, and store it in `field_mf1` and `field_mf2`*/
void BackwardTransform (amrex::MultiFab& field_mf1, int const field_index1,
- amrex::MultiFab& field_mf2, int const field_index2) {
- BL_PROFILE("SpectralSolverRZ::BackwardTransform");
- field_data.BackwardTransform(field_mf1, field_index1,
- field_mf2, field_index2);
- }
+ amrex::MultiFab& field_mf2, int const field_index2);
/* \brief Update the fields in spectral space, over one timestep */
- void pushSpectralFields () {
- BL_PROFILE("SpectralSolverRZ::pushSpectralFields");
- // Virtual function: the actual function used here depends
- // on the sub-class of `SpectralBaseAlgorithm` that was
- // initialized in the constructor of `SpectralSolverRZ`
- algorithm->pushSpectralFields(field_data);
- }
+ void pushSpectralFields ();
/* \brief Initialize K space filtering arrays */
void InitFilter (amrex::IntVect const & filter_npass_each_dir,
@@ -102,10 +83,8 @@ class SpectralSolverRZ
* \brief Public interface to call the member function ComputeSpectralDivE
* of the base class SpectralBaseAlgorithmRZ from objects of class SpectralSolverRZ
*/
- void ComputeSpectralDivE ( const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
- amrex::MultiFab& divE ) {
- algorithm->ComputeSpectralDivE( field_data, Efield, divE );
- }
+ void ComputeSpectralDivE (const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
+ amrex::MultiFab& divE);
/**
* \brief Public interface to call the virtual function \c CurrentCorrection,
@@ -117,10 +96,8 @@ class SpectralSolverRZ
* storing the three components of the current density
* \param[in] rho unique pointer to MultiFab storing the charge density
*/
- void CurrentCorrection ( std::array<std::unique_ptr<amrex::MultiFab>,3>& current,
- const std::unique_ptr<amrex::MultiFab>& rho ) {
- algorithm->CurrentCorrection( field_data, current, rho );
- }
+ void CurrentCorrection (std::array<std::unique_ptr<amrex::MultiFab>,3>& current,
+ const std::unique_ptr<amrex::MultiFab>& rho);
/**
* \brief Public interface to call the virtual function \c VayDeposition,
@@ -129,12 +106,9 @@ class SpectralSolverRZ
* unique pointer \c algorithm.
*
* \param[in,out] current Array of unique pointers to \c MultiFab storing
- * the three components of the current density
+ * the three components of the current density
*/
- void VayDeposition (std::array<std::unique_ptr<amrex::MultiFab>,3>& current)
- {
- algorithm->VayDeposition(field_data, current);
- }
+ void VayDeposition (std::array<std::unique_ptr<amrex::MultiFab>,3>& current);
private: