7#ifndef WARPX_SHAREDDEPOSITIONUTILS_H_
8#define WARPX_SHAREDDEPOSITIONUTILS_H_
28 const int nTiles = nCells / tilesize;
29 const int remainder = nCells % tilesize;
30 maxTilesize = tilesize + int(std::ceil((amrex::Real) remainder / nTiles));
40#if defined(AMREX_USE_HIP) || defined(AMREX_USE_CUDA)
46 using namespace amrex::literals;
50 for (
int icell = threadIdx.x; icell < bx.numPts(); icell += blockDim.x)
52 int k = icell / (len.x*len.y);
53 int j = (icell - k*(len.x*len.y)) / len.x;
54 int i = (icell - k*(len.x*len.y)) - j*len.x;
58 if (amrex::Math::abs(local(i, j, k)) > 0.0_rt) {
65#if (defined(AMREX_USE_HIP) || defined(AMREX_USE_CUDA)) && \
66 !(defined(WARPX_DIM_RCYLINDER) || defined(WARPX_DIM_RSPHERE))
67template <
int depos_order>
70 const amrex::ParticleReal *
const wp,
71 const amrex::ParticleReal *
const uxp,
72 const amrex::ParticleReal *
const uyp,
73 const amrex::ParticleReal *
const uzp,
77 const amrex::Real relative_time,
82 const int n_rz_azimuthal_modes,
83 const unsigned int ip,
84 const int zdir,
const int NODE,
const int CELL,
const int dir)
86 using namespace amrex::literals;
88#if !defined(WARPX_DIM_RZ)
94 const bool do_ionization = ion_lev;
96 const amrex::Real invvol = dinv.
x*dinv.
y*dinv.
z;
101 const amrex::Real gaminv = 1.0_rt/std::sqrt(1.0_rt + uxp[ip]*uxp[ip]*clightsq
102 + uyp[ip]*uyp[ip]*clightsq
103 + uzp[ip]*uzp[ip]*clightsq);
104 amrex::Real wq = q*wp[ip];
109 amrex::ParticleReal xp, yp, zp;
110 GetPosition(ip, xp, yp, zp);
112 const amrex::Real vx = uxp[ip]*gaminv;
113 const amrex::Real vy = uyp[ip]*gaminv;
114 const amrex::Real
vz = uzp[ip]*gaminv;
116#if defined(WARPX_DIM_RZ)
119 const amrex::Real xpmid = xp + relative_time*vx;
120 const amrex::Real ypmid = yp + relative_time*vy;
121 const amrex::Real rpmid = std::sqrt(xpmid*xpmid + ypmid*ypmid);
122 amrex::Real costheta;
123 amrex::Real sintheta;
125 costheta = xpmid/rpmid;
126 sintheta = ypmid/rpmid;
132 const amrex::Real wqx = wq*invvol*(+vx*costheta + vy*sintheta);
133 const amrex::Real wqy = wq*invvol*(-vx*sintheta + vy*costheta);
135 const amrex::Real wqx = wq*invvol*vx;
136 const amrex::Real wqy = wq*invvol*vy;
138 const amrex::Real wqz = wq*invvol*
vz;
140 amrex::Real pcurrent = 0.0;
143 }
else if (dir == 1) {
145 }
else if (dir == 2) {
151#if defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ) || defined(WARPX_DIM_3D)
155#if defined(WARPX_DIM_RZ)
157 const double xmid = (rpmid - xyzmin.
x)*dinv.
x;
159 const double xmid = ((xp - xyzmin.
x) + relative_time*vx)*dinv.
x;
166 double sx_node[depos_order + 1] = {0.};
167 double sx_cell[depos_order + 1] = {0.};
170 if (j_type[0] ==
NODE) {
171 j_node = compute_shape_factor(sx_node, xmid);
173 if (j_type[0] == CELL) {
174 j_cell = compute_shape_factor(sx_cell, xmid - 0.5);
177 amrex::Real sx_j[depos_order + 1] = {0._rt};
178 for (
int ix=0;
ix<=depos_order;
ix++)
180 sx_j[
ix] = ((j_type[0] ==
NODE) ? amrex::Real(sx_node[ix]) :
amrex::Real(sx_cell[
ix]));
183 int const j_j = ((j_type[0] ==
NODE) ? j_node : j_cell);
186#if defined(WARPX_DIM_3D)
189 const double ymid = ((yp - xyzmin.
y) + relative_time*vy)*dinv.
y;
190 double sy_node[depos_order + 1] = {0.};
191 double sy_cell[depos_order + 1] = {0.};
194 if (j_type[1] ==
NODE) {
195 k_node = compute_shape_factor(sy_node, ymid);
197 if (j_type[1] == CELL) {
198 k_cell = compute_shape_factor(sy_cell, ymid - 0.5);
200 amrex::Real sy_j[depos_order + 1] = {0._rt};
201 for (
int iy=0;
iy<=depos_order;
iy++)
203 sy_j[
iy] = ((j_type[1] ==
NODE) ? amrex::Real(sy_node[iy]) :
amrex::Real(sy_cell[
iy]));
205 int const k_j = ((j_type[1] ==
NODE) ? k_node : k_cell);
210 const double zmid = ((zp - xyzmin.
z) + relative_time*
vz)*dinv.
z;
211 double sz_node[depos_order + 1] = {0.};
212 double sz_cell[depos_order + 1] = {0.};
215 if (j_type[zdir] ==
NODE) {
216 l_node = compute_shape_factor(sz_node, zmid);
218 if (j_type[zdir] == CELL) {
219 l_cell = compute_shape_factor(sz_cell, zmid - 0.5);
221 amrex::Real sz_j[depos_order + 1] = {0._rt};
222 for (
int iz=0;
iz<=depos_order;
iz++)
224 sz_j[
iz] = ((j_type[zdir] ==
NODE) ? amrex::Real(sz_node[iz]) :
amrex::Real(sz_cell[
iz]));
226 int const l_j = ((j_type[zdir] ==
NODE) ? l_node : l_cell);
229#if defined(WARPX_DIM_1D_Z)
230 for (
int iz=0;
iz<=depos_order;
iz++){
232 &j_buff(lo.
x+l_j+iz, 0, 0, 0),
236#if defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
237 for (
int iz=0;
iz<=depos_order;
iz++){
238 for (
int ix=0;
ix<=depos_order;
ix++){
240 &j_buff(lo.
x+j_j+ix, lo.
y+l_j+iz, 0, 0),
241 sx_j[ix]*sz_j[iz]*pcurrent);
242#if defined(WARPX_DIM_RZ)
244 for (
int imode=1 ; imode < n_rz_azimuthal_modes ; imode++) {
253#elif defined(WARPX_DIM_3D)
254 for (
int iz=0;
iz<=depos_order;
iz++){
255 for (
int iy=0;
iy<=depos_order;
iy++){
256 for (
int ix=0;
ix<=depos_order;
ix++){
258 &j_buff(lo.
x+j_j+ix, lo.
y+k_j+iy, lo.
z+l_j+iz),
259 sx_j[ix]*sy_j[iy]*sz_j[iz]*pcurrent);
#define AMREX_FORCE_INLINE
#define AMREX_GPU_HOST_DEVICE
@ vz
Definition RigidInjectedParticleContainer.H:27
AMREX_FORCE_INLINE int getMaxTboxAlongDim(int nCells, int tilesize)
Definition SharedDepositionUtils.H:26
amrex::GpuComplex< amrex::Real > Complex
Definition WarpX_Complex.H:22
constexpr auto c
vacuum speed of light [m/s]
Definition constant.H:149
__host__ __device__ AMREX_FORCE_INLINE void AddNoRet(T *sum, T value) noexcept
__host__ __device__ void ignore_unused(const Ts &...)
__host__ __device__ Dim3 length(Array4< T > const &a) noexcept
__host__ __device__ Dim3 lbound(Array4< T > const &a) noexcept
Definition ShapeFactors.H:29
Functor that can be used to extract the positions of the macroparticles inside a ParallelFor kernel.
Definition GetAndSetPosition.H:75
__host__ __device__ constexpr T real() const noexcept
__host__ __device__ constexpr T imag() const noexcept