WarpX
SpectralKSpace.H
Go to the documentation of this file.
1 /* Copyright 2019 David Grote, Maxence Thevenet, Remi Lehe
2  *
3  *
4  * This file is part of WarpX.
5  *
6  * License: BSD-3-Clause-LBNL
7  */
8 #ifndef WARPX_SPECTRAL_K_SPACE_H_
9 #define WARPX_SPECTRAL_K_SPACE_H_
10 
11 #include "Utils/WarpX_Complex.H"
12 
13 #include <AMReX_BoxArray.H>
14 #include <AMReX_LayoutData.H>
15 
16 
17 // `KVectorComponent` and `SpectralShiftFactor` hold one 1D array
18 // ("DeviceVector") for each box ("LayoutData"). The arrays are
19 // only allocated if the corresponding box is owned by the local MPI rank.
20 using RealKVector = amrex::Gpu::DeviceVector<amrex::Real>;
21 using KVectorComponent = amrex::LayoutData< RealKVector >;
22 using SpectralShiftFactor = amrex::LayoutData<
23  amrex::Gpu::DeviceVector<Complex> >;
24 
25 // Indicate the type of correction "shift" factor to apply
26 // when the FFT is performed from/to a cell-centered grid in real space.
27 struct ShiftType {
29 };
30 
39 {
40  public:
41  amrex::BoxArray spectralspace_ba;
42  SpectralKSpace() : dx(amrex::RealVect::Zero) {}
43  SpectralKSpace( const amrex::BoxArray& realspace_ba,
44  const amrex::DistributionMapping& dm,
45  const amrex::RealVect realspace_dx );
46  KVectorComponent getKComponent(
47  const amrex::DistributionMapping& dm,
48  const amrex::BoxArray& realspace_ba,
49  const int i_dim, const bool only_positive_k ) const;
50  KVectorComponent getModifiedKComponent(
51  const amrex::DistributionMapping& dm, const int i_dim,
52  const int n_order, const bool nodal ) const;
53  SpectralShiftFactor getSpectralShiftFactor(
54  const amrex::DistributionMapping& dm, const int i_dim,
55  const int shift_type ) const;
56 
57  protected:
58  amrex::Array<KVectorComponent, AMREX_SPACEDIM> k_vec;
59  // 3D: k_vec is an Array of 3 components, corresponding to kx, ky, kz
60  // 2D: k_vec is an Array of 2 components, corresponding to kx, kz
61  amrex::RealVect dx;
62 };
63 
64 amrex::Vector<amrex::Real>
65 getFonbergStencilCoefficients( const int n_order, const bool nodal );
66 
67 #endif
amrex::Array< KVectorComponent, AMREX_SPACEDIM > k_vec
Definition: SpectralKSpace.H:58
Class that represents the spectral space.
Definition: SpectralKSpace.H:38
amrex::RealVect dx
Definition: SpectralKSpace.H:61
int dx
Definition: compute_domain.py:35
amrex::Vector< amrex::Real > getFonbergStencilCoefficients(const int n_order, const bool nodal)
Definition: SpectralKSpace.cpp:275
Definition: SpectralKSpace.H:27
amrex::BoxArray spectralspace_ba
Definition: SpectralKSpace.H:41
Definition: SpectralKSpace.H:28
SpectralKSpace()
Definition: SpectralKSpace.H:42
amrex::Gpu::DeviceVector< amrex::Real > RealKVector
Definition: SpectralKSpace.H:20
Definition: SpectralKSpace.H:28
amrex::LayoutData< amrex::Gpu::DeviceVector< Complex > > SpectralShiftFactor
Definition: SpectralKSpace.H:23
amrex::LayoutData< RealKVector > KVectorComponent
Definition: SpectralKSpace.H:21
Definition: PML.H:52