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 (int lev,
37  const amrex::BoxArray& realspace_ba,
38  const SpectralKSpaceRZ& k_space,
40  int n_field_required,
41  int n_modes);
42  SpectralFieldDataRZ () = default; // Default constructor
44 
53 
54 
55  void ForwardTransform (int lev, const amrex::MultiFab& mf, int field_index,
56  int i_comp=0);
57  void ForwardTransform (int lev, const amrex::MultiFab& mf_r, int field_index_r,
58  const amrex::MultiFab& mf_t, int field_index_t);
59  void BackwardTransform (int lev, amrex::MultiFab& mf, int field_index,
60  int i_comp=0);
61  void BackwardTransform (int lev, amrex::MultiFab& mf_r, int field_index_r,
62  amrex::MultiFab& mf_t, int field_index_t);
63 
64  void FABZForwardTransform (amrex::MFIter const & mfi, amrex::Box const & realspace_bx,
65  amrex::MultiFab const & tempHTransformedSplit,
66  int field_index, bool is_nodal_z);
67  void FABZBackwardTransform (amrex::MFIter const & mfi, amrex::Box const & realspace_bx,
68  int field_index,
69  amrex::MultiFab & tempHTransformedSplit,
70  bool is_nodal_z);
71 
79  void CopySpectralDataComp (int src_comp, int dest_comp)
80  {
81  // In spectral space fields of each mode are grouped together, so that the index
82  // of a field for a specific mode is given by field_index + mode*n_fields.
83  // That’s why, looping over all azimuthal modes, we need to take into account corresponding
84  // shift, given by imode * m_n_fields.
85  for (int imode=0 ; imode < n_rz_azimuthal_modes ; imode++)
86  {
87  const int shift_comp = imode * m_n_fields;
88  // The last two arguments represent the number of components and
89  // the number of ghost cells to perform this operation
90  Copy(fields, fields, src_comp + shift_comp, dest_comp + shift_comp, 1, 0);
91  }
92  }
93 
99  void ZeroOutDataComp(int icomp)
100  {
101  // In spectral space fields of each mode are grouped together, so that the index
102  // of a field for a specific mode is given by field_index + mode*n_fields.
103  // That’s why, looping over all azimuthal modes, we need to take into account corresponding
104  // shift, given by imode * m_n_fields.
105  for (int imode=0 ; imode < n_rz_azimuthal_modes ; imode++)
106  {
107  const int shift_comp = imode * m_n_fields;
108  // The last argument represents the number of components to perform this operation
109  fields.setVal(0., icomp + shift_comp, 1);
110  }
111  }
112 
119  void ScaleDataComp(int icomp, amrex::Real scale_factor)
120  {
121  // In spectral space fields of each mode are grouped together, so that the index
122  // of a field for a specific mode is given by field_index + mode*n_fields.
123  // That’s why, looping over all azimuthal modes, we need to take into account corresponding
124  // shift, given by imode * m_n_fields.
125  for (int imode=0 ; imode < n_rz_azimuthal_modes ; imode++)
126  {
127  const int shift_comp = imode * m_n_fields;
128  // The last argument represents the number of components to perform this operation
129  fields.mult(scale_factor, icomp + shift_comp, 1);
130  }
131  }
132 
133  void InitFilter (amrex::IntVect const & filter_npass_each_dir, bool compensation,
134  SpectralKSpaceRZ const & k_space);
135 
136  void ApplyFilter (int lev, int field_index);
137  void ApplyFilter (int lev, int field_index1,
138  int field_index2, int field_index3);
139 
140  // Returns an array that holds the kr for all of the modes
142  return multi_spectral_hankel_transformer[mfi].getKrArray();
143  }
144 
147 
151  int m_ncomps;
152 
153  private:
154 
157 
158  // tempHTransformed and tmpSpectralField store fields
159  // right before/after the z Fourier transform
160  SpectralField tempHTransformed; // contains Complexes
161  SpectralField tmpSpectralField; // contains Complexes
163  // Correcting "shift" factors when performing FFT from/to
164  // a cell-centered grid in real space, instead of a nodal grid
168 
169 };
170 
171 #endif // WARPX_SPECTRAL_FIELD_DATA_RZ_H_
Definition: SpectralFieldDataRZ.H:22
FFTplans backward_plan
Definition: SpectralFieldDataRZ.H:162
void ZeroOutDataComp(int icomp)
Set to zero the data on component icomp of fields.
Definition: SpectralFieldDataRZ.H:99
int m_n_fields
Definition: SpectralFieldDataRZ.H:156
void CopySpectralDataComp(int src_comp, int dest_comp)
Copy spectral data from component src_comp to component dest_comp of fields.
Definition: SpectralFieldDataRZ.H:79
MultiSpectralHankelTransformer multi_spectral_hankel_transformer
Definition: SpectralFieldDataRZ.H:166
SpectralFieldDataRZ()=default
SpectralFieldIndex m_spectral_index
Definition: SpectralFieldDataRZ.H:155
void BackwardTransform(int lev, amrex::MultiFab &mf, int field_index, int i_comp=0)
Definition: SpectralFieldDataRZ.cpp:568
SpectralFieldDataRZ & operator=(SpectralFieldDataRZ const &)=delete
SpectralFieldDataRZ & operator=(SpectralFieldDataRZ &&field_data)=default
SpectralField tempHTransformed
Definition: SpectralFieldDataRZ.H:160
void InitFilter(amrex::IntVect const &filter_npass_each_dir, bool compensation, SpectralKSpaceRZ const &k_space)
Definition: SpectralFieldDataRZ.cpp:748
HankelTransform::RealVector const & getKrArray(amrex::MFIter const &mfi) const
Definition: SpectralFieldDataRZ.H:141
void FABZBackwardTransform(amrex::MFIter const &mfi, amrex::Box const &realspace_bx, int field_index, amrex::MultiFab &tempHTransformedSplit, bool is_nodal_z)
Definition: SpectralFieldDataRZ.cpp:343
SpectralFieldDataRZ(SpectralFieldDataRZ const &)=delete
void ForwardTransform(int lev, const amrex::MultiFab &mf, int field_index, int i_comp=0)
Definition: SpectralFieldDataRZ.cpp:443
FFTplans forward_plan
Definition: SpectralFieldDataRZ.H:162
void FABZForwardTransform(amrex::MFIter const &mfi, amrex::Box const &realspace_bx, amrex::MultiFab const &tempHTransformedSplit, int field_index, bool is_nodal_z)
Definition: SpectralFieldDataRZ.cpp:237
void ScaleDataComp(int icomp, amrex::Real scale_factor)
Scale the data on component icomp of fields by a given scale factor.
Definition: SpectralFieldDataRZ.H:119
int m_ncomps
Number of MultiFab components, see WarpX::ncomps.
Definition: SpectralFieldDataRZ.H:151
SpectralField fields
fields stores fields in spectral space, as multicomponent FabArray
Definition: SpectralFieldDataRZ.H:146
int n_rz_azimuthal_modes
Number of modes for the RZ multi-mode version, see WarpX::n_rz_azimuthal_modes.
Definition: SpectralFieldDataRZ.H:149
SpectralShiftFactor zshift_FFTfromCell
Definition: SpectralFieldDataRZ.H:165
SpectralFieldDataRZ(SpectralFieldDataRZ &&)=default
~SpectralFieldDataRZ()
Definition: SpectralFieldDataRZ.cpp:204
SpectralField tmpSpectralField
Definition: SpectralFieldDataRZ.H:161
SpectralShiftFactor zshift_FFTtoCell
Definition: SpectralFieldDataRZ.H:165
void ApplyFilter(int lev, int field_index)
Definition: SpectralFieldDataRZ.cpp:765
BinomialFilter binomialfilter
Definition: SpectralFieldDataRZ.H:167
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)