WarpX
Loading...
Searching...
No Matches
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
23
24/* \brief Helper function to sample cosine of scattering angle from Nanbu's
25 * cumulative scattering distribution. See PhysRevE.55.4642 (1997).
26 * @param[in] s12 is the distribution variance
27 * @param[in] engine is a random number engine.
28 */
29
30template <typename T_PR>
33 T_PR const s12,
34 amrex::RandomEngine const& engine )
35{
36 using namespace amrex::literals;
37
38 T_PR r = amrex::Random(engine);
39
40 T_PR cosXs;
41 if ( s12 <= T_PR(0.1) )
42 {
43 while ( true )
44 {
45 cosXs = T_PR(1.0) + s12 * std::log(r);
46 // Avoid the bug when r is too small such that cosXs < -1
47 if ( cosXs >= T_PR(-1.0) ) { break; }
48 r = amrex::Random(engine);
49 }
50 }
51 else if ( s12 > T_PR(0.1) && s12 <= T_PR(3.0) )
52 {
53 T_PR const Ainv = static_cast<T_PR>(
54 0.0056958 + 0.9560202*s12 - 0.508139*s12*s12 +
55 0.47913906*s12*s12*s12 - 0.12788975*s12*s12*s12*s12 + 0.02389567*s12*s12*s12*s12*s12);
56 cosXs = Ainv * std::log( std::exp(T_PR(-1.0)/Ainv) +
57 T_PR(2.0) * r * std::sinh(T_PR(1.0)/Ainv) );
58 }
59 else if ( s12 > T_PR(3.0) && s12 <= T_PR(6.0) )
60 {
61 T_PR const A = T_PR(3.0) * std::exp(-s12);
62 cosXs = T_PR(1.0)/A * std::log( std::exp(-A) +
63 T_PR(2.0) * r * std::sinh(A) );
64 }
65 else
66 {
67 cosXs = T_PR(2.0) * r - T_PR(1.0);
68 }
69
70 return cosXs;
71}
72
73/* \brief Update particle velocities according to
74 * F. Perez et al., Phys.Plasmas.19.083104 (2012),
75 * which is based on Nanbu's method, PhysRevE.55.4642 (1997).
76 * @param[in] bmax is max(Debye length, minimal interparticle distance).
77 * @param[in] L is the Coulomb log. A fixed L will be used if L > 0,
78 * otherwise L will be calculated based on the algorithm.
79 * @param[in] n12 = max(w1,w2)*min(N1,N2)/dV is the effective density used for s12
80 * @param[in] sigma_max is the maximum cross section based on mfp = atomic spacing
81 * @param[in] scattering_model is the model used to obtain the scattering angle
82 * used for the normalized scattering length s12 (see Sec. II.C of Perez et al.)
83 * To see if there are nan or inf updated velocities,
84 * compile with USE_ASSERTION=TRUE.
85 *
86 * Updates and corrections to the original publication are documented in
87 * https://github.com/BLAST-WarpX/warpx/issues/429
88 * https://github.com/BLAST-WarpX/warpx/files/3799803/main.pdf
89 */
90
91template <typename T_PR, typename T_R>
94 T_PR& u1x, T_PR& u1y, T_PR& u1z, T_PR& u2x, T_PR& u2y, T_PR& u2z,
95 T_PR const q1, T_PR const m1, T_PR const w1,
96 T_PR const q2, T_PR const m2, T_PR const w2,
97 T_PR const n12, T_PR const sigma_max,
98 T_PR const L, T_PR const bmax,
99 CumulativeScatteringModel const scattering_model,
100 T_R const dt, amrex::RandomEngine const& engine )
101{
102
103 auto constexpr inv_c2 = PhysConst::inv_c2_v<T_PR>;
104
105 // Compute Lorentz factor gamma
106 T_PR const gb1sq = (u1x*u1x + u1y*u1y + u1z*u1z)*inv_c2;
107 T_PR const gb2sq = (u2x*u2x + u2y*u2y + u2z*u2z)*inv_c2;
108 T_PR const g1 = std::sqrt( T_PR(1.0) + gb1sq );
109 T_PR const g2 = std::sqrt( T_PR(1.0) + gb2sq );
110
111 T_PR const diffx = u1x-u2x;
112 T_PR const diffy = u1y-u2y;
113 T_PR const diffz = u1z-u2z;
114 T_PR const diffm = std::sqrt((diffx*diffx+diffy*diffy+diffz*diffz)*inv_c2);
115 T_PR const summm = std::sqrt(gb1sq) + std::sqrt(gb2sq);
116 // If g = u1 - u2 = 0, do not collide.
117 // Or if the relative difference is less than 1.0e-10.
118 if ( diffm < std::numeric_limits<T_PR>::min() || diffm/summm < 1.0e-10 ) { return; }
119
120 // Compute momenta
121 T_PR const p1x = u1x * m1;
122 T_PR const p1y = u1y * m1;
123 T_PR const p1z = u1z * m1;
124 T_PR const p2x = u2x * m2;
125 T_PR const p2y = u2y * m2;
126 T_PR const p2z = u2z * m2;
127
128 // Compute center-of-momentum (COM) velocity and gamma
129 T_PR const mass_g = m1 * g1 + m2 * g2;
130 T_PR const vcx = (p1x+p2x) / mass_g;
131 T_PR const vcy = (p1y+p2y) / mass_g;
132 T_PR const vcz = (p1z+p2z) / mass_g;
133 T_PR const vcms = vcx*vcx + vcy*vcy + vcz*vcz;
134 T_PR const gc = T_PR(1.0) / std::sqrt( T_PR(1.0) - vcms*inv_c2 );
135
136 // Compute vc dot v1 and v2
137 T_PR const vcDv1 = (vcx*u1x + vcy*u1y + vcz*u1z) / g1;
138 T_PR const vcDv2 = (vcx*u2x + vcy*u2y + vcz*u2z) / g2;
139
140 // Compute p1 star
141 // lorentz_transform_factor = ( (gc-1.0)/vcms*vcDv1 - gc )*m1*g1
142 // Rewrite by multiplying and dividing first term by (gc+1) to
143 // avoid loss of precision when gc is close to 1
144 T_PR const lorentz_transform_factor =
145 ( gc*gc*vcDv1*inv_c2/(T_PR(1.0) + gc) - gc )*m1*g1;
146 T_PR const p1sx = p1x + vcx*lorentz_transform_factor;
147 T_PR const p1sy = p1y + vcy*lorentz_transform_factor;
148 T_PR const p1sz = p1z + vcz*lorentz_transform_factor;
149 T_PR const p1sm = std::sqrt( p1sx*p1sx + p1sy*p1sy + p1sz*p1sz );
150
151 // Compute gamma star
152 T_PR const g1s = ( T_PR(1.0) - vcDv1*inv_c2 )*gc*g1;
153 T_PR const g2s = ( T_PR(1.0) - vcDv2*inv_c2 )*gc*g2;
154
155 // Compute relative velocity in center-of-momentum frame
156 // (not a Lorentz invariant quantity)
157 T_PR const muRst = g1s*m1*g2s*m2/(g1s*m1 + g2s*m2);
158 T_PR const vrelst = p1sm/muRst; // |v1s - v2s|
159
160 // Compute invariant relative velocity in center-of-momentum frame
161 // (Lorentz invariant quantity)
162 T_PR const denom = T_PR(1.0) + p1sm*p1sm/(m1*g1s*m2*g2s)*inv_c2; // (1.0 - v1s*v2s/c^2)
163 T_PR const vrelst_invar = vrelst/denom; // |v1s - v2s|/(1.0 - v1s*v2s/c^2)
164
165 // Compute s12
166 T_PR s12 = 0;
167 if (p1sm > std::numeric_limits<T_PR>::min()) {
168
169 // Writing b0 in a form that is directly analagous to the well-known non-relativistic form.
170 // See Eq. 3.3.2 in Principles of Plasma Discharges and Material Processing by
171 // M. A. Lieberman and A. J. Lichtenberg.
172 // Note that b0 on Eq. 22 of Perez POP 19 (2012) is bmin = b0/2,
173 // Note: there is a typo in Eq 22 of Perez, the last square is incorrect!
174 // See the SMILEI documentation: https://smileipic.github.io/Smilei/Understand/collisions.html
175 // and https://github.com/BLAST-WarpX/warpx/files/3799803/main.pdf from GitHub #429
176 T_PR const b0 = amrex::Math::abs(q1*q2) /
177 (T_PR(2.0)*MathConst::pi*PhysConst::epsilon_0*muRst*vrelst*vrelst_invar);
178
179 // Compute the Coulomb log lnLmd first
180 T_PR lnLmd;
181 if ( L > T_PR(0.0) ) { lnLmd = L; }
182 else
183 {
184 // Compute the minimum impact parameter from quantum: bqm = lDB/(4*pi) = hbar/(2*p1sm)
185 // See NRL formulary. Also see "An introduction to the physics of the Coulomb logarithm,
186 // with emphasis on quantum-mechanical effects", J. Plasma Phys. vol. 85 (2019). by J.A. Krommes.
187 // Note: The formula in Perez 2012 and in Lee and More 1984 uses h rather than
188 // hbar for bqm. If this is used, then the transition energy where bmin goes from classical
189 // to quantum is only 2.5 eV for electrons; compared to 100 eV when using hbar.
190 const T_PR bmin_qm = static_cast<T_PR>(PhysConst::hbar*0.5/p1sm);
191
192 // Set the minimum impact parameter
193 const T_PR bmin = amrex::max(bmin_qm, T_PR(0.5)*b0);
194
195 // Compute the Coulomb log lnLmd
196 lnLmd = amrex::max( T_PR(2.0),
197 T_PR(0.5)*std::log(T_PR(1.0) + bmax*bmax/(bmin*bmin)) );
198 }
199
200 // Compute s12 with sigma limited by sigma_max where mfp = atomic spacing
201 // See https://github.com/user-attachments/files/16555064/CoulombScattering_s12.pdf
202 // for a proof that this expression for s12 is the same as Eq. 9 of Perez 2012 when
203 // sigma_eff = pi*b0^2*lnLmd
204 const T_PR sigma_eff = amrex::min(T_PR(MathConst::pi)*b0*b0*lnLmd,sigma_max);
205 s12 = sigma_eff * n12 * dt * vrelst * g1s*g2s/(g1*g2);
206
207 }
208
209 // Only modify momenta if s12 is non-zero
210 if (s12 > std::numeric_limits<T_PR>::min()) {
211
212 T_PR cosXs = T_PR(0.0);
213 if (scattering_model == NANBU) {
214 // Obtain scattering angle by sampling from Nanbu's distibution
215 cosXs = SampleNanbu(s12, engine);
216 }
217 else if (scattering_model == BOBYLEV) {
218 // Obtain scattering angle from Bobylev's delta function distribution
219 cosXs = T_PR(1.0) - std::min(s12,T_PR(2.0));
220 }
221 T_PR const sinXs = std::sqrt(T_PR(1.0) - cosXs*cosXs);
222
223 // Get random azimuthal angle
224 T_PR const phis = amrex::Random(engine) * T_PR(2.0) * MathConst::pi;
225 T_PR const cosphis = std::cos(phis);
226 T_PR const sinphis = std::sin(phis);
227
228 // Compute post-collision momenta pfs in COM
229 T_PR p1fsx;
230 T_PR p1fsy;
231 T_PR p1fsz;
232 // p1sp is the p1s perpendicular
233 T_PR p1sp = std::sqrt( p1sx*p1sx + p1sy*p1sy );
234 // Make sure p1sp is not almost zero
235 if ( p1sp > std::numeric_limits<T_PR>::min() )
236 {
237 p1fsx = ( p1sx*p1sz/p1sp ) * sinXs*cosphis +
238 ( p1sy*p1sm/p1sp ) * sinXs*sinphis +
239 ( p1sx ) * cosXs;
240 p1fsy = ( p1sy*p1sz/p1sp ) * sinXs*cosphis +
241 (-p1sx*p1sm/p1sp ) * sinXs*sinphis +
242 ( p1sy ) * cosXs;
243 p1fsz = (-p1sp ) * sinXs*cosphis +
244 ( T_PR(0.0) ) * sinXs*sinphis +
245 ( p1sz ) * cosXs;
246 // Note a negative sign is different from
247 // Eq. (12) in Perez's paper,
248 // but they are the same due to the random nature of phis.
249 }
250 else
251 {
252 // If the previous p1sp is almost zero
253 // x->y y->z z->x
254 // This set is equivalent to the one in Nanbu's paper
255 p1sp = std::sqrt( p1sy*p1sy + p1sz*p1sz );
256 p1fsy = ( p1sy*p1sx/p1sp ) * sinXs*cosphis +
257 ( p1sz*p1sm/p1sp ) * sinXs*sinphis +
258 ( p1sy ) * cosXs;
259 p1fsz = ( p1sz*p1sx/p1sp ) * sinXs*cosphis +
260 (-p1sy*p1sm/p1sp ) * sinXs*sinphis +
261 ( p1sz ) * cosXs;
262 p1fsx = (-p1sp ) * sinXs*cosphis +
263 ( T_PR(0.0) ) * sinXs*sinphis +
264 ( p1sx ) * cosXs;
265 }
266
267 T_PR const p2fsx = -p1fsx;
268 T_PR const p2fsy = -p1fsy;
269 T_PR const p2fsz = -p1fsz;
270
271 // Transform from COM to lab frame
272 T_PR const vcDp1fs = vcx*p1fsx + vcy*p1fsy + vcz*p1fsz;
273 T_PR const vcDp2fs = vcx*p2fsx + vcy*p2fsy + vcz*p2fsz;
274 /* factor = (gc-1.0)/vcms; Rewrite to avoid subtraction losing precision when gc is close to 1 */
275 T_PR const factor = gc*gc*inv_c2/(gc+T_PR(1.0));
276 T_PR const factor1 = factor*vcDp1fs + m1*g1s*gc;
277 T_PR const factor2 = factor*vcDp2fs + m2*g2s*gc;
278 T_PR const p1fx = p1fsx + vcx * factor1;
279 T_PR const p1fy = p1fsy + vcy * factor1;
280 T_PR const p1fz = p1fsz + vcz * factor1;
281 T_PR const p2fx = p2fsx + vcx * factor2;
282 T_PR const p2fy = p2fsy + vcy * factor2;
283 T_PR const p2fz = p2fsz + vcz * factor2;
284
285 // Rejection method
286 T_PR r = amrex::Random(engine);
287 if ( w2 > r*amrex::max(w1, w2) )
288 {
289 u1x = p1fx / m1;
290 u1y = p1fy / m1;
291 u1z = p1fz / m1;
292 }
293 r = amrex::Random(engine);
294 if ( w1 > r*amrex::max(w1, w2) )
295 {
296 u2x = p2fx / m2;
297 u2y = p2fy / m2;
298 u2z = p2fz / m2;
299 }
300#ifndef AMREX_USE_DPCPP
301 AMREX_ASSERT(!std::isnan(u1x+u1y+u1z+u2x+u2y+u2z));
302 AMREX_ASSERT(!std::isinf(u1x+u1y+u1z+u2x+u2y+u2z));
303#endif
304
305 } // if s > std::numeric_limits<T_PR>::min()
306
307}
308
309#endif // WARPX_PARTICLES_COLLISION_UPDATE_MOMENTUM_PEREZ_ELASTIC_H_
#define AMREX_ASSERT(EX)
#define AMREX_INLINE
#define AMREX_GPU_HOST_DEVICE
CumulativeScatteringModel
Definition UpdateMomentumPerezElastic.H:19
@ NANBU
Definition UpdateMomentumPerezElastic.H:20
@ BOBYLEV
Definition UpdateMomentumPerezElastic.H:21
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 q1, T_PR const m1, T_PR const w1, T_PR const q2, T_PR const m2, T_PR const w2, T_PR const n12, T_PR const sigma_max, T_PR const L, T_PR const bmax, CumulativeScatteringModel const scattering_model, T_R const dt, amrex::RandomEngine const &engine)
Definition UpdateMomentumPerezElastic.H:93
AMREX_GPU_HOST_DEVICE AMREX_INLINE T_PR SampleNanbu(T_PR const s12, amrex::RandomEngine const &engine)
Definition UpdateMomentumPerezElastic.H:32
Real Random()
constexpr auto epsilon_0
vacuum permittivity: dielectric permittivity of vacuum [F/m]
Definition constant.H:163
constexpr auto inv_c2_v
inverse of the square of the vacuum speed of light [s^2/m^2] (variable template)
Definition constant.H:153
constexpr auto hbar
reduced Planck Constant = h / tau [J*s]
Definition constant.H:181
constexpr auto pi
ratio of a circle's circumference to its diameter
Definition constant.H:29
__host__ __device__ constexpr const T & min(const T &a, const T &b) noexcept
__host__ __device__ constexpr const T & max(const T &a, const T &b) noexcept