ARTEMIS
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_Particles.H>
15 
16 #include <AMReX_BaseFwd.H>
17 
18 namespace ParticleUtils {
19 
32  findParticlesInEachCell( int const lev, amrex::MFIter const& mfi,
34 
45  void getEnergy ( amrex::ParticleReal const u2, double const mass,
46  double& energy )
47  {
48  using std::sqrt;
49  using namespace amrex::literals;
50 
51  constexpr auto c2 = PhysConst::c * PhysConst::c;
52  energy = mass * u2 / (sqrt(1.0_rt + u2 / c2) + 1.0_rt) / PhysConst::q_e;
53  }
54 
67  void getCollisionEnergy ( amrex::ParticleReal const u2, double const m,
68  double const M, double& gamma, double& energy )
69  {
70  using std::sqrt;
71  using namespace amrex::literals;
72 
73  constexpr auto c2 = PhysConst::c * PhysConst::c;
74 
75  gamma = sqrt(1.0_rt + u2 / c2);
76  energy = (
77  2.0_rt * m * M * u2 / (gamma + 1.0_rt)
78  / (M + m + sqrt(m*m + M*M + 2.0_rt * m * M * gamma))
79  ) / PhysConst::q_e;
80  }
81 
92  void doLorentzTransform ( amrex::ParticleReal& ux, amrex::ParticleReal& uy,
93  amrex::ParticleReal& uz,
94  amrex::ParticleReal const Vx, amrex::ParticleReal const Vy,
95  amrex::ParticleReal const Vz )
96  {
97  // precompute repeatedly used quantities
98  constexpr auto c2 = PhysConst::c * PhysConst::c;
99  const auto V2 = (Vx*Vx + Vy*Vy + Vz*Vz);
100  const auto gamma_V = 1.0_prt / sqrt(1.0_prt - V2 / c2);
101  const auto gamma_u = sqrt(1.0_prt + (ux*ux + uy*uy + uz*uz) / c2);
102 
103  // copy velocity vector values
104  const auto vx = ux;
105  const auto vy = uy;
106  const auto vz = uz;
107 
108  ux = vx * (1.0_prt + (gamma_V - 1.0_prt) * Vx*Vx/V2)
109  + vy * (gamma_V - 1.0_prt) * Vx*Vy/V2
110  + vz * (gamma_V - 1.0_prt) * Vx*Vz/V2
111  - gamma_V * Vx * gamma_u;
112 
113  uy = vy * (1.0_prt + (gamma_V - 1.0_prt) * Vy*Vy/V2)
114  + vx * (gamma_V - 1.0_prt) * Vx*Vy/V2
115  + vz * (gamma_V - 1.0_prt) * Vy*Vz/V2
116  - gamma_V * Vy * gamma_u;
117 
118  uz = vz * (1.0_prt + (gamma_V - 1.0_prt) * Vz*Vz/V2)
119  + vx * (gamma_V - 1.0_prt) * Vx*Vz/V2
120  + vy * (gamma_V - 1.0_prt) * Vy*Vz/V2
121  - gamma_V * Vz * gamma_u;
122  }
123 
134  void getRandomVector ( amrex::ParticleReal& x, amrex::ParticleReal& y,
135  amrex::ParticleReal& z, amrex::RandomEngine const& engine )
136  {
137  using std::sqrt;
138  using std::cos;
139  using std::sin;
140  using namespace amrex::literals;
141 
142  auto const theta = amrex::Random(engine) * 2.0_prt * MathConst::pi;
143  z = 2.0_prt * amrex::Random(engine) - 1.0_prt;
144  auto const xy = sqrt(1_prt - z*z);
145  x = xy * cos(theta);
146  y = xy * sin(theta);
147  }
148 
149 
159  void RandomizeVelocity ( amrex::ParticleReal& ux, amrex::ParticleReal& uy,
160  amrex::ParticleReal& uz,
161  const amrex::ParticleReal vp,
162  amrex::RandomEngine const& engine )
163  {
164  amrex::ParticleReal x, y, z;
165  // generate random unit vector for the new velocity direction
166  getRandomVector(x, y, z, engine);
167 
168  // scale new vector to have the desired magnitude
169  ux = x * vp;
170  uy = y * vp;
171  uz = z * vp;
172  }
173 
174  /* \brief Determines whether the point is within the tilebox, inclusive of the boundaries.
175  * Note that this routine is needed since tilebox.contains excludes the boundaries.
176  * \param[in] tilebox The tilebox being checked
177  * \param[in] point The point being checked
178  * \result true if the point with within the boundary, otherwise false
179  */
181  bool containsInclusive (amrex::RealBox const& tilebox, amrex::XDim3 const point) {
182  const auto xlo = tilebox.lo();
183  const auto xhi = tilebox.hi();
184  return AMREX_D_TERM((xlo[0] <= point.x) && (point.x <= xhi[0]),
185  && (xlo[1] <= point.y) && (point.y <= xhi[1]),
186  && (xlo[2] <= point.z) && (point.z <= xhi[2]));
187  }
188 
189 }
190 
191 #endif // WARPX_PARTICLE_UTILS_H_
#define AMREX_INLINE
#define AMREX_GPU_HOST_DEVICE
#define AMREX_D_TERM(a, b, c)
AMREX_GPU_HOST_DEVICE const Real * lo() const &noexcept
AMREX_GPU_HOST_DEVICE const Real * hi() const &noexcept
def y
Definition: Excitation_Flag_Generator.py:76
def z
Definition: Excitation_Flag_Generator.py:77
Definition: ParticleUtils.cpp:25
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:67
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:92
AMREX_GPU_HOST_DEVICE AMREX_INLINE void getEnergy(amrex::ParticleReal const u2, double const mass, double &energy)
Return (relativistic) particle energy given velocity and mass. Note the use of double since this calc...
Definition: ParticleUtils.H:45
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:159
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:134
AMREX_GPU_HOST_DEVICE AMREX_INLINE bool containsInclusive(amrex::RealBox const &tilebox, amrex::XDim3 const point)
Definition: ParticleUtils.H:181
ParticleBins findParticlesInEachCell(int const lev, MFIter const &mfi, ParticleTileType const &ptile)
Find the particles and count the particles that are in each cell. More specifically this function ret...
Definition: ParticleUtils.cpp:37
static constexpr auto c
vacuum speed of light [m/s]
Definition: constant.H:44
static constexpr auto q_e
elementary charge [C]
Definition: constant.H:50
static constexpr amrex::Real pi
ratio of a circle's circumference to its diameter
Definition: constant.H:23
Real Random()
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE GpuComplex< T > sqrt(const GpuComplex< T > &a_z) noexcept
int gamma
Definition: stencil.py:446