7 #ifndef SHAREDDEPOSITIONUTILS_H_
8 #define SHAREDDEPOSITIONUTILS_H_
29 int nTiles = nCells / tilesize;
30 int remainder = nCells % tilesize;
31 maxTilesize = tilesize + int(std::ceil((amrex::Real) remainder / nTiles));
41 #if defined(AMREX_USE_HIP) || defined(AMREX_USE_CUDA)
49 for (
int icell = threadIdx.x; icell < bx.numPts(); icell += blockDim.x)
51 int k = icell / (len.x*len.y);
52 int j = (icell - k*(len.x*len.y)) / len.x;
53 int i = (icell - k*(len.x*len.y)) - j*len.x;
57 if (amrex::Math::abs(local(
i, j, k)) > 0.0_rt) {
64 #if defined(AMREX_USE_HIP) || defined(AMREX_USE_CUDA)
65 template <
int depos_order>
68 const amrex::ParticleReal *
const wp,
69 const amrex::ParticleReal *
const uxp,
70 const amrex::ParticleReal *
const uyp,
71 const amrex::ParticleReal *
const uzp,
72 const int *
const ion_lev,
75 const amrex::Real relative_time,
76 const std::array<amrex::Real,3>&
dx,
77 const std::array<amrex::Real,3>& xyzmin,
80 const int n_rz_azimuthal_modes,
81 const unsigned int ip,
82 const int zdir,
const int NODE,
const int CELL,
const int dir)
84 #if !defined(WARPX_DIM_RZ)
90 const bool do_ionization = ion_lev;
91 const amrex::Real dzi = 1.0_rt/
dx[2];
92 #if defined(WARPX_DIM_1D_Z)
93 const amrex::Real invvol = dzi;
95 #if defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
96 const amrex::Real dxi = 1.0_rt/
dx[0];
97 const amrex::Real invvol = dxi*dzi;
98 #elif defined(WARPX_DIM_3D)
99 const amrex::Real dxi = 1.0_rt/
dx[0];
100 const amrex::Real dyi = 1.0_rt/
dx[1];
101 const amrex::Real invvol = dxi*dyi*dzi;
104 #if (AMREX_SPACEDIM >= 2)
105 const amrex::Real xmin = xyzmin[0];
107 #if defined(WARPX_DIM_3D)
108 const amrex::Real ymin = xyzmin[1];
110 const amrex::Real zmin = xyzmin[2];
115 const amrex::Real gaminv = 1.0_rt/std::sqrt(1.0_rt + uxp[ip]*uxp[ip]*clightsq
116 + uyp[ip]*uyp[ip]*clightsq
117 + uzp[ip]*uzp[ip]*clightsq);
118 amrex::Real wq = q*wp[ip];
123 amrex::ParticleReal xp, yp, zp;
124 GetPosition(ip, xp, yp, zp);
126 const amrex::Real vx = uxp[ip]*gaminv;
127 const amrex::Real vy = uyp[ip]*gaminv;
128 const amrex::Real vz = uzp[ip]*gaminv;
130 #if defined(WARPX_DIM_RZ)
133 const amrex::Real xpmid = xp + relative_time*vx;
134 const amrex::Real ypmid = yp + relative_time*vy;
135 const amrex::Real rpmid = std::sqrt(xpmid*xpmid + ypmid*ypmid);
136 amrex::Real costheta;
137 amrex::Real sintheta;
139 costheta = xpmid/rpmid;
140 sintheta = ypmid/rpmid;
146 const amrex::Real wqx = wq*invvol*(+vx*costheta + vy*sintheta);
147 const amrex::Real wqy = wq*invvol*(-vx*sintheta + vy*costheta);
149 const amrex::Real wqx = wq*invvol*vx;
150 const amrex::Real wqy = wq*invvol*vy;
152 const amrex::Real wqz = wq*invvol*vz;
154 amrex::Real pcurrent = 0.0;
157 }
else if (dir == 1) {
159 }
else if (dir == 2) {
165 #if defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ) || defined(WARPX_DIM_3D)
169 #if defined(WARPX_DIM_RZ)
171 const double xmid = (rpmid - xmin)*dxi;
173 const double xmid = ((xp - xmin) + relative_time*vx)*dxi;
180 double sx_node[depos_order + 1] = {0.};
181 double sx_cell[depos_order + 1] = {0.};
184 if (j_type[0] == NODE) {
185 j_node = compute_shape_factor(sx_node, xmid);
187 if (j_type[0] == CELL) {
188 j_cell = compute_shape_factor(sx_cell, xmid - 0.5);
191 amrex::Real sx_j[depos_order + 1] = {0._rt};
192 for (
int ix=0;
ix<=depos_order;
ix++)
194 sx_j[
ix] = ((j_type[0] ==
NODE) ? amrex::Real(sx_node[ix]) : amrex::Real(sx_cell[ix]));
197 int const j_j = ((j_type[0] ==
NODE) ? j_node : j_cell);
200 #if defined(WARPX_DIM_3D)
203 const double ymid = ((yp - ymin) + relative_time*vy)*dyi;
204 double sy_node[depos_order + 1] = {0.};
205 double sy_cell[depos_order + 1] = {0.};
208 if (j_type[1] == NODE) {
209 k_node = compute_shape_factor(sy_node, ymid);
211 if (j_type[1] == CELL) {
212 k_cell = compute_shape_factor(sy_cell, ymid - 0.5);
214 amrex::Real sy_j[depos_order + 1] = {0._rt};
215 for (
int iy=0;
iy<=depos_order;
iy++)
217 sy_j[
iy] = ((j_type[1] ==
NODE) ? amrex::Real(sy_node[iy]) : amrex::Real(sy_cell[iy]));
219 int const k_j = ((j_type[1] ==
NODE) ? k_node : k_cell);
224 const double zmid = ((zp - zmin) + relative_time*vz)*dzi;
225 double sz_node[depos_order + 1] = {0.};
226 double sz_cell[depos_order + 1] = {0.};
229 if (j_type[zdir] == NODE) {
230 l_node = compute_shape_factor(sz_node, zmid);
232 if (j_type[zdir] == CELL) {
233 l_cell = compute_shape_factor(sz_cell, zmid - 0.5);
235 amrex::Real sz_j[depos_order + 1] = {0._rt};
236 for (
int iz=0;
iz<=depos_order;
iz++)
238 sz_j[
iz] = ((j_type[zdir] ==
NODE) ? amrex::Real(sz_node[iz]) : amrex::Real(sz_cell[iz]));
240 int const l_j = ((j_type[zdir] ==
NODE) ? l_node : l_cell);
243 #if defined(WARPX_DIM_1D_Z)
244 for (
int iz=0;
iz<=depos_order;
iz++){
246 &j_buff(lo.
x+l_j+iz, 0, 0, 0),
250 #if defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
251 for (
int iz=0;
iz<=depos_order;
iz++){
252 for (
int ix=0;
ix<=depos_order;
ix++){
254 &j_buff(lo.
x+j_j+ix, lo.
y+l_j+iz, 0, 0),
255 sx_j[ix]*sz_j[iz]*pcurrent);
256 #if defined(WARPX_DIM_RZ)
258 for (
int imode=1 ; imode < n_rz_azimuthal_modes ; imode++) {
267 #elif defined(WARPX_DIM_3D)
268 for (
int iz=0;
iz<=depos_order;
iz++){
269 for (
int iy=0;
iy<=depos_order;
iy++){
270 for (
int ix=0;
ix<=depos_order;
ix++){
272 &j_buff(lo.
x+j_j+ix, lo.
y+k_j+iy, lo.
z+l_j+iz),
273 sx_j[ix]*sy_j[iy]*sz_j[iz]*pcurrent);
#define AMREX_FORCE_INLINE
#define AMREX_GPU_HOST_DEVICE
AMREX_FORCE_INLINE int getMaxTboxAlongDim(int nCells, int tilesize)
Definition: SharedDepositionUtils.H:27
static constexpr auto c
vacuum speed of light [m/s]
Definition: constant.H:44
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void AddNoRet(T *sum, T value) noexcept
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Dim3 lbound(Array4< T > const &a) noexcept
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void ignore_unused(const Ts &...)
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Dim3 length(Array4< T > const &a) noexcept
i
Definition: check_interp_points_and_weights.py:174
int dx
Definition: stencil.py:436
Definition: ShapeFactors.H:27
Functor that can be used to extract the positions of the macroparticles inside a ParallelFor kernel.
Definition: GetAndSetPosition.H:53
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE constexpr T real() const noexcept
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE constexpr T imag() const noexcept