WarpX
ChargeDeposition.H
Go to the documentation of this file.
1 /* Copyright 2019 Axel Huebl, Andrew Myers, David Grote, Maxence Thevenet
2  * Weiqun Zhang
3  *
4  * This file is part of WarpX.
5  *
6  * License: BSD-3-Clause-LBNL
7  */
8 #ifndef CHARGEDEPOSITION_H_
9 #define CHARGEDEPOSITION_H_
10 
13 #include "Particles/ShapeFactors.H"
16 #ifdef WARPX_DIM_RZ
17 # include "Utils/WarpX_Complex.H"
18 #endif
19 
20 #include <AMReX.H>
21 
22 /* \brief Charge Deposition for thread thread_num
23  * \param GetPosition : A functor for returning the particle position.
24  * \param wp : Pointer to array of particle weights.
25  * \param ion_lev : Pointer to array of particle ionization level. This is
26  required to have the charge of each macroparticle
27  since q is a scalar. For non-ionizable species,
28  ion_lev is a null pointer.
29  * \param rho_fab : FArrayBox of charge density, either full array or tile.
30  * \param np_to_depose : Number of particles for which current is deposited.
31  * \param dx : 3D cell size
32  * \param xyzmin : Physical lower bounds of domain.
33  * \param lo : Index lower bounds of domain.
34  * \param q : species charge.
35  * \param n_rz_azimuthal_modes: Number of azimuthal modes when using RZ geometry.
36  * \param cost: Pointer to (load balancing) cost corresponding to box where present particles deposit current.
37  * \param load_balance_costs_update_algo: Selected method for updating load balance costs.
38  */
39 template <int depos_order>
41  const amrex::ParticleReal * const wp,
42  const int * const ion_lev,
43  amrex::FArrayBox& rho_fab,
44  const long np_to_depose,
45  const std::array<amrex::Real,3>& dx,
46  const std::array<amrex::Real, 3> xyzmin,
47  const amrex::Dim3 lo,
48  const amrex::Real q,
49  const int n_rz_azimuthal_modes,
50  amrex::Real* cost,
51  const long load_balance_costs_update_algo)
52 {
53  using namespace amrex;
54 
55 #if !defined(AMREX_USE_GPU)
56  amrex::ignore_unused(cost, load_balance_costs_update_algo);
57 #endif
58 
59  // Whether ion_lev is a null pointer (do_ionization=0) or a real pointer
60  // (do_ionization=1)
61  const bool do_ionization = ion_lev;
62  const amrex::Real dzi = 1.0_rt/dx[2];
63 #if defined(WARPX_DIM_1D_Z)
64  const amrex::Real invvol = dzi;
65 #endif
66 #if defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
67  const amrex::Real dxi = 1.0_rt/dx[0];
68  const amrex::Real invvol = dxi*dzi;
69 #elif defined(WARPX_DIM_3D)
70  const amrex::Real dxi = 1.0_rt/dx[0];
71  const amrex::Real dyi = 1.0_rt/dx[1];
72  const amrex::Real invvol = dxi*dyi*dzi;
73 #endif
74 
75 #if (AMREX_SPACEDIM >= 2)
76  const amrex::Real xmin = xyzmin[0];
77 #endif
78 #if defined(WARPX_DIM_3D)
79  const amrex::Real ymin = xyzmin[1];
80 #endif
81  const amrex::Real zmin = xyzmin[2];
82 
83  amrex::Array4<amrex::Real> const& rho_arr = rho_fab.array();
84  amrex::IntVect const rho_type = rho_fab.box().type();
85 
86  constexpr int NODE = amrex::IndexType::NODE;
87  constexpr int CELL = amrex::IndexType::CELL;
88 
89  // Loop over particles and deposit into rho_fab
90 #if defined(WARPX_USE_GPUCLOCK)
91  amrex::Real* cost_real = nullptr;
92  if( load_balance_costs_update_algo == LoadBalanceCostsUpdateAlgo::GpuClock) {
93  cost_real = (amrex::Real *) amrex::The_Managed_Arena()->alloc(sizeof(amrex::Real));
94  *cost_real = 0.;
95  }
96 #endif
98  np_to_depose,
99  [=] AMREX_GPU_DEVICE (long ip) {
100 #if defined(WARPX_USE_GPUCLOCK)
101  KernelTimer kernelTimer(cost && load_balance_costs_update_algo
103 #endif
104  // --- Get particle quantities
105  amrex::Real wq = q*wp[ip]*invvol;
106  if (do_ionization){
107  wq *= ion_lev[ip];
108  }
109 
110  amrex::ParticleReal xp, yp, zp;
111  GetPosition(ip, xp, yp, zp);
112 
113  // --- Compute shape factors
114  Compute_shape_factor< depos_order > const compute_shape_factor;
115 #if (AMREX_SPACEDIM >= 2)
116  // x direction
117  // Get particle position in grid coordinates
118 #if defined(WARPX_DIM_RZ)
119  const amrex::Real rp = std::sqrt(xp*xp + yp*yp);
120  amrex::Real costheta;
121  amrex::Real sintheta;
122  if (rp > 0.) {
123  costheta = xp/rp;
124  sintheta = yp/rp;
125  } else {
126  costheta = 1._rt;
127  sintheta = 0._rt;
128  }
129  const Complex xy0 = Complex{costheta, sintheta};
130  const amrex::Real x = (rp - xmin)*dxi;
131 #else
132  const amrex::Real x = (xp - xmin)*dxi;
133 #endif
134 
135  // Compute shape factor along x
136  // i: leftmost grid point that the particle touches
137  amrex::Real sx[depos_order + 1] = {0._rt};
138  int i = 0;
139  if (rho_type[0] == NODE) {
140  i = compute_shape_factor(sx, x);
141  } else if (rho_type[0] == CELL) {
142  i = compute_shape_factor(sx, x - 0.5_rt);
143  }
144 #endif //AMREX_SPACEDIM >= 2
145 #if defined(WARPX_DIM_3D)
146  // y direction
147  const amrex::Real y = (yp - ymin)*dyi;
148  amrex::Real sy[depos_order + 1] = {0._rt};
149  int j = 0;
150  if (rho_type[1] == NODE) {
151  j = compute_shape_factor(sy, y);
152  } else if (rho_type[1] == CELL) {
153  j = compute_shape_factor(sy, y - 0.5_rt);
154  }
155 #endif
156  // z direction
157  const amrex::Real z = (zp - zmin)*dzi;
158  amrex::Real sz[depos_order + 1] = {0._rt};
159  int k = 0;
160  if (rho_type[WARPX_ZINDEX] == NODE) {
161  k = compute_shape_factor(sz, z);
162  } else if (rho_type[WARPX_ZINDEX] == CELL) {
163  k = compute_shape_factor(sz, z - 0.5_rt);
164  }
165 
166  // Deposit charge into rho_arr
167 #if defined(WARPX_DIM_1D_Z)
168  for (int iz=0; iz<=depos_order; iz++){
170  &rho_arr(lo.x+k+iz, 0, 0, 0),
171  sz[iz]*wq);
172  }
173 #endif
174 #if defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
175  for (int iz=0; iz<=depos_order; iz++){
176  for (int ix=0; ix<=depos_order; ix++){
178  &rho_arr(lo.x+i+ix, lo.y+k+iz, 0, 0),
179  sx[ix]*sz[iz]*wq);
180 #if defined(WARPX_DIM_RZ)
181  Complex xy = xy0; // Throughout the following loop, xy takes the value e^{i m theta}
182  for (int imode=1 ; imode < n_rz_azimuthal_modes ; imode++) {
183  // The factor 2 on the weighting comes from the normalization of the modes
184  amrex::Gpu::Atomic::AddNoRet( &rho_arr(lo.x+i+ix, lo.y+k+iz, 0, 2*imode-1), 2._rt*sx[ix]*sz[iz]*wq*xy.real());
185  amrex::Gpu::Atomic::AddNoRet( &rho_arr(lo.x+i+ix, lo.y+k+iz, 0, 2*imode ), 2._rt*sx[ix]*sz[iz]*wq*xy.imag());
186  xy = xy*xy0;
187  }
188 #endif
189  }
190  }
191 #elif defined(WARPX_DIM_3D)
192  for (int iz=0; iz<=depos_order; iz++){
193  for (int iy=0; iy<=depos_order; iy++){
194  for (int ix=0; ix<=depos_order; ix++){
196  &rho_arr(lo.x+i+ix, lo.y+j+iy, lo.z+k+iz),
197  sx[ix]*sy[iy]*sz[iz]*wq);
198  }
199  }
200  }
201 #endif
202  }
203  );
204 #if defined(WARPX_USE_GPUCLOCK)
205  if(cost && load_balance_costs_update_algo == LoadBalanceCostsUpdateAlgo::GpuClock) {
207  *cost += *cost_real;
208  amrex::The_Managed_Arena()->free(cost_real);
209  }
210 #endif
211 
212 #ifndef WARPX_DIM_RZ
213  amrex::ignore_unused(n_rz_azimuthal_modes);
214 #endif
215 }
216 
217 #endif // CHARGEDEPOSITION_H_
const Box & box() const noexcept
AMREX_FORCE_INLINE Array4< Real const > array() const noexcept
void doChargeDepositionShapeN(const GetParticlePosition &GetPosition, const amrex::ParticleReal *const wp, const int *const ion_lev, amrex::FArrayBox &rho_fab, const long np_to_depose, const std::array< amrex::Real, 3 > &dx, const std::array< amrex::Real, 3 > xyzmin, const amrex::Dim3 lo, const amrex::Real q, const int n_rz_azimuthal_modes, amrex::Real *cost, const long load_balance_costs_update_algo)
Definition: ChargeDeposition.H:40
int dx
Definition: compute_domain.py:35
virtual void free(void *pt)=0
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void ignore_unused(const Ts &...)
Definition: ShapeFactors.H:26
AMREX_GPU_HOST_DEVICE constexpr AMREX_FORCE_INLINE T imag() const noexcept
std::enable_if_t< std::is_integral< T >::value > ParallelFor(TypeList< CTOs... >, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
#define AMREX_GPU_DEVICE
Definition: WarpXAlgorithmSelection.H:124
i
Definition: check_interp_points_and_weights.py:174
virtual void * alloc(std::size_t sz)=0
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void AddNoRet(T *sum, T value) noexcept
AMREX_GPU_HOST_DEVICE constexpr AMREX_FORCE_INLINE T real() const noexcept
Functor that can be used to extract the positions of the macroparticles inside a ParallelFor kernel...
Definition: GetAndSetPosition.H:52
void streamSynchronize() noexcept
Defines a timer object to be used on GPU; measures summed thread cycles.
Definition: KernelTimer.H:23
AMREX_GPU_HOST_DEVICE IntVect type() const noexcept
Arena * The_Managed_Arena()
NODE