7 #ifndef WARPX_PARTICLES_COLLISION_ELASTIC_COLLISION_PEREZ_H_
8 #define WARPX_PARTICLES_COLLISION_ELASTIC_COLLISION_PEREZ_H_
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,
55 bool const isSameSpecies)
72 if ( T1 <= T_PR(0.0) && L <= T_PR(0.0) )
77 if ( T2 <= T_PR(0.0) && L <= 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] ]; }
95 n1 = n1 / dV; n2 = n2 / dV;
97 int i1 = I1s;
int i2 = I2s;
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; }
107 if (isSameSpecies) n12 *= T_PR(2.0);
111 if ( T1t < T_PR(0.0) || T2t < T_PR(0.0) ) {
118 T_PR rmin = std::pow( T_PR(4.0) *
MathConst::pi / T_PR(3.0) *
122 #if (defined WARPX_DIM_RZ)
129 int i1 = I1s;
int i2 = I2s;
133 #if (defined WARPX_DIM_RZ)
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);
149 u1x[ I1[i1] ], u1y[ I1[i1] ], u1z[ I1[i1] ],
150 u2x[ I2[i2] ], u2y[ I2[i2] ], u2z[ I2[i2] ],
152 q1, m1, w1[ I1[i1] ], q2, m2, w2[ I2[i2] ],
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);
162 ++i1;
if ( i1 ==
static_cast<int>(I1e) ) { i1 = I1s; }
163 ++i2;
if ( i2 ==
static_cast<int>(I2e) ) { i2 = I2s; }
#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