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"
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)
29  using FFTplans = amrex::LayoutData<cufftHandle>;
30 #elif defined(AMREX_USE_HIP)
31  using FFTplans = amrex::LayoutData<rocfft_plan>;
32 #else
33  using FFTplans = amrex::LayoutData<fftw_plan>;
34 #endif
35  // Similarly, define the Hankel transformers and filter for each box.
36  using MultiSpectralHankelTransformer = amrex::LayoutData<SpectralHankelTransformer>;
37 
38  using BinomialFilter = amrex::LayoutData<SpectralBinomialFilter>;
39 
40  SpectralFieldDataRZ (const int lev,
41  const amrex::BoxArray& realspace_ba,
42  const SpectralKSpaceRZ& k_space,
43  const amrex::DistributionMapping& dm,
44  const int n_field_required,
45  const int n_modes);
46  SpectralFieldDataRZ () = default; // Default constructor
47  SpectralFieldDataRZ& operator=(SpectralFieldDataRZ&& field_data) = default;
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  // The last two arguments represent the number of components and
77  // the number of ghost cells to perform this operation
78  Copy(fields, fields, src_comp, dest_comp, n_rz_azimuthal_modes, 0);
79  }
80 
86  void ZeroOutDataComp(const int icomp)
87  {
88  // The last argument represents the number of components to perform this operation
89  fields.setVal(0., icomp, n_rz_azimuthal_modes);
90  }
91 
98  void ScaleDataComp(const int icomp, const amrex::Real scale_factor)
99  {
100  // The last argument represents the number of components to perform this operation
101  fields.mult(scale_factor, icomp, n_rz_azimuthal_modes);
102  }
103 
104  void InitFilter (amrex::IntVect const & filter_npass_each_dir, bool const compensation,
105  SpectralKSpaceRZ const & k_space);
106 
107  void ApplyFilter (const int lev, int const field_index);
108  void ApplyFilter (const int lev, int const field_index1,
109  int const field_index2, int const field_index3);
110 
111  // Returns an array that holds the kr for all of the modes
112  HankelTransform::RealVector const & getKrArray (amrex::MFIter const & mfi) const {
113  return multi_spectral_hankel_transformer[mfi].getKrArray();
114  }
115 
116  // `fields` stores fields in spectral space, as multicomponent FabArray
118 
120 
121  private:
122 
125 
126  // tempHTransformed and tmpSpectralField store fields
127  // right before/after the z Fourier transform
128  SpectralField tempHTransformed; // contains Complexes
129  SpectralField tmpSpectralField; // contains Complexes
131  // Correcting "shift" factors when performing FFT from/to
132  // a cell-centered grid in real space, instead of a nodal grid
136 
137 };
138 
139 #endif // WARPX_SPECTRAL_FIELD_DATA_RZ_H_
amrex::LayoutData< SpectralHankelTransformer > MultiSpectralHankelTransformer
Definition: SpectralFieldDataRZ.H:36
SpectralFieldDataRZ()=default
SpectralFieldDataRZ & operator=(SpectralFieldDataRZ &&field_data)=default
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
~SpectralFieldDataRZ()
Definition: SpectralFieldDataRZ.cpp:180
amrex::LayoutData< SpectralBinomialFilter > BinomialFilter
Definition: SpectralFieldDataRZ.H:38
SpectralField tempHTransformed
Definition: SpectralFieldDataRZ.H:128
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
HankelTransform::RealVector const & getKrArray(amrex::MFIter const &mfi) const
Definition: SpectralFieldDataRZ.H:112
SpectralShiftFactor zshift_FFTtoCell
Definition: SpectralFieldDataRZ.H:133
amrex::LayoutData< fftw_plan > FFTplans
Definition: SpectralFieldDataRZ.H:33
SpectralFieldIndex m_spectral_index
Definition: SpectralFieldDataRZ.H:123
FFTplans forward_plan
Definition: SpectralFieldDataRZ.H:130
amrex::Gpu::DeviceVector< amrex::Real > RealVector
Definition: HankelTransform.H:24
Definition: SpectralFieldDataRZ.H:20
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:208
BinomialFilter binomialfilter
Definition: SpectralFieldDataRZ.H:135
void ForwardTransform(const int lev, const amrex::MultiFab &mf, const int field_index, const int i_comp=0)
Definition: SpectralFieldDataRZ.cpp:398
SpectralField fields
Definition: SpectralFieldDataRZ.H:117
FFTplans backward_plan
Definition: SpectralFieldDataRZ.H:130
MultiSpectralHankelTransformer multi_spectral_hankel_transformer
Definition: SpectralFieldDataRZ.H:134
Definition: SpectralKSpaceRZ.H:18
void InitFilter(amrex::IntVect const &filter_npass_each_dir, bool const compensation, SpectralKSpaceRZ const &k_space)
Definition: SpectralFieldDataRZ.cpp:702
Definition: SpectralFieldData.H:32
int m_n_fields
Definition: SpectralFieldDataRZ.H:124
amrex::FabArray< amrex::BaseFab< Complex > > SpectralField
Definition: SpectralFieldData.H:30
void ApplyFilter(const int lev, int const field_index)
Definition: SpectralFieldDataRZ.cpp:719
int n_rz_azimuthal_modes
Definition: SpectralFieldDataRZ.H:119
SpectralShiftFactor zshift_FFTfromCell
Definition: SpectralFieldDataRZ.H:133
amrex::LayoutData< amrex::Gpu::DeviceVector< Complex > > SpectralShiftFactor
Definition: SpectralKSpace.H:32
void ZeroOutDataComp(const int icomp)
Set to zero the data on component icomp of fields.
Definition: SpectralFieldDataRZ.H:86
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:306
SpectralField tmpSpectralField
Definition: SpectralFieldDataRZ.H:129
void BackwardTransform(const int lev, amrex::MultiFab &mf, const int field_index, const int i_comp=0)
Definition: SpectralFieldDataRZ.cpp:522