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
33 {
34  public:
35 
72  SpectralSolver (int lev,
73  const amrex::BoxArray& realspace_ba,
75  int norder_x, int norder_y,
76  int norder_z, short grid_type,
77  const amrex::Vector<amrex::Real>& v_galilean,
78  const amrex::Vector<amrex::Real>& v_comoving,
80  amrex::Real dt,
81  bool pml,
82  bool periodic_single_box,
83  bool update_with_rho,
84  bool fft_do_time_averaging,
85  int psatd_solution_type,
86  int J_in_time,
87  int rho_in_time,
88  bool dive_cleaning,
89  bool divb_cleaning);
90 
100  void ForwardTransform (int lev,
101  const amrex::MultiFab& mf,
102  int field_index,
103  int i_comp = 0);
104 
109  void BackwardTransform( int lev,
110  amrex::MultiFab& mf,
111  int field_index,
112  const amrex::IntVect& fill_guards,
113  int i_comp=0 );
114 
118  void pushSpectralFields();
119 
124  void ComputeSpectralDivE ( int lev,
125  const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
126  amrex::MultiFab& divE ) {
127  algorithm->ComputeSpectralDivE( lev, field_data, Efield, divE );
128  }
129 
137  {
138  algorithm->CurrentCorrection(field_data);
139  }
140 
148  {
149  algorithm->VayDeposition(field_data);
150  }
151 
159  void CopySpectralDataComp (const int src_comp, const int dest_comp)
160  {
161  // The last two arguments represent the number of components and
162  // the number of ghost cells to perform this operation
163  Copy(field_data.fields, field_data.fields, src_comp, dest_comp, 1, 0);
164  }
165 
171  void ZeroOutDataComp (const int icomp)
172  {
173  // The last argument represents the number of components to perform this operation
174  field_data.fields.setVal(0., icomp, 1);
175  }
176 
184  void ScaleDataComp (const int icomp, const amrex::Real scale_factor)
185  {
186  // The last argument represents the number of components to perform this operation
187  field_data.fields.mult(scale_factor, icomp, 1);
188  }
189 
191 
192  protected:
193 
195 
196  private:
197 
198  void ReadParameters ();
199 
200  // Store field in spectral space and perform the Fourier transforms
202 
203  // Defines field update equation in spectral space and the associated coefficients.
204  // SpectralBaseAlgorithm is a base class; this pointer is meant to point
205  // to an instance of a sub-class defining a specific algorithm
206  std::unique_ptr<SpectralBaseAlgorithm> algorithm;
207 };
208 #endif // WARPX_USE_PSATD
209 #endif // WARPX_SPECTRAL_SOLVER_H_
Class that stores the fields in spectral space, and performs the Fourier transforms between real spac...
Definition: SpectralFieldData.H:143
SpectralField fields
Definition: SpectralFieldData.H:169
Definition: SpectralFieldData.H:34
Top-level class for the electromagnetic spectral solver.
Definition: SpectralSolver.H:33
void ReadParameters()
void ForwardTransform(int lev, const amrex::MultiFab &mf, int field_index, int i_comp=0)
Transform the component i_comp of the MultiFab mf to Fourier space, and store the result internally (...
Definition: SpectralSolver.cpp:118
SpectralFieldIndex m_spectral_index
Definition: SpectralSolver.H:190
void VayDeposition()
Public interface to call the virtual function VayDeposition, declared in the base class SpectralBaseA...
Definition: SpectralSolver.H:147
void ZeroOutDataComp(const int icomp)
Set to zero the data on component icomp of field_data.fields.
Definition: SpectralSolver.H:171
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:184
void pushSpectralFields()
Update the fields in spectral space, over one timestep.
Definition: SpectralSolver.cpp:139
void CurrentCorrection()
Public interface to call the virtual function CurrentCorrection, defined in the base class SpectralBa...
Definition: SpectralSolver.H:136
void ComputeSpectralDivE(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:124
SpectralFieldData field_data
Definition: SpectralSolver.H:201
amrex::IntVect m_fill_guards
Definition: SpectralSolver.H:194
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:159
void BackwardTransform(int lev, amrex::MultiFab &mf, int field_index, const amrex::IntVect &fill_guards, 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:128
std::unique_ptr< SpectralBaseAlgorithm > algorithm
Definition: SpectralSolver.H:206
SpectralSolver(int lev, const amrex::BoxArray &realspace_ba, const amrex::DistributionMapping &dm, int norder_x, int norder_y, int norder_z, short grid_type, const amrex::Vector< amrex::Real > &v_galilean, const amrex::Vector< amrex::Real > &v_comoving, amrex::RealVect dx, amrex::Real dt, bool pml, bool periodic_single_box, bool update_with_rho, bool fft_do_time_averaging, int psatd_solution_type, int J_in_time, int rho_in_time, bool dive_cleaning, bool divb_cleaning)
Constructor of the class SpectralSolver.
Definition: SpectralSolver.cpp:23
void setVal(value_type val)
void mult(value_type val, int comp, int num_comp, int nghost=0)
void Copy(FabArray< DFAB > &dst, FabArray< SFAB > const &src, int srccomp, int dstcomp, int numcomp, int nghost)
float dt
Definition: stencil.py:442
tuple dx
lab frame
Definition: stencil.py:429