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 /* \brief Top-level class for the electromagnetic spectral solver
16  *
17  * Stores the field in spectral space, and has member functions
18  * to Fourier-transform the fields between real space and spectral space
19  * and to update fields in spectral space over one time step.
20  */
22 {
23  public:
24  // Inline definition of the member functions of `SpectralSolverRZ`,
25  // except the constructor (see `SpectralSolverRZ.cpp`)
26  // The body of these functions is short, since the work is done in the
27  // underlying classes `SpectralFieldData` and `PsatdAlgorithm`
28 
29  // Constructor
30  SpectralSolverRZ (const int lev,
31  amrex::BoxArray const & realspace_ba,
32  amrex::DistributionMapping const & dm,
33  int const n_rz_azimuthal_modes,
34  int const norder_z, bool const nodal,
35  const amrex::Array<amrex::Real,3>& v_galilean,
36  amrex::RealVect const dx, amrex::Real const dt,
37  bool const update_with_rho,
38  const bool fft_do_time_averaging,
39  const bool J_linear_in_time,
40  const bool dive_cleaning,
41  const bool divb_cleaning);
42 
43  /* \brief Transform the component `i_comp` of MultiFab `field_mf`
44  * to spectral space, and store the corresponding result internally
45  * (in the spectral field specified by `field_index`) */
46  void ForwardTransform (const int lev, amrex::MultiFab const & field_mf, int const field_index,
47  int const i_comp=0);
48 
49  /* \brief Transform the two MultiFabs `field_mf1` and `field_mf2`
50  * to spectral space, and store the corresponding results internally
51  * (in the spectral field specified by `field_index1` and `field_index2`) */
52  void ForwardTransform (const int lev, amrex::MultiFab const & field_mf1, int const field_index1,
53  amrex::MultiFab const & field_mf2, int const field_index2);
54 
55  /* \brief Transform spectral field specified by `field_index` back to
56  * real space, and store it in the component `i_comp` of `field_mf` */
57  void BackwardTransform (const int lev, amrex::MultiFab& field_mf, int const field_index,
58  int const i_comp=0);
59 
60  /* \brief Transform spectral fields specified by `field_index1` and `field_index2`
61  * back to real space, and store it in `field_mf1` and `field_mf2`*/
62  void BackwardTransform (const int lev, amrex::MultiFab& field_mf1, int const field_index1,
63  amrex::MultiFab& field_mf2, int const field_index2);
64 
65  /* \brief Update the fields in spectral space, over one timestep */
66  void pushSpectralFields ();
67 
68  /* \brief Initialize K space filtering arrays */
69  void InitFilter (amrex::IntVect const & filter_npass_each_dir,
70  bool const compensation)
71  {
72  field_data.InitFilter(filter_npass_each_dir, compensation, k_space);
73  }
74 
75  /* \brief Apply K space filtering for a scalar */
76  void ApplyFilter (const int lev, int const field_index)
77  {
78  field_data.ApplyFilter(lev, field_index);
79  }
80 
81  /* \brief Apply K space filtering for a vector */
82  void ApplyFilter (const int lev, int const field_index1,
83  int const field_index2, int const field_index3)
84  {
85  field_data.ApplyFilter(lev, field_index1, field_index2, field_index3);
86  }
87 
92  void ComputeSpectralDivE (const int lev, const std::array<std::unique_ptr<amrex::MultiFab>,3>& Efield,
93  amrex::MultiFab& divE);
94 
105  void CurrentCorrection (const int lev, std::array<std::unique_ptr<amrex::MultiFab>,3>& current,
106  const std::unique_ptr<amrex::MultiFab>& rho);
107 
117  void VayDeposition (const int lev, std::array<std::unique_ptr<amrex::MultiFab>,3>& current);
118 
126  void CopySpectralDataComp (const int src_comp, const int dest_comp)
127  {
128  field_data.CopySpectralDataComp(src_comp, dest_comp);
129  }
130 
136  void ZeroOutDataComp (const int icomp)
137  {
139  }
140 
148  void ScaleDataComp (const int icomp, const amrex::Real scale_factor)
149  {
150  field_data.ScaleDataComp(icomp, scale_factor);
151  }
152 
154 
155  private:
156 
157  SpectralKSpaceRZ k_space; // Save the instance to initialize filtering
158  SpectralFieldDataRZ field_data; // Store field in spectral space
159  // and perform the Fourier transforms
160  std::unique_ptr<SpectralBaseAlgorithmRZ> algorithm;
161  // Defines field update equation in spectral space,
162  // and the associated coefficients.
163  // SpectralBaseAlgorithmRZ is a base class ; this pointer is meant
164  // to point an instance of a *sub-class* defining a specific algorithm
165 
166 };
167 
168 #endif // WARPX_SPECTRAL_SOLVER_RZ_H_
void ApplyFilter(const int lev, int const field_index)
Definition: SpectralSolverRZ.H:76
void InitFilter(amrex::IntVect const &filter_npass_each_dir, bool const compensation)
Definition: SpectralSolverRZ.H:69
void CopySpectralDataComp(const int src_comp, const int dest_comp)
Copy spectral data from component src_comp to component dest_comp of fields.
Definition: SpectralFieldDataRZ.H:74
int dx
Definition: compute_domain.py:35
void BackwardTransform(const int lev, amrex::MultiFab &field_mf, int const field_index, int const i_comp=0)
Definition: SpectralSolverRZ.cpp:99
void ScaleDataComp(const int icomp, const amrex::Real scale_factor)
Scale the data on component icomp of fields by a given scale factor.
Definition: SpectralFieldDataRZ.H:98
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: SpectralSolverRZ.cpp:133
void ForwardTransform(const int lev, amrex::MultiFab const &field_mf, int const field_index, int const i_comp=0)
Definition: SpectralSolverRZ.cpp:76
SpectralFieldIndex m_spectral_index
Definition: SpectralSolverRZ.H:153
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: SpectralSolverRZ.cpp:150
std::unique_ptr< SpectralBaseAlgorithmRZ > algorithm
Definition: SpectralSolverRZ.H:160
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:126
Definition: SpectralFieldDataRZ.H:20
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: SpectralSolverRZ.cpp:157
void pushSpectralFields()
Definition: SpectralSolverRZ.cpp:120
SpectralKSpaceRZ k_space
Definition: SpectralSolverRZ.H:157
void ZeroOutDataComp(const int icomp)
Set to zero the data on component icomp of field_data.fields.
Definition: SpectralSolverRZ.H:136
Definition: SpectralKSpaceRZ.H:18
void InitFilter(amrex::IntVect const &filter_npass_each_dir, bool const compensation, SpectralKSpaceRZ const &k_space)
Definition: SpectralFieldDataRZ.cpp:702
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:148
Definition: SpectralFieldData.H:32
void ApplyFilter(const int lev, int const field_index)
Definition: SpectralFieldDataRZ.cpp:719
void ZeroOutDataComp(const int icomp)
Set to zero the data on component icomp of fields.
Definition: SpectralFieldDataRZ.H:86
Definition: SpectralSolverRZ.H:21
SpectralFieldDataRZ field_data
Definition: SpectralSolverRZ.H:158
void ApplyFilter(const int lev, int const field_index1, int const field_index2, int const field_index3)
Definition: SpectralSolverRZ.H:82
SpectralSolverRZ(const int lev, 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, bool const update_with_rho, const bool fft_do_time_averaging, const bool J_linear_in_time, const bool dive_cleaning, const bool divb_cleaning)
Definition: SpectralSolverRZ.cpp:28