ARTEMIS
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 
10 #include "ComputeTemperature.H"
13 #include "Utils/WarpXConst.H"
14 
15 #include <AMReX_Random.H>
16 
17 
42 template <typename T_index, typename T_PR, typename T_R, typename SoaData_type>
45  T_index const I1s, T_index const I1e,
46  T_index const I2s, T_index const I2e,
47  T_index const* AMREX_RESTRICT I1,
48  T_index const* AMREX_RESTRICT I2,
49  SoaData_type soa_1, SoaData_type soa_2,
50  T_PR const q1, T_PR const q2,
51  T_PR const m1, T_PR const m2,
52  T_PR const T1, T_PR const T2,
53  T_R const dt, T_PR const L, T_R const dV,
54  amrex::RandomEngine const& engine,
55  bool const isSameSpecies)
56 {
57  int NI1 = I1e - I1s;
58  int NI2 = I2e - I2s;
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  // get local T1t and T2t
71  T_PR T1t; T_PR T2t;
72  if ( T1 <= T_PR(0.0) && L <= T_PR(0.0) )
73  {
74  T1t = ComputeTemperature(I1s,I1e,I1,u1x,u1y,u1z,m1);
75  }
76  else { T1t = T1; }
77  if ( T2 <= T_PR(0.0) && L <= T_PR(0.0) )
78  {
79  T2t = ComputeTemperature(I2s,I2e,I2,u2x,u2y,u2z,m2);
80  }
81  else { T2t = T2; }
82 
83  // local density
84  T_PR n1 = T_PR(0.0);
85  T_PR n2 = T_PR(0.0);
86  T_PR n12 = T_PR(0.0);
87  for (int i1=I1s; i1<static_cast<int>(I1e); ++i1) { n1 += w1[ I1[i1] ]; }
88  for (int i2=I2s; i2<static_cast<int>(I2e); ++i2) { n2 += w2[ I2[i2] ]; }
89  // Intra-species: the density is in fact the sum of the density of
90  // each sub-group of particles
91  if (isSameSpecies) {
92  n1 = n1 + n2;
93  n2 = n1;
94  }
95  n1 = n1 / dV; n2 = n2 / dV;
96  {
97  int i1 = I1s; int i2 = I2s;
98  for (int k = 0; k < amrex::max(NI1,NI2); ++k)
99  {
100  n12 += amrex::min( w1[ I1[i1] ], w2[ I2[i2] ] );
101  ++i1; if ( i1 == static_cast<int>(I1e) ) { i1 = I1s; }
102  ++i2; if ( i2 == static_cast<int>(I2e) ) { i2 = I2s; }
103  }
104  n12 = n12 / dV;
105  }
106  // Intra-species: Apply correction in collision rate
107  if (isSameSpecies) n12 *= T_PR(2.0);
108 
109  // compute Debye length lmdD
110  T_PR lmdD;
111  if ( T1t < T_PR(0.0) || T2t < T_PR(0.0) ) {
112  lmdD = T_PR(0.0);
113  }
114  else {
115  lmdD = T_PR(1.0)/std::sqrt( n1*q1*q1/(T1t*PhysConst::ep0) +
116  n2*q2*q2/(T2t*PhysConst::ep0) );
117  }
118  T_PR rmin = std::pow( T_PR(4.0) * MathConst::pi / T_PR(3.0) *
119  amrex::max(n1,n2), T_PR(-1.0/3.0) );
120  lmdD = amrex::max(lmdD, rmin);
121 
122 #if (defined WARPX_DIM_RZ)
123  T_PR * const AMREX_RESTRICT theta1 = soa_1.m_rdata[PIdx::theta];
124  T_PR * const AMREX_RESTRICT theta2 = soa_2.m_rdata[PIdx::theta];
125 #endif
126 
127  // call UpdateMomentumPerezElastic()
128  {
129  int i1 = I1s; int i2 = I2s;
130  for (int k = 0; k < amrex::max(NI1,NI2); ++k)
131  {
132 
133 #if (defined WARPX_DIM_RZ)
134  /* In RZ geometry, macroparticles can collide with other macroparticles
135  * in the same *cylindrical* cell. For this reason, collisions between macroparticles
136  * are actually not local in space. In this case, the underlying assumption is that
137  * particles within the same cylindrical cell represent a cylindrically-symmetry
138  * momentum distribution function. Therefore, here, we temporarily rotate the
139  * momentum of one of the macroparticles in agreement with this cylindrical symmetry.
140  * (This is technically only valid if we use only the m=0 azimuthal mode in the simulation;
141  * there is a corresponding assert statement at initialization.) */
142  T_PR const theta = theta2[I2[i2]]-theta1[I1[i1]];
143  T_PR const u1xbuf = u1x[I1[i1]];
144  u1x[I1[i1]] = u1xbuf*std::cos(theta) - u1y[I1[i1]]*std::sin(theta);
145  u1y[I1[i1]] = u1xbuf*std::sin(theta) + u1y[I1[i1]]*std::cos(theta);
146 #endif
147 
149  u1x[ I1[i1] ], u1y[ I1[i1] ], u1z[ I1[i1] ],
150  u2x[ I2[i2] ], u2y[ I2[i2] ], u2z[ I2[i2] ],
151  n1, n2, n12,
152  q1, m1, w1[ I1[i1] ], q2, m2, w2[ I2[i2] ],
153  dt, L, lmdD,
154  engine);
155 
156 #if (defined WARPX_DIM_RZ)
157  T_PR const u1xbuf_new = u1x[I1[i1]];
158  u1x[I1[i1]] = u1xbuf_new*std::cos(-theta) - u1y[I1[i1]]*std::sin(-theta);
159  u1y[I1[i1]] = u1xbuf_new*std::sin(-theta) + u1y[I1[i1]]*std::cos(-theta);
160 #endif
161 
162  ++i1; if ( i1 == static_cast<int>(I1e) ) { i1 = I1s; }
163  ++i2; if ( i2 == static_cast<int>(I2e) ) { i2 = I2s; }
164  }
165  }
166 
167 }
168 
169 #endif // WARPX_PARTICLES_COLLISION_ELASTIC_COLLISION_PEREZ_H_
#define AMREX_RESTRICT
#define AMREX_INLINE
#define AMREX_GPU_HOST_DEVICE
AMREX_GPU_HOST_DEVICE T_R ComputeTemperature(T_index const Is, T_index const Ie, T_index const *AMREX_RESTRICT I, T_R const *AMREX_RESTRICT ux, T_R const *AMREX_RESTRICT uy, T_R const *AMREX_RESTRICT uz, T_R const m)
Definition: ComputeTemperature.H:15
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 q1, T_PR const q2, T_PR const m1, T_PR const m2, T_PR const T1, T_PR const T2, T_R const dt, T_PR const L, T_R const dV, amrex::RandomEngine const &engine, bool const isSameSpecies)
Definition: ElasticCollisionPerez.H:44
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
static constexpr auto ep0
vacuum permittivity: dielectric permittivity of vacuum [F/m]
Definition: constant.H:46
static constexpr amrex::Real pi
ratio of a circle's circumference to its diameter
Definition: constant.H:23
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE constexpr const T & min(const T &a, const T &b) noexcept
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE constexpr const T & max(const T &a, const T &b) noexcept
int dt
Definition: stencil.py:440
@ theta
RZ needs all three position components.
Definition: NamedComponentParticleContainer.H:27
@ uz
Definition: NamedComponentParticleContainer.H:25
@ w
weight
Definition: NamedComponentParticleContainer.H:24
@ uy
Definition: NamedComponentParticleContainer.H:25
@ ux
Definition: NamedComponentParticleContainer.H:25