WarpX
SpectralFieldDataRZ.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_FIELD_DATA_RZ_H_
8 #define WARPX_SPECTRAL_FIELD_DATA_RZ_H_
9 
10 #include "SpectralBinomialFilter.H"
11 #include "SpectralFieldData.H"
13 #include "SpectralKSpaceRZ.H"
15 
16 #include <AMReX_MultiFab.H>
17 
18 /* \brief Class that stores the fields in spectral space, and performs the
19  * Fourier transforms between real space and spectral space
20  */
22 {
23 
24  public:
25 
26  // Define the FFTplans type, which holds one fft plan per box
27  // (plans are only initialized for the boxes that are owned by
28  // the local MPI rank)
30 
31  // Similarly, define the Hankel transformers and filter for each box.
33 
35 
36  SpectralFieldDataRZ (const int lev,
37  const amrex::BoxArray& realspace_ba,
38  const SpectralKSpaceRZ& k_space,
40  const int n_field_required,
41  const int n_modes);
42  SpectralFieldDataRZ () = default; // Default constructor
45 
46  void ForwardTransform (const int lev, const amrex::MultiFab& mf, const int field_index,
47  const int i_comp=0);
48  void ForwardTransform (const int lev, const amrex::MultiFab& mf_r, const int field_index_r,
49  const amrex::MultiFab& mf_t, const int field_index_t);
50  void BackwardTransform (const int lev, amrex::MultiFab& mf, const int field_index,
51  const int i_comp=0);
52  void BackwardTransform (const int lev, amrex::MultiFab& mf_r, const int field_index_r,
53  amrex::MultiFab& mf_t, const int field_index_t);
54 
55  void FABZForwardTransform (amrex::MFIter const & mfi, amrex::Box const & realspace_bx,
56  amrex::MultiFab const & tempHTransformedSplit,
57  int field_index, const bool is_nodal_z);
58  void FABZBackwardTransform (amrex::MFIter const & mfi, amrex::Box const & realspace_bx,
59  const int field_index,
60  amrex::MultiFab & tempHTransformedSplit,
61  const bool is_nodal_z);
62 
70  void CopySpectralDataComp (const int src_comp, const int dest_comp)
71  {
72  // In spectral space fields of each mode are grouped together, so that the index
73  // of a field for a specific mode is given by field_index + mode*n_fields.
74  // That’s why, looping over all azimuthal modes, we need to take into account corresponding
75  // shift, given by imode * m_n_fields.
76  for (int imode=0 ; imode < n_rz_azimuthal_modes ; imode++)
77  {
78  const int shift_comp = imode * m_n_fields;
79  // The last two arguments represent the number of components and
80  // the number of ghost cells to perform this operation
81  Copy(fields, fields, src_comp + shift_comp, dest_comp + shift_comp, 1, 0);
82  }
83  }
84 
90  void ZeroOutDataComp(const int icomp)
91  {
92  // In spectral space fields of each mode are grouped together, so that the index
93  // of a field for a specific mode is given by field_index + mode*n_fields.
94  // That’s why, looping over all azimuthal modes, we need to take into account corresponding
95  // shift, given by imode * m_n_fields.
96  for (int imode=0 ; imode < n_rz_azimuthal_modes ; imode++)
97  {
98  const int shift_comp = imode * m_n_fields;
99  // The last argument represents the number of components to perform this operation
100  fields.setVal(0., icomp + shift_comp, 1);
101  }
102  }
103 
110  void ScaleDataComp(const int icomp, const amrex::Real scale_factor)
111  {
112  // In spectral space fields of each mode are grouped together, so that the index
113  // of a field for a specific mode is given by field_index + mode*n_fields.
114  // That’s why, looping over all azimuthal modes, we need to take into account corresponding
115  // shift, given by imode * m_n_fields.
116  for (int imode=0 ; imode < n_rz_azimuthal_modes ; imode++)
117  {
118  const int shift_comp = imode * m_n_fields;
119  // The last argument represents the number of components to perform this operation
120  fields.mult(scale_factor, icomp + shift_comp, 1);
121  }
122  }
123 
124  void InitFilter (amrex::IntVect const & filter_npass_each_dir, bool const compensation,
125  SpectralKSpaceRZ const & k_space);
126 
127  void ApplyFilter (const int lev, int const field_index);
128  void ApplyFilter (const int lev, int const field_index1,
129  int const field_index2, int const field_index3);
130 
131  // Returns an array that holds the kr for all of the modes
133  return multi_spectral_hankel_transformer[mfi].getKrArray();
134  }
135 
138 
142  int m_ncomps;
143 
144  private:
145 
148 
149  // tempHTransformed and tmpSpectralField store fields
150  // right before/after the z Fourier transform
151  SpectralField tempHTransformed; // contains Complexes
152  SpectralField tmpSpectralField; // contains Complexes
154  // Correcting "shift" factors when performing FFT from/to
155  // a cell-centered grid in real space, instead of a nodal grid
159 
160 };
161 
162 #endif // WARPX_SPECTRAL_FIELD_DATA_RZ_H_
Definition: SpectralFieldDataRZ.H:22
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:110
void FABZForwardTransform(amrex::MFIter const &mfi, amrex::Box const &realspace_bx, amrex::MultiFab const &tempHTransformedSplit, int field_index, const bool is_nodal_z)
Definition: SpectralFieldDataRZ.cpp:238
FFTplans backward_plan
Definition: SpectralFieldDataRZ.H:153
int m_n_fields
Definition: SpectralFieldDataRZ.H:147
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:70
void FABZBackwardTransform(amrex::MFIter const &mfi, amrex::Box const &realspace_bx, const int field_index, amrex::MultiFab &tempHTransformedSplit, const bool is_nodal_z)
Definition: SpectralFieldDataRZ.cpp:344
MultiSpectralHankelTransformer multi_spectral_hankel_transformer
Definition: SpectralFieldDataRZ.H:157
SpectralFieldDataRZ()=default
SpectralFieldIndex m_spectral_index
Definition: SpectralFieldDataRZ.H:146
void BackwardTransform(const int lev, amrex::MultiFab &mf, const int field_index, const int i_comp=0)
Definition: SpectralFieldDataRZ.cpp:569
SpectralFieldDataRZ & operator=(SpectralFieldDataRZ &&field_data)=default
SpectralField tempHTransformed
Definition: SpectralFieldDataRZ.H:151
HankelTransform::RealVector const & getKrArray(amrex::MFIter const &mfi) const
Definition: SpectralFieldDataRZ.H:132
void ApplyFilter(const int lev, int const field_index)
Definition: SpectralFieldDataRZ.cpp:766
FFTplans forward_plan
Definition: SpectralFieldDataRZ.H:153
int m_ncomps
Number of MultiFab components, see WarpX::ncomps.
Definition: SpectralFieldDataRZ.H:142
void ForwardTransform(const int lev, const amrex::MultiFab &mf, const int field_index, const int i_comp=0)
Definition: SpectralFieldDataRZ.cpp:444
SpectralField fields
fields stores fields in spectral space, as multicomponent FabArray
Definition: SpectralFieldDataRZ.H:137
int n_rz_azimuthal_modes
Number of modes for the RZ multi-mode version, see WarpX::n_rz_azimuthal_modes.
Definition: SpectralFieldDataRZ.H:140
SpectralShiftFactor zshift_FFTfromCell
Definition: SpectralFieldDataRZ.H:156
void InitFilter(amrex::IntVect const &filter_npass_each_dir, bool const compensation, SpectralKSpaceRZ const &k_space)
Definition: SpectralFieldDataRZ.cpp:749
~SpectralFieldDataRZ()
Definition: SpectralFieldDataRZ.cpp:205
SpectralField tmpSpectralField
Definition: SpectralFieldDataRZ.H:152
SpectralShiftFactor zshift_FFTtoCell
Definition: SpectralFieldDataRZ.H:156
void ZeroOutDataComp(const int icomp)
Set to zero the data on component icomp of fields.
Definition: SpectralFieldDataRZ.H:90
BinomialFilter binomialfilter
Definition: SpectralFieldDataRZ.H:158
Definition: SpectralFieldData.H:33
Definition: SpectralKSpaceRZ.H:21
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)