WarpX
Loading...
Searching...
No Matches
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
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
38 void InitElementFinder (int lev, amrex::Real gamma_boost,
40 amrex::MFIter const& a_mfi,
41 AcceleratorLattice const& accelerator_lattice);
42
48 void AllocateIndices (AcceleratorLattice const& accelerator_lattice);
49
58 void UpdateIndices (int lev, amrex::MFIter const& a_mfi,
59 AcceleratorLattice const& accelerator_lattice,
60 const amrex::Vector<amrex::Real>& time);
61
62 /* Define the location and size of the index lookup table */
63 /* Use the type Real to be consistent with the way the main grid is defined */
64 int m_nz;
65 amrex::Real m_zmin;
66 amrex::Real m_dz;
67
68 /* Parameters needed for the Lorentz transforms into and out of the boosted frame */
69 /* The time for m_time is consistent with the main time variable */
70 amrex::ParticleReal m_gamma_boost;
71 amrex::ParticleReal m_uz_boost;
72 amrex::Real m_time;
73
83 WarpXParIter const& a_pti, int a_offset, AcceleratorLattice const& accelerator_lattice,
84 const amrex::Vector<amrex::Real>& dts) const;
85
86 /* The index lookup tables for each lattice element type */
89
100 amrex::Gpu::DeviceVector<int> & indices) const;
101};
102
108{
109
119 void
120 InitLatticeElementFinderDevice (WarpXParIter const& a_pti, int a_offset,
121 AcceleratorLattice const& accelerator_lattice,
122 LatticeElementFinder const & h_finder,
123 const amrex::Vector<amrex::Real>& dts);
124
125 /* Whether the class has been initialized */
126 bool m_initialized = false;
127
128 /* Size and location of the index lookup table */
129 amrex::Real m_zmin;
130 amrex::Real m_dz;
131 amrex::Real m_dt;
132
133 /* Parameters needed for the Lorentz transforms into and out of the boosted frame */
134 amrex::ParticleReal m_gamma_boost;
135 amrex::ParticleReal m_uz_boost;
136 amrex::Real m_time;
137
139 const amrex::ParticleReal* AMREX_RESTRICT m_ux = nullptr;
140 const amrex::ParticleReal* AMREX_RESTRICT m_uy = nullptr;
141 const amrex::ParticleReal* AMREX_RESTRICT m_uz = nullptr;
142
143 /* Device level instances for each lattice element type */
146
147 /* Device level index lookup tables for each element type */
148 int const* d_quad_indices_arr = nullptr;
149 int const* d_plasmalens_indices_arr = nullptr;
150
158 void operator () (const long i,
159 amrex::ParticleReal& field_Ex,
160 amrex::ParticleReal& field_Ey,
161 amrex::ParticleReal& field_Ez,
162 amrex::ParticleReal& field_Bx,
163 amrex::ParticleReal& field_By,
164 amrex::ParticleReal& field_Bz) const noexcept
165 {
166
167 using namespace amrex::literals;
168
169 amrex::ParticleReal x, y, z;
170 m_get_position(i, x, y, z);
171
172 // Find location of partice in the indices grid
173 // (which is in the boosted frame)
174 const int iz = static_cast<int>((z - m_zmin)/m_dz);
175
176 constexpr amrex::ParticleReal inv_c2 = 1._prt/(PhysConst::c*PhysConst::c);
177 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);
178 amrex::ParticleReal const vzp = m_uz[i]/gamma;
179
180 amrex::ParticleReal zpvdt = z + vzp*m_dt;
181
182 // The position passed to the get_field methods needs to be in the lab frame.
183 if (m_gamma_boost > 1._prt) {
185 zpvdt = m_gamma_boost*zpvdt + m_uz_boost*(m_time + m_dt);
186 }
187
188 amrex::ParticleReal Ex_sum = 0._prt;
189 amrex::ParticleReal Ey_sum = 0._prt;
190 const amrex::ParticleReal Ez_sum = 0._prt;
191 amrex::ParticleReal Bx_sum = 0._prt;
192 amrex::ParticleReal By_sum = 0._prt;
193 const amrex::ParticleReal Bz_sum = 0._prt;
194
195 if (d_quad.nelements > 0) {
196 if (d_quad_indices_arr[iz] > -1) {
197 const auto ielement = d_quad_indices_arr[iz];
198 amrex::ParticleReal Ex, Ey, Bx, By;
199 d_quad.get_field(ielement, x, y, z, zpvdt, Ex, Ey, Bx, By);
200 Ex_sum += Ex;
201 Ey_sum += Ey;
202 Bx_sum += Bx;
203 By_sum += By;
204 }
205 }
206
207 if (d_plasmalens.nelements > 0) {
208 if (d_plasmalens_indices_arr[iz] > -1) {
209 const auto ielement = d_plasmalens_indices_arr[iz];
210 amrex::ParticleReal Ex, Ey, Bx, By;
211 d_plasmalens.get_field(ielement, x, y, z, zpvdt, Ex, Ey, Bx, By);
212 Ex_sum += Ex;
213 Ey_sum += Ey;
214 Bx_sum += Bx;
215 By_sum += By;
216 }
217 }
218
219 if (m_gamma_boost > 1._prt) {
220 // The fields returned from get_field is in the lab frame
221 // Transform the fields to the boosted frame
222 const amrex::ParticleReal Ex_boost = m_gamma_boost*Ex_sum - m_uz_boost*By_sum;
223 const amrex::ParticleReal Ey_boost = m_gamma_boost*Ey_sum + m_uz_boost*Bx_sum;
224 const amrex::ParticleReal Bx_boost = m_gamma_boost*Bx_sum + m_uz_boost*Ey_sum*inv_c2;
225 const amrex::ParticleReal By_boost = m_gamma_boost*By_sum - m_uz_boost*Ex_sum*inv_c2;
226 Ex_sum = Ex_boost;
227 Ey_sum = Ey_boost;
228 Bx_sum = Bx_boost;
229 By_sum = By_boost;
230 }
231
232 field_Ex += Ex_sum;
233 field_Ey += Ey_sum;
234 field_Ez += Ez_sum;
235 field_Bx += Bx_sum;
236 field_By += By_sum;
237 field_Bz += Bz_sum;
238
239 }
240
241};
242
243#endif // WARPX_ACCELERATORLATTICE_LATTICEELEMENTS_LATTICEELEMENTFINDER_H_
#define AMREX_FORCE_INLINE
#define AMREX_RESTRICT
#define AMREX_GPU_HOST_DEVICE
Definition AcceleratorLattice.H:21
Definition WarpXParticleContainer.H:117
static constexpr auto c
vacuum speed of light [m/s]
Definition constant.H:44
PODVector< T, ArenaAllocator< T > > DeviceVector
Functor that can be used to extract the positions of the macroparticles inside a ParallelFor kernel.
Definition GetAndSetPosition.H:75
Definition HardEdgedPlasmaLens.H:64
Definition HardEdgedQuadrupole.H:64
The lattice element finder class that can be trivially copied to the device. This only has simple dat...
Definition LatticeElementFinder.H:108
const amrex::ParticleReal *AMREX_RESTRICT m_ux
Definition LatticeElementFinder.H:139
const amrex::ParticleReal *AMREX_RESTRICT m_uy
Definition LatticeElementFinder.H:140
int const * d_plasmalens_indices_arr
Definition LatticeElementFinder.H:149
amrex::Real m_time
Definition LatticeElementFinder.H:136
HardEdgedPlasmaLensDevice d_plasmalens
Definition LatticeElementFinder.H:145
amrex::Real m_dt
Definition LatticeElementFinder.H:131
bool m_initialized
Definition LatticeElementFinder.H:126
void InitLatticeElementFinderDevice(WarpXParIter const &a_pti, int a_offset, AcceleratorLattice const &accelerator_lattice, LatticeElementFinder const &h_finder, const amrex::Vector< amrex::Real > &dts)
Initialize the data needed to do the lookups.
Definition LatticeElementFinder.cpp:96
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void operator()(const long i, amrex::ParticleReal &field_Ex, amrex::ParticleReal &field_Ey, amrex::ParticleReal &field_Ez, amrex::ParticleReal &field_Bx, amrex::ParticleReal &field_By, amrex::ParticleReal &field_Bz) const noexcept
Gather the field for the particle from the lattice elements.
Definition LatticeElementFinder.H:158
amrex::Real m_dz
Definition LatticeElementFinder.H:130
int const * d_quad_indices_arr
Definition LatticeElementFinder.H:148
amrex::ParticleReal m_uz_boost
Definition LatticeElementFinder.H:135
const amrex::ParticleReal *AMREX_RESTRICT m_uz
Definition LatticeElementFinder.H:141
amrex::Real m_zmin
Definition LatticeElementFinder.H:129
HardEdgedQuadrupoleDevice d_quad
Definition LatticeElementFinder.H:144
GetParticlePosition< PIdx > m_get_position
Definition LatticeElementFinder.H:138
amrex::ParticleReal m_gamma_boost
Definition LatticeElementFinder.H:134
Definition LatticeElementFinder.H:27
void setup_lattice_indices(amrex::Gpu::DeviceVector< amrex::ParticleReal > const &zs, amrex::Gpu::DeviceVector< amrex::ParticleReal > const &ze, amrex::Gpu::DeviceVector< int > &indices) const
Fill in the index lookup tables This loops over the grid (in z) and finds the lattice element closest...
Definition LatticeElementFinder.cpp:133
amrex::Gpu::DeviceVector< int > d_quad_indices
Definition LatticeElementFinder.H:87
int m_nz
Definition LatticeElementFinder.H:64
amrex::Real m_zmin
Definition LatticeElementFinder.H:65
LatticeElementFinderDevice GetFinderDeviceInstance(WarpXParIter const &a_pti, int a_offset, AcceleratorLattice const &accelerator_lattice, const amrex::Vector< amrex::Real > &dts) const
Get the device level instance associated with this instance.
Definition LatticeElementFinder.cpp:86
amrex::Gpu::DeviceVector< int > d_plasmalens_indices
Definition LatticeElementFinder.H:88
amrex::ParticleReal m_gamma_boost
Definition LatticeElementFinder.H:70
amrex::Real m_time
Definition LatticeElementFinder.H:72
void UpdateIndices(int lev, amrex::MFIter const &a_mfi, AcceleratorLattice const &accelerator_lattice, const amrex::Vector< amrex::Real > &time)
Update the index lookup tables for each element type, filling in the values.
Definition LatticeElementFinder.cpp:61
amrex::Real m_dz
Definition LatticeElementFinder.H:66
amrex::ParticleReal m_uz_boost
Definition LatticeElementFinder.H:71
void InitElementFinder(int lev, amrex::Real gamma_boost, const amrex::Vector< amrex::Real > &time, 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:46