WarpX
Loading...
Searching...
No Matches
TemperatureDeposition.H
Go to the documentation of this file.
1/* Copyright 2025 The WarpX Community
2 *
3 * This file is part of WarpX.
4 *
5 * Authors: S. Eric Clark (Helion Energy)
6 *
7 * License: BSD-3-Clause-LBNL
8 */
9
10#ifndef WARPX_TEMPERATUREDEPOSITION_H_
11#define WARPX_TEMPERATUREDEPOSITION_H_
12
17#include "Utils/TextMsg.H"
19#include "Utils/WarpXConst.H"
20#ifdef WARPX_DIM_RZ
21# include "Utils/WarpX_Complex.H"
22#endif
23
25
26#include <AMReX.H>
27#include <AMReX_Arena.H>
28#include <AMReX_Array4.H>
29#include <AMReX_Dim3.H>
30#include <AMReX_REAL.H>
31
33
39 const amrex::IntVectND<3>& ixv,
40 const amrex::IntVectND<3>& iyv,
41 const amrex::IntVectND<3>& izv,
42 amrex::Real wpx_var,
43 amrex::Real wpy_var,
44 amrex::Real wpz_var,
45 const amrex::Array4<int> & nx_arr,
46 const amrex::Array4<int> & ny_arr,
47 const amrex::Array4<int> & nz_arr,
48 const amrex::Array4<amrex::Real> & wx_arr,
49 const amrex::Array4<amrex::Real> & wy_arr,
50 const amrex::Array4<amrex::Real> & wz_arr,
51 const amrex::Array4<amrex::Real> & w2x_arr,
52 const amrex::Array4<amrex::Real> & w2y_arr,
53 const amrex::Array4<amrex::Real> & w2z_arr,
54 const amrex::Array4<amrex::Real> & vxbar_arr,
55 const amrex::Array4<amrex::Real> & vybar_arr,
56 const amrex::Array4<amrex::Real> & vzbar_arr,
59)
60{
61 const auto vxr = static_cast<amrex::Real>(vx);
62 const auto vyr = static_cast<amrex::Real>(vy);
63 const auto vzr = static_cast<amrex::Real>(vz);
64
66 {
67 // Update sample count
69 &nx_arr(ixv[0], ixv[1], ixv[2], 0), 1);
71 &ny_arr(iyv[0], iyv[1], iyv[2], 0), 1);
73 &nz_arr(izv[0], izv[1], izv[2], 0), 1);
74
75 // Update w_sum
77 &wx_arr(ixv[0], ixv[1], ixv[2], 0), wpx_var);
79 &wy_arr(iyv[0], iyv[1], iyv[2], 0), wpy_var);
81 &wz_arr(izv[0], izv[1], izv[2], 0), wpz_var);
82
83 // Update wv_sum to be normalized to vbar
85 &vxbar_arr(ixv[0], ixv[1], ixv[2], 0), wpx_var*vxr);
87 &vybar_arr(iyv[0], iyv[1], iyv[2], 0), wpy_var*vyr);
89 &vzbar_arr(izv[0], izv[1], izv[2], 0), wpz_var*vzr);
90 }
91
93 {
94 // For a single pass run through, the vbar arrays will be zeroed, so this should work
95 // in either case. Will be non-demeaned in 1 pass
96 amrex::Real vxb = 0.;
97 amrex::Real vyb = 0.;
98 amrex::Real vzb = 0.;
99
100 // Normalization of accumulation arrays are not done until after all deposition, so normalize
101 // here for double pass algorithm
103 {
104 if (nx_arr(ixv[0], ixv[1], ixv[2], 0) > 0) {
105 vxb = vxbar_arr(ixv[0], ixv[1], ixv[2], 0)/wx_arr(ixv[0], ixv[1], ixv[2], 0);
106 }
107 if (ny_arr(iyv[0], iyv[1], iyv[2], 0) > 0) {
108 vyb = vybar_arr(iyv[0], iyv[1], iyv[2], 0)/wy_arr(iyv[0], iyv[1], iyv[2], 0);
109 }
110 if (nz_arr(izv[0], izv[1], izv[2], 0) > 0) {
111 vzb = vzbar_arr(izv[0], izv[1], izv[2], 0)/wz_arr(izv[0], izv[1], izv[2], 0);
112 }
113
114 }
115
116 const amrex::Real vxd = vxr - vxb;
117 const amrex::Real vyd = vyr - vyb;
118 const amrex::Real vzd = vzr - vzb;
119
120 // Update wv2_sum
122 &w2x_arr(ixv[0], ixv[1], ixv[2], 0), wpx_var*vxd*vxd);
124 &w2y_arr(iyv[0], iyv[1], iyv[2], 0), wpy_var*vyd*vyd);
126 &w2z_arr(izv[0], izv[1], izv[2], 0), wpz_var*vzd*vzd);
127 }
128}
129
152template <int depos_order>
155 [[maybe_unused]] const amrex::ParticleReal xp,
156 [[maybe_unused]] const amrex::ParticleReal yp,
157 [[maybe_unused]] const amrex::ParticleReal zp,
158 const amrex::ParticleReal wp,
159 const amrex::ParticleReal vx,
160 const amrex::ParticleReal vy,
162 const amrex::Array4<int> & nx_arr,
163 const amrex::Array4<int> & ny_arr,
164 const amrex::Array4<int> & nz_arr,
165 const amrex::Array4<amrex::Real> & wx_arr,
166 const amrex::Array4<amrex::Real> & wy_arr,
167 const amrex::Array4<amrex::Real> & wz_arr,
168 const amrex::Array4<amrex::Real> & w2x_arr,
169 const amrex::Array4<amrex::Real> & w2y_arr,
170 const amrex::Array4<amrex::Real> & w2z_arr,
171 const amrex::Array4<amrex::Real> & vxbar_arr,
172 const amrex::Array4<amrex::Real> & vybar_arr,
173 const amrex::Array4<amrex::Real> & vzbar_arr,
174 const amrex::IntVect & varx_type,
175 const amrex::IntVect & vary_type,
176 const amrex::IntVect & varz_type,
177 const TemperatureDepositionType type,
178 const TemperatureDepositionPass pass,
179 const amrex::Real relative_time,
180 const amrex::XDim3 & dinv,
181 const amrex::XDim3 & xyzmin,
182 const amrex::Dim3 lo,
183 [[maybe_unused]] const int n_rz_azimuthal_modes)
184{
185 using namespace amrex::literals;
186
187 constexpr int NODE = amrex::IndexType::NODE;
188 constexpr int CELL = amrex::IndexType::CELL;
189
190 // --- Compute shape factors
191 Compute_shape_factor< depos_order > const compute_shape_factor;
192#if !defined(WARPX_DIM_1D_Z)
193 // x direction
194 // Get particle position at dt=relative time from particle time
195 // Keep these double to avoid bug in single precision
196
197#if defined(WARPX_DIM_RZ) || defined(WARPX_DIM_RCYLINDER)
198 // In RZ, wpx is actually wpr, and wpy is wptheta
199 // Convert to cylindrical at the mid point
200 const amrex::Real xpmid = xp + relative_time*vx;
201 const amrex::Real ypmid = yp + relative_time*vy;
202 const amrex::Real rpmid = std::sqrt(xpmid*xpmid + ypmid*ypmid);
203
204 const double xmid = (rpmid - xyzmin.x)*dinv.x;
205#elif defined(WARPX_DIM_RSPHERE)
206 // Convert to spherical at the mid point
207 const amrex::Real xpmid = xp + relative_time*vx;
208 const amrex::Real ypmid = yp + relative_time*vy;
209 const amrex::Real zpmid = zp + relative_time*vz;
210 const amrex::Real rpmid = std::sqrt(xpmid*xpmid + ypmid*ypmid + zpmid*zpmid);
211
212 const double xmid = (rpmid - xyzmin.x)*dinv.x;
213#else
214 const double xmid = ((xp - xyzmin.x) + relative_time*vx)*dinv.x;
215#endif
216
217 // j_j[xyz] leftmost grid point in x that the particle touches for the centering of each current
218 // sx_j[xyz] shape factor along x for the centering of each current
219 // There are only two possible centerings, node or cell centered, so at most only two shape factor
220 // arrays will be needed.
221 // Keep these double to avoid bug in single precision
222 double sx_node[depos_order + 1] = {0.};
223 double sx_cell[depos_order + 1] = {0.};
224 int j_node = 0;
225 int j_cell = 0;
226 if (varx_type[0] == NODE || vary_type[0] == NODE || varz_type[0] == NODE) {
227 j_node = compute_shape_factor(sx_node, xmid);
228 }
229 if (varx_type[0] == CELL || vary_type[0] == CELL || varz_type[0] == CELL) {
230 j_cell = compute_shape_factor(sx_cell, xmid - 0.5);
231 }
232
233 amrex::Real sx_jx[depos_order + 1] = {0._rt};
234 amrex::Real sx_jy[depos_order + 1] = {0._rt};
235 amrex::Real sx_jz[depos_order + 1] = {0._rt};
236 for (int ix=0; ix<=depos_order; ix++)
237 {
238 sx_jx[ix] = ((varx_type[0] == NODE) ? amrex::Real(sx_node[ix]) : amrex::Real(sx_cell[ix]));
239 sx_jy[ix] = ((vary_type[0] == NODE) ? amrex::Real(sx_node[ix]) : amrex::Real(sx_cell[ix]));
240 sx_jz[ix] = ((varz_type[0] == NODE) ? amrex::Real(sx_node[ix]) : amrex::Real(sx_cell[ix]));
241 }
242
243 int const j_jx = ((varx_type[0] == NODE) ? j_node : j_cell);
244 int const j_jy = ((vary_type[0] == NODE) ? j_node : j_cell);
245 int const j_jz = ((varz_type[0] == NODE) ? j_node : j_cell);
246#endif
247
248#if defined(WARPX_DIM_3D)
249 // y direction
250 // Keep these double to avoid bug in single precision
251 const double ymid = ((yp - xyzmin.y) + relative_time*vy)*dinv.y;
252 double sy_node[depos_order + 1] = {0.};
253 double sy_cell[depos_order + 1] = {0.};
254 int k_node = 0;
255 int k_cell = 0;
256 if (varx_type[1] == NODE || vary_type[1] == NODE || varz_type[1] == NODE) {
257 k_node = compute_shape_factor(sy_node, ymid);
258 }
259 if (varx_type[1] == CELL || vary_type[1] == CELL || varz_type[1] == CELL) {
260 k_cell = compute_shape_factor(sy_cell, ymid - 0.5);
261 }
262 amrex::Real sy_jx[depos_order + 1] = {0._rt};
263 amrex::Real sy_jy[depos_order + 1] = {0._rt};
264 amrex::Real sy_jz[depos_order + 1] = {0._rt};
265 for (int iy=0; iy<=depos_order; iy++)
266 {
267 sy_jx[iy] = ((varx_type[1] == NODE) ? amrex::Real(sy_node[iy]) : amrex::Real(sy_cell[iy]));
268 sy_jy[iy] = ((vary_type[1] == NODE) ? amrex::Real(sy_node[iy]) : amrex::Real(sy_cell[iy]));
269 sy_jz[iy] = ((varz_type[1] == NODE) ? amrex::Real(sy_node[iy]) : amrex::Real(sy_cell[iy]));
270 }
271 int const k_jx = ((varx_type[1] == NODE) ? k_node : k_cell);
272 int const k_jy = ((vary_type[1] == NODE) ? k_node : k_cell);
273 int const k_jz = ((varz_type[1] == NODE) ? k_node : k_cell);
274#endif
275
276#if !defined(WARPX_DIM_RCYLINDER) && !defined(WARPX_DIM_RSPHERE)
277 // z direction
278 // Keep these double to avoid bug in single precision
279 constexpr int zdir = WARPX_ZINDEX;
280 const double zmid = ((zp - xyzmin.z) + relative_time*vz)*dinv.z;
281 double sz_node[depos_order + 1] = {0.};
282 double sz_cell[depos_order + 1] = {0.};
283 int l_node = 0;
284 int l_cell = 0;
285 if (varx_type[zdir] == NODE || vary_type[zdir] == NODE || varz_type[zdir] == NODE) {
286 l_node = compute_shape_factor(sz_node, zmid);
287 }
288 if (varx_type[zdir] == CELL || vary_type[zdir] == CELL || varz_type[zdir] == CELL) {
289 l_cell = compute_shape_factor(sz_cell, zmid - 0.5);
290 }
291 amrex::Real sz_jx[depos_order + 1] = {0._rt};
292 amrex::Real sz_jy[depos_order + 1] = {0._rt};
293 amrex::Real sz_jz[depos_order + 1] = {0._rt};
294 for (int iz=0; iz<=depos_order; iz++)
295 {
296 sz_jx[iz] = ((varx_type[zdir] == NODE) ? amrex::Real(sz_node[iz]) : amrex::Real(sz_cell[iz]));
297 sz_jy[iz] = ((vary_type[zdir] == NODE) ? amrex::Real(sz_node[iz]) : amrex::Real(sz_cell[iz]));
298 sz_jz[iz] = ((varz_type[zdir] == NODE) ? amrex::Real(sz_node[iz]) : amrex::Real(sz_cell[iz]));
299 }
300 int const l_jx = ((varx_type[zdir] == NODE) ? l_node : l_cell);
301 int const l_jy = ((vary_type[zdir] == NODE) ? l_node : l_cell);
302 int const l_jz = ((varz_type[zdir] == NODE) ? l_node : l_cell);
303#endif
304
305#if defined(WARPX_DIM_1D_Z)
306 for (int iz=0; iz<=depos_order; iz++){
307 const amrex::Real wpx_var = static_cast<amrex::Real>(wp)*sz_jx[iz];
308 const amrex::Real wpy_var = static_cast<amrex::Real>(wp)*sz_jy[iz];
309 const amrex::Real wpz_var = static_cast<amrex::Real>(wp)*sz_jz[iz];
310
311 const amrex::IntVectND<3> ixv{lo.z+l_jx+iz, 0, 0};
312 const amrex::IntVectND<3> iyv{lo.z+l_jy+iz, 0, 0};
313 const amrex::IntVectND<3> izv{lo.z+l_jz+iz, 0, 0};
314
316 vx, vy, vz,
317 ixv, iyv, izv,
318 wpx_var, wpy_var, wpz_var,
319 nx_arr, ny_arr, nz_arr,
320 wx_arr, wy_arr, wz_arr,
321 w2x_arr, w2y_arr, w2z_arr,
322 vxbar_arr, vybar_arr, vzbar_arr,
323 type, pass
324 );
325 }
326#elif defined(WARPX_DIM_RCYLINDER) || defined(WARPX_DIM_RSPHERE)
327 for (int ix=0; ix<=depos_order; ix++){
328 const amrex::Real wpx_var = wp*sx_jx[ix];
329 const amrex::Real wpy_var = wp*sx_jy[ix];
330 const amrex::Real wpz_var = wp*sx_jz[ix];
331
332 amrex::IntVectND<3> ixv{lo.x+j_jx+ix, 0, 0};
333 amrex::IntVectND<3> iyv{lo.x+j_jy+ix, 0, 0};
334 amrex::IntVectND<3> izv{lo.x+j_jz+ix, 0, 0};
335
337 vx, vy, vz,
338 ixv, iyv, izv,
339 wpx_var, wpy_var, wpz_var,
340 nx_arr, ny_arr, nz_arr,
341 wx_arr, wy_arr, wz_arr,
342 w2x_arr, w2y_arr, w2z_arr,
343 vxbar_arr, vybar_arr, vzbar_arr,
344 type, pass
345 );
346 }
347#elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
348 for (int iz=0; iz<=depos_order; iz++){
349 for (int ix=0; ix<=depos_order; ix++){
350 const amrex::Real wpx_var = wp*sx_jx[ix]*sz_jx[iz];
351 const amrex::Real wpy_var = wp*sx_jy[ix]*sz_jy[iz];
352 const amrex::Real wpz_var = wp*sx_jz[ix]*sz_jz[iz];
353
354 const amrex::IntVectND<3> ixv{lo.x+j_jx+ix, lo.z+l_jx+iz, 0};
355 const amrex::IntVectND<3> iyv{lo.x+j_jy+ix, lo.z+l_jy+iz, 0};
356 const amrex::IntVectND<3> izv{lo.x+j_jz+ix, lo.z+l_jz+iz, 0};
357
359 vx, vy, vz,
360 ixv, iyv, izv,
361 wpx_var, wpy_var, wpz_var,
362 nx_arr, ny_arr, nz_arr,
363 wx_arr, wy_arr, wz_arr,
364 w2x_arr, w2y_arr, w2z_arr,
365 vxbar_arr, vybar_arr, vzbar_arr,
366 type, pass
367 );
368 }
369 }
370#elif defined(WARPX_DIM_3D)
371 for (int iz=0; iz<=depos_order; iz++){
372 for (int iy=0; iy<=depos_order; iy++){
373 for (int ix=0; ix<=depos_order; ix++){
374 const amrex::Real wpx_var = wp*sx_jx[ix]*sy_jx[iy]*sz_jx[iz];
375 const amrex::Real wpy_var = wp*sx_jy[ix]*sy_jy[iy]*sz_jy[iz];
376 const amrex::Real wpz_var = wp*sx_jz[ix]*sy_jz[iy]*sz_jz[iz];
377
378 const amrex::IntVectND<3> ixv{lo.x+j_jx+ix, lo.y+k_jx+iy, lo.z+l_jx+iz};
379 const amrex::IntVectND<3> iyv{lo.x+j_jy+ix, lo.y+k_jy+iy, lo.z+l_jy+iz};
380 const amrex::IntVectND<3> izv{lo.x+j_jz+ix, lo.y+k_jz+iy, lo.z+l_jz+iz};
381
383 vx, vy, vz,
384 ixv, iyv, izv,
385 wpx_var, wpy_var, wpz_var,
386 nx_arr, ny_arr, nz_arr,
387 wx_arr, wy_arr, wz_arr,
388 w2x_arr, w2y_arr, w2z_arr,
389 vxbar_arr, vybar_arr, vzbar_arr,
390 type, pass
391 );
392 }
393 }
394 }
395#endif
396}
397
425template <int depos_order>
427 const amrex::ParticleReal * wp,
428 const amrex::ParticleReal * uxp,
429 const amrex::ParticleReal * uyp,
430 const amrex::ParticleReal * uzp,
431 amrex::FArrayBox& varx_fab,
432 amrex::FArrayBox& vary_fab,
433 amrex::FArrayBox& varz_fab,
434 amrex::IArrayBox& nx_iab,
435 amrex::IArrayBox& ny_iab,
436 amrex::IArrayBox& nz_iab,
437 amrex::FArrayBox& wx_fab,
438 amrex::FArrayBox& wy_fab,
439 amrex::FArrayBox& wz_fab,
440 amrex::FArrayBox& w2x_fab,
441 amrex::FArrayBox& w2y_fab,
442 amrex::FArrayBox& w2z_fab,
443 amrex::FArrayBox& vxbar_fab,
444 amrex::FArrayBox& vybar_fab,
445 amrex::FArrayBox& vzbar_fab,
446 const TemperatureDepositionType type,
447 const TemperatureDepositionPass pass,
448 const long np_to_deposit,
449 const amrex::Real relative_time,
450 const amrex::XDim3 & dinv,
451 const amrex::XDim3 & xyzmin,
452 const amrex::Dim3 lo,
453 [[maybe_unused]]const int n_rz_azimuthal_modes)
454{
455 using namespace amrex::literals;
456
457#if defined(WARPX_DIM_RZ)
459 n_rz_azimuthal_modes == 1,
460 "Azimuthal Fourier decomposition for temperature deposition is not implemented."
461 " Only mode=0 supported in RZ."
462 );
463#endif
464
465 const amrex::Real inv_c2 = PhysConst::inv_c2;
466
467 const amrex::Array4<int> & nx_arr = nx_iab.array();
468 const amrex::Array4<int> & ny_arr = ny_iab.array();
469 const amrex::Array4<int> & nz_arr = nz_iab.array();
470 const amrex::Array4<amrex::Real> & wx_arr = wx_fab.array();
471 const amrex::Array4<amrex::Real> & wy_arr = wy_fab.array();
472 const amrex::Array4<amrex::Real> & wz_arr = wz_fab.array();
473 const amrex::Array4<amrex::Real> & w2x_arr = w2x_fab.array();
474 const amrex::Array4<amrex::Real> & w2y_arr = w2y_fab.array();
475 const amrex::Array4<amrex::Real> & w2z_arr = w2z_fab.array();
476 const amrex::Array4<amrex::Real> & vxbar_arr = vxbar_fab.array();
477 const amrex::Array4<amrex::Real> & vybar_arr = vybar_fab.array();
478 const amrex::Array4<amrex::Real> & vzbar_arr = vzbar_fab.array();
479 const amrex::IntVect & varx_type = varx_fab.box().type();
480 const amrex::IntVect & vary_type = vary_fab.box().type();
481 const amrex::IntVect & varz_type = varz_fab.box().type();
482
483 // Loop over particles and deposit into the deposition buffers
485 np_to_deposit,
486 [=] AMREX_GPU_DEVICE (long ip) {
487 amrex::ParticleReal xp, yp, zp;
488 GetPosition(ip, xp, yp, zp);
489
490 // --- Get particle quantities
491 const amrex::ParticleReal gaminv = 1.0_rt/std::sqrt(1.0_rt + uxp[ip]*uxp[ip]*inv_c2
492 + uyp[ip]*uyp[ip]*inv_c2
493 + uzp[ip]*uzp[ip]*inv_c2);
494 const amrex::ParticleReal vx = uxp[ip]*gaminv;
495 const amrex::ParticleReal vy = uyp[ip]*gaminv;
496 const amrex::ParticleReal vz = uzp[ip]*gaminv;
497
499 xp, yp, zp, wp[ip], vx, vy, vz,
500 nx_arr, ny_arr, nz_arr,
501 wx_arr, wy_arr, wz_arr,
502 w2x_arr, w2y_arr, w2z_arr,
503 vxbar_arr, vybar_arr, vzbar_arr,
504 varx_type, vary_type, varz_type,
505 type, pass,
506 relative_time, dinv, xyzmin,
507 lo, n_rz_azimuthal_modes);
508 }
509 );
510
512}
513
514} // namespace warpx::particles::deposition
515#endif // WARPX_TEMPERATUREDEPOSITION_H_
#define AMREX_INLINE
#define AMREX_GPU_DEVICE
#define AMREX_GPU_HOST_DEVICE
@ vz
Definition RigidInjectedParticleContainer.H:27
#define WARPX_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition TextMsg.H:13
const Box & box() const noexcept
Array4< T const > array() const noexcept
__host__ __device__ IntVectND< dim > type() const noexcept
amrex_real Real
amrex_particle_real ParticleReal
ArrayND< T, 4, true > Array4
constexpr auto inv_c2
inverse of the square of the vacuum speed of light [s^2/m^2]
Definition constant.H:220
__host__ __device__ AMREX_FORCE_INLINE void AddNoRet(T *sum, T value) noexcept
void streamSynchronize() noexcept
std::enable_if_t< std::is_integral_v< T > > ParallelFor(TypeList< CTOs... > ctos, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
IntVectND< 3 > IntVect
Definition TemperatureDeposition.H:32
AMREX_GPU_HOST_DEVICE AMREX_INLINE void varianceDepositionSubKernel(amrex::ParticleReal vx, amrex::ParticleReal vy, amrex::ParticleReal vz, const amrex::IntVectND< 3 > &ixv, const amrex::IntVectND< 3 > &iyv, const amrex::IntVectND< 3 > &izv, amrex::Real wpx_var, amrex::Real wpy_var, amrex::Real wpz_var, const amrex::Array4< int > &nx_arr, const amrex::Array4< int > &ny_arr, const amrex::Array4< int > &nz_arr, const amrex::Array4< amrex::Real > &wx_arr, const amrex::Array4< amrex::Real > &wy_arr, const amrex::Array4< amrex::Real > &wz_arr, const amrex::Array4< amrex::Real > &w2x_arr, const amrex::Array4< amrex::Real > &w2y_arr, const amrex::Array4< amrex::Real > &w2z_arr, const amrex::Array4< amrex::Real > &vxbar_arr, const amrex::Array4< amrex::Real > &vybar_arr, const amrex::Array4< amrex::Real > &vzbar_arr, const TemperatureDepositionType type, const TemperatureDepositionPass pass)
Definition TemperatureDeposition.H:35
TemperatureDepositionType
Definition TemperatureDepositionTypes.H:20
@ SINGLE_PASS
Definition TemperatureDepositionTypes.H:20
@ DOUBLE_PASS
Definition TemperatureDepositionTypes.H:20
AMREX_GPU_HOST_DEVICE AMREX_INLINE void doVarianceDepositionShapeNKernel(const amrex::ParticleReal xp, const amrex::ParticleReal yp, const amrex::ParticleReal zp, const amrex::ParticleReal wp, const amrex::ParticleReal vx, const amrex::ParticleReal vy, const amrex::ParticleReal vz, const amrex::Array4< int > &nx_arr, const amrex::Array4< int > &ny_arr, const amrex::Array4< int > &nz_arr, const amrex::Array4< amrex::Real > &wx_arr, const amrex::Array4< amrex::Real > &wy_arr, const amrex::Array4< amrex::Real > &wz_arr, const amrex::Array4< amrex::Real > &w2x_arr, const amrex::Array4< amrex::Real > &w2y_arr, const amrex::Array4< amrex::Real > &w2z_arr, const amrex::Array4< amrex::Real > &vxbar_arr, const amrex::Array4< amrex::Real > &vybar_arr, const amrex::Array4< amrex::Real > &vzbar_arr, const amrex::IntVect &varx_type, const amrex::IntVect &vary_type, const amrex::IntVect &varz_type, const TemperatureDepositionType type, const TemperatureDepositionPass pass, const amrex::Real relative_time, const amrex::XDim3 &dinv, const amrex::XDim3 &xyzmin, const amrex::Dim3 lo, const int n_rz_azimuthal_modes)
Kernel for the direct current deposition for thread thread_num.
Definition TemperatureDeposition.H:154
TemperatureDepositionPass
Definition TemperatureDepositionTypes.H:25
@ SECOND
Definition TemperatureDepositionTypes.H:25
@ FIRST
Definition TemperatureDepositionTypes.H:25
void doVarianceDepositionShapeN(const GetParticlePosition< PIdx > &GetPosition, const amrex::ParticleReal *wp, const amrex::ParticleReal *uxp, const amrex::ParticleReal *uyp, const amrex::ParticleReal *uzp, amrex::FArrayBox &varx_fab, amrex::FArrayBox &vary_fab, amrex::FArrayBox &varz_fab, amrex::IArrayBox &nx_iab, amrex::IArrayBox &ny_iab, amrex::IArrayBox &nz_iab, amrex::FArrayBox &wx_fab, amrex::FArrayBox &wy_fab, amrex::FArrayBox &wz_fab, amrex::FArrayBox &w2x_fab, amrex::FArrayBox &w2y_fab, amrex::FArrayBox &w2z_fab, amrex::FArrayBox &vxbar_fab, amrex::FArrayBox &vybar_fab, amrex::FArrayBox &vzbar_fab, const TemperatureDepositionType type, const TemperatureDepositionPass pass, const long np_to_deposit, const amrex::Real relative_time, const amrex::XDim3 &dinv, const amrex::XDim3 &xyzmin, const amrex::Dim3 lo, const int n_rz_azimuthal_modes)
Temperature Deposition Algorithm from Bell (1979) https://dl.acm.org/doi/pdf/10.1145/359146....
Definition TemperatureDeposition.H:426
Definition ShapeFactors.H:29
Functor that can be used to extract the positions of the macroparticles inside a ParallelFor kernel.
Definition GetAndSetPosition.H:75