tback to cell-center f_pf and F_pf - sphere - GPU-based 3D discrete element met… | |
git clone git://src.adamsgaard.dk/sphere | |
Log | |
Files | |
Refs | |
LICENSE | |
--- | |
commit 03a1e5cfa91cabac2788716303dc592da5d261f1 | |
parent bcc75acd3cd0e3e2b9e32510e870232204eb3cc5 | |
Author: Anders Damsgaard <[email protected]> | |
Date: Sat, 7 Jun 2014 11:36:43 +0200 | |
back to cell-center f_pf and F_pf | |
Diffstat: | |
M src/device.cu | 8 ++++++++ | |
M src/navierstokes.cuh | 8 ++++---- | |
M src/sphere.h | 3 +++ | |
3 files changed, 15 insertions(+), 4 deletions(-) | |
--- | |
diff --git a/src/device.cu b/src/device.cu | |
t@@ -929,6 +929,11 @@ __host__ void DEM::startTime() | |
cudaThreadSynchronize(); | |
checkForCudaErrorsIter("Post setNSghostNodes(dev_ns_p)", iter); | |
+ setNSghostNodes<Float><<<dimGridFluid, dimBlockFluid>>>( | |
+ dev_ns_phi, ns.bc_bot, ns.bc_top); | |
+ cudaThreadSynchronize(); | |
+ checkForCudaErrorsIter("Post setNSghostNodes(dev_ns_p)", iter); | |
+ | |
// Per particle, find the fluid-particle interaction force f_pf | |
// and apply it to the particle | |
findInteractionForce<<<dimGrid, dimBlock>>>( | |
t@@ -959,6 +964,9 @@ __host__ void DEM::startTime() | |
dev_cellEnd, | |
dev_ns_f_pf, | |
dev_ns_F_pf); | |
+ //dev_ns_F_pf_x, | |
+ //dev_ns_F_pf_y, | |
+ //dev_ns_F_pf_z); | |
cudaThreadSynchronize(); | |
checkForCudaErrorsIter("Post applyInteractionForceToFluid", | |
iter); | |
diff --git a/src/navierstokes.cuh b/src/navierstokes.cuh | |
t@@ -2953,10 +2953,10 @@ __global__ void findInteractionForce( | |
const Float v_x = dev_ns_v_x[vidx(i_x,i_y,i_z)]; | |
const Float v_x_p = dev_ns_v_x[vidx(i_x+1,i_y,i_z)]; | |
- const Float v_y = dev_ns_v_x[vidx(i_x,i_y,i_z)]; | |
- const Float v_y_p = dev_ns_v_x[vidx(i_x,i_y+1,i_z)]; | |
- const Float v_z = dev_ns_v_x[vidx(i_x,i_y,i_z)]; | |
- const Float v_z_p = dev_ns_v_x[vidx(i_x,i_y,i_z+1)]; | |
+ const Float v_y = dev_ns_v_y[vidx(i_x,i_y,i_z)]; | |
+ const Float v_y_p = dev_ns_v_y[vidx(i_x,i_y+1,i_z)]; | |
+ const Float v_z = dev_ns_v_z[vidx(i_x,i_y,i_z)]; | |
+ const Float v_z_p = dev_ns_v_z[vidx(i_x,i_y,i_z+1)]; | |
const Float div_tau_x = dev_ns_div_tau_x[vidx(i_x,i_y,i_z)]; | |
const Float div_tau_x_p = dev_ns_div_tau_x[vidx(i_x+1,i_y,i_z)]; | |
diff --git a/src/sphere.h b/src/sphere.h | |
t@@ -182,6 +182,9 @@ class DEM { | |
Float3* dev_ns_vp_avg; // Average particle velocity in cell | |
Float* dev_ns_d_avg; // Average particle diameter in cell | |
Float3* dev_ns_F_pf; // Interaction force on fluid | |
+ //Float* dev_ns_F_pf_x; // Interaction force on fluid | |
+ //Float* dev_ns_F_pf_y; // Interaction force on fluid | |
+ //Float* dev_ns_F_pf_z; // Interaction force on fluid | |
Float* dev_ns_phi; // Cell porosity | |
Float* dev_ns_dphi; // Cell porosity change | |
//Float3* dev_ns_div_phi_v_v; // Divegence used in velocity prediction |