tfix arguments - sphere - GPU-based 3D discrete element method algorithm with o… | |
git clone git://src.adamsgaard.dk/sphere | |
Log | |
Files | |
Refs | |
LICENSE | |
--- | |
commit e1e3c1d366af936d429548770cd277edc14d2e10 | |
parent d009927b87ac9bdd75a225ef9af5a7433930b2bc | |
Author: Anders Damsgaard <[email protected]> | |
Date: Tue, 3 Jun 2014 16:21:38 +0200 | |
fix arguments | |
Diffstat: | |
M src/device.cu | 14 +++++++------- | |
M src/navierstokes.cuh | 24 ++++++++++++------------ | |
2 files changed, 19 insertions(+), 19 deletions(-) | |
--- | |
diff --git a/src/device.cu b/src/device.cu | |
t@@ -1199,17 +1199,17 @@ __host__ void DEM::startTime() | |
startTimer(&kernel_tic); | |
findNSforcing<<<dimGridFluid, dimBlockFluid>>>( | |
dev_ns_epsilon, | |
- dev_ns_f1, | |
- dev_ns_f2, | |
- dev_ns_f, | |
dev_ns_phi, | |
dev_ns_dphi, | |
dev_ns_v_p, | |
- dev_ns_v_x, | |
- dev_ns_v_y, | |
- dev_ns_v_z, | |
+ dev_ns_v_p_x, | |
+ dev_ns_v_p_y, | |
+ dev_ns_v_p_z, | |
nijac, | |
- ns.ndem); | |
+ ns.ndem, | |
+ dev_ns_f1, | |
+ dev_ns_f2, | |
+ dev_ns_f); | |
cudaThreadSynchronize(); | |
if (PROFILING == 1) | |
stopTimer(&kernel_tic, &kernel_toc, &kernel_elapsed, | |
diff --git a/src/navierstokes.cuh b/src/navierstokes.cuh | |
t@@ -2242,18 +2242,18 @@ __global__ void findPredNSvelocities( | |
// At each iteration, the value of the forcing function is found as: | |
// f = f1 - f2 dot grad(epsilon) | |
__global__ void findNSforcing( | |
- Float* dev_ns_epsilon, | |
- Float* dev_ns_f1, | |
- Float3* dev_ns_f2, | |
- Float* dev_ns_f, | |
- Float* dev_ns_phi, | |
- Float* dev_ns_dphi, | |
- Float3* dev_ns_v_p, | |
- Float* dev_ns_v_x, | |
- Float* dev_ns_v_y, | |
- Float* dev_ns_v_z, | |
- unsigned int nijac, | |
- unsigned int ndem) | |
+ Float* dev_ns_epsilon, // in | |
+ Float* dev_ns_phi, // in | |
+ Float* dev_ns_dphi, // in | |
+ Float3* dev_ns_v_p, // in | |
+ Float* dev_ns_v_p_x, // in | |
+ Float* dev_ns_v_p_y, // in | |
+ Float* dev_ns_v_p_z, // in | |
+ unsigned int nijac, // in | |
+ unsigned int ndem, // in | |
+ Float* dev_ns_f1, // out | |
+ Float3* dev_ns_f2, // out | |
+ Float* dev_ns_f) // out | |
{ | |
// 3D thread index | |
const unsigned int x = blockDim.x * blockIdx.x + threadIdx.x; |