WarpX
UpdateMomentumPerezElastic.H
Go to the documentation of this file.
1 /* Copyright 2019 Yinjian Zhao
2  *
3  * This file is part of WarpX.
4  *
5  * License: BSD-3-Clause-LBNL
6  */
7 #ifndef WARPX_PARTICLES_COLLISION_UPDATE_MOMENTUM_PEREZ_ELASTIC_H_
8 #define WARPX_PARTICLES_COLLISION_UPDATE_MOMENTUM_PEREZ_ELASTIC_H_
9 
10 #include "Utils/WarpXConst.H"
11 
12 #include <AMReX_Math.H>
13 #include <AMReX_Random.H>
14 
15 #include <cmath> // isnan() isinf()
16 #include <limits> // numeric_limits<float>::min()
17 
18 /* \brief Update particle velocities according to
19  * F. Perez et al., Phys.Plasmas.19.083104 (2012),
20  * which is based on Nanbu's method, PhysRevE.55.4642 (1997).
21  * @param[in] bmax is max(Debye length, minimal interparticle distance).
22  * @param[in] L is the Coulomb log. A fixed L will be used if L > 0,
23  * otherwise L will be calculated based on the algorithm.
24  * @param[in] n12 = max(w1,w2)*min(N1,N2)/dV is the effective density used for s12
25  * To see if there are nan or inf updated velocities,
26  * compile with USE_ASSERTION=TRUE.
27  *
28  * Updates and corrections to the original publication are documented in
29  * https://github.com/ECP-WarpX/WarpX/issues/429
30  * https://github.com/ECP-WarpX/WarpX/files/3799803/main.pdf
31  */
32 
33 template <typename T_PR, typename T_R>
36  T_PR& u1x, T_PR& u1y, T_PR& u1z, T_PR& u2x, T_PR& u2y, T_PR& u2z,
37  T_PR const n1, T_PR const n2, T_PR const n12,
38  T_PR const q1, T_PR const m1, T_PR const w1,
39  T_PR const q2, T_PR const m2, T_PR const w2,
40  T_R const dt, T_PR const L, T_PR const bmax,
41  amrex::RandomEngine const& engine)
42 {
43 
44  T_PR const diffx = amrex::Math::abs(u1x-u2x);
45  T_PR const diffy = amrex::Math::abs(u1y-u2y);
46  T_PR const diffz = amrex::Math::abs(u1z-u2z);
47  T_PR const diffm = std::sqrt(diffx*diffx+diffy*diffy+diffz*diffz);
48  T_PR const summm = std::sqrt(u1x*u1x+u1y*u1y+u1z*u1z) + std::sqrt(u2x*u2x+u2y*u2y+u2z*u2z);
49  // If g = u1 - u2 = 0, do not collide.
50  // Or if the relative difference is less than 1.0e-10.
51  if ( diffm < std::numeric_limits<T_PR>::min() || diffm/summm < 1.0e-10 ) { return; }
52 
53  T_PR constexpr inv_c2 = T_PR(1.0) / ( PhysConst::c * PhysConst::c );
54 
55  // Compute Lorentz factor gamma
56  T_PR const g1 = std::sqrt( T_PR(1.0) + (u1x*u1x+u1y*u1y+u1z*u1z)*inv_c2 );
57  T_PR const g2 = std::sqrt( T_PR(1.0) + (u2x*u2x+u2y*u2y+u2z*u2z)*inv_c2 );
58 
59  // Compute momenta
60  T_PR const p1x = u1x * m1;
61  T_PR const p1y = u1y * m1;
62  T_PR const p1z = u1z * m1;
63  T_PR const p2x = u2x * m2;
64  T_PR const p2y = u2y * m2;
65  T_PR const p2z = u2z * m2;
66 
67  // Compute center-of-mass (COM) velocity and gamma
68  T_PR const mass_g = m1 * g1 + m2 * g2;
69  T_PR const vcx = (p1x+p2x) / mass_g;
70  T_PR const vcy = (p1y+p2y) / mass_g;
71  T_PR const vcz = (p1z+p2z) / mass_g;
72  T_PR const vcms = vcx*vcx + vcy*vcy + vcz*vcz;
73  T_PR const gc = T_PR(1.0) / std::sqrt( T_PR(1.0) - vcms*inv_c2 );
74 
75  // Compute vc dot v1 and v2
76  T_PR const vcDv1 = (vcx*u1x + vcy*u1y + vcz*u1z) / g1;
77  T_PR const vcDv2 = (vcx*u2x + vcy*u2y + vcz*u2z) / g2;
78 
79  // Compute p1 star
80  T_PR p1sx;
81  T_PR p1sy;
82  T_PR p1sz;
83  if ( vcms > std::numeric_limits<T_PR>::min() )
84  {
85  /* lorentz_transform_factor = ( (gc-1.0)/vcms*vcDv1 - gc )*m1*g1;
86  * Rewrite to avoid loss of precision from subtracting similar
87  * numbers when gc is close to 1 */
88  T_PR const lorentz_transform_factor =
89  ( (gc*gc*vcms*inv_c2/(T_PR(1.0) + gc))/vcms*vcDv1 - gc )*m1*g1;
90  p1sx = p1x + vcx*lorentz_transform_factor;
91  p1sy = p1y + vcy*lorentz_transform_factor;
92  p1sz = p1z + vcz*lorentz_transform_factor;
93  }
94  else // If vcms = 0, don't do Lorentz-transform.
95  {
96  p1sx = p1x;
97  p1sy = p1y;
98  p1sz = p1z;
99  }
100  T_PR const p1sm = std::sqrt( p1sx*p1sx + p1sy*p1sy + p1sz*p1sz );
101 
102  // Compute gamma star
103  T_PR const g1s = ( T_PR(1.0) - vcDv1*inv_c2 )*gc*g1;
104  T_PR const g2s = ( T_PR(1.0) - vcDv2*inv_c2 )*gc*g2;
105 
106  // Compute s
107  T_PR s = 0;
108  if (p1sm > std::numeric_limits<T_PR>::min()) {
109 
110  // s is non-zero (i.e. particles scatter) only if the relative
111  // motion between particles is not negligible (p1sm non-zero)
112 
113  // Compute the Coulomb log lnLmd first
114  T_PR lnLmd;
115  if ( L > T_PR(0.0) ) { lnLmd = L; }
116  else
117  {
118  // Compute b0 according to eq (22) from Perez et al., Phys.Plasmas.19.083104 (2012)
119  // Note: there is a typo in the equation, the last square is incorrect!
120  // See the SMILEI documentation: https://smileipic.github.io/Smilei/Understand/collisions.html
121  // and https://github.com/ECP-WarpX/WarpX/files/3799803/main.pdf from GitHub #429
122  T_PR const b0 = amrex::Math::abs(q1*q2) * inv_c2 /
123  (T_PR(4.0)*MathConst::pi*PhysConst::ep0) * gc/mass_g *
124  ( m1*g1s*m2*g2s/(p1sm*p1sm*inv_c2) + T_PR(1.0) );
125 
126  // Compute the minimal impact parameter
127  constexpr T_PR hbar_pi = static_cast<T_PR>(PhysConst::hbar*MathConst::pi);
128  const T_PR bmin = amrex::max(hbar_pi/p1sm, b0);
129 
130  // Compute the Coulomb log lnLmd
131  lnLmd = amrex::max( T_PR(2.0),
132  T_PR(0.5)*std::log(T_PR(1.0) + bmax*bmax/(bmin*bmin)) );
133  }
134 
135  // Compute s
136  const auto tts = m1*g1s*m2*g2s/(inv_c2*p1sm*p1sm) + T_PR(1.0);
137  const auto tts2 = tts*tts;
138  s = n12 * dt*lnLmd*q1*q1*q2*q2 /
139  ( T_PR(4.0) * MathConst::pi * PhysConst::ep0 * PhysConst::ep0 *
140  m1*g1*m2*g2/(inv_c2*inv_c2) ) * gc*p1sm/mass_g * tts2;
141 
142  // Compute s'
143  const auto cbrt_n1 = std::cbrt(n1);
144  const auto cbrt_n2 = std::cbrt(n2);
145  const auto coeff = static_cast<T_PR>(
146  std::pow(4.0*MathConst::pi/3.0,1.0/3.0));
147  T_PR const vrel = mass_g*p1sm/(m1*g1s*m2*g2s*gc);
148  T_PR const sp = coeff * n12 * dt * vrel * (m1+m2) /
149  amrex::max( m1*cbrt_n1*cbrt_n1,
150  m2*cbrt_n2*cbrt_n2);
151 
152  // Determine s
153  s = amrex::min(s,sp);
154  }
155 
156  // Only modify momenta if is s is non-zero
157  if (s > std::numeric_limits<T_PR>::min()) {
158 
159  // Get random numbers
160  T_PR r = amrex::Random(engine);
161 
162  // Compute scattering angle
163  T_PR cosXs;
164  T_PR sinXs;
165  if ( s <= T_PR(0.1) )
166  {
167  while ( true )
168  {
169  cosXs = T_PR(1.0) + s * std::log(r);
170  // Avoid the bug when r is too small such that cosXs < -1
171  if ( cosXs >= T_PR(-1.0) ) { break; }
172  r = amrex::Random(engine);
173  }
174  }
175  else if ( s > T_PR(0.1) && s <= T_PR(3.0) )
176  {
177  T_PR const Ainv = static_cast<T_PR>(
178  0.0056958 + 0.9560202*s - 0.508139*s*s +
179  0.47913906*s*s*s - 0.12788975*s*s*s*s + 0.02389567*s*s*s*s*s);
180  cosXs = Ainv * std::log( std::exp(T_PR(-1.0)/Ainv) +
181  T_PR(2.0) * r * std::sinh(T_PR(1.0)/Ainv) );
182  }
183  else if ( s > T_PR(3.0) && s <= T_PR(6.0) )
184  {
185  T_PR const A = T_PR(3.0) * std::exp(-s);
186  cosXs = T_PR(1.0)/A * std::log( std::exp(-A) +
187  T_PR(2.0) * r * std::sinh(A) );
188  }
189  else
190  {
191  cosXs = T_PR(2.0) * r - T_PR(1.0);
192  }
193  sinXs = std::sqrt(T_PR(1.0) - cosXs*cosXs);
194 
195  // Get random azimuthal angle
196  T_PR const phis = amrex::Random(engine) * T_PR(2.0) * MathConst::pi;
197  T_PR const cosphis = std::cos(phis);
198  T_PR const sinphis = std::sin(phis);
199 
200  // Compute post-collision momenta pfs in COM
201  T_PR p1fsx;
202  T_PR p1fsy;
203  T_PR p1fsz;
204  // p1sp is the p1s perpendicular
205  T_PR p1sp = std::sqrt( p1sx*p1sx + p1sy*p1sy );
206  // Make sure p1sp is not almost zero
207  if ( p1sp > std::numeric_limits<T_PR>::min() )
208  {
209  p1fsx = ( p1sx*p1sz/p1sp ) * sinXs*cosphis +
210  ( p1sy*p1sm/p1sp ) * sinXs*sinphis +
211  ( p1sx ) * cosXs;
212  p1fsy = ( p1sy*p1sz/p1sp ) * sinXs*cosphis +
213  (-p1sx*p1sm/p1sp ) * sinXs*sinphis +
214  ( p1sy ) * cosXs;
215  p1fsz = (-p1sp ) * sinXs*cosphis +
216  ( T_PR(0.0) ) * sinXs*sinphis +
217  ( p1sz ) * cosXs;
218  // Note a negative sign is different from
219  // Eq. (12) in Perez's paper,
220  // but they are the same due to the random nature of phis.
221  }
222  else
223  {
224  // If the previous p1sp is almost zero
225  // x->y y->z z->x
226  // This set is equivalent to the one in Nanbu's paper
227  p1sp = std::sqrt( p1sy*p1sy + p1sz*p1sz );
228  p1fsy = ( p1sy*p1sx/p1sp ) * sinXs*cosphis +
229  ( p1sz*p1sm/p1sp ) * sinXs*sinphis +
230  ( p1sy ) * cosXs;
231  p1fsz = ( p1sz*p1sx/p1sp ) * sinXs*cosphis +
232  (-p1sy*p1sm/p1sp ) * sinXs*sinphis +
233  ( p1sz ) * cosXs;
234  p1fsx = (-p1sp ) * sinXs*cosphis +
235  ( T_PR(0.0) ) * sinXs*sinphis +
236  ( p1sx ) * cosXs;
237  }
238 
239  T_PR const p2fsx = -p1fsx;
240  T_PR const p2fsy = -p1fsy;
241  T_PR const p2fsz = -p1fsz;
242 
243  // Transform from COM to lab frame
244  T_PR p1fx; T_PR p2fx;
245  T_PR p1fy; T_PR p2fy;
246  T_PR p1fz; T_PR p2fz;
247  if ( vcms > std::numeric_limits<T_PR>::min() )
248  {
249  T_PR const vcDp1fs = vcx*p1fsx + vcy*p1fsy + vcz*p1fsz;
250  T_PR const vcDp2fs = vcx*p2fsx + vcy*p2fsy + vcz*p2fsz;
251  /* factor = (gc-1.0)/vcms; Rewrite to avoid subtraction losing precision when gc is close to 1 */
252  T_PR const factor = gc*gc*inv_c2/(gc+T_PR(1.0));
253  T_PR const factor1 = factor*vcDp1fs + m1*g1s*gc;
254  T_PR const factor2 = factor*vcDp2fs + m2*g2s*gc;
255  p1fx = p1fsx + vcx * factor1;
256  p1fy = p1fsy + vcy * factor1;
257  p1fz = p1fsz + vcz * factor1;
258  p2fx = p2fsx + vcx * factor2;
259  p2fy = p2fsy + vcy * factor2;
260  p2fz = p2fsz + vcz * factor2;
261  }
262  else // If vcms = 0, don't do Lorentz-transform.
263  {
264  p1fx = p1fsx;
265  p1fy = p1fsy;
266  p1fz = p1fsz;
267  p2fx = p2fsx;
268  p2fy = p2fsy;
269  p2fz = p2fsz;
270  }
271 
272  // Rejection method
273  r = amrex::Random(engine);
274  if ( w2 > r*amrex::max(w1, w2) )
275  {
276  u1x = p1fx / m1;
277  u1y = p1fy / m1;
278  u1z = p1fz / m1;
279  }
280  r = amrex::Random(engine);
281  if ( w1 > r*amrex::max(w1, w2) )
282  {
283  u2x = p2fx / m2;
284  u2y = p2fy / m2;
285  u2z = p2fz / m2;
286  }
287 #ifndef AMREX_USE_DPCPP
288  AMREX_ASSERT(!std::isnan(u1x+u1y+u1z+u2x+u2y+u2z));
289  AMREX_ASSERT(!std::isinf(u1x+u1y+u1z+u2x+u2y+u2z));
290 #endif
291 
292  } // if s > std::numeric_limits<T_PR>::min()
293 
294 }
295 
296 #endif // WARPX_PARTICLES_COLLISION_UPDATE_MOMENTUM_PEREZ_ELASTIC_H_
#define AMREX_ASSERT(EX)
#define AMREX_INLINE
#define AMREX_GPU_HOST_DEVICE
AMREX_GPU_HOST_DEVICE AMREX_INLINE void UpdateMomentumPerezElastic(T_PR &u1x, T_PR &u1y, T_PR &u1z, T_PR &u2x, T_PR &u2y, T_PR &u2z, T_PR const n1, T_PR const n2, T_PR const n12, T_PR const q1, T_PR const m1, T_PR const w1, T_PR const q2, T_PR const m2, T_PR const w2, T_R const dt, T_PR const L, T_PR const bmax, amrex::RandomEngine const &engine)
Definition: UpdateMomentumPerezElastic.H:35
static constexpr auto c
vacuum speed of light [m/s]
Definition: constant.H:44
static constexpr auto ep0
vacuum permittivity: dielectric permittivity of vacuum [F/m]
Definition: constant.H:46
static constexpr auto hbar
reduced Planck Constant = h / tau [J*s]
Definition: constant.H:59
Real Random()
AMREX_GPU_HOST_DEVICE constexpr AMREX_FORCE_INLINE const T & max(const T &a, const T &b) noexcept
AMREX_GPU_HOST_DEVICE constexpr AMREX_FORCE_INLINE const T & min(const T &a, const T &b) noexcept
float dt
Definition: stencil.py:442
sp
Definition: stencil.py:27