WarpX
SpectralSolverRZ.H
Go to the documentation of this file.
1 /* Copyright 2019-2020 David Grote
2  *
3  * This file is part of WarpX.
4  *
5  * License: BSD-3-Clause-LBNL
6  */
7 #ifndef WARPX_SPECTRAL_SOLVER_RZ_H_
8 #define WARPX_SPECTRAL_SOLVER_RZ_H_
9 
11 #include "SpectralFieldDataRZ.H"
12 
13 /* \brief Top-level class for the electromagnetic spectral solver
14  *
15  * Stores the field in spectral space, and has member functions
16  * to Fourier-transform the fields between real space and spectral space
17  * and to update fields in spectral space over one time step.
18  */
20 {
21  public:
22  // Inline definition of the member functions of `SpectralSolverRZ`,
23  // except the constructor (see `SpectralSolverRZ.cpp`)
24  // The body of these functions is short, since the work is done in the
25  // underlying classes `SpectralFieldData` and `PsatdAlgorithm`
26 
27  // Constructor
28  SpectralSolverRZ (amrex::BoxArray const & realspace_ba,
29  amrex::DistributionMapping const & dm,
30  int const n_rz_azimuthal_modes,
31  int const norder_z, bool const nodal,
32  const amrex::Array<amrex::Real,3>& v_galilean,
33  amrex::RealVect const dx, amrex::Real const dt,
34  int const lev);
35 
36  /* \brief Transform the component `i_comp` of MultiFab `field_mf`
37  * to spectral space, and store the corresponding result internally
38  * (in the spectral field specified by `field_index`) */
39  void ForwardTransform (amrex::MultiFab const & field_mf, int const field_index,
40  int const i_comp=0);
41 
42  /* \brief Transform the two MultiFabs `field_mf1` and `field_mf2`
43  * to spectral space, and store the corresponding results internally
44  * (in the spectral field specified by `field_index1` and `field_index2`) */
45  void ForwardTransform (amrex::MultiFab const & field_mf1, int const field_index1,
46  amrex::MultiFab const & field_mf2, int const field_index2);
47 
48  /* \brief Transform spectral field specified by `field_index` back to
49  * real space, and store it in the component `i_comp` of `field_mf` */
50  void BackwardTransform (amrex::MultiFab& field_mf, int const field_index,
51  int const i_comp=0);
52 
53  /* \brief Transform spectral fields specified by `field_index1` and `field_index2`
54  * back to real space, and store it in `field_mf1` and `field_mf2`*/
55  void BackwardTransform (amrex::MultiFab& field_mf1, int const field_index1,
56  amrex::MultiFab& field_mf2, int const field_index2);
57 
58  /* \brief Update the fields in spectral space, over one timestep */
59  void pushSpectralFields ();
60 
61  /* \brief Initialize K space filtering arrays */
62  void InitFilter (amrex::IntVect const & filter_npass_each_dir,
63  bool const compensation)
64  {
65  field_data.InitFilter(filter_npass_each_dir, compensation, k_space);
66  }
67 
68  /* \brief Apply K space filtering for a scalar */
69  void ApplyFilter (int const field_index)
70  {
71  field_data.ApplyFilter(field_index);
72  }
73 
74  /* \brief Apply K space filtering for a vector */
75  void ApplyFilter (int const field_index1, int const field_index2, int const field_index3)
76  {
77  field_data.ApplyFilter(field_index1, field_index2, field_index3);
78  }
79 
84  void ComputeSpectralDivE (const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
85  amrex::MultiFab& divE);
86 
97  void CurrentCorrection (std::array<std::unique_ptr<amrex::MultiFab>,3>& current,
98  const std::unique_ptr<amrex::MultiFab>& rho);
99 
109  void VayDeposition (std::array<std::unique_ptr<amrex::MultiFab>,3>& current);
110 
111  private:
112 
113  SpectralKSpaceRZ k_space; // Save the instance to initialize filtering
114  SpectralFieldDataRZ field_data; // Store field in spectral space
115  // and perform the Fourier transforms
116  std::unique_ptr<SpectralBaseAlgorithmRZ> algorithm;
117  // Defines field update equation in spectral space,
118  // and the associated coefficients.
119  // SpectralBaseAlgorithmRZ is a base class ; this pointer is meant
120  // to point an instance of a *sub-class* defining a specific algorithm
121 
122 };
123 
124 #endif // WARPX_SPECTRAL_SOLVER_RZ_H_
void InitFilter(amrex::IntVect const &filter_npass_each_dir, bool const compensation)
Definition: SpectralSolverRZ.H:62
void ApplyFilter(int const field_index)
Definition: SpectralFieldDataRZ.cpp:479
int dx
Definition: compute_domain.py:35
void ApplyFilter(int const field_index1, int const field_index2, int const field_index3)
Definition: SpectralSolverRZ.H:75
void ApplyFilter(int const field_index)
Definition: SpectralSolverRZ.H:69
std::unique_ptr< SpectralBaseAlgorithmRZ > algorithm
Definition: SpectralSolverRZ.H:116
Definition: SpectralFieldDataRZ.H:19
void BackwardTransform(amrex::MultiFab &field_mf, int const field_index, int const i_comp=0)
Definition: SpectralSolverRZ.cpp:86
void pushSpectralFields()
Definition: SpectralSolverRZ.cpp:104
SpectralKSpaceRZ k_space
Definition: SpectralSolverRZ.H:113
Definition: SpectralKSpaceRZ.H:18
void InitFilter(amrex::IntVect const &filter_npass_each_dir, bool const compensation, SpectralKSpaceRZ const &k_space)
Definition: SpectralFieldDataRZ.cpp:462
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: SpectralSolverRZ.cpp:133
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: SpectralSolverRZ.cpp:117
SpectralSolverRZ(amrex::BoxArray const &realspace_ba, amrex::DistributionMapping const &dm, int const n_rz_azimuthal_modes, int const norder_z, bool const nodal, const amrex::Array< amrex::Real, 3 > &v_galilean, amrex::RealVect const dx, amrex::Real const dt, int const lev)
Definition: SpectralSolverRZ.cpp:28
void ForwardTransform(amrex::MultiFab const &field_mf, int const field_index, int const i_comp=0)
Definition: SpectralSolverRZ.cpp:66
Definition: SpectralSolverRZ.H:19
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: SpectralSolverRZ.cpp:139
SpectralFieldDataRZ field_data
Definition: SpectralSolverRZ.H:114