aboutsummaryrefslogtreecommitdiff
path: root/Source/FieldSolver
diff options
context:
space:
mode:
Diffstat (limited to 'Source/FieldSolver')
-rw-r--r--Source/FieldSolver/SpectralSolver/SpectralAlgorithms/Make.package1
-rw-r--r--Source/FieldSolver/SpectralSolver/SpectralAlgorithms/PsatdAlgorithm.H20
-rw-r--r--Source/FieldSolver/SpectralSolver/SpectralAlgorithms/PsatdAlgorithm.cpp11
-rw-r--r--Source/FieldSolver/SpectralSolver/SpectralAlgorithms/SpectralBaseAlgorithm.H44
-rw-r--r--Source/FieldSolver/SpectralSolver/SpectralSolver.H13
5 files changed, 62 insertions, 27 deletions
diff --git a/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/Make.package b/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/Make.package
index 1caf69397..c62c21f44 100644
--- a/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/Make.package
+++ b/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/Make.package
@@ -1,3 +1,4 @@
+CEXE_headers += SpectralBaseAlgorithm.H
CEXE_headers += PsatdAlgorithm.H
CEXE_sources += PsatdAlgorithm.cpp
diff --git a/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/PsatdAlgorithm.H b/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/PsatdAlgorithm.H
index acefcc466..0487e5226 100644
--- a/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/PsatdAlgorithm.H
+++ b/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/PsatdAlgorithm.H
@@ -1,31 +1,23 @@
#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
* and stores the coefficients of the corresponding update equation.
*/
-class PsatdAlgorithm
+class PsatdAlgorithm : public SpectralBaseAlgorithm
{
- using SpectralCoefficients = amrex::FabArray< amrex::BaseFab <amrex::Real> >;
-
public:
PsatdAlgorithm(const SpectralKSpace& spectral_kspace,
const amrex::DistributionMapping& dm,
const int norder_x, const int norder_y,
- const int norder_z, const bool nodal, const amrex::Real dt);
- PsatdAlgorithm() = default; // Default constructor
- PsatdAlgorithm& operator=(PsatdAlgorithm&& algorithm) = default;
- void pushSpectralFields(SpectralFieldData& f) const;
+ const int norder_z, const bool nodal,
+ const amrex::Real dt);
+ // Redefine update equation from base class
+ virtual void pushSpectralFields(SpectralFieldData& f) const override final;
private:
- // Modified finite-order vectors
- KVectorComponent modified_kx_vec, modified_kz_vec;
-#if (AMREX_SPACEDIM==3)
- KVectorComponent modified_ky_vec;
-#endif
SpectralCoefficients C_coef, S_ck_coef, X1_coef, X2_coef, X3_coef;
};
diff --git a/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/PsatdAlgorithm.cpp b/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/PsatdAlgorithm.cpp
index ada7506c3..37892d35a 100644
--- a/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/PsatdAlgorithm.cpp
+++ b/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/PsatdAlgorithm.cpp
@@ -9,14 +9,9 @@ PsatdAlgorithm::PsatdAlgorithm(const SpectralKSpace& spectral_kspace,
const DistributionMapping& dm,
const int norder_x, const int norder_y,
const int norder_z, const bool nodal, const Real dt)
-// 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)),
- modified_kz_vec(spectral_kspace.getModifiedKComponent(dm,2,norder_z,nodal))
-#else
- modified_kz_vec(spectral_kspace.getModifiedKComponent(dm,1,norder_z,nodal))
-#endif
+ // Initialize members of base class
+ : SpectralBaseAlgorithm( spectral_kspace, dm,
+ norder_x, norder_y, norder_z, nodal )
{
const BoxArray& ba = spectral_kspace.spectralspace_ba;
diff --git a/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/SpectralBaseAlgorithm.H b/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/SpectralBaseAlgorithm.H
new file mode 100644
index 000000000..18d26e0c8
--- /dev/null
+++ b/Source/FieldSolver/SpectralSolver/SpectralAlgorithms/SpectralBaseAlgorithm.H
@@ -0,0 +1,44 @@
+#ifndef WARPX_SPECTRAL_BASE_ALGORITHM_H_
+#define WARPX_SPECTRAL_BASE_ALGORITHM_H_
+
+#include <SpectralKSpace.H>
+#include <SpectralFieldData.H>
+
+/* \brief Class that updates the field in spectral space
+ * and stores the coefficients of the corresponding update equation.
+ * TODO: Mention base class
+ */
+class SpectralBaseAlgorithm
+{
+ public:
+ // Member function that updates the fields in spectral space ;
+ // meant to be overridden in subclasses
+ virtual void pushSpectralFields(SpectralFieldData& f) const = 0;
+
+ protected: // Meant to be used in the subclasses
+
+ using SpectralCoefficients = amrex::FabArray< amrex::BaseFab <amrex::Real> >;
+
+ // Constructor
+ SpectralBaseAlgorithm(const SpectralKSpace& spectral_kspace,
+ const amrex::DistributionMapping& dm,
+ const int norder_x, const int norder_y,
+ const int norder_z, const bool nodal)
+ // 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)),
+ modified_kz_vec(spectral_kspace.getModifiedKComponent(dm,2,norder_z,nodal))
+#else
+ modified_kz_vec(spectral_kspace.getModifiedKComponent(dm,1,norder_z,nodal))
+#endif
+ {};
+
+ // Modified finite-order vectors
+ KVectorComponent modified_kx_vec, modified_kz_vec;
+#if (AMREX_SPACEDIM==3)
+ KVectorComponent modified_ky_vec;
+#endif
+};
+
+#endif // WARPX_SPECTRAL_BASE_ALGORITHM_H_
diff --git a/Source/FieldSolver/SpectralSolver/SpectralSolver.H b/Source/FieldSolver/SpectralSolver/SpectralSolver.H
index 7444452af..5266025b6 100644
--- a/Source/FieldSolver/SpectralSolver/SpectralSolver.H
+++ b/Source/FieldSolver/SpectralSolver/SpectralSolver.H
@@ -31,8 +31,8 @@ class SpectralSolver
// as well as the value of the corresponding k coordinates)
const SpectralKSpace k_space= SpectralKSpace(realspace_ba, dm, dx);
// - Initialize the algorithm (coefficients) over this space
- algorithm = PsatdAlgorithm( k_space, dm, norder_x, norder_y,
- norder_z, nodal, dt );
+ algorithm = std::unique_ptr<PsatdAlgorithm>( new PsatdAlgorithm(
+ k_space, dm, norder_x, norder_y, norder_z, nodal, dt ) );
// - Initialize arrays for fields in Fourier space + FFT plans
field_data = SpectralFieldData( realspace_ba, k_space, dm );
};
@@ -59,14 +59,17 @@ class SpectralSolver
/* \brief Update the fields in spectral space, over one timestep */
void pushSpectralFields(){
BL_PROFILE("SpectralSolver::pushSpectralFields");
- algorithm.pushSpectralFields( field_data );
+ algorithm->pushSpectralFields( field_data );
};
private:
SpectralFieldData field_data; // Store field in spectral space
// and perform the Fourier transforms
- PsatdAlgorithm algorithm; // Contains Psatd coefficients
- // and field update equation
+ std::unique_ptr<SpectralBaseAlgorithm> algorithm;
+ // Defines field update equation in spectral space,
+ // and the associated coefficients.
+ // SpectralBaseAlgorithm is a base class ; this pointer is meant
+ // to point an instance of a *sub-class* defining a specific algorithm
};
#endif // WARPX_SPECTRAL_SOLVER_H_