WarpX
Loading...
Searching...
No Matches
PhotonCreationFunc.H
Go to the documentation of this file.
1/* Copyright 2024 David Grote
2 *
3 * This file is part of WarpX.
4 *
5 * License: BSD-3-Clause-LBNL
6 */
7
8#ifndef WARPX_PHOTON_CREATION_FUNC_H_
9#define WARPX_PHOTON_CREATION_FUNC_H_
10
12
16#include "Utils/ParticleUtils.H"
17
18#include <AMReX_DenseBins.H>
19#include <AMReX_GpuAtomic.H>
20#include <AMReX_GpuDevice.H>
21#include <AMReX_GpuContainers.H>
22#include <AMReX_INT.H>
23#include <AMReX_Random.H>
24#include <AMReX_REAL.H>
25#include <AMReX_Vector.H>
26
32{
33 // Define shortcuts for frequently-used type names
36 using ParticleTileDataType = typename ParticleTileType::ParticleTileDataType;
39 using SoaData_type = typename WarpXParticleContainer::ParticleTileType::ParticleTileDataType;
40
41public:
45 PhotonCreationFunc () = default;
46
53 PhotonCreationFunc (const std::string& collision_name, MultiParticleContainer const * mypc);
54
90 const index_type& n_total_pairs,
91 ParticleTileType& ptile1, ParticleTileType& ptile2,
93 ParticleTileType** AMREX_RESTRICT tile_products,
94 const amrex::ParticleReal& m1, const amrex::ParticleReal& m2,
95 const amrex::Vector<amrex::ParticleReal>& /*products_mass*/,
96 const index_type* AMREX_RESTRICT p_mask,
97 const amrex::Vector<index_type>& products_np,
98 const SmartCopy* AMREX_RESTRICT copy_species1,
99 const SmartCopy* /*copy_species2*/,
100 const index_type* AMREX_RESTRICT p_pair_indices_1,
101 const index_type* AMREX_RESTRICT p_pair_indices_2,
102 const amrex::ParticleReal* AMREX_RESTRICT p_pair_reaction_weight,
103 const amrex::ParticleReal* AMREX_RESTRICT p_product_data
104 ) const
105 {
106 using namespace amrex::literals;
107
108 if (n_total_pairs == 0) { return amrex::Vector<int>(m_num_product_species, 0); }
109
110 // Compute offset array and allocate memory for the produced species
111 amrex::Gpu::DeviceVector<index_type> offsets(n_total_pairs);
112 const index_type* AMREX_RESTRICT p_offsets = offsets.dataPtr();
113
114 int total;
115 if (m_create_photons) {
116 total = amrex::Scan::ExclusiveSum(n_total_pairs, p_mask, offsets.data());
117 } else {
118 total = 0;
119 }
120
122 for (int i = 0; i < m_num_product_species; i++)
123 {
124 const index_type num_added = total;
125 num_added_vec[i] = static_cast<int>(num_added);
126 tile_products[i]->resize(products_np[i] + num_added);
127 }
128
129 const auto soa_1 = ptile1.getParticleTileData();
130 const auto soa_2 = ptile2.getParticleTileData();
131
132 // Create necessary GPU vectors, that will be used in the kernel below
133 amrex::Vector<SoaData_type> soa_products;
134 for (int i = 0; i < m_num_product_species; i++)
135 {
136 soa_products.push_back(tile_products[i]->getParticleTileData());
137 }
138#ifdef AMREX_USE_GPU
141
143 soa_products.end(),
144 device_soa_products.begin());
146 products_np.end(),
147 device_products_np.begin());
148
150 SoaData_type* AMREX_RESTRICT soa_products_data = device_soa_products.data();
151 const index_type* AMREX_RESTRICT products_np_data = device_products_np.data();
152#else
153 SoaData_type* AMREX_RESTRICT soa_products_data = soa_products.data();
154 const index_type* AMREX_RESTRICT products_np_data = products_np.data();
155#endif
156
157 bool const create_photons = m_create_photons;
158
159 amrex::ParallelForRNG(n_total_pairs,
160 [=] AMREX_GPU_DEVICE (int i, amrex::RandomEngine const& engine) noexcept
161 {
162 if (p_mask[i])
163 {
164 constexpr auto c2 = PhysConst::c * PhysConst::c;
165
166 int const i1 = static_cast<int>(p_pair_indices_1[i]);
167 int const i2 = static_cast<int>(p_pair_indices_2[i]);
168
169 // Move to the ion rest frame and calculate the momentum and energy of the three particles
170 // This is done here to avoid having to save the photon momentum components
171 amrex::ParticleReal & u1x = soa_1.m_rdata[PIdx::ux][i1];
172 amrex::ParticleReal & u1y = soa_1.m_rdata[PIdx::uy][i1];
173 amrex::ParticleReal & u1z = soa_1.m_rdata[PIdx::uz][i1];
174 amrex::ParticleReal const & w1 = soa_1.m_rdata[PIdx::w][i1];
175 amrex::ParticleReal & u2x = soa_2.m_rdata[PIdx::ux][i2];
176 amrex::ParticleReal & u2y = soa_2.m_rdata[PIdx::uy][i2];
177 amrex::ParticleReal & u2z = soa_2.m_rdata[PIdx::uz][i2];
178 amrex::ParticleReal const & w2 = soa_2.m_rdata[PIdx::w][i2];
179
180 amrex::ParticleReal u1x_rel = u1x;
181 amrex::ParticleReal u1y_rel = u1y;
182 amrex::ParticleReal u1z_rel = u1z;
183 ParticleUtils::doLorentzTransformWithU(u1x_rel, u1y_rel, u1z_rel, u2x, u2y, u2z);
184
185 amrex::ParticleReal const u1sq_rel = (u1x_rel*u1x_rel + u1y_rel*u1y_rel + u1z_rel*u1z_rel);
186 amrex::ParticleReal const gamma1_rel = std::sqrt(1.0_prt + u1sq_rel/c2);
187
188 amrex::ParticleReal const Ephoton = p_product_data[i];
189 amrex::ParticleReal const wp = p_pair_reaction_weight[i];
190
191 // Calculate electron and ion energy and momentum using conservation
192 amrex::ParticleReal const u1_rel = std::sqrt(u1sq_rel);
193 amrex::ParticleReal const A = u1_rel/PhysConst::c - Ephoton*wp/(w1*m1*c2);
194 amrex::ParticleReal const B = gamma1_rel + w2*m2/(w1*m1) - Ephoton*wp/(w1*m1*c2);
195
196 // Expanding A**2 - B**2 cancels some terms and reduces the round off error
197 amrex::ParticleReal const AAmBB = -1.0_prt - 2.0_prt*gamma1_rel*w2*m2/(w1*m1) - w2*m2/(w1*m1)*w2*m2/(w1*m1)
198 +2.0_prt*(gamma1_rel + w2*m2/(w1*m1) - u1_rel/PhysConst::c)*Ephoton*wp/(w1*m1*c2);
199 amrex::ParticleReal const D = w2*m2*w2*m2/(w1*m1*w1*m1) + AAmBB - 1.0_prt;
200
201 amrex::ParticleReal const a = 4.0_prt*(-AAmBB);
202 amrex::ParticleReal const b = 4.0_prt*A*D;
203 amrex::ParticleReal const c = 4.0_prt*B*B - D*D;
204
205 // c could also be expanded but it doesn't seem to improve the round off
206 // u = u1_rel/PhysConst::c
207 // E = Ephoton*wp/(w1*m1*c2)
208 // w = w2*m2/(w1*m1)
209 // g = gamma1_rel
210 // cc1 = w*w*(- 4.0_prt*u*u + 8.0_prt*E*g - 4.0_prt*E*E)
211 // cc2 = w*(+ 8.0_prt*E*u*u - 8.0_prt*E*u*g + 8.0_prt*E + 8.0_prt*E*E*u - 8.0_prt*E*E*g)
212 // cc3 = (- 8.0_prt*E*u + 4.0_prt*u*u - 8.0_prt*E*E*u*u + 8.0_prt*E*E*u*g)
213 // c = cc1 + cc2 + cc3
214
215 // Make sure that the root term is always non-negative.
216 amrex::ParticleReal const bac = std::max(b*b - 4.0_prt*a*c, 0._prt);
217
218 // Use the "plus" solution since it keeps the electon motion forward
219 amrex::ParticleReal const u1_rel_after = (-b + std::sqrt(bac))/(2.0_prt*a)*PhysConst::c;
220
221 // The terms inside the sqrt can also be expanded but doesn't seem to improve the round off error
222 // u1_rel_after = (-b + 4.0_prt*std::sqrt(B*B*D*D + 4.0_prt*A*A*B*B - 4.0_prt*B*B*B*B))/(2.0_prt*a)*PhysConst::c
223
224 amrex::ParticleReal const u2_rel_after = (A*PhysConst::c - u1_rel_after)*w1*m1/(w2*m2);
225
226 u1x_rel *= u1_rel_after/u1_rel;
227 u1y_rel *= u1_rel_after/u1_rel;
228 u1z_rel *= u1_rel_after/u1_rel;
229 amrex::ParticleReal u2x_rel = u2_rel_after*u1x_rel/u1_rel_after;
230 amrex::ParticleReal u2y_rel = u2_rel_after*u1y_rel/u1_rel_after;
231 amrex::ParticleReal u2z_rel = u2_rel_after*u1z_rel/u1_rel_after;
232 amrex::ParticleReal p3x_rel = (Ephoton/PhysConst::c)*u1x_rel/u1_rel_after;
233 amrex::ParticleReal p3y_rel = (Ephoton/PhysConst::c)*u1y_rel/u1_rel_after;
234 amrex::ParticleReal p3z_rel = (Ephoton/PhysConst::c)*u1z_rel/u1_rel_after;
235
236 // Lorentz transform electron and ion back to lab frame
237 ParticleUtils::doLorentzTransformWithU(u1x_rel, u1y_rel, u1z_rel, -u2x, -u2y, -u2z);
238 ParticleUtils::doLorentzTransformWithU(u2x_rel, u2y_rel, u2z_rel, -u2x, -u2y, -u2z);
239 ParticleUtils::doLorentzTransformWithP(p3x_rel, p3y_rel, p3z_rel, 0.0_prt, -u2x, -u2y, -u2z);
240 u1x = u1x_rel;
241 u1y = u1y_rel;
242 u1z = u1z_rel;
243 u2x = u2x_rel;
244 u2y = u2y_rel;
245 u2z = u2z_rel;
246
247 if (create_photons) {
248
249 // Create product particle at position of particle 1
250 const auto product_index = products_np_data[0] + p_offsets[i];
251 copy_species1[0](soa_products_data[0], soa_1, i1,
252 static_cast<int>(product_index), engine);
253
254 // Set the weight of the new particle to p_pair_reaction_weight[i]
255 soa_products_data[0].m_rdata[PIdx::w][product_index] = wp;
256
257 amrex::ParticleReal & upx = soa_products_data[0].m_rdata[PIdx::ux][product_index];
258 amrex::ParticleReal & upy = soa_products_data[0].m_rdata[PIdx::uy][product_index];
259 amrex::ParticleReal & upz = soa_products_data[0].m_rdata[PIdx::uz][product_index];
260
261 // Lorentz transform to lab frame (from ion rest frame)
262 upx = p3x_rel/PhysConst::m_e;
263 upy = p3y_rel/PhysConst::m_e;
264 upz = p3z_rel/PhysConst::m_e;
265 }
266 }
267 });
268
269 // Initialize the user runtime components
270 for (int i = 0; i < m_num_product_species; i++)
271 {
272 const auto start_index = int(products_np[i]);
273 const auto stop_index = int(products_np[i] + num_added_vec[i]);
275 0, 0,
276 pc_products[i]->getUserRealAttribs(), pc_products[i]->getUserIntAttribs(),
277 pc_products[i]->GetRealSoANames(), pc_products[i]->GetIntSoANames(),
278 pc_products[i]->getUserRealAttribParser(),
279 pc_products[i]->getUserIntAttribParser(),
280#ifdef WARPX_QED
281 false, // do not initialize QED quantities, since they were initialized
282 // when calling the SmartCopy functors
283 pc_products[i]->get_breit_wheeler_engine_ptr(),
284 pc_products[i]->get_quantum_sync_engine_ptr(),
285#endif
286 pc_products[i]->getIonizationInitialLevel(),
287 start_index, stop_index);
288 }
289
291
292 return num_added_vec;
293 }
294
295private:
296 // How many different type of species the collision produces
298
299 // Whether photons are created
301
303
304};
305
306#endif // WARPX_PHOTON_CREATION_FUNC_H_
#define AMREX_RESTRICT
#define AMREX_INLINE
#define AMREX_GPU_DEVICE
CollisionType
Definition BinaryCollisionUtils.H:17
Definition MultiParticleContainer.H:68
bool m_create_photons
Definition PhotonCreationFunc.H:300
PhotonCreationFunc()=default
Default constructor of the PhotonCreationFunc class.
typename WarpXParticleContainer::ParticleType ParticleType
Definition PhotonCreationFunc.H:34
amrex::DenseBins< ParticleTileDataType > ParticleBins
Definition PhotonCreationFunc.H:37
CollisionType m_collision_type
Definition PhotonCreationFunc.H:302
AMREX_INLINE amrex::Vector< int > operator()(const index_type &n_total_pairs, ParticleTileType &ptile1, ParticleTileType &ptile2, const amrex::Vector< WarpXParticleContainer * > &pc_products, ParticleTileType **AMREX_RESTRICT tile_products, const amrex::ParticleReal &m1, const amrex::ParticleReal &m2, const amrex::Vector< amrex::ParticleReal > &, const index_type *AMREX_RESTRICT p_mask, const amrex::Vector< index_type > &products_np, const SmartCopy *AMREX_RESTRICT copy_species1, const SmartCopy *, const index_type *AMREX_RESTRICT p_pair_indices_1, const index_type *AMREX_RESTRICT p_pair_indices_2, const amrex::ParticleReal *AMREX_RESTRICT p_pair_reaction_weight, const amrex::ParticleReal *AMREX_RESTRICT p_product_data) const
operator() of the PhotonCreationFunc class. It creates new photon particles from binary collisions....
Definition PhotonCreationFunc.H:89
int m_num_product_species
Definition PhotonCreationFunc.H:297
typename WarpXParticleContainer::ParticleTileType ParticleTileType
Definition PhotonCreationFunc.H:35
typename WarpXParticleContainer::ParticleTileType::ParticleTileDataType SoaData_type
Definition PhotonCreationFunc.H:39
typename ParticleBins::index_type index_type
Definition PhotonCreationFunc.H:38
typename ParticleTileType::ParticleTileDataType ParticleTileDataType
Definition PhotonCreationFunc.H:36
iterator begin() noexcept
T * dataPtr() noexcept
T * data() noexcept
void DefaultInitializeRuntimeAttributes(PTile &ptile, const int n_external_attr_real, const int n_external_attr_int, const std::vector< std::string > &user_real_attribs, const std::vector< std::string > &user_int_attribs, const std::vector< std::string > &particle_comps, const std::vector< std::string > &particle_icomps, const std::vector< amrex::Parser * > &user_real_attrib_parser, const std::vector< amrex::Parser * > &user_int_attrib_parser, const bool do_qed_comps, BreitWheelerEngine *p_bw_engine, QuantumSynchrotronEngine *p_qs_engine, const int ionization_initial_level, int start, int stop)
Default initialize runtime attributes in a tile. This routine does not initialize the first n_externa...
Definition DefaultInitialization.H:118
AMREX_GPU_HOST_DEVICE AMREX_INLINE void doLorentzTransformWithP(amrex::ParticleReal &px, amrex::ParticleReal &py, amrex::ParticleReal &pz, amrex::ParticleReal mass, amrex::ParticleReal const Ux, amrex::ParticleReal const Uy, amrex::ParticleReal const Uz)
Perform a Lorentz transformation of the given momentum to a frame moving with gamma*velocity (Ux,...
Definition ParticleUtils.H:159
AMREX_GPU_HOST_DEVICE AMREX_INLINE void doLorentzTransformWithU(amrex::ParticleReal &ux, amrex::ParticleReal &uy, amrex::ParticleReal &uz, amrex::ParticleReal const Ux, amrex::ParticleReal const Uy, amrex::ParticleReal const Uz)
Perform a Lorentz transformation of the given velocity to a frame moving with gamma*velocity (Ux,...
Definition ParticleUtils.H:122
static constexpr auto c
vacuum speed of light [m/s]
Definition constant.H:44
static constexpr auto m_e
electron mass [kg]
Definition constant.H:54
void synchronize() noexcept
void copyAsync(HostToDevice, InIter begin, InIter end, OutIter result) noexcept
static constexpr HostToDevice hostToDevice
void streamSynchronize() noexcept
PODVector< T, ArenaAllocator< T > > DeviceVector
T ExclusiveSum(N n, T const *in, T *out, RetSum a_ret_sum=retSum)
AMREX_ATTRIBUTE_FLATTEN_FOR void ParallelForRNG(T n, L const &f) noexcept
@ uz
Definition WarpXParticleContainer.H:70
@ w
weight
Definition WarpXParticleContainer.H:69
@ uy
Definition WarpXParticleContainer.H:70
@ ux
Definition WarpXParticleContainer.H:70
This is a functor for performing a "smart copy" that works in both host and device code.
Definition SmartCopy.H:34