WarpX
Loading...
Searching...
No Matches
LinearBreitWheelerCollisionFunc.H
Go to the documentation of this file.
1/* Copyright 2023 Arianna Formenti
2 *
3 * This file is part of WarpX.
4 *
5 * License: BSD-3-Clause-LBNL
6 */
7
8#ifndef LINEAR_BREIT_WHEELER_COLLISION_FUNC_H_
9#define LINEAR_BREIT_WHEELER_COLLISION_FUNC_H_
10
12
18#include "Utils/TextMsg.H"
19
20#include <AMReX_Algorithm.H>
21#include <AMReX_DenseBins.H>
22#include <AMReX_ParmParse.H>
23#include <AMReX_Random.H>
24#include <AMReX_REAL.H>
25#include <AMReX_Vector.H>
26
36{
37 // Define shortcuts for frequently-used type names
40 using ParticleTileDataType = ParticleTileType::ParticleTileDataType;
43 using SoaData_type = WarpXParticleContainer::ParticleTileType::ParticleTileDataType;
44
45public:
50
58 const std::string& collision_name,
59 MultiParticleContainer const * const /*mypc*/,
60 const bool isSameSpecies) :
61 m_isSameSpecies{isSameSpecies}
62 {
63 using namespace amrex::literals;
64
65#ifdef AMREX_SINGLE_PRECISION_PARTICLES
66 WARPX_ABORT_WITH_MESSAGE("Linear Breit-Wheeler collision module does not currently work with single precision");
67#endif
68
69 const amrex::ParmParse pp_collision_name(collision_name);
70
71 constexpr auto default_event_multiplier = 1.0_prt;
72 m_event_multiplier = default_event_multiplier;
74 pp_collision_name, "event_multiplier", m_event_multiplier);
75
76 constexpr auto default_probability_threshold = 0.02_prt;
77 m_probability_threshold = default_probability_threshold;
79 pp_collision_name, "probability_threshold", m_probability_threshold);
80
81 constexpr auto default_probability_target_value = 0.002_prt;
82 m_probability_target_value = default_probability_target_value;
84 pp_collision_name, "probability_target_value",
86
87 m_exe.m_event_multiplier = m_event_multiplier;
88 m_exe.m_probability_threshold = m_probability_threshold;
89 m_exe.m_probability_target_value = m_probability_target_value;
90 m_exe.m_isSameSpecies = m_isSameSpecies;
91 }
92
93 struct Executor {
137 index_type const I1s, index_type const I1e,
138 index_type const I2s, index_type const I2e,
139 index_type const* AMREX_RESTRICT I1,
140 index_type const* AMREX_RESTRICT I2,
141 const SoaData_type& soa_1, const SoaData_type& soa_2,
142 GetParticlePosition<PIdx> /*get_position_1*/, GetParticlePosition<PIdx> /*get_position_2*/,
143 amrex::ParticleReal const /*n1*/, amrex::ParticleReal const /*n2*/,
144 amrex::ParticleReal const /*T1*/, amrex::ParticleReal const /*T2*/,
145 amrex::Real const /*global_lamdb*/,
146 amrex::ParticleReal const /*q1*/, amrex::ParticleReal const /*q2*/,
147 amrex::ParticleReal const /*m1*/, amrex::ParticleReal const /*m2*/,
148 amrex::Real const dt, amrex::Real const dV, index_type coll_idx,
149 index_type const cell_start_pair, index_type* AMREX_RESTRICT p_mask,
150 index_type* AMREX_RESTRICT p_pair_indices_1, index_type* AMREX_RESTRICT p_pair_indices_2,
151 amrex::ParticleReal* AMREX_RESTRICT p_pair_reaction_weight,
152 amrex::ParticleReal* /*p_product_data*/,
153 amrex::RandomEngine const& engine) const
154 {
155 amrex::ParticleReal * const AMREX_RESTRICT w1 = soa_1.m_rdata[PIdx::w];
156 amrex::ParticleReal * const AMREX_RESTRICT u1x = soa_1.m_rdata[PIdx::ux];
157 amrex::ParticleReal * const AMREX_RESTRICT u1y = soa_1.m_rdata[PIdx::uy];
158 amrex::ParticleReal * const AMREX_RESTRICT u1z = soa_1.m_rdata[PIdx::uz];
159 uint64_t* AMREX_RESTRICT idcpu1 = soa_1.m_idcpu;
160
161 amrex::ParticleReal * const AMREX_RESTRICT w2 = soa_2.m_rdata[PIdx::w];
162 amrex::ParticleReal * const AMREX_RESTRICT u2x = soa_2.m_rdata[PIdx::ux];
163 amrex::ParticleReal * const AMREX_RESTRICT u2y = soa_2.m_rdata[PIdx::uy];
164 amrex::ParticleReal * const AMREX_RESTRICT u2z = soa_2.m_rdata[PIdx::uz];
165 uint64_t* AMREX_RESTRICT idcpu2 = soa_2.m_idcpu;
166
167 // Number of macroparticles of each species
168 const index_type NI1 = I1e - I1s;
169 const index_type NI2 = I2e - I2s;
170 const index_type max_N = amrex::max(NI1,NI2);
171 const index_type min_N = amrex::min(NI1,NI2);
172
173 index_type pair_index = cell_start_pair + coll_idx;
174
175 // multiplier ratio to take into account unsampled pairs
176 const auto multiplier_ratio = static_cast<int>(
177 m_isSameSpecies ? min_N + max_N - 1 : min_N);
178
179#if defined(WARPX_DIM_RZ) || defined(WARPX_DIM_RCYLINDER)
180 amrex::ParticleReal * const AMREX_RESTRICT theta1 = soa_1.m_rdata[PIdx::theta];
181 amrex::ParticleReal * const AMREX_RESTRICT theta2 = soa_2.m_rdata[PIdx::theta];
182#endif
183 index_type i1 = I1s + coll_idx;
184 index_type i2 = I2s + coll_idx;
185 // we will start from collision number = coll_idx and then add
186 // stride (smaller set size) until we do all collisions (larger set size)
187 for (index_type k = coll_idx; k < max_N; k += min_N)
188 {
189 // do not check for collision if a particle's weight was
190 // reduced to zero from a previous collision
191 if (idcpu1[ I1[i1] ]==amrex::ParticleIdCpus::Invalid ||
192 idcpu2[ I2[i2] ]==amrex::ParticleIdCpus::Invalid ) {
193 p_mask[pair_index] = false;
194 }
195 else {
196#if defined(WARPX_DIM_RZ) || defined(WARPX_DIM_RCYLINDER)
197 /* In RZ and RCYLINDER geometry, macroparticles can collide with other macroparticles
198 * in the same *cylindrical* cell. For this reason, collisions between macroparticles
199 * are actually not local in space. In this case, the underlying assumption is that
200 * particles within the same cylindrical cell represent a cylindrically-symmetry
201 * momentum distribution function. Therefore, here, we temporarily rotate the
202 * momentum of one of the macroparticles in agreement with this cylindrical symmetry.
203 * (This is technically only valid if we use only the m=0 azimuthal mode in the simulation;
204 * there is a corresponding assert statement at initialization.) */
205 amrex::ParticleReal const theta = theta2[I2[i2]]-theta1[I1[i1]];
206 amrex::ParticleReal const u1xbuf = u1x[I1[i1]];
207 u1x[I1[i1]] = u1xbuf*std::cos(theta) - u1y[I1[i1]]*std::sin(theta);
208 u1y[I1[i1]] = u1xbuf*std::sin(theta) + u1y[I1[i1]]*std::cos(theta);
209#endif
211 u1x[ I1[i1] ], u1y[ I1[i1] ], u1z[ I1[i1] ],
212 u2x[ I2[i2] ], u2y[ I2[i2] ], u2z[ I2[i2] ],
213 w1[ I1[i1] ], w2[ I2[i2] ],
214 dt, dV, static_cast<int>(pair_index), p_mask, p_pair_reaction_weight,
215 m_event_multiplier, multiplier_ratio,
218 engine);
219#if defined(WARPX_DIM_RZ) || defined(WARPX_DIM_RCYLINDER)
220 amrex::ParticleReal const u1xbuf_new = u1x[I1[i1]];
221 u1x[I1[i1]] = u1xbuf_new*std::cos(-theta) - u1y[I1[i1]]*std::sin(-theta);
222 u1y[I1[i1]] = u1xbuf_new*std::sin(-theta) + u1y[I1[i1]]*std::cos(-theta);
223#endif
224 // Remove pair reaction weight from the colliding particles' weights
225 if (p_mask[pair_index]) {
227 w1[ I1[i1] ], idcpu1[ I1[i1] ], p_pair_reaction_weight[pair_index]);
229 w2[ I2[i2] ], idcpu2[ I2[i2] ], p_pair_reaction_weight[pair_index]);
230 }
231 }
232
233 p_pair_indices_1[pair_index] = I1[i1];
234 p_pair_indices_2[pair_index] = I2[i2];
235 if (max_N == NI1) { i1 += min_N; }
236 if (max_N == NI2) { i2 += min_N; }
237 pair_index += min_N;
238 }
239 }
240
248 };
249
250 [[nodiscard]] Executor const& executor () const { return m_exe; }
251
252 bool use_global_debye_length() { return false; }
253
254private:
255 // Factor used to increase the number of collisions by decreasing the weight of the
256 // produced particles
258 // If the event multiplier is too high and results in a collision probability that approaches
259 // 1, there is a risk of underestimating the total yield of electrons and positrons.
260 // In these cases, we reduce the event multiplier used in a given collision.
261 // m_probability_threshold is the probability above which we reduce the event multiplier.
262 // m_probability_target_value is the target probability used to determine by how much
263 // the event multiplier should be reduced.
267
269};
270
271#endif // LINEAR_BREIT_WHEELER_COLLISION_FUNC_H_
#define AMREX_RESTRICT
#define AMREX_INLINE
#define AMREX_GPU_HOST_DEVICE
AMREX_GPU_HOST_DEVICE AMREX_INLINE void SingleLinearBreitWheelerCollisionEvent(const amrex::ParticleReal &u1x, const amrex::ParticleReal &u1y, const amrex::ParticleReal &u1z, const amrex::ParticleReal &u2x, const amrex::ParticleReal &u2y, const amrex::ParticleReal &u2z, amrex::ParticleReal w1, amrex::ParticleReal w2, const amrex::Real &dt, const amrex::ParticleReal &dV, const int &pair_index, index_type *AMREX_RESTRICT p_mask, amrex::ParticleReal *AMREX_RESTRICT p_pair_reaction_weight, const amrex::ParticleReal &event_multiplier, const int &multiplier_ratio, const amrex::ParticleReal &probability_threshold, const amrex::ParticleReal &probability_target_value, const amrex::RandomEngine &engine)
This function computes whether the collision between two photons results in a pair-producing (electro...
Definition SingleLinearBreitWheelerCollisionEvent.H:56
#define WARPX_ABORT_WITH_MESSAGE(MSG)
Definition TextMsg.H:15
WarpXParticleContainer::ParticleTileType ParticleTileType
Definition LinearBreitWheelerCollisionFunc.H:39
amrex::ParticleReal m_probability_threshold
Definition LinearBreitWheelerCollisionFunc.H:264
bool use_global_debye_length()
Definition LinearBreitWheelerCollisionFunc.H:252
bool m_isSameSpecies
Definition LinearBreitWheelerCollisionFunc.H:266
ParticleBins::index_type index_type
Definition LinearBreitWheelerCollisionFunc.H:42
Executor const & executor() const
Definition LinearBreitWheelerCollisionFunc.H:250
amrex::DenseBins< ParticleTileDataType > ParticleBins
Definition LinearBreitWheelerCollisionFunc.H:41
Executor m_exe
Definition LinearBreitWheelerCollisionFunc.H:268
LinearBreitWheelerCollisionFunc()=default
Default constructor of the LinearBreitWheelerCollisionFunc class.
WarpXParticleContainer::ParticleTileType::ParticleTileDataType SoaData_type
Definition LinearBreitWheelerCollisionFunc.H:43
amrex::ParticleReal m_probability_target_value
Definition LinearBreitWheelerCollisionFunc.H:265
ParticleTileType::ParticleTileDataType ParticleTileDataType
Definition LinearBreitWheelerCollisionFunc.H:40
LinearBreitWheelerCollisionFunc(const std::string &collision_name, MultiParticleContainer const *const, const bool isSameSpecies)
Constructor of the LinearBreitWheelerCollisionFunc class.
Definition LinearBreitWheelerCollisionFunc.H:57
WarpXParticleContainer::ParticleType ParticleType
Definition LinearBreitWheelerCollisionFunc.H:38
amrex::ParticleReal m_event_multiplier
Definition LinearBreitWheelerCollisionFunc.H:257
Definition MultiParticleContainer.H:69
std::conditional_t< is_rtsoa_pc, ParticleTileRT< typename ParticleType::RealType, typename ParticleType::IntType >, ParticleTile< ParticleType, NArrayReal, NArrayInt, Allocator > > ParticleTileType
amrex_real Real
amrex_particle_real ParticleReal
AMREX_GPU_HOST_DEVICE AMREX_INLINE void remove_weight_from_colliding_particle(amrex::ParticleReal &weight, uint64_t &idcpu, const amrex::ParticleReal reaction_weight)
Subtract given weight from particle and set its ID to invalid if the weight reaches zero.
Definition BinaryCollisionUtils.H:193
constexpr std::uint64_t Invalid
__host__ __device__ constexpr const T & min(const T &a, const T &b) noexcept
__host__ __device__ constexpr const T & max(const T &a, const T &b) noexcept
int queryWithParser(const amrex::ParmParse &a_pp, char const *const str, T &val)
Definition ParserUtils.H:102
Functor that can be used to extract the positions of the macroparticles inside a ParallelFor kernel.
Definition GetAndSetPosition.H:75
Definition LinearBreitWheelerCollisionFunc.H:93
amrex::ParticleReal m_event_multiplier
Definition LinearBreitWheelerCollisionFunc.H:241
amrex::ParticleReal m_probability_threshold
Definition LinearBreitWheelerCollisionFunc.H:242
AMREX_GPU_HOST_DEVICE AMREX_INLINE void operator()(index_type const I1s, index_type const I1e, index_type const I2s, index_type const I2e, index_type const *AMREX_RESTRICT I1, index_type const *AMREX_RESTRICT I2, const SoaData_type &soa_1, const SoaData_type &soa_2, GetParticlePosition< PIdx >, GetParticlePosition< PIdx >, amrex::ParticleReal const, amrex::ParticleReal const, amrex::ParticleReal const, amrex::ParticleReal const, amrex::Real const, amrex::ParticleReal const, amrex::ParticleReal const, amrex::ParticleReal const, amrex::ParticleReal const, amrex::Real const dt, amrex::Real const dV, index_type coll_idx, index_type const cell_start_pair, index_type *AMREX_RESTRICT p_mask, index_type *AMREX_RESTRICT p_pair_indices_1, index_type *AMREX_RESTRICT p_pair_indices_2, amrex::ParticleReal *AMREX_RESTRICT p_pair_reaction_weight, amrex::ParticleReal *, amrex::RandomEngine const &engine) const
Executor of the LinearBreitWheelerCollisionFunc class. Performs two-photon collisions at the cell lev...
Definition LinearBreitWheelerCollisionFunc.H:136
bool m_computeSpeciesDensities
Definition LinearBreitWheelerCollisionFunc.H:245
bool m_computeSpeciesTemperatures
Definition LinearBreitWheelerCollisionFunc.H:246
amrex::ParticleReal m_probability_target_value
Definition LinearBreitWheelerCollisionFunc.H:243
bool m_isSameSpecies
Definition LinearBreitWheelerCollisionFunc.H:244
bool m_need_product_data
Definition LinearBreitWheelerCollisionFunc.H:247
@ uz
Definition WarpXParticleContainer.H:70
@ w
Definition WarpXParticleContainer.H:70
@ uy
Definition WarpXParticleContainer.H:70
@ ux
Definition WarpXParticleContainer.H:70