ARTEMIS
LatticeElementFinder.H
Go to the documentation of this file.
1 /* Copyright 2022 David Grote
2  *
3  * This file is part of WarpX.
4  *
5  * License: BSD-3-Clause-LBNL
6  */
7 #ifndef WARPX_ACCELERATORLATTICE_LATTICEELEMENTS_LATTICEELEMENTFINDER_H_
8 #define WARPX_ACCELERATORLATTICE_LATTICEELEMENTS_LATTICEELEMENTFINDER_H_
9 
14 
15 #include <AMReX_REAL.H>
16 #include <AMReX_GpuContainers.H>
17 
18 class AcceleratorLattice;
20 
21 // Instances of the LatticeElementFinder class are saved in the AcceleratorLattice class
22 // as the objects in a LayoutData.
23 // The LatticeElementFinder handles the lookup needed to find the lattice elements at
24 // particle locations.
25 
27 {
28 
36  void InitElementFinder (int const lev, amrex::MFIter const& a_mfi,
37  AcceleratorLattice const& accelerator_lattice);
38 
44  void AllocateIndices (AcceleratorLattice const& accelerator_lattice);
45 
53  void UpdateIndices (int const lev, amrex::MFIter const& a_mfi,
54  AcceleratorLattice const& accelerator_lattice);
55 
56  /* Define the location and size of the index lookup table */
57  /* Use the type Real to be consistent with the way the main grid is defined */
58  int m_nz;
59  amrex::Real m_zmin;
60  amrex::Real m_dz;
61 
62  /* Parameters needed for the Lorentz transforms into and out of the boosted frame */
63  /* The time for m_time is consistent with the main time variable */
64  amrex::ParticleReal m_gamma_boost;
65  amrex::ParticleReal m_uz_boost;
66  amrex::Real m_time;
67 
75  LatticeElementFinderDevice GetFinderDeviceInstance (WarpXParIter const& a_pti, int const a_offset,
76  AcceleratorLattice const& accelerator_lattice);
77 
78  /* The index lookup tables for each lattice element type */
81 
93  {
94 
95  using namespace amrex::literals;
96 
97  int nelements = static_cast<int>(zs.size());
98  amrex::ParticleReal const * zs_arr = zs.data();
99  amrex::ParticleReal const * ze_arr = ze.data();
100  int * indices_arr = indices.data();
101 
102  amrex::Real const zmin = m_zmin;
103  amrex::Real const dz = m_dz;
104 
105  amrex::ParticleReal const gamma_boost = m_gamma_boost;
106  amrex::ParticleReal const uz_boost = m_uz_boost;
107  amrex::Real const time = m_time;
108 
110  [=] AMREX_GPU_DEVICE (int iz) {
111 
112  // Get the location of the grid node
113  amrex::Real z_node = zmin + iz*dz;
114 
115  if (gamma_boost > 1._prt) {
116  // Transform to lab frame
117  z_node = gamma_boost*z_node + uz_boost*time;
118  }
119 
120  // Find the index to the element that is closest to the grid cell.
121  // For now, this assumes that there is no overlap among elements of the same type.
122  for (int ie = 0 ; ie < nelements ; ie++) {
123  // Find the mid points between element ie and the ones before and after it.
124  // The first and last element need special handling.
125  amrex::ParticleReal zcenter_left, zcenter_right;
126  if (ie == 0) {
127  zcenter_left = std::numeric_limits<amrex::ParticleReal>::lowest();
128  } else {
129  zcenter_left = 0.5_prt*(ze_arr[ie-1] + zs_arr[ie]);
130  }
131  if (ie < nelements - 1) {
132  zcenter_right = 0.5_prt*(ze_arr[ie] + zs_arr[ie+1]);
133  } else {
134  zcenter_right = std::numeric_limits<amrex::ParticleReal>::max();
135  }
136 
137  if (zcenter_left <= z_node && z_node < zcenter_right) {
138  indices_arr[iz] = ie;
139  }
140 
141  }
142  }
143  );
144  }
145 
146 };
147 
153 {
154 
163  void
164  InitLatticeElementFinderDevice (WarpXParIter const& a_pti, int const a_offset,
165  AcceleratorLattice const& accelerator_lattice,
166  LatticeElementFinder const & h_finder);
167 
168  /* Size and location of the index lookup table */
169  amrex::Real m_zmin;
170  amrex::Real m_dz;
171  amrex::Real m_dt;
172 
173  /* Parameters needed for the Lorentz transforms into and out of the boosted frame */
174  amrex::ParticleReal m_gamma_boost;
175  amrex::ParticleReal m_uz_boost;
176  amrex::Real m_time;
177 
179  const amrex::ParticleReal* AMREX_RESTRICT m_ux = nullptr;
180  const amrex::ParticleReal* AMREX_RESTRICT m_uy = nullptr;
181  const amrex::ParticleReal* AMREX_RESTRICT m_uz = nullptr;
182 
183  /* Device level instances for each lattice element type */
186 
187  /* Device level index lookup tables for each element type */
188  int const* d_quad_indices_arr = nullptr;
189  int const* d_plasmalens_indices_arr = nullptr;
190 
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
205  {
206 
207  using namespace amrex::literals;
208 
209  amrex::ParticleReal x, y, z;
210  m_get_position(i, x, y, z);
211 
212  // Find location of partice in the indices grid
213  // (which is in the boosted frame)
214  const int iz = static_cast<int>((z - m_zmin)/m_dz);
215 
216  constexpr amrex::ParticleReal inv_c2 = 1._prt/(PhysConst::c*PhysConst::c);
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;
219 
220  amrex::ParticleReal zpvdt = z + vzp*m_dt;
221 
222  // The position passed to the get_field methods needs to be in the lab frame.
223  if (m_gamma_boost > 1._prt) {
225  zpvdt = m_gamma_boost*zpvdt + m_uz_boost*(m_time + m_dt);
226  }
227 
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;
234 
235  if (d_quad.nelements > 0) {
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);
240  Ex_sum += Ex;
241  Ey_sum += Ey;
242  Bx_sum += Bx;
243  By_sum += By;
244  }
245  }
246 
247  if (d_plasmalens.nelements > 0) {
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);
252  Ex_sum += Ex;
253  Ey_sum += Ey;
254  Bx_sum += Bx;
255  By_sum += By;
256  }
257  }
258 
259  if (m_gamma_boost > 1._prt) {
260  // The fields returned from get_field is in the lab frame
261  // Transform the fields to the boosted frame
262  const amrex::ParticleReal Ex_boost = m_gamma_boost*Ex_sum - m_uz_boost*By_sum;
263  const amrex::ParticleReal Ey_boost = m_gamma_boost*Ey_sum + m_uz_boost*Bx_sum;
264  const amrex::ParticleReal Bx_boost = m_gamma_boost*Bx_sum + m_uz_boost*Ey_sum*inv_c2;
265  const amrex::ParticleReal By_boost = m_gamma_boost*By_sum - m_uz_boost*Ex_sum*inv_c2;
266  Ex_sum = Ex_boost;
267  Ey_sum = Ey_boost;
268  Bx_sum = Bx_boost;
269  By_sum = By_boost;
270  }
271 
272  field_Ex += Ex_sum;
273  field_Ey += Ey_sum;
274  field_Ez += Ez_sum;
275  field_Bx += Bx_sum;
276  field_By += By_sum;
277  field_Bz += Bz_sum;
278 
279  }
280 
281 };
282 
283 #endif // WARPX_ACCELERATORLATTICE_LATTICEELEMENTS_LATTICEELEMENTFINDER_H_
#define AMREX_FORCE_INLINE
#define AMREX_RESTRICT
#define AMREX_GPU_DEVICE
#define AMREX_GPU_HOST_DEVICE
Definition: AcceleratorLattice.H:21
Definition: WarpXParticleContainer.H:52
size_type size() const noexcept
T * data() 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