WarpX
Loading...
Searching...
No Matches
ParticleUtils.H
Go to the documentation of this file.
1/* Copyright 2019-2020 Neil Zaim, Yinjian Zhao
2 *
3 * This file is part of WarpX.
4 *
5 * License: BSD-3-Clause-LBNL
6 */
7#ifndef WARPX_PARTICLE_UTILS_H_
8#define WARPX_PARTICLE_UTILS_H_
9
11#include "Utils/WarpXConst.H"
12
13#include <AMReX_DenseBins.H>
14#include <AMReX_Geometry.H>
15#include <AMReX_Particles.H>
16
17#include <AMReX_BaseFwd.H>
18
19namespace ParticleUtils {
20
21 // Define shortcuts for frequently-used type names
23 using ParticleTileDataType = typename ParticleTileType::ParticleTileDataType;
24
38 amrex::MFIter const & mfi,
40
53 void getCollisionEnergy ( amrex::ParticleReal const u2, double const m,
54 double const M, double& gamma, double& energy )
55 {
56 using std::sqrt;
57 using namespace amrex::literals;
58
59 gamma = sqrt(1.0_rt + u2 * PhysConst::inv_c2);
60 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))
64 }
65
79 amrex::ParticleReal const Vz )
80 {
81 using namespace amrex::literals;
82
83 // precompute repeatedly used quantities
84 const auto V2 = (Vx*Vx + Vy*Vy + Vz*Vz);
85 const auto gamma_V = 1.0_prt / std::sqrt(1.0_prt - V2 * PhysConst::inv_c2);
86 const auto gamma_u = std::sqrt(1.0_prt + (ux*ux + uy*uy + uz*uz)* PhysConst::inv_c2);
87
88 // copy velocity vector values
89 const auto vx = ux;
90 const auto vy = uy;
91 const auto vz = uz;
92
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;
97
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;
102
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;
107 }
108
121 amrex::ParticleReal const Ux, amrex::ParticleReal const Uy,
122 amrex::ParticleReal const Uz )
123 {
124 using namespace amrex::literals;
125
126 constexpr auto c2 = PhysConst::c2;
127
128 amrex::ParticleReal const Usq = Ux*Ux + Uy*Uy + Uz*Uz;
129 amrex::ParticleReal const U = std::sqrt(Usq);
130 if (U > 0._prt) {
131 amrex::ParticleReal const gamma = std::sqrt(1._prt + Usq/c2);
132 // A nice way of calculating (gamma - 1) when it is small
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;
142 }
143
144 }
145
158 amrex::ParticleReal const Ux, amrex::ParticleReal const Uy,
159 amrex::ParticleReal const Uz )
160 {
161 using namespace amrex::literals;
162
163 constexpr auto c2 = PhysConst::c2;
164
165 amrex::ParticleReal const Usq = Ux*Ux + Uy*Uy + Uz*Uz;
166 amrex::ParticleReal const U = std::sqrt(Usq);
167 if (U > 0._prt) {
168 amrex::ParticleReal const gamma = std::sqrt(1._prt + Usq/c2);
169 // A nice way of calculating (gamma - 1) when it is small
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;
176 px = px + gammam1*pdotn*Ux_hat - Ux*P0/PhysConst::c;
177 py = py + gammam1*pdotn*Uy_hat - Uy*P0/PhysConst::c;
178 pz = pz + gammam1*pdotn*Uz_hat - Uz*P0/PhysConst::c;
179 }
180
181 }
182
195 amrex::ParticleReal& z, amrex::RandomEngine const& engine )
196 {
197 using namespace amrex::literals;
198
199 auto const phi = 2.0_prt * MathConst::pi * amrex::Random(engine);
200 auto const cos_theta = 2.0_prt * amrex::Random(engine) - 1.0_prt; // opposite sign of that in Higginson's paper
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);
204 z = cos_theta;
205 }
206
207
219 const amrex::ParticleReal vp,
220 amrex::RandomEngine const& engine )
221 {
222 amrex::ParticleReal x, y, z;
223 // generate random unit vector for the new velocity direction
224 getRandomVector(x, y, z, engine);
225
226 // scale new vector to have the desired magnitude
227 ux = x * vp;
228 uy = y * vp;
229 uz = z * vp;
230 }
231
232 /* \brief Determines whether the point is within the tilebox, inclusive of the boundaries.
233 * Note that this routine is needed since tilebox.contains excludes the boundaries.
234 * \param[in] tilebox The tilebox being checked
235 * \param[in] point The point being checked
236 * \result true if the point with within the boundary, otherwise false
237 */
239 bool containsInclusive (amrex::RealBox const& tilebox, amrex::XDim3 const point) {
240 const auto *const xlo = tilebox.lo();
241 const auto *const xhi = tilebox.hi();
242 return AMREX_D_TERM((xlo[0] <= point.x) && (point.x <= xhi[0]),
243 && (xlo[1] <= point.y) && (point.y <= xhi[1]),
244 && (xlo[2] <= point.z) && (point.z <= xhi[2]));
245 }
246
247 /* \brief Crops the position at the specified boundary if do_cropping is true, 1D version
248 * \param[in] x0 The initial particle position
249 * \param[in/out] x1 The final particle position, with the returned value being cropped
250 * \param[in] xmin The lower boundary of the domain in grid units
251 * \param[in] xmax The upper boundary of the domain in grid units
252 * \param[in] do_cropping Whether to crop at each of the boundaries
253 */
255 void crop_at_boundary(double & x1,
256 double xmin, double xmax, amrex::GpuArray<bool,2> const& do_cropping) {
257 if (do_cropping[0] && x1 < xmin) {
258 x1 = xmin;
259 }
260 else if (do_cropping[1] && x1 > xmax) {
261 x1 = xmax;
262 }
263 }
264
265 /* \brief Crops the position at the specified boundary if do_cropping is true, 2D version
266 * \param[in] x0, z0 The initial particle position
267 * \param[in/out] x1, z1 The final particle position, with the returned value being cropped
268 * \param[in] xmin The lower boundary of the domain in grid units
269 * \param[in] xmax The upper boundary of the domain in grid units
270 * \param[in] do_cropping Whether to crop at each of the boundaries
271 */
273 void crop_at_boundary(double x0, double z0,
274 double & x1, double & z1,
275 double xmin, double xmax, amrex::GpuArray<bool,2> const& do_cropping) {
276 if (do_cropping[0] && x1 < xmin) {
277 double const fx = ((x0 - x1) != 0. ? (x0 - xmin)/(x0 - x1) : 0.);
278 x1 = xmin;
279 z1 = z0 + fx*(z1 - z0);
280 }
281 else if (do_cropping[1] && x1 > xmax) {
282 double const fx = ((x0 - x1) != 0. ? (xmax - x0)/(x1 - x0) : 0.);
283 x1 = xmax;
284 z1 = z0 + fx*(z1 - z0);
285 }
286 }
287
288 /* \brief Crops the position at the specified boundary if do_cropping is true, 3D version
289 * \param[in] x0, y0, z0 The initial particle position
290 * \param[in/out] x1, y1, z1 The final particle position, with the returned value being cropped
291 * \param[in] xmin The lower boundary of the domain in grid units
292 * \param[in] xmax The upper boundary of the domain in grid units
293 * \param[in] do_cropping Whether to crop at each of the boundaries
294 */
296 void crop_at_boundary(double x0, double y0, double z0,
297 double & x1, double & y1, double & z1,
298 double xmin, double xmax, amrex::GpuArray<bool,2> const& do_cropping) {
299 if (do_cropping[0] && x1 < xmin) {
300 double const fx = ((x0 - x1) != 0. ? (x0 - xmin)/(x0 - x1) : 0.);
301 x1 = xmin;
302 y1 = y0 + fx*(y1 - y0);
303 z1 = z0 + fx*(z1 - z0);
304 }
305 else if (do_cropping[1] && x1 > xmax) {
306 double const fx = ((x0 - x1) != 0. ? (xmax - x0)/(x1 - x0) : 0.);
307 x1 = xmax;
308 y1 = y0 + fx*(y1 - y0);
309 z1 = z0 + fx*(z1 - z0);
310 }
311 }
312
313 /* \brief Check if the particle position is out of bounds at a boundary
314 * where orbit cropping is being performed.
315 * \param[in] xp, yp, zp The particle position
316 * \param[in] dinv 3D cell size inverse
317 * \param[in] xyzmin Physical lower bounds of domain
318 * \param[in] domain_double The domain boundary indicies cast to double values
319 * \param[in] do_cropping Whether to crop particle orbits at domain boundaries
320 * return whether position is out of bounds
321 */
323 bool is_out_of_bounds ([[maybe_unused]] double xp,
324 [[maybe_unused]] double yp,
325 [[maybe_unused]] double zp,
326 const amrex::XDim3 & dinv,
327 const amrex::XDim3 & xyzmin,
328 const amrex::GpuArray<amrex::GpuArray<double,2>, AMREX_SPACEDIM> & domain_double,
329 const amrex::GpuArray<amrex::GpuArray<bool,2>, AMREX_SPACEDIM> & do_cropping)
330 {
331
332 // Store logical particle position (in grid units) in temporary array
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;
348#endif
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;
352#endif
353
354 // Check if particle's position is out of bounds
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])) {
358 return true;
359 }
360 }
361 return false;
362
363 }
364
365 /* \brief Add or subtract a small percent of energy from a
366 * pair of particles such that momentum is not disturbed.
367 *
368 * This is done by treating it as an inelastic scattering event
369 * with potential U = deltaE and zero scattering angle
370 * deltaE > 0 ==> need to take energy away
371 * deltaE < 0 ==> need to add energy
372 *
373 * \param[in/out] uxp, uyp, uzp momenta of the particles
374 * \param[in] w weight of the particles
375 * \param[in] indices indices of particles, listed in cell order
376 * \param[in] cell_start start index of particles in the cell
377 * \param[in] cell_stop start index of particles in the next cell
378 * \param[in] mass mass of the particles
379 * \param[in] energy_fraction fraction of pair's relative energy to change
380 * \param[in] energy_fraction_max max fracton of total relative energy that can be changed
381 * \param[in/out] deltaE amount of energy to be distributed in this species
382 */
383 template <typename T_index>
385 bool ModifyEnergyPairwise (amrex::ParticleReal * const AMREX_RESTRICT uxp, //NOLINT(readability-non-const-parameter)
386 amrex::ParticleReal * const AMREX_RESTRICT uyp, //NOLINT(readability-non-const-parameter)
387 amrex::ParticleReal * const AMREX_RESTRICT uzp, //NOLINT(readability-non-const-parameter)
388 const amrex::ParticleReal * const AMREX_RESTRICT w,
389 T_index const * const AMREX_RESTRICT indices,
390 T_index const cell_start,
391 T_index const cell_stop,
392 amrex::ParticleReal const mass,
393 amrex::ParticleReal const energy_fraction,
394 amrex::ParticleReal const energy_fraction_max,
395 amrex::ParticleReal & deltaE )
396 {
397 using namespace amrex::literals;
398
399 constexpr auto c2 = PhysConst::c2;
400
401 int const numCell = (cell_stop - cell_start);
402
403 // Adjust particle's momentum to absorb deltaE.
404 // This loops over the particles, two at a time, distributing
405 // the energy until all of it has been distributed or failing if not.
406 // If deltaE < 0, it should never fail since energy can always be added
407 // to the particles.
408 // if deltaE > 0, it will sometimes fail because energy cannot always be
409 // removed from the particles without changing the momentum.
410 // The loop is setup to switch between even first and odd first pairs
411 // each loop over the particles.
412 int loop_count = 0;
413 int const max_loop_count = 10;
414 amrex::ParticleReal Erel_cumm = 0.0_prt;
415 amrex::ParticleReal fmult_fact = 1.0_prt;
416 amrex::ParticleReal energy_frac_eff = 0.0_prt;
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) {
421
422 T_index i2 = i1 + 1;
423 if (i2 == cell_stop) {
424 // If i1 is the last particle, set the other particle, i2, to be the first particle
425 i2 = cell_start;
426 }
427
428 int const sign = (deltaE < 0.0 ? -1 : 1);
429
430 const amrex::ParticleReal ux1 = uxp[ indices[i1] ];
431 const amrex::ParticleReal uy1 = uyp[ indices[i1] ];
432 const amrex::ParticleReal uz1 = uzp[ indices[i1] ];
433 const amrex::ParticleReal ux2 = uxp[ indices[i2] ];
434 const amrex::ParticleReal uy2 = uyp[ indices[i2] ];
435 const amrex::ParticleReal uz2 = uzp[ indices[i2] ];
436 const amrex::ParticleReal wpmp1 = mass*w[ indices[i1] ];
437 const amrex::ParticleReal wpmp2 = mass*w[ indices[i2] ];
438
439 // For the calculations below, since it involves taking differences
440 // between similar numbers, which can lose precision, to achieve machine
441 // level accuracy on the energy conservation long double precision would
442 // be required. However, long double is implemented inconsistently across
443 // platforms, so standard precision is used here with the accordant loss
444 // of precision. For practical purposes however, the loss of precision is negligible.
445
446 const amrex::ParticleReal ux = ux1 - ux2;
447 const amrex::ParticleReal uy = uy1 - uy2;
448 const amrex::ParticleReal uz = uz1 - uz2;
449 const amrex::ParticleReal dbetasq = (ux*ux + uy*uy + uz*uz)/c2;
450 const amrex::ParticleReal gbsq1 = (ux1*ux1 + uy1*uy1 + uz1*uz1)/c2;
451 const amrex::ParticleReal gbsq2 = (ux2*ux2 + uy2*uy2 + uz2*uz2)/c2;
452 const amrex::ParticleReal gamma1 = std::sqrt(1.0_prt + gbsq1);
453 const amrex::ParticleReal gamma2 = std::sqrt(1.0_prt + gbsq2);
454 const amrex::ParticleReal E1 = wpmp1*gamma1*c2;
455 const amrex::ParticleReal E2 = wpmp2*gamma2*c2;
456 const amrex::ParticleReal Etot = E1 + E2;
457 const amrex::ParticleReal pxtot = wpmp1*ux1 + wpmp2*ux2;
458 const amrex::ParticleReal pytot = wpmp1*uy1 + wpmp2*uy2;
459 const amrex::ParticleReal pztot = wpmp1*uz1 + wpmp2*uz2;
460 const amrex::ParticleReal Ecm = std::sqrt(Etot*Etot - (pxtot*pxtot + pytot*pytot + pztot*pztot)*c2);
461 const amrex::ParticleReal Erel = Ecm - wpmp1*c2 - wpmp2*c2;
462 if (Erel <= 0.0) { continue; }
463
464 // Set the amount of energy change for this pair based on the
465 // relative energy in the center of momentum.
466 // If deltaE > 0, no more that Erel can be removed from this pair.
467 // If deltaE < 0, limit the energy added by Erel at first, but
468 // after the first loop, add whatever is left (using the adjusted fmult_fact).
469 amrex::ParticleReal deltaEpair = static_cast<amrex::ParticleReal>(sign*energy_fraction*fmult_fact)*Erel;
470 if (std::abs(deltaEpair) >= std::abs(deltaE)) {
471 deltaEpair = deltaE;
472 deltaE = 0.0_prt;
473 deltaE_is_zero = true;
474 } else {
475 deltaE -= deltaEpair;
476 }
477
478 const amrex::ParticleReal A = Etot - deltaEpair;
479 const amrex::ParticleReal D = A*A + E2*E2 - E1*E1;
480 const auto p2dotu = static_cast<amrex::ParticleReal>(wpmp2*(ux2*ux + uy2*uy + uz2*uz));
481 const amrex::ParticleReal ptdotu = (pxtot*ux + pytot*uy + pztot*uz);
482
483 // Compute coefficients for quadratic equation for alpha.
484 const amrex::ParticleReal a = A*A*dbetasq - ptdotu*ptdotu;
485 const amrex::ParticleReal b = D*ptdotu - 2.0_prt*A*A*p2dotu;
486 const amrex::ParticleReal c = A*A*E2*E2 - D*D/4.0_prt;
487
488 const amrex::ParticleReal root = b*b - 4.0_prt*a*c;
489 if (root < 0.0 || a == 0.0) { continue; }
490 const amrex::ParticleReal alpha = (-b + std::sqrt(root))/(2.0_prt*a);
491
492 // Update particle velocities.
493 const amrex::ParticleReal ratio1 = alpha/(wpmp1*c2);
494 const amrex::ParticleReal ratio2 = alpha/(wpmp2*c2);
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;
501
502 Erel_cumm += Erel - deltaEpair;
503
504 if (i1 < cell_stop - 2) {
505 // If not done with the loop, increment to the next pair of particles
506 i1 = i2 + 1;
507 } else {
508 loop_count++;
509 // On next time through the loop, use a different set of pairs.
510 // What index to start with depends on whether the number of particles
511 // is even or odd.
512 if (i1 == cell_stop - 2) {
513 if ((numCell % 2) == 0) { i1 = cell_start + 1; }
514 if ((numCell % 2) != 0) { i1 = cell_start; }
515 } else {
516 if ((numCell % 2) == 0) { i1 = cell_start; }
517 if ((numCell % 2) != 0) { i1 = cell_start + 1; }
518 }
519 // Compute the energy fraction to correct with one more loop.
520 // This is the remaining energy correction divided by the cummulated Erel,
521 // the amount of "free" energy available.
522 // If deltaE < 0, this check always passes since the amount of energy that
523 // can be added is unlimited.
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;
527 break;
528 } else if (deltaE < 0.) {
529 // If after the first time around, there is energy left to distribute,
530 // increase the multiplier to distribute the rest faster.
531 fmult_fact = std::max(1._prt, energy_frac_eff/energy_fraction);
532 }
533 }
534
535 }
536
537 return correction_failed;
538 }
539}
540
541#endif // WARPX_PARTICLE_UTILS_H_
#define AMREX_FORCE_INLINE
#define AMREX_RESTRICT
#define AMREX_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
Real Random()
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