trevert to previous logic - sphere - GPU-based 3D discrete element method algor… | |
git clone git://src.adamsgaard.dk/sphere | |
Log | |
Files | |
Refs | |
LICENSE | |
--- | |
commit b0e5f3b5ebb17b6639a2464a9deea444dee67e50 | |
parent c24a98931ae68aaabf77055f855b50ef0c118be1 | |
Author: Anders Damsgaard Christensen <[email protected]> | |
Date: Wed, 15 Jun 2016 15:41:19 -0700 | |
revert to previous logic | |
Diffstat: | |
M src/contactsearch.cuh | 80 -----------------------------… | |
M src/darcy.cuh | 12 ++++++++++-- | |
2 files changed, 10 insertions(+), 82 deletions(-) | |
--- | |
diff --git a/src/contactsearch.cuh b/src/contactsearch.cuh | |
t@@ -80,86 +80,6 @@ __device__ int findDistMod(int3* targetCell, Float3* distmo… | |
return 0; | |
} | |
-/* Copy of findDistMod, but modified for mirror BC at frictionless boundaries | |
- during porosity esimation */ | |
-__device__ int findDistModPorosity(int3* targetCell, Float3* distmod) | |
-{ | |
- // Check whether x- and y boundaries are to be treated as periodic | |
- // 1: x- and y boundaries periodic | |
- // 2: x boundaries periodic | |
- if (devC_grid.periodic == 1) { | |
- | |
- // Periodic x-boundary | |
- if (targetCell->x < 0) { | |
- //targetCell->x = devC_grid.num[0] - 1; | |
- targetCell->x += devC_grid.num[0]; | |
- *distmod += MAKE_FLOAT3(devC_grid.L[0], 0.0, 0.0); | |
- } | |
- if (targetCell->x >= devC_grid.num[0]) { | |
- //targetCell->x = 0; | |
- targetCell->x -= devC_grid.num[0]; | |
- *distmod -= MAKE_FLOAT3(devC_grid.L[0], 0.0, 0.0); | |
- } | |
- | |
- // Periodic y-boundary | |
- if (targetCell->y < 0) { | |
- //targetCell->y = devC_grid.num[1] - 1; | |
- targetCell->y += devC_grid.num[1]; | |
- *distmod += MAKE_FLOAT3(0.0, devC_grid.L[1], 0.0); | |
- } | |
- if (targetCell->y >= devC_grid.num[1]) { | |
- //targetCell->y = 0; | |
- targetCell->y -= devC_grid.num[1]; | |
- *distmod -= MAKE_FLOAT3(0.0, devC_grid.L[1], 0.0); | |
- } | |
- | |
- | |
- // Only x-boundaries are periodic | |
- } else if (devC_grid.periodic == 2) { | |
- | |
- // Periodic x-boundary | |
- if (targetCell->x < 0) { | |
- //targetCell->x = devC_grid.num[0] - 1; | |
- targetCell->x += devC_grid.num[0]; | |
- *distmod += MAKE_FLOAT3(devC_grid.L[0], 0.0, 0.0); | |
- } | |
- if (targetCell->x >= devC_grid.num[0]) { | |
- //targetCell->x = 0; | |
- targetCell->x -= devC_grid.num[0]; | |
- *distmod -= MAKE_FLOAT3(devC_grid.L[0], 0.0, 0.0); | |
- } | |
- | |
- // Hande out-of grid cases on y-axis by mirroring the grid | |
- if (targetCell->y < 0) { | |
- targetCell->y += 2; | |
- //*distmod += MAKE_FLOAT3(devC_grid.L[0], 0.0, 0.0); | |
- } | |
- if (targetCell->y >= devC_grid.num[1]) { | |
- targetCell->y -= 2; | |
- //*distmod += MAKE_FLOAT3(devC_grid.L[0], 0.0, 0.0); | |
- } | |
- | |
- | |
- // No periodic boundaries | |
- } else { | |
- | |
- // Don't modify targetCell or distmod. | |
- | |
- // Hande out-of grid cases on x- and y-axes | |
- if (targetCell->x < 0 || targetCell->x >= devC_grid.num[0]) | |
- return -1; | |
- if (targetCell->y < 0 || targetCell->y >= devC_grid.num[1]) | |
- return -1; | |
- } | |
- | |
- // Handle out-of-grid cases on z-axis | |
- if (targetCell->z < 0 || targetCell->z >= devC_grid.num[2]) | |
- return -1; | |
- else | |
- return 0; | |
-} | |
- | |
- | |
// Find overlaps between particle no. 'idx' and particles in cell 'gridpos'. | |
// Contacts are processed as soon as they are identified. | |
diff --git a/src/darcy.cuh b/src/darcy.cuh | |
t@@ -463,10 +463,18 @@ __global__ void findDarcyPorositiesLinear( | |
if (targetCell.z == -1) | |
targetCell.z = 1; | |
+ // Mirror particle grid at frictionless boundaries | |
+ /*if (devC_grid.periodic == 2) { | |
+ if (targetCell.y == -1) | |
+ targetCell.y = 1; | |
+ else if (targetCell.y == devC_grid.num[1]) | |
+ targetCell.y = devC_grid.num[1] - 2; | |
+ }*/ | |
+ | |
// Get distance modifier for interparticle | |
// vector, if it crosses a periodic boundary | |
distmod = MAKE_FLOAT3(0.0, 0.0, 0.0); | |
- if (findDistModPorosity(&targetCell, &distmod) != -1) { | |
+ if (findDistMod(&targetCell, &distmod) != -1) { | |
// Calculate linear cell ID | |
cellID = targetCell.x | |
t@@ -767,7 +775,7 @@ __global__ void findDarcyPorosities( | |
// Get distance modifier for interparticle | |
// vector, if it crosses a periodic boundary | |
distmod = MAKE_FLOAT3(0.0, 0.0, 0.0); | |
- if (findDistModPorosity(&targetCell, &distmod) != -1) { | |
+ if (findDistMod(&targetCell, &distmod) != -1) { | |
// Calculate linear cell ID | |
cellID = targetCell.x |