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 amrex::ParticleReal zcenter_left, zcenter_right;
127 zcenter_left = std::numeric_limits<amrex::ParticleReal>::lowest();
129 zcenter_left = 0.5_prt*(ze_arr[ie-1] + zs_arr[ie]);
132 zcenter_right = 0.5_prt*(ze_arr[ie] + zs_arr[ie+1]);
134 zcenter_right = std::numeric_limits<amrex::ParticleReal>::max();
137 if (zcenter_left <= z_node && z_node < zcenter_right) {
138 indices_arr[iz] = ie;
164 InitLatticeElementFinderDevice (
WarpXParIter const& a_pti,
int const a_offset,
188 int const* d_quad_indices_arr =
nullptr;
189 int const* d_plasmalens_indices_arr =
nullptr;
198 void operator () (
const long i,
199 amrex::ParticleReal& field_Ex,
200 amrex::ParticleReal& field_Ey,
201 amrex::ParticleReal& field_Ez,
202 amrex::ParticleReal& field_Bx,
203 amrex::ParticleReal& field_By,
204 amrex::ParticleReal& field_Bz)
const noexcept
207 using namespace amrex::literals;
209 amrex::ParticleReal
x,
y,
z;
210 m_get_position(
i,
x,
y,
z);
214 const int iz =
static_cast<int>((
z -
m_zmin)/
m_dz);
217 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);
218 amrex::ParticleReal
const vzp = m_uz[
i]/
gamma;
220 amrex::ParticleReal zpvdt =
z + vzp*m_dt;
228 amrex::ParticleReal Ex_sum = 0._prt;
229 amrex::ParticleReal Ey_sum = 0._prt;
230 amrex::ParticleReal Ez_sum = 0._prt;
231 amrex::ParticleReal Bx_sum = 0._prt;
232 amrex::ParticleReal By_sum = 0._prt;
233 amrex::ParticleReal Bz_sum = 0._prt;
236 if (d_quad_indices_arr[iz] > -1) {
237 int ielement = d_quad_indices_arr[iz];
238 amrex::ParticleReal Ex, Ey, Bx, By;
239 d_quad.
get_field(ielement,
x,
y,
z, zpvdt, Ex, Ey, Bx, By);
248 if (d_plasmalens_indices_arr[iz] > -1) {
249 int ielement = d_plasmalens_indices_arr[iz];
250 amrex::ParticleReal Ex, Ey, Bx, By;
251 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
def y
Definition: Excitation_Flag_Generator.py:76
def x
Formats datastring to remove "+" at the end of the string #####.
Definition: Excitation_Flag_Generator.py:75
def z
Definition: Excitation_Flag_Generator.py:77
static constexpr auto c
vacuum speed of light [m/s]
Definition: constant.H:44
std::enable_if_t< std::is_integral_v< T > > ParallelFor(TypeList< CTOs... > 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 gamma
Definition: stencil.py:446
int dz
Definition: stencil.py:438
constexpr int nelements
Definition: IonizationEnergiesTable.H:116
Functor that can be used to extract the positions of the macroparticles inside a ParallelFor kernel.
Definition: GetAndSetPosition.H:53
Definition: HardEdgedPlasmaLens.H:67
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:93
int nelements
Definition: HardEdgedPlasmaLens.H:76
Definition: HardEdgedQuadrupole.H:67
int nelements
Definition: HardEdgedQuadrupole.H:76
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:93
The lattice element finder class that can be trivially copied to the device. This only has simple dat...
Definition: LatticeElementFinder.H:153
amrex::Real m_time
Definition: LatticeElementFinder.H:176
HardEdgedPlasmaLensDevice d_plasmalens
Definition: LatticeElementFinder.H:185
amrex::Real m_dt
Definition: LatticeElementFinder.H:171
amrex::Real m_dz
Definition: LatticeElementFinder.H:170
amrex::ParticleReal m_uz_boost
Definition: LatticeElementFinder.H:175
amrex::Real m_zmin
Definition: LatticeElementFinder.H:169
HardEdgedQuadrupoleDevice d_quad
Definition: LatticeElementFinder.H:184
GetParticlePosition m_get_position
Definition: LatticeElementFinder.H:178
amrex::ParticleReal m_gamma_boost
Definition: LatticeElementFinder.H:174
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
void UpdateIndices(int const 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::ParticleReal m_gamma_boost
Definition: LatticeElementFinder.H:64
amrex::Real m_time
Definition: LatticeElementFinder.H:66
amrex::Real m_dz
Definition: LatticeElementFinder.H:60
amrex::ParticleReal m_uz_boost
Definition: LatticeElementFinder.H:65
void AllocateIndices(AcceleratorLattice const &accelerator_lattice)
Allocate the index lookup tables for each element type.
Definition: LatticeElementFinder.cpp:39
void InitElementFinder(int const lev, amrex::MFIter const &a_mfi, AcceleratorLattice const &accelerator_lattice)
Initialize the element finder at the level and grid.
Definition: LatticeElementFinder.cpp:18
LatticeElementFinderDevice GetFinderDeviceInstance(WarpXParIter const &a_pti, int const a_offset, AcceleratorLattice const &accelerator_lattice)
Get the device level instance associated with this instance.
Definition: LatticeElementFinder.cpp:80