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 
10 #include "SpectralSolverRZ_fwd.H"
11 
13 #include "SpectralFieldDataRZ.H"
14 
15 #include <ablastr/utils/Enums.H>
16 
17 
18 /* \brief Top-level class for the electromagnetic spectral solver
19  *
20  * Stores the field in spectral space, and has member functions
21  * to Fourier-transform the fields between real space and spectral space
22  * and to update fields in spectral space over one time step.
23  */
25 {
26  public:
27  // Inline definition of the member functions of `SpectralSolverRZ`,
28  // except the constructor (see `SpectralSolverRZ.cpp`)
29  // The body of these functions is short, since the work is done in the
30  // underlying classes `SpectralFieldData` and `PsatdAlgorithm`
31 
32  // Constructor
33  SpectralSolverRZ (int lev,
34  amrex::BoxArray const & realspace_ba,
35  amrex::DistributionMapping const & dm,
36  int n_rz_azimuthal_modes,
37  int norder_z,
39  const amrex::Vector<amrex::Real>& v_galilean,
40  amrex::RealVect dx, amrex::Real dt,
41  bool with_pml,
42  bool update_with_rho,
43  bool fft_do_time_averaging,
44  int J_in_time,
45  int rho_in_time,
46  bool dive_cleaning,
47  bool divb_cleaning);
48 
49  /* \brief Transform the component `i_comp` of MultiFab `field_mf`
50  * to spectral space, and store the corresponding result internally
51  * (in the spectral field specified by `field_index`) */
52  void ForwardTransform (int lev, amrex::MultiFab const & field_mf, int field_index,
53  int i_comp=0);
54 
55  /* \brief Transform the two MultiFabs `field_mf1` and `field_mf2`
56  * to spectral space, and store the corresponding results internally
57  * (in the spectral field specified by `field_index1` and `field_index2`) */
58  void ForwardTransform (int lev, amrex::MultiFab const & field_mf1, int field_index1,
59  amrex::MultiFab const & field_mf2, int field_index2);
60 
61  /* \brief Transform spectral field specified by `field_index` back to
62  * real space, and store it in the component `i_comp` of `field_mf` */
63  void BackwardTransform (int lev, amrex::MultiFab& field_mf, int field_index,
64  int i_comp=0);
65 
66  /* \brief Transform spectral fields specified by `field_index1` and `field_index2`
67  * back to real space, and store it in `field_mf1` and `field_mf2`*/
68  void BackwardTransform (int lev, amrex::MultiFab& field_mf1, int field_index1,
69  amrex::MultiFab& field_mf2, int field_index2);
70 
71  /* \brief Update the fields in spectral space, over one timestep */
72  void pushSpectralFields (bool doing_pml=false);
73 
74  /* \brief Initialize K space filtering arrays */
75  void InitFilter (amrex::IntVect const & filter_npass_each_dir,
76  bool const compensation)
77  {
78  field_data.InitFilter(filter_npass_each_dir, compensation, k_space);
79  }
80 
81  /* \brief Apply K space filtering for a scalar */
82  void ApplyFilter (const int lev, int const field_index)
83  {
84  field_data.ApplyFilter(lev, field_index);
85  }
86 
87  /* \brief Apply K space filtering for a vector */
88  void ApplyFilter (const int lev, int const field_index1,
89  int const field_index2, int const field_index3)
90  {
91  field_data.ApplyFilter(lev, field_index1, field_index2, field_index3);
92  }
93 
98  void ComputeSpectralDivE (int lev, const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
99  amrex::MultiFab& divE);
100 
107  void CurrentCorrection ();
108 
115  void VayDeposition ();
116 
124  void CopySpectralDataComp (const int src_comp, const int dest_comp)
125  {
126  field_data.CopySpectralDataComp(src_comp, dest_comp);
127  }
128 
134  void ZeroOutDataComp (const int icomp)
135  {
137  }
138 
146  void ScaleDataComp (const int icomp, const amrex::Real scale_factor)
147  {
148  field_data.ScaleDataComp(icomp, scale_factor);
149  }
150 
152 
153  private:
154 
155  SpectralKSpaceRZ k_space; // Save the instance to initialize filtering
156  SpectralFieldDataRZ field_data; // Store field in spectral space
157  // and perform the Fourier transforms
158  std::unique_ptr<SpectralBaseAlgorithmRZ> algorithm;
159  std::unique_ptr<SpectralBaseAlgorithmRZ> PML_algorithm;
160  // Defines field update equation in spectral space,
161  // and the associated coefficients.
162  // SpectralBaseAlgorithmRZ is a base class ; this pointer is meant
163  // to point an instance of a *sub-class* defining a specific algorithm
164 
165 };
166 
167 #endif // WARPX_SPECTRAL_SOLVER_RZ_H_
Definition: SpectralFieldDataRZ.H:23
void ZeroOutDataComp(int icomp)
Set to zero the data on component icomp of fields.
Definition: SpectralFieldDataRZ.H:100
void CopySpectralDataComp(int src_comp, int dest_comp)
Copy spectral data from component src_comp to component dest_comp of fields.
Definition: SpectralFieldDataRZ.H:80
void InitFilter(amrex::IntVect const &filter_npass_each_dir, bool compensation, SpectralKSpaceRZ const &k_space)
Definition: SpectralFieldDataRZ.cpp:756
void ScaleDataComp(int icomp, amrex::Real scale_factor)
Scale the data on component icomp of fields by a given scale factor.
Definition: SpectralFieldDataRZ.H:120
void ApplyFilter(int lev, int field_index)
Definition: SpectralFieldDataRZ.cpp:773
Definition: SpectralFieldData.H:34
Definition: SpectralKSpaceRZ.H:21
Definition: SpectralSolverRZ.H:25
SpectralFieldIndex m_spectral_index
Definition: SpectralSolverRZ.H:151
void ApplyFilter(const int lev, int const field_index1, int const field_index2, int const field_index3)
Definition: SpectralSolverRZ.H:88
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: SpectralSolverRZ.H:124
void ZeroOutDataComp(const int icomp)
Set to zero the data on component icomp of field_data.fields.
Definition: SpectralSolverRZ.H:134
void ApplyFilter(const int lev, int const field_index)
Definition: SpectralSolverRZ.H:82
void InitFilter(amrex::IntVect const &filter_npass_each_dir, bool const compensation)
Definition: SpectralSolverRZ.H:75
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: SpectralSolverRZ.H:146
SpectralKSpaceRZ k_space
Definition: SpectralSolverRZ.H:155
SpectralFieldDataRZ field_data
Definition: SpectralSolverRZ.H:156
void VayDeposition()
Public interface to call the virtual function VayDeposition, declared in the base class SpectralBaseA...
Definition: SpectralSolverRZ.cpp:163
void BackwardTransform(int lev, amrex::MultiFab &field_mf, int field_index, int i_comp=0)
Definition: SpectralSolverRZ.cpp:106
std::unique_ptr< SpectralBaseAlgorithmRZ > algorithm
Definition: SpectralSolverRZ.H:158
void pushSpectralFields(bool doing_pml=false)
Definition: SpectralSolverRZ.cpp:127
void ForwardTransform(int lev, amrex::MultiFab const &field_mf, int field_index, int i_comp=0)
Definition: SpectralSolverRZ.cpp:83
std::unique_ptr< SpectralBaseAlgorithmRZ > PML_algorithm
Definition: SpectralSolverRZ.H:159
void CurrentCorrection()
Public interface to call the virtual function CurrentCorrection, defined in the base class SpectralBa...
Definition: SpectralSolverRZ.cpp:157
SpectralSolverRZ(int lev, amrex::BoxArray const &realspace_ba, amrex::DistributionMapping const &dm, int n_rz_azimuthal_modes, int norder_z, ablastr::utils::enums::GridType grid_type, const amrex::Vector< amrex::Real > &v_galilean, amrex::RealVect dx, amrex::Real dt, bool with_pml, bool update_with_rho, bool fft_do_time_averaging, int J_in_time, int rho_in_time, bool dive_cleaning, bool divb_cleaning)
Definition: SpectralSolverRZ.cpp:28
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: SpectralSolverRZ.cpp:144
GridType
Definition: Enums.H:17
float dt
Definition: stencil.py:442
tuple dx
lab frame
Definition: stencil.py:429