7#ifndef WARPX_PARTICLE_UTILS_H_
8#define WARPX_PARTICLE_UTILS_H_
54 double const M,
double& gamma,
double& energy )
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))
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;
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;
200 auto const cos_theta = 2.0_prt *
amrex::Random(engine) - 1.0_prt;
201 auto const sin_theta = std::sqrt(1.0_prt - cos_theta*cos_theta);
202 x = sin_theta * std::cos(phi);
203 y = sin_theta * std::sin(phi);
240 const auto *
const xlo = tilebox.
lo();
241 const auto *
const xhi = tilebox.
hi();
243 && (xlo[1] <= point.
y) && (point.
y <= xhi[1]),
244 && (xlo[2] <= point.
z) && (point.
z <= xhi[2]));
257 if (do_cropping[0] && x1 < xmin) {
260 else if (do_cropping[1] && x1 > xmax) {
274 double & x1,
double & z1,
276 if (do_cropping[0] && x1 < xmin) {
277 double const fx = ((x0 - x1) != 0. ? (x0 - xmin)/(x0 - x1) : 0.);
279 z1 = z0 + fx*(z1 - z0);
281 else if (do_cropping[1] && x1 > xmax) {
282 double const fx = ((x0 - x1) != 0. ? (xmax - x0)/(x1 - x0) : 0.);
284 z1 = z0 + fx*(z1 - z0);
297 double & x1,
double & y1,
double & z1,
299 if (do_cropping[0] && x1 < xmin) {
300 double const fx = ((x0 - x1) != 0. ? (x0 - xmin)/(x0 - x1) : 0.);
302 y1 = y0 + fx*(y1 - y0);
303 z1 = z0 + fx*(z1 - z0);
305 else if (do_cropping[1] && x1 > xmax) {
306 double const fx = ((x0 - x1) != 0. ? (xmax - x0)/(x1 - x0) : 0.);
308 y1 = y0 + fx*(y1 - y0);
309 z1 = z0 + fx*(z1 - z0);
324 [[maybe_unused]]
double yp,
325 [[maybe_unused]]
double zp,
334#if defined(WARPX_DIM_3D)
335 xp_grid[0] = (xp - xyzmin.
x)*dinv.
x;
336 xp_grid[1] = (yp - xyzmin.
y)*dinv.
y;
337 xp_grid[2] = (zp - xyzmin.
z)*dinv.
z;
338#elif defined(WARPX_DIM_XZ)
339 xp_grid[0] = (xp - xyzmin.
x)*dinv.
x;
340 xp_grid[1] = (zp - xyzmin.
z)*dinv.
z;
341#elif defined(WARPX_DIM_1D_Z)
342 xp_grid[0] = (zp - xyzmin.
z)*dinv.
z;
343#elif defined(WARPX_DIM_RZ) || defined(WARPX_DIM_RCYLINDER)
344 const double rp = std::sqrt(xp*xp + yp*yp);
345 xp_grid[0] = (rp - xyzmin.
x)*dinv.
x;
346#if defined(WARPX_DIM_RZ)
347 xp_grid[1] = (zp - xyzmin.
z)*dinv.
z;
349#elif defined(WARPX_DIM_RSPHERE)
350 const double rp = std::sqrt(xp*xp + yp*yp + zp*zp);
351 xp_grid[0] = (rp - xyzmin.
x)*dinv.
x;
355 for (
int dim = 0; dim < AMREX_SPACEDIM; dim++) {
356 if ( (do_cropping[dim][0] && xp_grid[dim] <= domain_double[dim][0]) ||
357 (do_cropping[dim][1] && xp_grid[dim] >= domain_double[dim][1])) {
383 template <
typename T_index>
390 T_index
const cell_start,
391 T_index
const cell_stop,
401 int const numCell = (cell_stop - cell_start);
413 int const max_loop_count = 10;
417 T_index i1 = cell_start;
418 bool correction_failed =
false;
419 bool deltaE_is_zero =
false;
420 while (!deltaE_is_zero && !correction_failed && loop_count <= max_loop_count) {
423 if (i2 == cell_stop) {
428 int const sign = (deltaE < 0.0 ? -1 : 1);
460 const amrex::ParticleReal Ecm = std::sqrt(Etot*Etot - (pxtot*pxtot + pytot*pytot + pztot*pztot)*c2);
462 if (Erel <= 0.0) {
continue; }
470 if (std::abs(deltaEpair) >= std::abs(deltaE)) {
473 deltaE_is_zero =
true;
475 deltaE -= deltaEpair;
489 if (root < 0.0 || a == 0.0) {
continue; }
495 uxp[ indices[i1] ] += ratio1*ux;
496 uyp[ indices[i1] ] += ratio1*uy;
497 uzp[ indices[i1] ] += ratio1*uz;
498 uxp[ indices[i2] ] -= ratio2*ux;
499 uyp[ indices[i2] ] -= ratio2*uy;
500 uzp[ indices[i2] ] -= ratio2*uz;
502 Erel_cumm += Erel - deltaEpair;
504 if (i1 < cell_stop - 2) {
512 if (i1 == cell_stop - 2) {
513 if ((numCell % 2) == 0) { i1 = cell_start + 1; }
514 if ((numCell % 2) != 0) { i1 = cell_start; }
516 if ((numCell % 2) == 0) { i1 = cell_start; }
517 if ((numCell % 2) != 0) { i1 = cell_start + 1; }
524 energy_frac_eff = std::abs(deltaE)/Erel_cumm;
525 if (energy_frac_eff > energy_fraction_max || loop_count > max_loop_count) {
526 correction_failed =
true;
528 }
else if (deltaE < 0.) {
531 fmult_fact = std::max(1._prt, energy_frac_eff/energy_fraction);
537 return correction_failed;
#define AMREX_FORCE_INLINE
#define AMREX_GPU_HOST_DEVICE
#define AMREX_D_TERM(a, b, c)
@ vz
Definition RigidInjectedParticleContainer.H:27
std::conditional_t< is_rtsoa_pc, ParticleTileRT< typename ParticleType::RealType, typename ParticleType::IntType >, ParticleTile< ParticleType, NArrayReal, NArrayInt, Allocator > > ParticleTileType
__host__ __device__ const Real * lo() const &noexcept
__host__ __device__ const Real * hi() const &noexcept
amrex_particle_real ParticleReal
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:217
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_FORCE_INLINE bool is_out_of_bounds(double xp, double yp, double zp, const amrex::XDim3 &dinv, const amrex::XDim3 &xyzmin, const amrex::GpuArray< amrex::GpuArray< double, 2 >, 3 > &domain_double, const amrex::GpuArray< amrex::GpuArray< bool, 2 >, 3 > &do_cropping)
Definition ParticleUtils.H:323
AMREX_GPU_HOST_DEVICE AMREX_INLINE void getRandomVector(amrex::ParticleReal &x, amrex::ParticleReal &y, amrex::ParticleReal &z, amrex::RandomEngine const &engine)
Generate a random unit vector with isotropic distribution, following the method described in equation...
Definition ParticleUtils.H:194
AMREX_GPU_HOST_DEVICE AMREX_INLINE bool containsInclusive(amrex::RealBox const &tilebox, amrex::XDim3 const point)
Definition ParticleUtils.H:239
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void crop_at_boundary(double &x1, double xmin, double xmax, amrex::GpuArray< bool, 2 > const &do_cropping)
Definition ParticleUtils.H:255
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, const 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:385
typename ParticleTileType::ParticleTileDataType ParticleTileDataType
Definition ParticleUtils.H:23
constexpr auto c
vacuum speed of light [m/s]
Definition constant.H:153
constexpr auto c2
square of the vacuum speed of light [m^2/s^2]
Definition constant.H:214
constexpr auto inv_c2
inverse of the square of the vacuum speed of light [s^2/m^2]
Definition constant.H:220
constexpr auto q_e
elementary charge [C]
Definition constant.H:162
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