7#ifndef WARPX_PARTICLE_UTILS_H_
8#define WARPX_PARTICLE_UTILS_H_
54 double const M,
double& gamma,
double& energy )
57 using namespace amrex::literals;
61 2.0_rt * m * M * u2 / (gamma + 1.0_rt)
62 / (M + m +
sqrt(m*m + M*M + 2.0_rt * m * M * gamma))
77 amrex::ParticleReal& uz,
78 amrex::ParticleReal
const Vx, amrex::ParticleReal
const Vy,
79 amrex::ParticleReal
const Vz )
81 using namespace amrex::literals;
84 const auto V2 = (Vx*Vx + Vy*Vy + Vz*Vz);
86 const auto gamma_u = std::sqrt(1.0_prt + (ux*ux + uy*uy + uz*uz)*
PhysConst::inv_c2);
93 ux = vx * (1.0_prt + (gamma_V - 1.0_prt) * Vx*Vx/V2)
94 + vy * (gamma_V - 1.0_prt) * Vx*Vy/V2
95 +
vz * (gamma_V - 1.0_prt) * Vx*Vz/V2
96 - gamma_V * Vx * gamma_u;
98 uy = vy * (1.0_prt + (gamma_V - 1.0_prt) * Vy*Vy/V2)
99 + vx * (gamma_V - 1.0_prt) * Vx*Vy/V2
100 +
vz * (gamma_V - 1.0_prt) * Vy*Vz/V2
101 - gamma_V * Vy * gamma_u;
103 uz =
vz * (1.0_prt + (gamma_V - 1.0_prt) * Vz*Vz/V2)
104 + vx * (gamma_V - 1.0_prt) * Vx*Vz/V2
105 + vy * (gamma_V - 1.0_prt) * Vy*Vz/V2
106 - gamma_V * Vz * gamma_u;
120 amrex::ParticleReal& uz,
121 amrex::ParticleReal
const Ux, amrex::ParticleReal
const Uy,
122 amrex::ParticleReal
const Uz )
124 using namespace amrex::literals;
128 amrex::ParticleReal
const Usq = Ux*Ux + Uy*Uy + Uz*Uz;
129 amrex::ParticleReal
const U = std::sqrt(Usq);
131 amrex::ParticleReal
const gamma = std::sqrt(1._prt + Usq/c2);
133 amrex::ParticleReal
const gammam1 = Usq/c2/(gamma + 1.0_prt);
134 amrex::ParticleReal
const Ux_hat = Ux/U;
135 amrex::ParticleReal
const Uy_hat = Uy/U;
136 amrex::ParticleReal
const Uz_hat = Uz/U;
137 amrex::ParticleReal
const P0 = std::sqrt(1._prt + (ux*ux + uy*uy + uz*uz)/c2);
138 amrex::ParticleReal
const udotn = ux*Ux_hat + uy*Uy_hat + uz*Uz_hat;
139 ux = ux + gammam1*udotn*Ux_hat - Ux*P0;
140 uy = uy + gammam1*udotn*Uy_hat - Uy*P0;
141 uz = uz + gammam1*udotn*Uz_hat - Uz*P0;
157 amrex::ParticleReal& pz, amrex::ParticleReal mass,
158 amrex::ParticleReal
const Ux, amrex::ParticleReal
const Uy,
159 amrex::ParticleReal
const Uz )
161 using namespace amrex::literals;
165 amrex::ParticleReal
const Usq = Ux*Ux + Uy*Uy + Uz*Uz;
166 amrex::ParticleReal
const U = std::sqrt(Usq);
168 amrex::ParticleReal
const gamma = std::sqrt(1._prt + Usq/c2);
170 amrex::ParticleReal
const gammam1 = Usq/c2/(gamma + 1.0_prt);
171 amrex::ParticleReal
const Ux_hat = Ux/U;
172 amrex::ParticleReal
const Uy_hat = Uy/U;
173 amrex::ParticleReal
const Uz_hat = Uz/U;
174 amrex::ParticleReal
const P0 = std::sqrt(px*px + py*py + pz*pz + mass*mass*c2);
175 amrex::ParticleReal
const pdotn = px*Ux_hat + py*Uy_hat + pz*Uz_hat;
199 using namespace amrex::literals;
203 auto const xy =
sqrt(1_prt - z*z);
219 amrex::ParticleReal& uz,
220 const amrex::ParticleReal vp,
223 amrex::ParticleReal x, y, z;
241 const auto *
const xlo = tilebox.
lo();
242 const auto *
const xhi = tilebox.
hi();
244 && (xlo[1] <= point.
y) && (point.
y <= xhi[1]),
245 && (xlo[2] <= point.
z) && (point.
z <= xhi[2]));
267 template <
typename T_index>
274 T_index
const cell_start,
275 T_index
const cell_stop,
276 amrex::ParticleReal
const mass,
277 amrex::ParticleReal
const energy_fraction,
278 amrex::ParticleReal
const energy_fraction_max,
279 amrex::ParticleReal & deltaE )
281 using namespace amrex::literals;
285 int const numCell = (cell_stop - cell_start);
297 int const max_loop_count = 10;
298 amrex::ParticleReal Erel_cumm = 0.0_prt;
299 amrex::ParticleReal fmult_fact = 1.0_prt;
300 amrex::ParticleReal energy_frac_eff = 0.0_prt;
301 T_index i1 = cell_start;
302 bool correction_failed =
false;
303 bool deltaE_is_zero =
false;
304 while (!deltaE_is_zero && !correction_failed && loop_count <= max_loop_count) {
307 if (i2 == cell_stop) {
312 int const sign = (deltaE < 0.0 ? -1 : 1);
314 const amrex::ParticleReal ux1 = uxp[ indices[i1] ];
315 const amrex::ParticleReal uy1 = uyp[ indices[i1] ];
316 const amrex::ParticleReal uz1 = uzp[ indices[i1] ];
317 const amrex::ParticleReal ux2 = uxp[ indices[i2] ];
318 const amrex::ParticleReal uy2 = uyp[ indices[i2] ];
319 const amrex::ParticleReal uz2 = uzp[ indices[i2] ];
320 const amrex::ParticleReal wpmp1 = mass*w[ indices[i1] ];
321 const amrex::ParticleReal wpmp2 = mass*w[ indices[i2] ];
330 const amrex::ParticleReal ux = ux1 - ux2;
331 const amrex::ParticleReal uy = uy1 - uy2;
332 const amrex::ParticleReal uz = uz1 - uz2;
333 const amrex::ParticleReal dbetasq = (ux*ux + uy*uy + uz*uz)/c2;
334 const amrex::ParticleReal gbsq1 = (ux1*ux1 + uy1*uy1 + uz1*uz1)/c2;
335 const amrex::ParticleReal gbsq2 = (ux2*ux2 + uy2*uy2 + uz2*uz2)/c2;
336 const amrex::ParticleReal gamma1 = std::sqrt(1.0_prt + gbsq1);
337 const amrex::ParticleReal gamma2 = std::sqrt(1.0_prt + gbsq2);
338 const amrex::ParticleReal E1 = wpmp1*gamma1*c2;
339 const amrex::ParticleReal E2 = wpmp2*gamma2*c2;
340 const amrex::ParticleReal Etot = E1 + E2;
341 const amrex::ParticleReal pxtot = wpmp1*ux1 + wpmp2*ux2;
342 const amrex::ParticleReal pytot = wpmp1*uy1 + wpmp2*uy2;
343 const amrex::ParticleReal pztot = wpmp1*uz1 + wpmp2*uz2;
344 const amrex::ParticleReal Ecm = std::sqrt(Etot*Etot - (pxtot*pxtot + pytot*pytot + pztot*pztot)*c2);
345 const amrex::ParticleReal Erel = Ecm - wpmp1*c2 - wpmp2*c2;
346 if (Erel <= 0.0) {
continue; }
353 amrex::ParticleReal deltaEpair =
static_cast<amrex::ParticleReal
>(sign*energy_fraction*fmult_fact)*Erel;
354 if (std::abs(deltaEpair) >= std::abs(deltaE)) {
357 deltaE_is_zero =
true;
359 deltaE -= deltaEpair;
362 const amrex::ParticleReal A = Etot - deltaEpair;
363 const amrex::ParticleReal D = A*A + E2*E2 - E1*E1;
364 const amrex::ParticleReal p2dotu =
static_cast<amrex::ParticleReal
>(wpmp2*(ux2*ux + uy2*uy + uz2*uz));
365 const amrex::ParticleReal ptdotu = (pxtot*ux + pytot*uy + pztot*uz);
368 const amrex::ParticleReal a = A*A*dbetasq - ptdotu*ptdotu;
369 const amrex::ParticleReal b = D*ptdotu - 2.0_prt*A*A*p2dotu;
370 const amrex::ParticleReal c = A*A*E2*E2 - D*D/4.0_prt;
372 const amrex::ParticleReal root = b*b - 4.0_prt*a*c;
373 if (root < 0.0 || a == 0.0) {
continue; }
374 const amrex::ParticleReal
alpha = (-b + std::sqrt(root))/(2.0_prt*a);
377 const amrex::ParticleReal ratio1 =
alpha/(wpmp1*c2);
378 const amrex::ParticleReal ratio2 =
alpha/(wpmp2*c2);
379 uxp[ indices[i1] ] += ratio1*ux;
380 uyp[ indices[i1] ] += ratio1*uy;
381 uzp[ indices[i1] ] += ratio1*uz;
382 uxp[ indices[i2] ] -= ratio2*ux;
383 uyp[ indices[i2] ] -= ratio2*uy;
384 uzp[ indices[i2] ] -= ratio2*uz;
386 Erel_cumm += Erel - deltaEpair;
388 if (i1 < cell_stop - 2) {
396 if (i1 == cell_stop - 2) {
397 if ((numCell % 2) == 0) { i1 = cell_start + 1; }
398 if ((numCell % 2) != 0) { i1 = cell_start; }
400 if ((numCell % 2) == 0) { i1 = cell_start; }
401 if ((numCell % 2) != 0) { i1 = cell_start + 1; }
408 energy_frac_eff = std::abs(deltaE)/Erel_cumm;
409 if (energy_frac_eff > energy_fraction_max || loop_count > max_loop_count) {
410 correction_failed =
true;
412 }
else if (deltaE < 0.) {
415 fmult_fact = std::max(1._prt, energy_frac_eff/energy_fraction);
421 return correction_failed;
#define AMREX_GPU_HOST_DEVICE
#define AMREX_D_TERM(a, b, c)
@ vz
Definition RigidInjectedParticleContainer.H:27
@ alpha
Definition SpeciesPhysicalProperties.H:18
ParticleTile< ParticleType, NArrayReal, NArrayInt, Allocator > ParticleTileType
__host__ __device__ const Real * lo() const &noexcept
__host__ __device__ const Real * hi() const &noexcept
Definition ParticleUtils.cpp:24
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:156
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:119
AMREX_GPU_HOST_DEVICE AMREX_INLINE void getCollisionEnergy(amrex::ParticleReal const u2, double const m, double const M, double &gamma, double &energy)
Return (relativistic) collision energy assuming the target (with mass M) is stationary and the projec...
Definition ParticleUtils.H:53
AMREX_GPU_HOST_DEVICE AMREX_INLINE void doLorentzTransform(amrex::ParticleReal &ux, amrex::ParticleReal &uy, amrex::ParticleReal &uz, amrex::ParticleReal const Vx, amrex::ParticleReal const Vy, amrex::ParticleReal const Vz)
Perform a Lorentz transformation of the given velocity to a frame moving with velocity (Vx,...
Definition ParticleUtils.H:76
typename WarpXParticleContainer::ParticleTileType ParticleTileType
Definition ParticleUtils.H:22
AMREX_GPU_HOST_DEVICE AMREX_INLINE void RandomizeVelocity(amrex::ParticleReal &ux, amrex::ParticleReal &uy, amrex::ParticleReal &uz, const amrex::ParticleReal vp, amrex::RandomEngine const &engine)
Function to perform scattering of a particle that results in a random velocity vector with given magn...
Definition ParticleUtils.H:218
amrex::DenseBins< ParticleTileDataType > findParticlesInEachCell(amrex::Geometry const &geom_lev, amrex::MFIter const &mfi, ParticleTileType &ptile)
Find the particles and count the particles that are in each cell. More specifically this function ret...
Definition ParticleUtils.cpp:35
AMREX_GPU_HOST_DEVICE AMREX_INLINE bool ModifyEnergyPairwise(amrex::ParticleReal *const AMREX_RESTRICT uxp, amrex::ParticleReal *const AMREX_RESTRICT uyp, amrex::ParticleReal *const AMREX_RESTRICT uzp, amrex::ParticleReal *const AMREX_RESTRICT w, T_index const *const AMREX_RESTRICT indices, T_index const cell_start, T_index const cell_stop, amrex::ParticleReal const mass, amrex::ParticleReal const energy_fraction, amrex::ParticleReal const energy_fraction_max, amrex::ParticleReal &deltaE)
Definition ParticleUtils.H:269
AMREX_GPU_HOST_DEVICE AMREX_INLINE void getRandomVector(amrex::ParticleReal &x, amrex::ParticleReal &y, amrex::ParticleReal &z, amrex::RandomEngine const &engine)
Generate random unit vector in 3 dimensions https://mathworld.wolfram.com/SpherePointPicking....
Definition ParticleUtils.H:193
AMREX_GPU_HOST_DEVICE AMREX_INLINE bool containsInclusive(amrex::RealBox const &tilebox, amrex::XDim3 const point)
Definition ParticleUtils.H:240
typename ParticleTileType::ParticleTileDataType ParticleTileDataType
Definition ParticleUtils.H:23
constexpr auto c
vacuum speed of light [m/s]
Definition constant.H:149
constexpr auto inv_c2
inverse of the square of the vacuum speed of light [s^2/m^2]
Definition constant.H:216
constexpr auto q_e
elementary charge [C]
Definition constant.H:158
constexpr auto pi
ratio of a circle's circumference to its diameter
Definition constant.H:29
__host__ __device__ GpuComplex< T > sqrt(const GpuComplex< T > &a_z) noexcept