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