7#ifndef WARPX_PARTICLES_COLLISION_ELASTIC_COLLISION_PEREZ_H_
8#define WARPX_PARTICLES_COLLISION_ELASTIC_COLLISION_PEREZ_H_
42template <
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,
49 SoaData_type soa_1, SoaData_type soa_2,
50 T_PR
const n1, T_PR
const n2,
51 T_PR
const T1, T_PR
const T2,
52 T_PR
const q1, T_PR
const q2,
53 T_PR
const m1, T_PR
const m2,
54 T_R
const dt, T_R
const global_lamdb, T_PR
const L, T_R
const dV,
56 bool const isSameSpecies, T_index coll_idx)
58 const T_index NI1 = I1e - I1s;
59 const T_index NI2 = I2e - I2s;
74 T_PR lmdD = T_PR(-1.0);
75 if ( L <= T_PR(0.0) ) {
76 if (global_lamdb <= T_PR(0.0)) {
86 const auto rmin =
static_cast<T_PR
>( 1.0/std::cbrt(4.0*
MathConst::pi/3.0*maxn) );
92 const T_PR sigma_max = T_PR(1.0) / (maxn * rmin);
94#if (defined WARPX_DIM_RZ)
101 T_index i1 = I1s + coll_idx;
102 T_index i2 = I2s + coll_idx;
106 for (T_index k = coll_idx; k < max_N; k += min_N)
108#if (defined WARPX_DIM_RZ)
117 T_PR
const theta = theta2[I2[i2]]-theta1[I1[i1]];
118 T_PR
const u1xbuf = u1x[I1[i1]];
119 u1x[I1[i1]] = u1xbuf*std::cos(theta) - u1y[I1[i1]]*std::sin(theta);
120 u1y[I1[i1]] = u1xbuf*std::sin(theta) + u1y[I1[i1]]*std::cos(theta);
131 const T_PR wpmax =
amrex::max(w1[ I1[i1] ],w2[ I2[i2] ]);
132 if (isSameSpecies) { n12 = wpmax*
static_cast<T_PR
>(min_N+max_N-1)/dV; }
133 else { n12 = wpmax*
static_cast<T_PR
>(min_N)/dV; }
136 u1x[ I1[i1] ], u1y[ I1[i1] ], u1z[ I1[i1] ],
137 u2x[ I2[i2] ], u2y[ I2[i2] ], u2z[ I2[i2] ],
138 q1, m1, w1[ I1[i1] ], q2, m2, w2[ I2[i2] ],
139 n12, sigma_max, L, bmax, dt, engine );
141#if (defined WARPX_DIM_RZ)
142 T_PR
const u1xbuf_new = u1x[I1[i1]];
143 u1x[I1[i1]] = u1xbuf_new*std::cos(-theta) - u1y[I1[i1]]*std::sin(-theta);
144 u1y[I1[i1]] = u1xbuf_new*std::sin(-theta) + u1y[I1[i1]]*std::cos(-theta);
#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_R const global_lamdb, T_PR const L, T_R const dV, amrex::RandomEngine const &engine, bool const isSameSpecies, T_index coll_idx)
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 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, T_R const dt, amrex::RandomEngine const &engine)
Definition UpdateMomentumPerezElastic.H:37
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
__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
@ theta
RZ needs all three position components.
Definition WarpXParticleContainer.H:72
@ uz
Definition WarpXParticleContainer.H:70
@ w
weight
Definition WarpXParticleContainer.H:69
@ uy
Definition WarpXParticleContainer.H:70
@ ux
Definition WarpXParticleContainer.H:70