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 
10 #include "SpectralSolver_fwd.H"
11 
13 #include "SpectralFieldData.H"
14 
15 #include <AMReX_Array.H>
16 #include <AMReX_REAL.H>
17 #include <AMReX_RealVect.H>
18 
19 #include <AMReX_BaseFwd.H>
20 
21 #include <array>
22 #include <memory>
23 
24 #ifdef WARPX_USE_PSATD
25 
33 {
34  public:
35 
69  SpectralSolver (const int lev,
70  const amrex::BoxArray& realspace_ba,
71  const amrex::DistributionMapping& dm,
72  const int norder_x, const int norder_y,
73  const int norder_z, const bool nodal,
74  const amrex::IntVect& fill_guards,
75  const amrex::Array<amrex::Real,3>& v_galilean,
76  const amrex::Array<amrex::Real,3>& v_comoving,
77  const amrex::RealVect dx,
78  const amrex::Real dt,
79  const bool pml,
80  const bool periodic_single_box,
81  const bool update_with_rho,
82  const bool fft_do_time_averaging,
83  const bool J_linear_in_time,
84  const bool dive_cleaning,
85  const bool divb_cleaning);
86 
91  void ForwardTransform( const int lev,
92  const amrex::MultiFab& mf,
93  const int field_index,
94  const int i_comp=0 );
95 
100  void BackwardTransform( const int lev,
101  amrex::MultiFab& mf,
102  const int field_index,
103  const int i_comp=0 );
104 
108  void pushSpectralFields();
109 
114  void ComputeSpectralDivE ( const int lev,
115  const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
116  amrex::MultiFab& divE ) {
117  algorithm->ComputeSpectralDivE( lev, field_data, Efield, divE );
118  }
119 
130  void CurrentCorrection ( const int lev,
131  std::array<std::unique_ptr<amrex::MultiFab>,3>& current,
132  const std::unique_ptr<amrex::MultiFab>& rho ) {
133  algorithm->CurrentCorrection( lev, field_data, current, rho );
134  }
135 
145  void VayDeposition (const int lev, std::array<std::unique_ptr<amrex::MultiFab>,3>& current)
146  {
147  algorithm->VayDeposition(lev, field_data, current);
148  }
149 
157  void CopySpectralDataComp (const int src_comp, const int dest_comp)
158  {
159  // The last two arguments represent the number of components and
160  // the number of ghost cells to perform this operation
161  Copy(field_data.fields, field_data.fields, src_comp, dest_comp, 1, 0);
162  }
163 
169  void ZeroOutDataComp (const int icomp)
170  {
171  // The last argument represents the number of components to perform this operation
172  field_data.fields.setVal(0., icomp, 1);
173  }
174 
182  void ScaleDataComp (const int icomp, const amrex::Real scale_factor)
183  {
184  // The last argument represents the number of components to perform this operation
185  field_data.fields.mult(scale_factor, icomp, 1);
186  }
187 
189 
190  protected:
191 
192  amrex::IntVect m_fill_guards;
193 
194  private:
195 
196  void ReadParameters ();
197 
198  // Store field in spectral space and perform the Fourier transforms
200 
201  // Defines field update equation in spectral space and the associated coefficients.
202  // SpectralBaseAlgorithm is a base class; this pointer is meant to point
203  // to an instance of a sub-class defining a specific algorithm
204  std::unique_ptr<SpectralBaseAlgorithm> algorithm;
205 };
206 #endif // WARPX_USE_PSATD
207 #endif // WARPX_SPECTRAL_SOLVER_H_
SpectralFieldData field_data
Definition: SpectralSolver.H:199
amrex::IntVect m_fill_guards
Definition: SpectralSolver.H:192
int dx
Definition: compute_domain.py:35
void CopySpectralDataComp(const int src_comp, const int dest_comp)
Copy spectral data from component src_comp to component dest_comp of field_data.fields.
Definition: SpectralSolver.H:157
void ForwardTransform(const int lev, 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:79
void pushSpectralFields()
Update the fields in spectral space, over one timestep.
Definition: SpectralSolver.cpp:99
void VayDeposition(const int lev, 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:145
void ZeroOutDataComp(const int icomp)
Set to zero the data on component icomp of field_data.fields.
Definition: SpectralSolver.H:169
Class that stores the fields in spectral space, and performs the Fourier transforms between real spac...
Definition: SpectralFieldData.H:105
SpectralFieldIndex m_spectral_index
Definition: SpectralSolver.H:188
void ReadParameters()
Definition: SpectralFieldData.H:32
std::unique_ptr< SpectralBaseAlgorithm > algorithm
Definition: SpectralSolver.H:204
void CurrentCorrection(const int lev, 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:130
SpectralSolver(const int lev, 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::IntVect &fill_guards, const amrex::Array< amrex::Real, 3 > &v_galilean, const amrex::Array< amrex::Real, 3 > &v_comoving, const amrex::RealVect dx, const amrex::Real dt, const bool pml, const bool periodic_single_box, const bool update_with_rho, const bool fft_do_time_averaging, const bool J_linear_in_time, const bool dive_cleaning, const bool divb_cleaning)
Constructor of the class SpectralSolver.
Definition: SpectralSolver.cpp:20
void ScaleDataComp(const int icomp, const amrex::Real scale_factor)
Scale the data on component icomp of field_data.fields by a given scale factor.
Definition: SpectralSolver.H:182
SpectralField fields
Definition: SpectralFieldData.H:133
void BackwardTransform(const int lev, 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:89
void ComputeSpectralDivE(const int lev, 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:114
Top-level class for the electromagnetic spectral solver.
Definition: SpectralSolver.H:32