WarpX
SpectralSolver.H
Go to the documentation of this file.
1 /* Copyright 2019-2020 Maxence Thevenet, Remi Lehe, Edoardo Zoni
2  *
3  * This file is part of WarpX.
4  *
5  * License: BSD-3-Clause-LBNL
6  */
7 #ifndef WARPX_SPECTRAL_SOLVER_H_
8 #define WARPX_SPECTRAL_SOLVER_H_
9 
11 #include "SpectralFieldData.H"
12 
13 
14 #ifdef WARPX_USE_PSATD
15 
23 {
24  public:
25  // Inline definition of the member functions of `SpectralSolver`,
26  // except the constructor (see `SpectralSolver.cpp`)
27  // The body of these functions is short, since the work is done in the
28  // underlying classes `SpectralFieldData` and `PsatdAlgorithm`
29 
30  // Constructor
31  SpectralSolver( const amrex::BoxArray& realspace_ba,
32  const amrex::DistributionMapping& dm,
33  const int norder_x, const int norder_y,
34  const int norder_z, const bool nodal,
35  const amrex::Array<amrex::Real,3>& v_galilean,
36  const amrex::RealVect dx, const amrex::Real dt,
37  const bool pml=false,
38  const bool periodic_single_box=false,
39  const bool update_with_rho=false,
40  const bool fft_do_time_averaging=false);
41 
46  void ForwardTransform( const amrex::MultiFab& mf,
47  const int field_index,
48  const int i_comp=0 );
49 
54  void BackwardTransform( amrex::MultiFab& mf,
55  const int field_index,
56  const int i_comp=0 );
57 
61  void pushSpectralFields();
62 
67  void ComputeSpectralDivE ( const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
68  amrex::MultiFab& divE ) {
69  algorithm->ComputeSpectralDivE( field_data, Efield, divE );
70  }
71 
82  void CurrentCorrection ( std::array<std::unique_ptr<amrex::MultiFab>,3>& current,
83  const std::unique_ptr<amrex::MultiFab>& rho ) {
84  algorithm->CurrentCorrection( field_data, current, rho );
85  }
86 
96  void VayDeposition (std::array<std::unique_ptr<amrex::MultiFab>,3>& current)
97  {
98  algorithm->VayDeposition(field_data, current);
99  }
100 
101  private:
102  void ReadParameters ();
103 
104  // Store field in spectral space and perform the Fourier transforms
106 
107  // Defines field update equation in spectral space and the associated coefficients.
108  // SpectralBaseAlgorithm is a base class; this pointer is meant to point
109  // to an instance of a sub-class defining a specific algorithm
110  std::unique_ptr<SpectralBaseAlgorithm> algorithm;
111 };
112 #endif // WARPX_USE_PSATD
113 #endif // WARPX_SPECTRAL_SOLVER_H_
SpectralFieldData field_data
Definition: SpectralSolver.H:105
int dx
Definition: compute_domain.py:35
SpectralSolver(const amrex::BoxArray &realspace_ba, const amrex::DistributionMapping &dm, const int norder_x, const int norder_y, const int norder_z, const bool nodal, const amrex::Array< amrex::Real, 3 > &v_galilean, const amrex::RealVect dx, const amrex::Real dt, const bool pml=false, const bool periodic_single_box=false, const bool update_with_rho=false, const bool fft_do_time_averaging=false)
Definition: SpectralSolver.cpp:34
void VayDeposition(std::array< std::unique_ptr< amrex::MultiFab >, 3 > &current)
Public interface to call the virtual function VayDeposition, declared in the base class SpectralBaseA...
Definition: SpectralSolver.H:96
void CurrentCorrection(std::array< std::unique_ptr< amrex::MultiFab >, 3 > &current, const std::unique_ptr< amrex::MultiFab > &rho)
Public interface to call the virtual function CurrentCorrection, defined in the base class SpectralBa...
Definition: SpectralSolver.H:82
void ComputeSpectralDivE(const std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Efield, amrex::MultiFab &divE)
Public interface to call the member function ComputeSpectralDivE of the base class SpectralBaseAlgori...
Definition: SpectralSolver.H:67
void pushSpectralFields()
Update the fields in spectral space, over one timestep.
Definition: SpectralSolver.cpp:102
void ForwardTransform(const amrex::MultiFab &mf, const int field_index, const int i_comp=0)
Transform the component i_comp of MultiFab mf to spectral space, and store the corresponding result i...
Definition: SpectralSolver.cpp:84
Class that stores the fields in spectral space, and performs the Fourier transforms between real spac...
Definition: SpectralFieldData.H:46
void ReadParameters()
void BackwardTransform(amrex::MultiFab &mf, const int field_index, const int i_comp=0)
Transform spectral field specified by field_index back to real space, and store it in the component i...
Definition: SpectralSolver.cpp:93
std::unique_ptr< SpectralBaseAlgorithm > algorithm
Definition: SpectralSolver.H:110
Top-level class for the electromagnetic spectral solver.
Definition: SpectralSolver.H:22