aboutsummaryrefslogtreecommitdiff
path: root/Source/FieldSolver/SpectralSolver
diff options
context:
space:
mode:
authorGravatar Remi Lehe <remi.lehe@normalesup.org> 2019-04-19 10:43:18 -0700
committerGravatar Remi Lehe <remi.lehe@normalesup.org> 2019-04-23 12:43:53 -0700
commit551e934fdee50f2321076b0dd1882a74cc92fb30 (patch)
tree4cf460bd13bc6061cb51da902df15e527858b645 /Source/FieldSolver/SpectralSolver
parent84cf3b88ad82fa25183f9281136fe154eaa00992 (diff)
downloadWarpX-551e934fdee50f2321076b0dd1882a74cc92fb30.tar.gz
WarpX-551e934fdee50f2321076b0dd1882a74cc92fb30.tar.zst
WarpX-551e934fdee50f2321076b0dd1882a74cc92fb30.zip
Link spectral solver to the rest of the code
Diffstat (limited to 'Source/FieldSolver/SpectralSolver')
-rw-r--r--Source/FieldSolver/SpectralSolver/PsatdAlgorithm.H10
-rw-r--r--Source/FieldSolver/SpectralSolver/SpectralFieldData.H2
-rw-r--r--Source/FieldSolver/SpectralSolver/SpectralKSpace.H7
-rw-r--r--Source/FieldSolver/SpectralSolver/SpectralKSpace.cpp5
-rw-r--r--Source/FieldSolver/SpectralSolver/SpectralSolver.H40
5 files changed, 35 insertions, 29 deletions
diff --git a/Source/FieldSolver/SpectralSolver/PsatdAlgorithm.H b/Source/FieldSolver/SpectralSolver/PsatdAlgorithm.H
index 9fbdc7073..375fa48bb 100644
--- a/Source/FieldSolver/SpectralSolver/PsatdAlgorithm.H
+++ b/Source/FieldSolver/SpectralSolver/PsatdAlgorithm.H
@@ -11,10 +11,12 @@ class PsatdAlgorithm
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 amrex::Real dt);
+ PsatdAlgorithm( const SpectralKSpace& spectral_kspace,
+ const amrex::DistributionMapping& dm,
+ const int norder_x, const int norder_y,
+ const int norder_z, const amrex::Real dt);
+ PsatdAlgorithm() = default; // Default constructor
+ PsatdAlgorithm& operator=(PsatdAlgorithm&& algorithm) = default; // Default move assignment
void pushSpectralFields( SpectralFieldData& f ) const;
private:
diff --git a/Source/FieldSolver/SpectralSolver/SpectralFieldData.H b/Source/FieldSolver/SpectralSolver/SpectralFieldData.H
index 9f5b85a2b..a74bfe22b 100644
--- a/Source/FieldSolver/SpectralSolver/SpectralFieldData.H
+++ b/Source/FieldSolver/SpectralSolver/SpectralFieldData.H
@@ -30,6 +30,8 @@ class SpectralFieldData
SpectralFieldData( const amrex::BoxArray& realspace_ba,
const SpectralKSpace& k_space,
const amrex::DistributionMapping& dm );
+ 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 );
diff --git a/Source/FieldSolver/SpectralSolver/SpectralKSpace.H b/Source/FieldSolver/SpectralSolver/SpectralKSpace.H
index f61cffe14..a790be94f 100644
--- a/Source/FieldSolver/SpectralSolver/SpectralKSpace.H
+++ b/Source/FieldSolver/SpectralSolver/SpectralKSpace.H
@@ -13,15 +13,16 @@ class SpectralKSpace
{
public:
SpectralKSpace( const amrex::BoxArray& realspace_ba,
- const amrex::DistributionMapping& dm, const amrex::Real* dx );
+ const amrex::DistributionMapping& dm,
+ const amrex::Array<amrex::Real,3> dx );
amrex::BoxArray spectralspace_ba;
SpectralKVector kx_vec, ky_vec, kz_vec;
- const amrex::Real* dx;
+ amrex::Array<amrex::Real,3> dx;
};
void
AllocateAndFillKvector( amrex::Gpu::ManagedVector<amrex::Real>& k,
- const amrex::Box& bx, const amrex::Real* dx, const int i_dim );
+ const amrex::Box& bx, const amrex::Array<amrex::Real,3> dx, const int i_dim );
void
ComputeModifiedKVector( amrex::Gpu::ManagedVector<amrex::Real>& modified_k,
diff --git a/Source/FieldSolver/SpectralSolver/SpectralKSpace.cpp b/Source/FieldSolver/SpectralSolver/SpectralKSpace.cpp
index 0bd57c7ea..2b1e7ee33 100644
--- a/Source/FieldSolver/SpectralSolver/SpectralKSpace.cpp
+++ b/Source/FieldSolver/SpectralSolver/SpectralKSpace.cpp
@@ -6,7 +6,7 @@ using namespace Gpu;
SpectralKSpace::SpectralKSpace( const BoxArray& realspace_ba,
const DistributionMapping& dm,
- const Real* realspace_dx )
+ const Array<Real,3> realspace_dx )
{
// Create the box array that corresponds to spectral space
BoxList spectral_bl; // Create empty box list
@@ -37,7 +37,8 @@ SpectralKSpace::SpectralKSpace( const BoxArray& realspace_ba,
}
void
-AllocateAndFillKvector( ManagedVector<Real>& k, const Box& bx, const Real* dx, const int i_dim )
+AllocateAndFillKvector( ManagedVector<Real>& k, const Box& bx,
+ const Array<Real,3> dx, const int i_dim )
{
// Alllocate k to the right size
int N = bx.length( i_dim );
diff --git a/Source/FieldSolver/SpectralSolver/SpectralSolver.H b/Source/FieldSolver/SpectralSolver/SpectralSolver.H
index cde1fccd2..521e558ba 100644
--- a/Source/FieldSolver/SpectralSolver/SpectralSolver.H
+++ b/Source/FieldSolver/SpectralSolver/SpectralSolver.H
@@ -1,24 +1,39 @@
#ifndef WARPX_SPECTRAL_SOLVER_H_
#define WARPX_SPECTRAL_SOLVER_H_
+#include <SpectralKSpace.H>
+#include <PsatdAlgorithm.H>
+#include <SpectralFieldData.H>
+
/* \brief
* TODO
*/
class SpectralSolver
{
private:
- SpectralKSpace k_space; // Contains size of each box in spectral space,
- // and corresponding values of the k vectors
SpectralFieldData field_data; // Store field in spectral space
// and perform the Fourier transforms
PsatdAlgorithm algorithm; // Contains Psatd coefficients
// and field update equation
public:
- SpectralSolver();
- void ForwardTransform( const MultiFab& mf, const int field_index ){
+ SpectralSolver( const amrex::BoxArray& realspace_ba,
+ const amrex::DistributionMapping& dm,
+ const int norder_x, const int norder_y, const int norder_z,
+ const amrex::Array<amrex::Real,3> dx, const amrex::Real dt )
+ {
+ // Initialize all structures using the same distribution mapping dm
+ // - Initialize k space (Contains size of each box in spectral space,
+ // and corresponding values of the k vectors)
+ const SpectralKSpace k_space = SpectralKSpace( realspace_ba, dm, dx );
+ // - Initialize algorithm (coefficients) over this space
+ algorithm = PsatdAlgorithm( k_space, dm, norder_x, norder_y, norder_z, dt );
+ // - 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 BackwardTransform( MultiFab& mf, const int field_index ){
+ void BackwardTransform( amrex::MultiFab& mf, const int field_index ){
field_data.BackwardTransform( mf, field_index );
};
void pushSpectralFields(){
@@ -26,19 +41,4 @@ class SpectralSolver
};
};
-SpectralSolver::SpectralSolver( const SpectralKSpace& realspace_ba,
- const DistributionMapping& dm,
- const int norder_x, const int norder_y,
- const int norder_z, const Real* dx,
- const Real dt )
-{
- // Initialize all structures using the same distribution mapping dm
- // - Initialize k space values
- k_space = SpectralKSpace( realspace_ba, dm, dx );
- // - Initialize algorithm (coefficients) over this space
- algorithm = PsatdAlgorithm( k_space, dm, norder_x, norder_y, norder_z );
- // - Initialize arrays for fields in Fourier space + FFT plans
- field_data = SpectralFieldData( realspace_ba, k_space, dm );
-};
-
#endif // WARPX_SPECTRAL_SOLVER_H_