1#include "../../include/integrator/device_explicit_euler.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) {
13 particles->
vx[bodyIndex + offset] += dt * (particles->
ax[bodyIndex + offset] + particles->
g_ax[bodyIndex + offset]);
15 particles->
vy[bodyIndex + offset] += dt * (particles->
ay[bodyIndex + offset] + particles->
g_ay[bodyIndex + offset]);
17 particles->
vz[bodyIndex + offset] += dt * (particles->
az[bodyIndex + offset] + particles->
g_az[bodyIndex + offset]);
27 particles->
x[bodyIndex + offset] += dt * particles->
vx[bodyIndex + offset];
29 particles->
y[bodyIndex + offset] += dt * particles->
vy[bodyIndex + offset];
31 particles->
z[bodyIndex + offset] += dt * particles->
vz[bodyIndex + offset];
36 particles->
rho[bodyIndex + offset] = particles->
rho[bodyIndex + offset] + dt * particles->
drhodt[bodyIndex + offset];
40 particles->
e[bodyIndex + offset] += dt * particles->
dedt[bodyIndex + offset];
41 if (particles->
e[bodyIndex + offset] < 1e-50) {
42 particles->
e[bodyIndex + offset] = 1e-50;
48 particles->
sml[bodyIndex + offset] = particles->
sml[bodyIndex + offset] + dt * particles->dsmldt[bodyIndex + offset];
Execution policy/instruction for CUDA kernel execution.
Particle(s) class based on SoA (Structur of Arrays).
real * e
(pointer to) internal energy (array)
real * x
(pointer to) x position (array)
real * rho
(pointer to) density (array)
real * drhodt
(pointer to) time derivative of density (array)
real * g_ay
(pointer to) y gravitational acceleration (array)
real * ay
(pointer to) y acceleration (array)
real * az
(pointer to) z acceleration (array)
real * y
(pointer to) y position (array)
real * ax
(pointer to) x acceleration (array)
real * g_ax
(pointer to) x gravitational acceleration (array)
real * sml
(pointer to) smoothing length (array)
real * dedt
(pointer to) time derivative of internal energy (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 update(Particles *particles, integer n, real dt)
Wrapper for ExplicitEulerNS::Kernel::update().
__global__ void update(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.