WarpX
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 lev, amrex::MFIter const& a_mfi,
37  AcceleratorLattice const& accelerator_lattice);
38 
44  void AllocateIndices (AcceleratorLattice const& accelerator_lattice);
45 
53  void UpdateIndices (int 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 
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  const auto 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  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;
131  }
132 
133  }
134  }
135  );
136  }
137 
138 };
139 
145 {
146 
155  void
156  InitLatticeElementFinderDevice (WarpXParIter const& a_pti, int a_offset,
157  AcceleratorLattice const& accelerator_lattice,
158  LatticeElementFinder const & h_finder);
159 
160  /* Size and location of the index lookup table */
161  amrex::Real m_zmin;
162  amrex::Real m_dz;
163  amrex::Real m_dt;
164 
165  /* Parameters needed for the Lorentz transforms into and out of the boosted frame */
166  amrex::ParticleReal m_gamma_boost;
167  amrex::ParticleReal m_uz_boost;
168  amrex::Real m_time;
169 
171  const amrex::ParticleReal* AMREX_RESTRICT m_ux = nullptr;
172  const amrex::ParticleReal* AMREX_RESTRICT m_uy = nullptr;
173  const amrex::ParticleReal* AMREX_RESTRICT m_uz = nullptr;
174 
175  /* Device level instances for each lattice element type */
178 
179  /* Device level index lookup tables for each element type */
180  int const* d_quad_indices_arr = nullptr;
181  int const* d_plasmalens_indices_arr = nullptr;
182 
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
197  {
198 
199  using namespace amrex::literals;
200 
201  amrex::ParticleReal x, y, z;
202  m_get_position(i, x, y, z);
203 
204  // Find location of partice in the indices grid
205  // (which is in the boosted frame)
206  const int iz = static_cast<int>((z - m_zmin)/m_dz);
207 
208  constexpr amrex::ParticleReal inv_c2 = 1._prt/(PhysConst::c*PhysConst::c);
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;
211 
212  amrex::ParticleReal zpvdt = z + vzp*m_dt;
213 
214  // The position passed to the get_field methods needs to be in the lab frame.
215  if (m_gamma_boost > 1._prt) {
217  zpvdt = m_gamma_boost*zpvdt + m_uz_boost*(m_time + m_dt);
218  }
219 
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;
226 
227  if (d_quad.nelements > 0) {
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);
232  Ex_sum += Ex;
233  Ey_sum += Ey;
234  Bx_sum += Bx;
235  By_sum += By;
236  }
237  }
238 
239  if (d_plasmalens.nelements > 0) {
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);
244  Ex_sum += Ex;
245  Ey_sum += Ey;
246  Bx_sum += Bx;
247  By_sum += By;
248  }
249  }
250 
251  if (m_gamma_boost > 1._prt) {
252  // The fields returned from get_field is in the lab frame
253  // Transform the fields to the boosted frame
254  const amrex::ParticleReal Ex_boost = m_gamma_boost*Ex_sum - m_uz_boost*By_sum;
255  const amrex::ParticleReal Ey_boost = m_gamma_boost*Ey_sum + m_uz_boost*Bx_sum;
256  const amrex::ParticleReal Bx_boost = m_gamma_boost*Bx_sum + m_uz_boost*Ey_sum*inv_c2;
257  const amrex::ParticleReal By_boost = m_gamma_boost*By_sum - m_uz_boost*Ex_sum*inv_c2;
258  Ex_sum = Ex_boost;
259  Ey_sum = Ey_boost;
260  Bx_sum = Bx_boost;
261  By_sum = By_boost;
262  }
263 
264  field_Ex += Ex_sum;
265  field_Ey += Ey_sum;
266  field_Ez += Ez_sum;
267  field_Bx += Bx_sum;
268  field_By += By_sum;
269  field_Bz += Bz_sum;
270 
271  }
272 
273 };
274 
275 #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
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