WarpX
ElasticCollisionPerez.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_ELASTIC_COLLISION_PEREZ_H_
8 #define WARPX_PARTICLES_COLLISION_ELASTIC_COLLISION_PEREZ_H_
9 
12 #include "Utils/WarpXConst.H"
13 
14 #include <AMReX_Random.H>
15 
16 
39 template <typename T_index, typename T_PR, typename T_R, typename SoaData_type>
42  T_index const I1s, T_index const I1e,
43  T_index const I2s, T_index const I2e,
44  T_index const* AMREX_RESTRICT I1,
45  T_index const* AMREX_RESTRICT I2,
46  SoaData_type soa_1, SoaData_type soa_2,
47  T_PR const n1, T_PR const n2,
48  T_PR const T1, T_PR const T2,
49  T_PR const q1, T_PR const q2,
50  T_PR const m1, T_PR const m2,
51  T_R const dt, T_PR const L, T_R const dV,
52  amrex::RandomEngine const& engine,
53  bool const isSameSpecies, T_index coll_idx)
54 {
55  const T_index NI1 = I1e - I1s;
56  const T_index NI2 = I2e - I2s;
57  const T_index max_N = amrex::max(NI1,NI2);
58  const T_index min_N = amrex::min(NI1,NI2);
59 
60  T_PR * const AMREX_RESTRICT w1 = soa_1.m_rdata[PIdx::w];
61  T_PR * const AMREX_RESTRICT u1x = soa_1.m_rdata[PIdx::ux];
62  T_PR * const AMREX_RESTRICT u1y = soa_1.m_rdata[PIdx::uy];
63  T_PR * const AMREX_RESTRICT u1z = soa_1.m_rdata[PIdx::uz];
64 
65  T_PR * const AMREX_RESTRICT w2 = soa_2.m_rdata[PIdx::w];
66  T_PR * const AMREX_RESTRICT u2x = soa_2.m_rdata[PIdx::ux];
67  T_PR * const AMREX_RESTRICT u2y = soa_2.m_rdata[PIdx::uy];
68  T_PR * const AMREX_RESTRICT u2z = soa_2.m_rdata[PIdx::uz];
69 
70  // compute Debye length lmdD (if not using a fixed L = Coulomb log)
71  T_PR lmdD = T_PR(-1.0);
72  if ( L <= T_PR(0.0) ) {
73  lmdD = T_PR(1.0)/std::sqrt( n1*q1*q1/(T1*PhysConst::ep0) +
74  n2*q2*q2/(T2*PhysConst::ep0) );
75  }
76 
77  // compute atomic spacing
78  const T_PR maxn = amrex::max(n1,n2);
79  const auto rmin = static_cast<T_PR>( 1.0/std::cbrt(4.0*MathConst::pi/3.0*maxn) );
80 
81  // bmax (screening length) cannot be smaller than atomic spacing
82  const T_PR bmax = amrex::max(lmdD, rmin);
83 
84 #if (defined WARPX_DIM_RZ)
85  T_PR * const AMREX_RESTRICT theta1 = soa_1.m_rdata[PIdx::theta];
86  T_PR * const AMREX_RESTRICT theta2 = soa_2.m_rdata[PIdx::theta];
87 #endif
88 
89  // call UpdateMomentumPerezElastic()
90  {
91  T_index i1 = I1s + coll_idx;
92  T_index i2 = I2s + coll_idx;
93 
94  // we will start from collision number = coll_idx and then add
95  // stride (smaller set size) until we do all collisions (larger set size)
96  for (T_index k = coll_idx; k < max_N; k += min_N)
97  {
98 #if (defined WARPX_DIM_RZ)
99  /* In RZ geometry, macroparticles can collide with other macroparticles
100  * in the same *cylindrical* cell. For this reason, collisions between macroparticles
101  * are actually not local in space. In this case, the underlying assumption is that
102  * particles within the same cylindrical cell represent a cylindrically-symmetry
103  * momentum distribution function. Therefore, here, we temporarily rotate the
104  * momentum of one of the macroparticles in agreement with this cylindrical symmetry.
105  * (This is technically only valid if we use only the m=0 azimuthal mode in the simulation;
106  * there is a corresponding assert statement at initialization.) */
107  T_PR const theta = theta2[I2[i2]]-theta1[I1[i1]];
108  T_PR const u1xbuf = u1x[I1[i1]];
109  u1x[I1[i1]] = u1xbuf*std::cos(theta) - u1y[I1[i1]]*std::sin(theta);
110  u1y[I1[i1]] = u1xbuf*std::sin(theta) + u1y[I1[i1]]*std::cos(theta);
111 #endif
112 
113  // Compute the effective density n12 used to compute the normalized
114  // scattering path s12 in UpdateMomentumPerezElastic().
115  // s12 is defined such that the expected value of the change in particle
116  // velocity is equal to that from the full NxN pairing method, as described
117  // here https://arxiv.org/submit/5758216/view. This method is a direct extension
118  // of the original method by Takizuka and Abe JCP 25 (1977) to weighted particles.
119  T_PR n12;
120  const T_PR wpmax = amrex::max(w1[ I1[i1] ],w2[ I2[i2] ]);
121  if (isSameSpecies) { n12 = wpmax*static_cast<T_PR>(min_N+max_N-1)/dV; }
122  else { n12 = wpmax*static_cast<T_PR>(min_N)/dV; }
123 
125  u1x[ I1[i1] ], u1y[ I1[i1] ], u1z[ I1[i1] ],
126  u2x[ I2[i2] ], u2y[ I2[i2] ], u2z[ I2[i2] ],
127  n1, n2, n12,
128  q1, m1, w1[ I1[i1] ], q2, m2, w2[ I2[i2] ],
129  dt, L, bmax,
130  engine);
131 
132 #if (defined WARPX_DIM_RZ)
133  T_PR const u1xbuf_new = u1x[I1[i1]];
134  u1x[I1[i1]] = u1xbuf_new*std::cos(-theta) - u1y[I1[i1]]*std::sin(-theta);
135  u1y[I1[i1]] = u1xbuf_new*std::sin(-theta) + u1y[I1[i1]]*std::cos(-theta);
136 #endif
137 
138  if (max_N == NI1) {
139  i1 += min_N;
140  }
141  if (max_N == NI2) {
142  i2 += min_N;
143  }
144  }
145  }
146 
147 }
148 
149 #endif // WARPX_PARTICLES_COLLISION_ELASTIC_COLLISION_PEREZ_H_
#define AMREX_RESTRICT
#define AMREX_INLINE
#define AMREX_GPU_HOST_DEVICE
AMREX_GPU_HOST_DEVICE AMREX_INLINE void ElasticCollisionPerez(T_index const I1s, T_index const I1e, T_index const I2s, T_index const I2e, T_index const *AMREX_RESTRICT I1, T_index const *AMREX_RESTRICT I2, SoaData_type soa_1, SoaData_type soa_2, T_PR const n1, T_PR const n2, T_PR const T1, T_PR const T2, T_PR const q1, T_PR const q2, T_PR const m1, T_PR const m2, T_R const dt, T_PR const L, T_R const dV, amrex::RandomEngine const &engine, bool const isSameSpecies, T_index coll_idx)
Definition: ElasticCollisionPerez.H:41
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 ep0
vacuum permittivity: dielectric permittivity of vacuum [F/m]
Definition: constant.H:46
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
@ theta
RZ needs all three position components.
Definition: NamedComponentParticleContainer.H:36
@ uz
Definition: NamedComponentParticleContainer.H:34
@ w
weight
Definition: NamedComponentParticleContainer.H:33
@ uy
Definition: NamedComponentParticleContainer.H:34
@ ux
Definition: NamedComponentParticleContainer.H:34