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