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 "SpectralKSpaceRZ.H"
11 #include "SpectralFieldData.H"
13 #include "SpectralBinomialFilter.H"
14 #include <AMReX_MultiFab.H>
15 
16 /* \brief Class that stores the fields in spectral space, and performs the
17  * Fourier transforms between real space and spectral space
18  */
20 {
21 
22  public:
23 
24  // Define the FFTplans type, which holds one fft plan per box
25  // (plans are only initialized for the boxes that are owned by
26  // the local MPI rank)
27 #ifdef AMREX_USE_GPU
28  using FFTplans = amrex::LayoutData<cufftHandle>;
29 #else
30  using FFTplans = amrex::LayoutData<fftw_plan>;
31 #endif
32  // Similarly, define the Hankel transformers and filter for each box.
33  using MultiSpectralHankelTransformer = amrex::LayoutData<SpectralHankelTransformer>;
34 
35  using BinomialFilter = amrex::LayoutData<SpectralBinomialFilter>;
36 
37  SpectralFieldDataRZ (const amrex::BoxArray& realspace_ba,
38  const SpectralKSpaceRZ& k_space,
39  const amrex::DistributionMapping& dm,
40  const int n_field_required,
41  const int n_modes,
42  const int lev);
43  SpectralFieldDataRZ () = default; // Default constructor
44  SpectralFieldDataRZ& operator=(SpectralFieldDataRZ&& field_data) = default;
46 
47  void ForwardTransform (const amrex::MultiFab& mf, const int field_index,
48  const int i_comp=0);
49  void ForwardTransform (const amrex::MultiFab& mf_r, const int field_index_r,
50  const amrex::MultiFab& mf_t, const int field_index_t);
51  void BackwardTransform (amrex::MultiFab& mf, const int field_index,
52  const int i_comp=0);
53  void BackwardTransform (amrex::MultiFab& mf_r, const int field_index_r,
54  amrex::MultiFab& mf_t, const int field_index_t);
55 
56  void FABZForwardTransform (amrex::MFIter const & mfi, amrex::Box const & realspace_bx,
57  amrex::MultiFab const & tempHTransformedSplit,
58  int field_index, const bool is_nodal_z);
59  void FABZBackwardTransform (amrex::MFIter const & mfi, amrex::Box const & realspace_bx,
60  const int field_index,
61  amrex::MultiFab & tempHTransformedSplit,
62  const bool is_nodal_z);
63 
64  void InitFilter (amrex::IntVect const & filter_npass_each_dir, bool const compensation,
65  SpectralKSpaceRZ const & k_space);
66 
67  void ApplyFilter (int const field_index);
68  void ApplyFilter (int const field_index1, int const field_index2, int const field_index3);
69 
70  // Returns an array that holds the kr for all of the modes
71  HankelTransform::RealVector const & getKrArray (amrex::MFIter const & mfi) const {
72  return multi_spectral_hankel_transformer[mfi].getKrArray();
73  }
74 
75  // `fields` stores fields in spectral space, as multicomponent FabArray
77 
79 
80  private:
81 
82  // tempHTransformed and tmpSpectralField store fields
83  // right before/after the z Fourier transform
84  SpectralField tempHTransformed; // contains Complexes
85  SpectralField tmpSpectralField; // contains Complexes
87  // Correcting "shift" factors when performing FFT from/to
88  // a cell-centered grid in real space, instead of a nodal grid
92 
93 };
94 
95 #endif // WARPX_SPECTRAL_FIELD_DATA_RZ_H_
amrex::LayoutData< SpectralHankelTransformer > MultiSpectralHankelTransformer
Definition: SpectralFieldDataRZ.H:33
SpectralFieldDataRZ()=default
SpectralFieldDataRZ & operator=(SpectralFieldDataRZ &&field_data)=default
void ApplyFilter(int const field_index)
Definition: SpectralFieldDataRZ.cpp:479
~SpectralFieldDataRZ()
Definition: SpectralFieldDataRZ.cpp:128
amrex::LayoutData< SpectralBinomialFilter > BinomialFilter
Definition: SpectralFieldDataRZ.H:35
SpectralField tempHTransformed
Definition: SpectralFieldDataRZ.H:84
HankelTransform::RealVector const & getKrArray(amrex::MFIter const &mfi) const
Definition: SpectralFieldDataRZ.H:71
SpectralShiftFactor zshift_FFTtoCell
Definition: SpectralFieldDataRZ.H:89
amrex::LayoutData< fftw_plan > FFTplans
Definition: SpectralFieldDataRZ.H:30
FFTplans forward_plan
Definition: SpectralFieldDataRZ.H:86
amrex::Gpu::DeviceVector< amrex::Real > RealVector
Definition: HankelTransform.H:24
Definition: SpectralFieldDataRZ.H:19
void BackwardTransform(amrex::MultiFab &mf, const int field_index, const int i_comp=0)
Definition: SpectralFieldDataRZ.cpp:376
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:153
BinomialFilter binomialfilter
Definition: SpectralFieldDataRZ.H:91
SpectralField fields
Definition: SpectralFieldDataRZ.H:76
FFTplans backward_plan
Definition: SpectralFieldDataRZ.H:86
MultiSpectralHankelTransformer multi_spectral_hankel_transformer
Definition: SpectralFieldDataRZ.H:90
Definition: SpectralKSpaceRZ.H:18
void InitFilter(amrex::IntVect const &filter_npass_each_dir, bool const compensation, SpectralKSpaceRZ const &k_space)
Definition: SpectralFieldDataRZ.cpp:462
amrex::FabArray< amrex::BaseFab< Complex > > SpectralField
Definition: SpectralFieldData.H:20
int n_rz_azimuthal_modes
Definition: SpectralFieldDataRZ.H:78
SpectralShiftFactor zshift_FFTfromCell
Definition: SpectralFieldDataRZ.H:89
amrex::LayoutData< amrex::Gpu::DeviceVector< Complex > > SpectralShiftFactor
Definition: SpectralKSpace.H:23
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:226
void ForwardTransform(const amrex::MultiFab &mf, const int field_index, const int i_comp=0)
Definition: SpectralFieldDataRZ.cpp:293
SpectralField tmpSpectralField
Definition: SpectralFieldDataRZ.H:85