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] LmdD 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  * To see if there are nan or inf updated velocities,
25  * compile with USE_ASSERTION=TRUE.
26 */
27 
28 template <typename T_PR, typename T_R>
31  T_PR& u1x, T_PR& u1y, T_PR& u1z, T_PR& u2x, T_PR& u2y, T_PR& u2z,
32  T_PR const n1, T_PR const n2, T_PR const n12,
33  T_PR const q1, T_PR const m1, T_PR const w1,
34  T_PR const q2, T_PR const m2, T_PR const w2,
35  T_R const dt, T_PR const L, T_PR const lmdD,
36  amrex::RandomEngine const& engine)
37 {
38 
39  T_PR const diffx = amrex::Math::abs(u1x-u2x);
40  T_PR const diffy = amrex::Math::abs(u1y-u2y);
41  T_PR const diffz = amrex::Math::abs(u1z-u2z);
42  T_PR const diffm = std::sqrt(diffx*diffx+diffy*diffy+diffz*diffz);
43  T_PR const summm = std::sqrt(u1x*u1x+u1y*u1y+u1z*u1z) + std::sqrt(u2x*u2x+u2y*u2y+u2z*u2z);
44  // If g = u1 - u2 = 0, do not collide.
45  // Or if the relative difference is less than 1.0e-10.
46  if ( diffm < std::numeric_limits<T_PR>::min() || diffm/summm < 1.0e-10 ) { return; }
47 
48  T_PR constexpr inv_c2 = T_PR(1.0) / ( PhysConst::c * PhysConst::c );
49 
50  // Compute Lorentz factor gamma
51  T_PR const g1 = std::sqrt( T_PR(1.0) + (u1x*u1x+u1y*u1y+u1z*u1z)*inv_c2 );
52  T_PR const g2 = std::sqrt( T_PR(1.0) + (u2x*u2x+u2y*u2y+u2z*u2z)*inv_c2 );
53 
54  // Compute momenta
55  T_PR const p1x = u1x * m1;
56  T_PR const p1y = u1y * m1;
57  T_PR const p1z = u1z * m1;
58  T_PR const p2x = u2x * m2;
59  T_PR const p2y = u2y * m2;
60  T_PR const p2z = u2z * m2;
61 
62  // Compute center-of-mass (COM) velocity and gamma
63  T_PR const mass_g = m1 * g1 + m2 * g2;
64  T_PR const vcx = (p1x+p2x) / mass_g;
65  T_PR const vcy = (p1y+p2y) / mass_g;
66  T_PR const vcz = (p1z+p2z) / mass_g;
67  T_PR const vcms = vcx*vcx + vcy*vcy + vcz*vcz;
68  T_PR const gc = T_PR(1.0) / std::sqrt( T_PR(1.0) - vcms*inv_c2 );
69 
70  // Compute vc dot v1 and v2
71  T_PR const vcDv1 = (vcx*u1x + vcy*u1y + vcz*u1z) / g1;
72  T_PR const vcDv2 = (vcx*u2x + vcy*u2y + vcz*u2z) / g2;
73 
74  // Compute p1 star
75  T_PR p1sx;
76  T_PR p1sy;
77  T_PR p1sz;
78  if ( vcms > std::numeric_limits<T_PR>::min() )
79  {
80  T_PR const lorentz_tansform_factor =
81  ( (gc-T_PR(1.0))/vcms*vcDv1 - gc )*m1*g1;
82  p1sx = p1x + vcx*lorentz_tansform_factor;
83  p1sy = p1y + vcy*lorentz_tansform_factor;
84  p1sz = p1z + vcz*lorentz_tansform_factor;
85  }
86  else // If vcms = 0, don't do Lorentz-transform.
87  {
88  p1sx = p1x;
89  p1sy = p1y;
90  p1sz = p1z;
91  }
92  T_PR const p1sm = std::sqrt( p1sx*p1sx + p1sy*p1sy + p1sz*p1sz );
93 
94  // Compute gamma star
95  T_PR const g1s = ( T_PR(1.0) - vcDv1*inv_c2 )*gc*g1;
96  T_PR const g2s = ( T_PR(1.0) - vcDv2*inv_c2 )*gc*g2;
97 
98  // Compute the Coulomb log lnLmd
99  T_PR lnLmd;
100  if ( L > T_PR(0.0) ) { lnLmd = L; }
101  else
102  {
103  // Compute b0
104  T_PR const b0 = amrex::Math::abs(q1*q2) * inv_c2 /
105  (T_PR(4.0)*MathConst::pi*PhysConst::ep0) * gc/mass_g *
106  ( m1*g1s*m2*g2s/(p1sm*p1sm*inv_c2) + T_PR(1.0) );
107 
108  // Compute the minimal impact parameter
109  constexpr T_PR hbar_pi = static_cast<T_PR>(PhysConst::hbar*MathConst::pi);
110  const T_PR bmin = amrex::max(hbar_pi/p1sm, b0);
111 
112  // Compute the Coulomb log lnLmd
113  lnLmd = amrex::max( T_PR(2.0),
114  T_PR(0.5)*std::log(T_PR(1.0)+lmdD*lmdD/(bmin*bmin)) );
115  }
116 
117  // Compute s
118  const auto tts = m1*g1s*m2*g2s/(inv_c2*p1sm*p1sm) + T_PR(1.0);
119  const auto tts2 = tts*tts;
120  T_PR s = n1*n2/n12 * dt*lnLmd*q1*q1*q2*q2 /
121  ( T_PR(4.0) * MathConst::pi * PhysConst::ep0 * PhysConst::ep0 *
122  m1*g1*m2*g2/(inv_c2*inv_c2) ) * gc*p1sm/mass_g * tts2;
123 
124  // Compute s'
125  const auto cbrt_n1 = std::cbrt(n1);
126  const auto cbrt_n2 = std::cbrt(n2);
127  const auto coeff = static_cast<T_PR>(
128  std::pow(4.0*MathConst::pi/3.0,1.0/3.0));
129  T_PR const vrel = mass_g*p1sm/(m1*g1s*m2*g2s*gc);
130  T_PR const sp = coeff * n1*n2/n12 * dt * vrel * (m1+m2) /
131  amrex::max( m1*cbrt_n1*cbrt_n1,
132  m2*cbrt_n2*cbrt_n2);
133 
134  // Determine s
135  s = amrex::min(s,sp);
136 
137  // Get random numbers
138  T_PR r = amrex::Random(engine);
139 
140  // Compute scattering angle
141  T_PR cosXs;
142  T_PR sinXs;
143  if ( s <= T_PR(0.1) )
144  {
145  while ( true )
146  {
147  cosXs = T_PR(1.0) + s * std::log(r);
148  // Avoid the bug when r is too small such that cosXs < -1
149  if ( cosXs >= T_PR(-1.0) ) { break; }
150  r = amrex::Random(engine);
151  }
152  }
153  else if ( s > T_PR(0.1) && s <= T_PR(3.0) )
154  {
155  T_PR const Ainv = static_cast<T_PR>(
156  0.0056958 + 0.9560202*s - 0.508139*s*s +
157  0.47913906*s*s*s - 0.12788975*s*s*s*s + 0.02389567*s*s*s*s*s);
158  cosXs = Ainv * std::log( std::exp(T_PR(-1.0)/Ainv) +
159  T_PR(2.0) * r * std::sinh(T_PR(1.0)/Ainv) );
160  }
161  else if ( s > T_PR(3.0) && s <= T_PR(6.0) )
162  {
163  T_PR const A = T_PR(3.0) * std::exp(-s);
164  cosXs = T_PR(1.0)/A * std::log( std::exp(-A) +
165  T_PR(2.0) * r * std::sinh(A) );
166  }
167  else
168  {
169  cosXs = T_PR(2.0) * r - T_PR(1.0);
170  }
171  sinXs = std::sqrt(T_PR(1.0) - cosXs*cosXs);
172 
173  // Get random azimuthal angle
174  T_PR const phis = amrex::Random(engine) * T_PR(2.0) * MathConst::pi;
175  T_PR const cosphis = std::cos(phis);
176  T_PR const sinphis = std::sin(phis);
177 
178  // Compute post-collision momenta pfs in COM
179  T_PR p1fsx;
180  T_PR p1fsy;
181  T_PR p1fsz;
182  // p1sp is the p1s perpendicular
183  T_PR p1sp = std::sqrt( p1sx*p1sx + p1sy*p1sy );
184  // Make sure p1sp is not almost zero
185  if ( p1sp > std::numeric_limits<T_PR>::min() )
186  {
187  p1fsx = ( p1sx*p1sz/p1sp ) * sinXs*cosphis +
188  ( p1sy*p1sm/p1sp ) * sinXs*sinphis +
189  ( p1sx ) * cosXs;
190  p1fsy = ( p1sy*p1sz/p1sp ) * sinXs*cosphis +
191  (-p1sx*p1sm/p1sp ) * sinXs*sinphis +
192  ( p1sy ) * cosXs;
193  p1fsz = (-p1sp ) * sinXs*cosphis +
194  ( T_PR(0.0) ) * sinXs*sinphis +
195  ( p1sz ) * cosXs;
196  // Note a negative sign is different from
197  // Eq. (12) in Perez's paper,
198  // but they are the same due to the random nature of phis.
199  }
200  else
201  {
202  // If the previous p1sp is almost zero
203  // x->y y->z z->x
204  // This set is equivalent to the one in Nanbu's paper
205  p1sp = std::sqrt( p1sy*p1sy + p1sz*p1sz );
206  p1fsy = ( p1sy*p1sx/p1sp ) * sinXs*cosphis +
207  ( p1sz*p1sm/p1sp ) * sinXs*sinphis +
208  ( p1sy ) * cosXs;
209  p1fsz = ( p1sz*p1sx/p1sp ) * sinXs*cosphis +
210  (-p1sy*p1sm/p1sp ) * sinXs*sinphis +
211  ( p1sz ) * cosXs;
212  p1fsx = (-p1sp ) * sinXs*cosphis +
213  ( T_PR(0.0) ) * sinXs*sinphis +
214  ( p1sx ) * cosXs;
215  }
216 
217  T_PR const p2fsx = -p1fsx;
218  T_PR const p2fsy = -p1fsy;
219  T_PR const p2fsz = -p1fsz;
220 
221  // Transform from COM to lab frame
222  T_PR p1fx; T_PR p2fx;
223  T_PR p1fy; T_PR p2fy;
224  T_PR p1fz; T_PR p2fz;
225  if ( vcms > std::numeric_limits<T_PR>::min() )
226  {
227  T_PR const vcDp1fs = vcx*p1fsx + vcy*p1fsy + vcz*p1fsz;
228  T_PR const vcDp2fs = vcx*p2fsx + vcy*p2fsy + vcz*p2fsz;
229  T_PR const factor = (gc-T_PR(1.0))/vcms;
230  T_PR const factor1 = factor*vcDp1fs + m1*g1s*gc;
231  T_PR const factor2 = factor*vcDp2fs + m2*g2s*gc;
232  p1fx = p1fsx + vcx * factor1;
233  p1fy = p1fsy + vcy * factor1;
234  p1fz = p1fsz + vcz * factor1;
235  p2fx = p2fsx + vcx * factor2;
236  p2fy = p2fsy + vcy * factor2;
237  p2fz = p2fsz + vcz * factor2;
238  }
239  else // If vcms = 0, don't do Lorentz-transform.
240  {
241  p1fx = p1fsx;
242  p1fy = p1fsy;
243  p1fz = p1fsz;
244  p2fx = p2fsx;
245  p2fy = p2fsy;
246  p2fz = p2fsz;
247  }
248 
249  // Rejection method
250  r = amrex::Random(engine);
251  if ( w2 > r*amrex::max(w1, w2) )
252  {
253  u1x = p1fx / m1;
254  u1y = p1fy / m1;
255  u1z = p1fz / m1;
256 #ifndef AMREX_USE_DPCPP
257  AMREX_ASSERT(!std::isnan(u1x+u1y+u1z+u2x+u2y+u2z));
258  AMREX_ASSERT(!std::isinf(u1x+u1y+u1z+u2x+u2y+u2z));
259 #endif
260  }
261  r = amrex::Random(engine);
262  if ( w1 > r*amrex::max(w1, w2) )
263  {
264  u2x = p2fx / m2;
265  u2y = p2fy / m2;
266  u2z = p2fz / m2;
267 #ifndef AMREX_USE_DPCPP
268  AMREX_ASSERT(!std::isnan(u1x+u1y+u1z+u2x+u2y+u2z));
269  AMREX_ASSERT(!std::isinf(u1x+u1y+u1z+u2x+u2y+u2z));
270 #endif
271  }
272 
273 }
274 
275 #endif // WARPX_PARTICLES_COLLISION_UPDATE_MOMENTUM_PEREZ_ELASTIC_H_
AMREX_GPU_HOST_DEVICE constexpr const AMREX_FORCE_INLINE T & max(const T &a, const T &b) noexcept
Real Random()
static constexpr auto c
vacuum speed of light [m/s]
Definition: constant.H:44
int dt
Definition: Stencil.py:468
#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 lmdD, amrex::RandomEngine const &engine)
Definition: UpdateMomentumPerezElastic.H:30
#define AMREX_INLINE
static constexpr auto hbar
reduced Planck Constant = h / tau [J*s]
Definition: constant.H:59
AMREX_GPU_HOST_DEVICE constexpr const AMREX_FORCE_INLINE T & min(const T &a, const T &b) noexcept
static constexpr auto ep0
vacuum permittivity: dielectric permittivity of vacuum [F/m]
Definition: constant.H:46
s
Definition: plot_results.py:104
static constexpr amrex::Real pi
ratio of a circle&#39;s circumference to its diameter
Definition: constant.H:23
#define AMREX_ASSERT(EX)