7 #ifndef WARPX_ACCELERATORLATTICE_LATTICEELEMENTS_LATTICEELEMENTFINDER_H_
8 #define WARPX_ACCELERATORLATTICE_LATTICEELEMENTS_LATTICEELEMENTFINDER_H_
95 using namespace amrex::literals;
98 amrex::ParticleReal
const * zs_arr = zs.
data();
99 amrex::ParticleReal
const * ze_arr = ze.
data();
100 int * indices_arr = indices.
data();
102 amrex::Real
const zmin =
m_zmin;
103 amrex::Real
const dz =
m_dz;
106 amrex::ParticleReal
const uz_boost =
m_uz_boost;
107 amrex::Real
const time =
m_time;
113 amrex::Real z_node = zmin + iz*
dz;
122 for (
int ie = 0 ; ie <
nelements ; ie++) {
125 const amrex::ParticleReal zcenter_left = (ie == 0)?
126 (std::numeric_limits<amrex::ParticleReal>::lowest()) : (0.5_prt*(ze_arr[ie-1] + zs_arr[ie]));
127 const amrex::ParticleReal zcenter_right = (ie <
nelements - 1)?
128 (0.5_prt*(ze_arr[ie] + zs_arr[ie+1])) : (std::numeric_limits<amrex::ParticleReal>::max());
129 if (zcenter_left <= z_node && z_node < zcenter_right) {
130 indices_arr[iz] = ie;
156 InitLatticeElementFinderDevice (
WarpXParIter const& a_pti,
int a_offset,
180 int const* d_quad_indices_arr =
nullptr;
181 int const* d_plasmalens_indices_arr =
nullptr;
190 void operator () (
const long i,
191 amrex::ParticleReal& field_Ex,
192 amrex::ParticleReal& field_Ey,
193 amrex::ParticleReal& field_Ez,
194 amrex::ParticleReal& field_Bx,
195 amrex::ParticleReal& field_By,
196 amrex::ParticleReal& field_Bz)
const noexcept
199 using namespace amrex::literals;
201 amrex::ParticleReal x, y, z;
202 m_get_position(
i, x, y, z);
206 const int iz =
static_cast<int>((z -
m_zmin)/
m_dz);
209 amrex::ParticleReal
const gamma = std::sqrt(1._prt + (m_ux[
i]*m_ux[
i] + m_uy[
i]*m_uy[
i] + m_uz[
i]*m_uz[
i])*inv_c2);
210 amrex::ParticleReal
const vzp = m_uz[
i]/
gamma;
212 amrex::ParticleReal zpvdt = z + vzp*m_dt;
220 amrex::ParticleReal Ex_sum = 0._prt;
221 amrex::ParticleReal Ey_sum = 0._prt;
222 const amrex::ParticleReal Ez_sum = 0._prt;
223 amrex::ParticleReal Bx_sum = 0._prt;
224 amrex::ParticleReal By_sum = 0._prt;
225 const amrex::ParticleReal Bz_sum = 0._prt;
228 if (d_quad_indices_arr[iz] > -1) {
229 const auto ielement = d_quad_indices_arr[iz];
230 amrex::ParticleReal Ex, Ey, Bx, By;
231 d_quad.
get_field(ielement, x, y, z, zpvdt, Ex, Ey, Bx, By);
240 if (d_plasmalens_indices_arr[iz] > -1) {
241 const auto ielement = d_plasmalens_indices_arr[iz];
242 amrex::ParticleReal Ex, Ey, Bx, By;
243 d_plasmalens.
get_field(ielement, x, y, z, zpvdt, Ex, Ey, Bx, By);
#define AMREX_FORCE_INLINE
#define AMREX_GPU_HOST_DEVICE
Definition: AcceleratorLattice.H:21
Definition: WarpXParticleContainer.H:52
size_type size() const noexcept
static constexpr auto c
vacuum speed of light [m/s]
Definition: constant.H:44
std::enable_if_t< std::is_integral< T >::value > ParallelFor(TypeList< CTOs... >, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
i
Definition: check_interp_points_and_weights.py:174
int gamma_boost
Definition: compute_domain.py:41
int dz
Definition: compute_domain.py:36
int gamma
boosted frame
Definition: stencil.py:431
constexpr int nelements
Definition: IonizationEnergiesTable.H:116
Definition: HardEdgedPlasmaLens.H:64
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void get_field(const int ielement, const amrex::ParticleReal x, const amrex::ParticleReal y, const amrex::ParticleReal z, const amrex::ParticleReal zpvdt, amrex::ParticleReal &Ex, amrex::ParticleReal &Ey, amrex::ParticleReal &Bx, amrex::ParticleReal &By) const
Fetch the field of the specified element at the given location.
Definition: HardEdgedPlasmaLens.H:90
int nelements
Definition: HardEdgedPlasmaLens.H:73
Definition: HardEdgedQuadrupole.H:64
int nelements
Definition: HardEdgedQuadrupole.H:73
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void get_field(const int ielement, const amrex::ParticleReal x, const amrex::ParticleReal y, const amrex::ParticleReal z, const amrex::ParticleReal zpvdt, amrex::ParticleReal &Ex, amrex::ParticleReal &Ey, amrex::ParticleReal &Bx, amrex::ParticleReal &By) const
Fetch the field of the specified element at the given location.
Definition: HardEdgedQuadrupole.H:90
The lattice element finder class that can be trivially copied to the device. This only has simple dat...
Definition: LatticeElementFinder.H:145
amrex::Real m_time
Definition: LatticeElementFinder.H:168
HardEdgedPlasmaLensDevice d_plasmalens
Definition: LatticeElementFinder.H:177
amrex::Real m_dt
Definition: LatticeElementFinder.H:163
amrex::Real m_dz
Definition: LatticeElementFinder.H:162
amrex::ParticleReal m_uz_boost
Definition: LatticeElementFinder.H:167
amrex::Real m_zmin
Definition: LatticeElementFinder.H:161
HardEdgedQuadrupoleDevice d_quad
Definition: LatticeElementFinder.H:176
GetParticlePosition< PIdx > m_get_position
Definition: LatticeElementFinder.H:170
amrex::ParticleReal m_gamma_boost
Definition: LatticeElementFinder.H:166
Definition: LatticeElementFinder.H:27
amrex::Gpu::DeviceVector< int > d_quad_indices
Definition: LatticeElementFinder.H:79
int m_nz
Definition: LatticeElementFinder.H:58
void setup_lattice_indices(amrex::Gpu::DeviceVector< amrex::ParticleReal > const &zs, amrex::Gpu::DeviceVector< amrex::ParticleReal > const &ze, amrex::Gpu::DeviceVector< int > &indices)
Fill in the index lookup tables This loops over the grid (in z) and finds the lattice element closest...
Definition: LatticeElementFinder.H:90
amrex::Real m_zmin
Definition: LatticeElementFinder.H:59
amrex::Gpu::DeviceVector< int > d_plasmalens_indices
Definition: LatticeElementFinder.H:80
amrex::ParticleReal m_gamma_boost
Definition: LatticeElementFinder.H:64
amrex::Real m_time
Definition: LatticeElementFinder.H:66
void UpdateIndices(int lev, amrex::MFIter const &a_mfi, AcceleratorLattice const &accelerator_lattice)
Update the index lookup tables for each element type, filling in the values.
Definition: LatticeElementFinder.cpp:54
amrex::Real m_dz
Definition: LatticeElementFinder.H:60
amrex::ParticleReal m_uz_boost
Definition: LatticeElementFinder.H:65
void InitElementFinder(int lev, amrex::MFIter const &a_mfi, AcceleratorLattice const &accelerator_lattice)
Initialize the element finder at the level and grid.
Definition: LatticeElementFinder.cpp:18
void AllocateIndices(AcceleratorLattice const &accelerator_lattice)
Allocate the index lookup tables for each element type.
Definition: LatticeElementFinder.cpp:39
LatticeElementFinderDevice GetFinderDeviceInstance(WarpXParIter const &a_pti, int a_offset, AcceleratorLattice const &accelerator_lattice)
Get the device level instance associated with this instance.
Definition: LatticeElementFinder.cpp:80