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 
70  SpectralSolver (const int lev,
71  const amrex::BoxArray& realspace_ba,
73  const int norder_x, const int norder_y,
74  const int norder_z, const bool nodal,
75  const amrex::Vector<amrex::Real>& v_galilean,
76  const amrex::Vector<amrex::Real>& 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 int J_in_time,
84  const int rho_in_time,
85  const bool dive_cleaning,
86  const bool divb_cleaning);
87 
97  void ForwardTransform (const int lev,
98  const amrex::MultiFab& mf,
99  const int field_index,
100  const int i_comp = 0);
101 
106  void BackwardTransform( const int lev,
107  amrex::MultiFab& mf,
108  const int field_index,
109  const amrex::IntVect& fill_guards,
110  const int i_comp=0 );
111 
115  void pushSpectralFields();
116 
121  void ComputeSpectralDivE ( const int lev,
122  const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
123  amrex::MultiFab& divE ) {
124  algorithm->ComputeSpectralDivE( lev, field_data, Efield, divE );
125  }
126 
134  {
135  algorithm->CurrentCorrection(field_data);
136  }
137 
149  {
150  algorithm->VayDeposition(field_data);
151  }
152 
160  void CopySpectralDataComp (const int src_comp, const int dest_comp)
161  {
162  // The last two arguments represent the number of components and
163  // the number of ghost cells to perform this operation
164  Copy(field_data.fields, field_data.fields, src_comp, dest_comp, 1, 0);
165  }
166 
172  void ZeroOutDataComp (const int icomp)
173  {
174  // The last argument represents the number of components to perform this operation
175  field_data.fields.setVal(0., icomp, 1);
176  }
177 
185  void ScaleDataComp (const int icomp, const amrex::Real scale_factor)
186  {
187  // The last argument represents the number of components to perform this operation
188  field_data.fields.mult(scale_factor, icomp, 1);
189  }
190 
192 
193  protected:
194 
196 
197  private:
198 
199  void ReadParameters ();
200 
201  // Store field in spectral space and perform the Fourier transforms
203 
204  // Defines field update equation in spectral space and the associated coefficients.
205  // SpectralBaseAlgorithm is a base class; this pointer is meant to point
206  // to an instance of a sub-class defining a specific algorithm
207  std::unique_ptr<SpectralBaseAlgorithm> algorithm;
208 };
209 #endif // WARPX_USE_PSATD
210 #endif // WARPX_SPECTRAL_SOLVER_H_
void setVal(value_type val)
SpectralFieldData field_data
Definition: SpectralSolver.H:202
amrex::IntVect m_fill_guards
Definition: SpectralSolver.H:195
int dx
Definition: compute_domain.py:35
void CurrentCorrection()
Public interface to call the virtual function CurrentCorrection, defined in the base class SpectralBa...
Definition: SpectralSolver.H:133
int dt
Definition: Stencil.py:468
void BackwardTransform(const int lev, amrex::MultiFab &mf, const int field_index, const amrex::IntVect &fill_guards, 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:101
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:160
void ForwardTransform(const int lev, const amrex::MultiFab &mf, const int field_index, const int i_comp=0)
Transform the component i_comp of the MultiFab mf to Fourier space, and store the result internally (...
Definition: SpectralSolver.cpp:91
void pushSpectralFields()
Update the fields in spectral space, over one timestep.
Definition: SpectralSolver.cpp:112
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::Vector< amrex::Real > &v_galilean, const amrex::Vector< amrex::Real > &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 int J_in_time, const int rho_in_time, const bool dive_cleaning, const bool divb_cleaning)
Constructor of the class SpectralSolver.
Definition: SpectralSolver.cpp:21
void Copy(FabArray< FAB > &dst, FabArray< FAB > const &src, int srccomp, int dstcomp, int numcomp, int nghost)
void ZeroOutDataComp(const int icomp)
Set to zero the data on component icomp of field_data.fields.
Definition: SpectralSolver.H:172
Class that stores the fields in spectral space, and performs the Fourier transforms between real spac...
Definition: SpectralFieldData.H:117
void mult(value_type val, int comp, int num_comp, int nghost=0)
SpectralFieldIndex m_spectral_index
Definition: SpectralSolver.H:191
void ReadParameters()
Definition: SpectralFieldData.H:32
std::unique_ptr< SpectralBaseAlgorithm > algorithm
Definition: SpectralSolver.H:207
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:185
void VayDeposition()
Public interface to call the virtual function VayDeposition, declared in the base class SpectralBaseA...
Definition: SpectralSolver.H:148
SpectralField fields
Definition: SpectralFieldData.H:139
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:121
Top-level class for the electromagnetic spectral solver.
Definition: SpectralSolver.H:32