1#include "../../include/integrator/device_leapfrog.cuh"
2#include "../../include/cuda_utils/cuda_launcher.cuh"
6 integer bodyIndex = threadIdx.x + blockIdx.x * blockDim.x;
7 integer stride = blockDim.x * gridDim.x;
10 while (bodyIndex + offset < n) {
19 particles->
x[bodyIndex + offset] += dt * (particles->
vx[bodyIndex + offset] + 0.5 * dt * particles->
g_ax[bodyIndex + offset]);
20 particles->
g_ax_old[bodyIndex + offset] = particles->
g_ax[bodyIndex + offset];
22 particles->
y[bodyIndex + offset] += dt * (particles->
vy[bodyIndex + offset] + 0.5 * dt * particles->
g_ay[bodyIndex + offset]);
23 particles->
g_ay_old[bodyIndex + offset] = particles->
g_ay[bodyIndex + offset];
25 particles->
z[bodyIndex + offset] += dt * (particles->
vz[bodyIndex + offset] + 0.5 * dt * particles->
g_az[bodyIndex + offset]);
26 particles->
g_az_old[bodyIndex + offset] = particles->
g_az[bodyIndex + offset];
36 integer bodyIndex = threadIdx.x + blockIdx.x * blockDim.x;
37 integer stride = blockDim.x * gridDim.x;
40 while (bodyIndex + offset < n) {
48 particles->
vx[bodyIndex + offset] += 0.5 * dt * (particles->
g_ax[bodyIndex + offset] + particles->
g_ax_old[bodyIndex + offset]);
50 particles->
vy[bodyIndex + offset] += 0.5 * dt * (particles->
g_ay[bodyIndex + offset] + particles->
g_ay_old[bodyIndex + offset]);
52 particles->
vz[bodyIndex + offset] += 0.5 * dt * (particles->
g_az[bodyIndex + offset] + particles->
g_az_old[bodyIndex + offset]);
Execution policy/instruction for CUDA kernel execution.
Particle(s) class based on SoA (Structur of Arrays).
real * x
(pointer to) x position (array)
real * g_ay
(pointer to) y gravitational acceleration (array)
real * y
(pointer to) y position (array)
real * g_ax
(pointer to) x gravitational acceleration (array)
real * z
(pointer to) z position (array)
real * vz
(pointer to) z velocity (array)
real * g_az
(pointer to) z gravitational acceleration (array)
real * vx
(pointer to) x velocity (array)
real * vy
(pointer to) y velocity (array)
real updateX(Particles *particles, integer n, real dt)
Wrapper for LeapfrogNS::Kernel::update().
real updateV(Particles *particles, integer n, real dt)
__global__ void updateV(Particles *particles, integer n, real dt)
__global__ void updateX(Particles *particles, integer n, real dt)
Update/move/advance particles.
real launch(bool timeKernel, const ExecutionPolicy &policy, void(*f)(Arguments...), Arguments... args)
CUDA execution wrapper function.