WarpX
Loading...
Searching...
No Matches
VirtualPhotonCreation.H
Go to the documentation of this file.
1/* Copyright 2025 Arianna Formenti
2 *
3 * This file is part of WarpX.
4 *
5 * License: BSD-3-Clause-LBNL
6 */
7
8#ifndef WARPX_VIRTUAL_PHOTONS_H_
9#define WARPX_VIRTUAL_PHOTONS_H_
10
15
16#include "Utils/TextMsg.H"
18
19#include <AMReX_INT.H>
20#include <AMReX_REAL.H>
21#include <AMReX_Particle.H>
22
23#include <cmath>
24
25
27
28using namespace amrex::literals;
29using SoaData_type = typename WarpXParticleContainer::ParticleTileType::ParticleTileDataType;
30
31
54
55 WARPX_PROFILE("collision::binarycollision::virtualphotons::GenerateVirtualPhotons()");
56
57 // Loop through the species
58 for (int i_s = 0; i_s < mypc->nSpecies(); ++i_s) {
59
60 auto& primary = mypc->GetParticleContainer(i_s);
61
62 if(!primary.has_virtual_photons()){
63 continue;
64 }
65
66 // Get the virtual photon species corresponding to this primary species
67 int vphotons_index = primary.getVirtualPhotonSpeciesIndex();
68 auto& vphotons = mypc->GetParticleContainer(vphotons_index);
69 const amrex::ParmParse pp_species_name(mypc->GetSpeciesNames()[vphotons_index]);
70
71 // Minimum allowed energy of the virtual photons
72 amrex::Real vphoton_min_energy = 0.0_rt;
73 utils::parser::getWithParser(pp_species_name, "qed_virtual_photons_min_energy", vphoton_min_energy);
74
75 // Sampling factor (a.k.a. multiplier):
76 // the number of virtual photons generated is multiplied by this factor,
77 // the weight of each virtual photon is divided by this factor
78 amrex::Real sampling_factor = 0.0_rt;
79 utils::parser::getWithParser(pp_species_name, "qed_virtual_photons_multiplier", sampling_factor);
80
81 amrex::Real const alpha_over_pi = PhysConst::alpha / MathConst::pi;
82 amrex::Real const inv_c2 = 1._rt / (PhysConst::c * PhysConst::c);
83 amrex::Real const mass = primary.getMass();
84
85 int const nlevs = std::max(0, primary.finestLevel()+1);
86 for (int lev = 0; lev < nlevs; ++lev) {
87#ifdef AMREX_USE_OMP
88#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
89#endif
90 for (amrex::MFIter mfi = primary.MakeMFIter(lev); mfi.isValid(); ++mfi)
91 {
92 // Notation: _vp means virtual photon
93 // Primary particles (leptons) in the current tile
94 ParticleTileType& ptile = primary.ParticlesAt(lev, mfi);
95 const auto soa = ptile.getParticleTileData();
96
97 // Number of primary particles in the current tile
98 amrex::Long const num = ptile.numParticles();
99
100 // Vector that will contain the number of virtual photons for each primary particle
102 auto* num_vp_data = num_vp.dataPtr();
103
104 // First pass: compute the number of virtual photons for each primary particle
105 // and fill the corresponding vector
107 [=] AMREX_GPU_DEVICE (amrex::Long i, amrex::RandomEngine const& engine) noexcept
108 {
109 amrex::ParticleReal ux = soa.m_rdata[PIdx::ux][i]; // u=v*gamma=p/m_e
110 amrex::ParticleReal uy = soa.m_rdata[PIdx::uy][i];
111 amrex::ParticleReal uz = soa.m_rdata[PIdx::uz][i];
112
113 // Formula 99.16 in Berestetskii et al., Quantum Electrodynamics
114 // integrated over the photon energies from vphoton_min_energy to the energy of the primary particle
115 // A similar formula is 15.58 in Jackson's, Classical Electrodynamics
116 // but neglect longitudinal field, assume relativistic velocities, and integrate in energy
117 amrex::ParticleReal gamma = std::sqrt( 1.0_rt + (ux*ux + uy*uy + uz*uz) * inv_c2 );
118 // Minimum fractional (w.r.t. the primary) photon energy
119 amrex::ParticleReal y_min = vphoton_min_energy * inv_c2 / (gamma * mass);
120 amrex::ParticleReal lny = std::log( y_min );
121 // Number of virtual photons per primary particle
122 amrex::Real r_photons = alpha_over_pi * lny * lny * sampling_factor;
123
124 // `n_photons` must be an integer, but must average to `r_photons` over many realizations
125 // This is achieved by adding a random number between 0 and 1, and taking the integer part.
126 amrex::Long n_photons = (amrex::Long)( r_photons + amrex::Random(engine) );
127
128 num_vp_data[i] = n_photons;
129 });
130
131 // Compute the offsets vector as the cumulative sum of the elements of num_vp excluding the current element,
132 // i.e., offsets[i] = sum_{j=0}^{i-1} num_vp[j],
133 // and return the total number of virtual photons to be generated in the current tile
134 // (which is the last element of the offsets vector)
136 const amrex::Long total_num_vp = amrex::Scan::ExclusiveSum(num_vp.size(), num_vp.data(), offsets_vp.data());
137 auto *const offset_vp_data = offsets_vp.dataPtr();
138
139 // Now we can allocate and build the virtual photon species in the current tile
140 // Note that this operation will overwrite any virtual photons that were previously generated by mypc
141 // namely the ones that were created in the previous time step.
142 ParticleTileType& ptile_vp = vphotons.ParticlesAt(lev, mfi);
143 ptile_vp.resize(total_num_vp);
144
145 // Get the starting particle ID on CPU and reserve IDs for all virtual photons
146 // This must be done on CPU because NextID() is not thread-safe and cannot be called from GPU
147 amrex::Long pid;
148#ifdef AMREX_USE_OMP
149#pragma omp critical (virtual_photon_nextid)
150#endif
151 {
152 pid = ParticleTileType::ParticleType::NextID();
153 ParticleTileType::ParticleType::NextID(pid + total_num_vp);
154 }
155
156 const int cpuid = amrex::ParallelDescriptor::MyProc();
157
158 // SoA that will contain the virtual photons data
159 auto &soa_vp = ptile_vp.GetStructOfArrays();
160
161 // Array with the PIDs of the virtual photons
162 uint64_t * AMREX_RESTRICT pid_vp = soa_vp.GetIdCPUData().data();
163
164 // Pointers to the arrays that will contain the particle attributes of the virtual photons
166 for (int ia = 0; ia < PIdx::nattribs; ++ia) {
167 pa_vp[ia] = soa_vp.GetRealData(ia).data();
168 }
169
170 // Capture the starting PID for use in the GPU kernel
171 const amrex::Long pid_start = pid;
172
173 // Second pass: populate the virtual photon species
175 [=] AMREX_GPU_DEVICE (amrex::Long i, amrex::RandomEngine const& engine) noexcept
176 {
177 // Primary particle
178 amrex::ParticleReal ux_primary = soa.m_rdata[PIdx::ux][i];
179 amrex::ParticleReal uy_primary = soa.m_rdata[PIdx::uy][i];
180 amrex::ParticleReal uz_primary = soa.m_rdata[PIdx::uz][i];
181 amrex::ParticleReal u_primary = std::sqrt(ux_primary*ux_primary + uy_primary*uy_primary + uz_primary*uz_primary);
182 amrex::ParticleReal nx = ux_primary / u_primary; // normalized ux
183 amrex::ParticleReal ny = uy_primary / u_primary; // normalized uy
184 amrex::ParticleReal nz = uz_primary / u_primary; // normalized uz
185 amrex::ParticleReal gamma_primary = std::sqrt( 1.0_rt + (ux_primary*ux_primary + uy_primary*uy_primary + uz_primary*uz_primary)*inv_c2 );
186
187#if defined (WARPX_DIM_3D)
188 amrex::ParticleReal x = soa.m_rdata[PIdx::x][i];
189 amrex::ParticleReal y = soa.m_rdata[PIdx::y][i];
190 amrex::ParticleReal z = soa.m_rdata[PIdx::z][i];
191#elif defined (WARPX_DIM_XZ)
192 amrex::ParticleReal x = soa.m_rdata[PIdx::x][i];
193 amrex::ParticleReal z = soa.m_rdata[PIdx::z][i];
194#elif defined (WARPX_DIM_RZ)
195 amrex::ParticleReal x = soa.m_rdata[PIdx::x][i];
196 amrex::ParticleReal z = soa.m_rdata[PIdx::z][i];
197 amrex::ParticleReal theta = soa.m_rdata[PIdx::theta][i];
198#elif defined (WARPX_DIM_1D_Z)
199 amrex::ParticleReal z = soa.m_rdata[PIdx::z][i];
200#elif defined (WARPX_DIM_RCYLINDER)
201 amrex::ParticleReal x = soa.m_rdata[PIdx::x][i];
202 amrex::ParticleReal theta = soa.m_rdata[PIdx::theta][i];
203#elif defined(WARPX_DIM_RSPHERE)
204 amrex::ParticleReal x = soa.m_rdata[PIdx::x][i];
205 amrex::ParticleReal theta = soa.m_rdata[PIdx::theta][i];
206 amrex::ParticleReal phi = soa.m_rdata[PIdx::phi][i];
207#endif
208 amrex::ParticleReal w = soa.m_rdata[PIdx::w][i];
209
210 // TODO: add a runtime attribute to the virtual photon species
211 // that containes the pid of the parent particle = soa.m_idcpu[i]
212 // This will allow to update the parent lepton if needed
213
214 // Minimum fractional (wrt primary particle) photon energy
215 amrex::Real y_min = vphoton_min_energy / (mass * gamma_primary * PhysConst::c * PhysConst::c);
216 amrex::Real umin = 0._rt;
217 amrex::Real umax = std::log(y_min) * std::log(y_min);
218
219 for (int j = 0; j < num_vp_data[i]; j++)
220 {
221 // Sample frac_energy from a probability distribution function
222 // that is proportional to log(frac_energy)/frac_energy
223 // (formula 99.16 in Berestetskii et al.)
224 // using the method of the inverse cumulative distributionfunction
225
226 // Draw a random number between umin and umax
227 amrex::ParticleReal rnd = (umax - umin) * amrex::Random(engine) + umin ;
228 // Fractional energy of the photon, often denoted as y (or x)
229 amrex::ParticleReal frac_energy = std::exp( - std::sqrt(rnd) );
230 // Energy of the virtual photon
231 amrex::ParticleReal vphoton_energy = frac_energy * gamma_primary * PhysConst::c;
232
233 // Photon index for the current primary
234 amrex::Long ip = offset_vp_data[i] + j;
235 pa_vp[PIdx::ux][ip] = vphoton_energy * nx; // will be multiplied by m_e before dumping the outputs
236 pa_vp[PIdx::uy][ip] = vphoton_energy * ny; // will be multiplied by m_e before dumping the outputs
237 pa_vp[PIdx::uz][ip] = vphoton_energy * nz; // will be multiplied by m_e before dumping the outputs
238
239 // TODO: add beam size effect - displace the photon position
240#if defined (WARPX_DIM_3D)
241 pa_vp[PIdx::x][ip] = x;
242 pa_vp[PIdx::y][ip] = y;
243 pa_vp[PIdx::z][ip] = z;
244#elif defined (WARPX_DIM_XZ)
245 pa_vp[PIdx::x][ip] = x;
246 pa_vp[PIdx::z][ip] = z;
247#elif defined (WARPX_DIM_RZ)
248 pa_vp[PIdx::x][ip] = x;
249 pa_vp[PIdx::z][ip] = z;
250 pa_vp[PIdx::theta][ip] = theta;
251#elif defined (WARPX_DIM_1D_Z)
252 pa_vp[PIdx::z][ip] = z;
253#elif defined (WARPX_DIM_RCYLINDER)
254 pa_vp[PIdx::x][ip] = x;
255 pa_vp[PIdx::theta][ip] = theta;
256#elif defined(WARPX_DIM_RSPHERE)
257 pa_vp[PIdx::x][ip] = x;
258 pa_vp[PIdx::theta][ip] = theta;
259 pa_vp[PIdx::phi][ip] = phi;
260#endif
261 pa_vp[PIdx::w][ip] = w / sampling_factor;
262 pid_vp[ip] = amrex::SetParticleIDandCPU(pid_start + ip, cpuid);
263 }
264 });
265 } // mfi
266 } // lev
267 } // species
268} // function
269} // close namespace
270#endif // WARPX_VIRTUAL_PHOTONS_H_
#define AMREX_RESTRICT
#define AMREX_GPU_DEVICE
#define WARPX_PROFILE(fname)
Definition WarpXProfilerWrapper.H:13
Definition MultiParticleContainer.H:69
int nSpecies() const
Definition MultiParticleContainer.H:267
WarpXParticleContainer & GetParticleContainer(int index) const
Definition MultiParticleContainer.H:83
std::vector< std::string > GetSpeciesNames() const
Definition MultiParticleContainer.H:320
virtual int getVirtualPhotonSpeciesIndex() const
Definition WarpXParticleContainer.H:537
size_type size() const noexcept
T * dataPtr() noexcept
T * data() noexcept
amrex_real Real
amrex_particle_real ParticleReal
amrex_long Long
T ExclusiveSum(N n, T const *in, T *out, RetSum a_ret_sum=retSum)
PODVector< T, ArenaAllocator< T > > DeviceVector
Real Random()
constexpr auto alpha
fine-structure constant = mu0/(4*pi)*q_e*q_e*c/hbar [dimensionless]
Definition constant.H:177
constexpr auto c
vacuum speed of light [m/s]
Definition constant.H:153
constexpr auto pi
ratio of a circle's circumference to its diameter
Definition constant.H:29
__host__ __device__ std::uint64_t SetParticleIDandCPU(Long id, int cpu) noexcept
AMREX_ATTRIBUTE_FLATTEN_FOR void ParallelForRNG(T n, L const &f) noexcept
Definition VirtualPhotonCreation.H:26
void GenerateVirtualPhotons(MultiParticleContainer *mypc)
Definition VirtualPhotonCreation.H:53
typename WarpXParticleContainer::ParticleTileType::ParticleTileDataType SoaData_type
Definition VirtualPhotonCreation.H:29
void getWithParser(const amrex::ParmParse &a_pp, char const *const str, T &val)
Definition ParserUtils.H:166
@ nattribs
Definition WarpXParticleContainer.H:70
@ x
Definition WarpXParticleContainer.H:70
@ uz
Definition WarpXParticleContainer.H:70
@ w
Definition WarpXParticleContainer.H:70
@ uy
Definition WarpXParticleContainer.H:70
@ z
Definition WarpXParticleContainer.H:70
@ ux
Definition WarpXParticleContainer.H:70
__host__ __device__ const T * data() const noexcept
ParticleTileDataType getParticleTileData()
void resize(std::size_t count, GrowthStrategy strategy=GrowthStrategy::Poisson)