WarpX
Loading...
Searching...
No Matches
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
40
41template <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 n1, T_PR const n2,
50 T_PR const T1, T_PR const T2,
51 T_PR const q1, T_PR const q2,
52 T_PR const m1, T_PR const m2,
53 T_R const dt, T_R const global_lamdb, T_PR const L, T_R const dV,
54 amrex::RandomEngine const& engine,
55 bool const isSameSpecies, T_index coll_idx)
56{
57 const T_index NI1 = I1e - I1s;
58 const T_index NI2 = I2e - I2s;
59 const T_index max_N = amrex::max(NI1,NI2);
60 const T_index min_N = amrex::min(NI1,NI2);
61
62 T_PR * const AMREX_RESTRICT w1 = soa_1.m_rdata[PIdx::w];
63 T_PR * const AMREX_RESTRICT u1x = soa_1.m_rdata[PIdx::ux];
64 T_PR * const AMREX_RESTRICT u1y = soa_1.m_rdata[PIdx::uy];
65 T_PR * const AMREX_RESTRICT u1z = soa_1.m_rdata[PIdx::uz];
66
67 T_PR * const AMREX_RESTRICT w2 = soa_2.m_rdata[PIdx::w];
68 T_PR * const AMREX_RESTRICT u2x = soa_2.m_rdata[PIdx::ux];
69 T_PR * const AMREX_RESTRICT u2y = soa_2.m_rdata[PIdx::uy];
70 T_PR * const AMREX_RESTRICT u2z = soa_2.m_rdata[PIdx::uz];
71
72 // compute Debye length lmdD (if not using a fixed L = Coulomb log)
73 T_PR lmdD = T_PR(-1.0);
74 if ( L <= T_PR(0.0) ) {
75 if (global_lamdb <= T_PR(0.0)) {
76 lmdD = T_PR(1.0)/std::sqrt( n1*q1*q1/(T1*PhysConst::epsilon_0) +
77 n2*q2*q2/(T2*PhysConst::epsilon_0) );
78 } else {
79 lmdD = global_lamdb;
80 }
81 }
82
83 // compute atomic spacing
84 const T_PR maxn = amrex::max(n1,n2);
85 const auto rmin = static_cast<T_PR>( 1.0/std::cbrt(4.0*MathConst::pi/3.0*maxn) );
86
87 // bmax (screening length) cannot be smaller than atomic spacing
88 const T_PR bmax = amrex::max(lmdD, rmin);
89
90 // set max cross section based on mfp = atomic spacing
91 const T_PR sigma_max = T_PR(1.0) / (maxn * rmin);
92
93 // call UpdateMomentumPerezElastic()
94 {
95 T_index i1 = I1s + coll_idx;
96 T_index i2 = I2s + coll_idx;
97
98 // we will start from collision number = coll_idx and then add
99 // stride (smaller set size) until we do all collisions (larger set size)
100 for (T_index k = coll_idx; k < max_N; k += min_N)
101 {
102
103 // Compute the effective density n12 used to compute the normalized
104 // transport mean free path (s12) in UpdateMomentumPerezElastic().
105 // s12 is defined such that the expected value of the change in particle
106 // velocity is equal to that from the full NxN pairing method for each particle, as described in
107 // Justin Ray Angus, Yichen Fu, Vasily Geyko, Dave Grote, David Larson, JCP 531, 113927 (2025).
108 // This method is a direct extension of the original binary-pairing method for equally weighted
109 // particles by Takizuka and Abe JCP 25 (1977).
110 T_PR n12;
111 const T_PR wpmax = amrex::max(w1[ I1[i1] ],w2[ I2[i2] ]);
112 if (isSameSpecies) { n12 = wpmax*static_cast<T_PR>(min_N+max_N-1)/dV; }
113 else { n12 = wpmax*static_cast<T_PR>(min_N)/dV; }
114
116 u1x[ I1[i1] ], u1y[ I1[i1] ], u1z[ I1[i1] ],
117 u2x[ I2[i2] ], u2y[ I2[i2] ], u2z[ I2[i2] ],
118 q1, m1, w1[ I1[i1] ], q2, m2, w2[ I2[i2] ],
119 n12, sigma_max, L, bmax, dt, engine );
120
121 if (max_N == NI1) {
122 i1 += min_N;
123 }
124 if (max_N == NI2) {
125 i2 += min_N;
126 }
127 }
128 }
129
130}
131
132#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_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:43
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
constexpr auto epsilon_0
vacuum permittivity: dielectric permittivity of vacuum [F/m]
Definition constant.H:156
constexpr auto pi
ratio of a circle's circumference to its diameter
Definition constant.H:29
__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
@ uz
Definition WarpXParticleContainer.H:70
@ w
Definition WarpXParticleContainer.H:70
@ uy
Definition WarpXParticleContainer.H:70
@ ux
Definition WarpXParticleContainer.H:70