tfixed out-of-grid errors in porosity estimation routine on non-periodic grids - sphere - GPU-based 3D discrete element method algorithm with optional fluid coupling
(HTM) git clone git://src.adamsgaard.dk/sphere
(DIR) Log
(DIR) Files
(DIR) Refs
(DIR) LICENSE
---
(DIR) commit fb03b34bb82f5d8cfac1c1f5d0c87961afe57dca
(DIR) parent 6ce8fd65468c3a399a4b131a3bce345a4cc3fdac
(HTM) Author: Anders Damsgaard <anders.damsgaard@geo.au.dk>
Date: Fri, 15 Aug 2014 19:51:47 +0200
fixed out-of-grid errors in porosity estimation routine on non-periodic grids
Diffstat:
M CMakeLists.txt | 2 +-
M src/contactsearch.cuh | 6 +++---
M src/navierstokes.cuh | 4 +++-
M tests/fluid_particle_interaction.py | 2 ++
4 files changed, 9 insertions(+), 5 deletions(-)
---
(DIR) diff --git a/CMakeLists.txt b/CMakeLists.txt
t@@ -31,7 +31,7 @@ enable_testing()
# Set build type = Debug
#set(CMAKE_BUILD_TYPE Debug)
#if (CUDA_FOUND)
-# set(CUDA_NVCC_FLAGS -g;-G)
+# set(CUDA_NVCC_FLAGS ${CUDA_NVCC_FLAGS};-g -G)
#endif()
# Set build type = Release
(DIR) diff --git a/src/contactsearch.cuh b/src/contactsearch.cuh
t@@ -57,7 +57,7 @@ __device__ int findDistMod(int3* targetCell, Float3* distmod)
}
// Hande out-of grid cases on y-axis
- if (targetCell->y < 0 || targetCell->y == devC_grid.num[1])
+ if (targetCell->y < 0 || targetCell->y >= devC_grid.num[1])
return -1;
t@@ -65,9 +65,9 @@ __device__ int findDistMod(int3* targetCell, Float3* distmod)
} else {
// Hande out-of grid cases on x- and y-axes
- if (targetCell->x < 0 || targetCell->x == devC_grid.num[0])
+ if (targetCell->x < 0 || targetCell->x >= devC_grid.num[0])
return -1;
- if (targetCell->y < 0 || targetCell->y == devC_grid.num[1])
+ if (targetCell->y < 0 || targetCell->y >= devC_grid.num[1])
return -1;
}
(DIR) diff --git a/src/navierstokes.cuh b/src/navierstokes.cuh
t@@ -1184,15 +1184,17 @@ __global__ void findPorositiesVelocitiesDiametersSpherical(
cellID = targetCell.x
+ targetCell.y * devC_grid.num[0]
+ (devC_grid.num[0] * devC_grid.num[1])
- * targetCell.z;
+ * targetCell.z;
// Lowest particle index in cell
+ __syncthreads();
startIdx = dev_cellStart[cellID];
// Make sure cell is not empty
if (startIdx != 0xffffffff) {
// Highest particle index in cell
+ __syncthreads();
endIdx = dev_cellEnd[cellID];
// Iterate over cell particles
(DIR) diff --git a/tests/fluid_particle_interaction.py b/tests/fluid_particle_interaction.py
t@@ -20,6 +20,8 @@ sim.addParticle([0.5, 0.5, 0.5], 0.05)
sim.initTemporal(total=0.01, file_dt=0.001)
sim.run(verbose=False)
+#sim.run(dry=True)
+#sim.run(cudamemcheck=True)
#sim.writeVTKall()
sim.readlast()