WarpX
Loading...
Searching...
No Matches
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
11#include "SpectralFieldData.H"
13#include "SpectralKSpaceRZ.H"
14
16
17#include <AMReX_MultiFab.H>
18
19/* \brief Class that stores the fields in spectral space, and performs the
20 * Fourier transforms between real space and spectral space
21 */
23{
24
25 public:
26
27 // Define the FFTplans type, which holds one fft plan per box
28 // (plans are only initialized for the boxes that are owned by
29 // the local MPI rank)
31
32 // Similarly, define the Hankel transformers and filter for each box.
34
36
37 SpectralFieldDataRZ (int lev,
38 const amrex::BoxArray& realspace_ba,
39 const SpectralKSpaceRZ& k_space,
41 int n_field_required,
42 int n_modes);
43 SpectralFieldDataRZ () = default; // Default constructor
45
54
55
56 void ForwardTransform (int lev, const amrex::MultiFab& mf, int field_index,
57 int i_comp=0);
58 void ForwardTransform (int lev, const amrex::MultiFab& mf_r, int field_index_r,
59 const amrex::MultiFab& mf_t, int field_index_t);
60 void BackwardTransform (int lev, amrex::MultiFab& mf, int field_index,
61 int i_comp=0);
62 void BackwardTransform (int lev, amrex::MultiFab& mf_r, int field_index_r,
63 amrex::MultiFab& mf_t, int field_index_t);
64
65 void FABZForwardTransform (amrex::MFIter const & mfi, amrex::Box const & realspace_bx,
66 amrex::MultiFab const & tempHTransformedSplit,
67 int field_index, bool is_nodal_z);
68 void FABZBackwardTransform (amrex::MFIter const & mfi, amrex::Box const & realspace_bx,
69 int field_index,
70 amrex::MultiFab & tempHTransformedSplit,
71 bool is_nodal_z);
72
80 void CopySpectralDataComp (int src_comp, int dest_comp);
81
87 void ZeroOutDataComp(int icomp);
88
95 void ScaleDataComp(int icomp, amrex::Real scale_factor);
96
97 void InitFilter (amrex::IntVect const & filter_npass_each_dir, bool compensation,
98 SpectralKSpaceRZ const & k_space);
99
100 void ApplyFilter (int lev, int field_index);
101 void ApplyFilter (int lev, int field_index1,
102 int field_index2, int field_index3);
103
104 // Returns an array that holds the kr for all of the modes
106 return multi_spectral_hankel_transformer[mfi].getKrArray();
107 }
108
111
116
117 private:
118
121
122 // tempHTransformed and tmpSpectralField store fields
123 // right before/after the z Fourier transform
124 SpectralField tempHTransformed; // contains Complexes
125 SpectralField tmpSpectralField; // contains Complexes
127 // Correcting "shift" factors when performing FFT from/to
128 // a cell-centered grid in real space, instead of a nodal grid
132
133};
134
135#endif // WARPX_SPECTRAL_FIELD_DATA_RZ_H_
amrex::FabArray< amrex::BaseFab< Complex > > SpectralField
Definition SpectralFieldData.H:32
amrex::LayoutData< amrex::Gpu::DeviceVector< Complex > > SpectralShiftFactor
Definition SpectralKSpace.H:34
amrex::Gpu::DeviceVector< amrex::Real > RealVector
Definition HankelTransform.H:32
HankelTransform::RealVector const & getKrArray(amrex::MFIter const &mfi) const
Definition SpectralFieldDataRZ.H:105
SpectralFieldDataRZ(int lev, const amrex::BoxArray &realspace_ba, const SpectralKSpaceRZ &k_space, const amrex::DistributionMapping &dm, int n_field_required, int n_modes)
Definition SpectralFieldDataRZ.cpp:28
FFTplans backward_plan
Definition SpectralFieldDataRZ.H:126
void ZeroOutDataComp(int icomp)
Set to zero the data on component icomp of fields.
Definition SpectralFieldDataRZ.cpp:770
int m_n_fields
Definition SpectralFieldDataRZ.H:120
void CopySpectralDataComp(int src_comp, int dest_comp)
Copy spectral data from component src_comp to component dest_comp of fields.
Definition SpectralFieldDataRZ.cpp:755
MultiSpectralHankelTransformer multi_spectral_hankel_transformer
Definition SpectralFieldDataRZ.H:130
amrex::LayoutData< SpectralHankelTransformer > MultiSpectralHankelTransformer
Definition SpectralFieldDataRZ.H:33
SpectralFieldDataRZ()=default
SpectralFieldIndex m_spectral_index
Definition SpectralFieldDataRZ.H:119
void BackwardTransform(int lev, amrex::MultiFab &mf, int field_index, int i_comp=0)
Definition SpectralFieldDataRZ.cpp:576
SpectralField tempHTransformed
Definition SpectralFieldDataRZ.H:124
void InitFilter(amrex::IntVect const &filter_npass_each_dir, bool compensation, SpectralKSpaceRZ const &k_space)
Definition SpectralFieldDataRZ.cpp:800
SpectralFieldDataRZ & operator=(SpectralFieldDataRZ const &)=delete
void FABZBackwardTransform(amrex::MFIter const &mfi, amrex::Box const &realspace_bx, int field_index, amrex::MultiFab &tempHTransformedSplit, bool is_nodal_z)
Definition SpectralFieldDataRZ.cpp:349
SpectralFieldDataRZ(SpectralFieldDataRZ const &)=delete
void ForwardTransform(int lev, const amrex::MultiFab &mf, int field_index, int i_comp=0)
Definition SpectralFieldDataRZ.cpp:451
FFTplans forward_plan
Definition SpectralFieldDataRZ.H:126
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:241
void ScaleDataComp(int icomp, amrex::Real scale_factor)
Scale the data on component icomp of fields by a given scale factor.
Definition SpectralFieldDataRZ.cpp:784
int m_ncomps
Number of MultiFab components, see WarpX::ncomps.
Definition SpectralFieldDataRZ.H:115
SpectralField fields
fields stores fields in spectral space, as multicomponent FabArray
Definition SpectralFieldDataRZ.H:110
int n_rz_azimuthal_modes
Number of modes for the RZ multi-mode version, see WarpX::n_rz_azimuthal_modes.
Definition SpectralFieldDataRZ.H:113
SpectralShiftFactor zshift_FFTfromCell
Definition SpectralFieldDataRZ.H:129
amrex::LayoutData< ablastr::math::anyfft::VendorFFTPlan > FFTplans
Definition SpectralFieldDataRZ.H:30
SpectralFieldDataRZ(SpectralFieldDataRZ &&)=default
~SpectralFieldDataRZ()
Definition SpectralFieldDataRZ.cpp:208
SpectralField tmpSpectralField
Definition SpectralFieldDataRZ.H:125
SpectralFieldDataRZ & operator=(SpectralFieldDataRZ &&field_data)=default
SpectralShiftFactor zshift_FFTtoCell
Definition SpectralFieldDataRZ.H:129
amrex::LayoutData< SpectralBinomialFilter > BinomialFilter
Definition SpectralFieldDataRZ.H:35
void ApplyFilter(int lev, int field_index)
Definition SpectralFieldDataRZ.cpp:817
BinomialFilter binomialfilter
Definition SpectralFieldDataRZ.H:131
Definition SpectralFieldData.H:35
Definition SpectralKSpaceRZ.H:21
BoxND< 3 > Box
IntVectND< 3 > IntVect