Introduction
Introduction Statistics Contact Development Disclaimer Help
tRemoved NC error reporting, an error exception system is needed - sphere - GPU…
git clone git://src.adamsgaard.dk/sphere
Log
Files
Refs
LICENSE
---
commit 72e8a2458f7da3e8d00b2476c966715afa57fef7
parent 377088874c49108d71af8048d568abb87baebd83
Author: Anders Damsgaard <[email protected]>
Date: Wed, 10 Oct 2012 11:02:59 +0200
Removed NC error reporting, an error exception system is needed
Diffstat:
M src/contactsearch.cuh | 25 +++++++++----------------
M src/datatypes.h | 4 ++--
M src/device.cu | 17 ++++++++---------
3 files changed, 19 insertions(+), 27 deletions(-)
---
diff --git a/src/contactsearch.cuh b/src/contactsearch.cuh
t@@ -179,7 +179,7 @@ __device__ void findAndProcessContactsInCell(int3 targetCe…
// Used for shearmodel=2, where bookkeeping of contact history is necessary.
// Kernel executed on device, and callable from device only.
// Function is called from topology().
-__device__ int findContactsInCell(int3 targetCell,
+__device__ void findContactsInCell(int3 targetCell,
unsigned int idx_a,
Float4 x_a, Float radius_a,
Float4* dev_x_sorted,
t@@ -223,6 +223,7 @@ __device__ int findContactsInCell(int3 targetCell,
for (unsigned int idx_b = startIdx; idx_b<endIdx; ++idx_b) {
if (idx_b != idx_a) { // Do not collide particle with itself
+
// Fetch position and radius of particle B.
Float4 x_b = dev_x_sorted[idx_b];
Float radius_b = dev_radius_sorted[idx_b];
t@@ -275,11 +276,12 @@ __device__ int findContactsInCell(int3 targetCell,
// Increment contact counter
++*nc;
-
+
// Check if the number of contacts of particle A
// exceeds the max. number of contacts per particle
- if (*nc >= NC)
- return 1;
+ if (*nc > devC_nc)
+ return; // I would like to throw an error, but it doesn't seem pos…
+
}
// Write the inter-particle position vector correction and radius of p…
t@@ -300,9 +302,6 @@ __device__ int findContactsInCell(int3 targetCell,
} // Do not collide particle with itself end
} // Iterate over cell particles end
} // Check wether cell is empty end
-
- // Return without errors
- return 0;
} // End of findContactsInCell(...)
t@@ -310,7 +309,7 @@ __device__ int findContactsInCell(int3 targetCell,
// Search for neighbors to particle 'idx' inside the 27 closest cells,
// and save the contact pairs in global memory.
// Function is called from mainGPU loop.
-__global__ int topology(unsigned int* dev_cellStart,
+__global__ void topology(unsigned int* dev_cellStart,
unsigned int* dev_cellEnd, // Input: Particles in…
unsigned int* dev_gridParticleIndex, // Input: Unsort…
Float4* dev_x_sorted, Float* dev_radius_sorted,
t@@ -326,7 +325,6 @@ __global__ int topology(unsigned int* dev_cellStart,
// Count the number of contacts in this time step
int nc = 0;
- int status = 0;
// Grid position of host and neighbor cells in uniform, cubic grid
int3 gridPos;
t@@ -343,20 +341,15 @@ __global__ int topology(unsigned int* dev_cellStart,
for (int y_dim=-1; y_dim<2; ++y_dim) { // y-axis
for (int x_dim=-1; x_dim<2; ++x_dim) { // x-axis
targetPos = gridPos + make_int3(x_dim, y_dim, z_dim);
- if (findContactsInCell(targetPos, idx_a, x_a, radius_a,
+ findContactsInCell(targetPos, idx_a, x_a, radius_a,
dev_x_sorted, dev_radius_sorted,
dev_cellStart, dev_cellEnd,
dev_gridParticleIndex,
- &nc, dev_contacts, dev_distmod) == 1)
- status = 1;
+ &nc, dev_contacts, dev_distmod);
}
}
}
}
-
- // Returns 0 if no particles have more contacts than allowed,
- // and 1 otherwise
- return status;
} // End of topology(...)
diff --git a/src/datatypes.h b/src/datatypes.h
t@@ -60,8 +60,8 @@ const unsigned int ND = 3;
const Float VERS = 0.25;
// Max. number of contacts per particle
-//#define NC 16
-#define NC 32
+//const char NC = 16;
+const char NC = 32;
///////////////////////////
diff --git a/src/device.cu b/src/device.cu
t@@ -484,16 +484,15 @@ __host__ void gpuMain(Float4* host_x,
// For each particle: Search contacts in neighbor cells
if (PROFILING == 1)
startTimer(&kernel_tic);
- if(topology<<<dimGrid, dimBlock>>>(dev_cellStart,
- dev_cellEnd,
- dev_gridParticleIndex,
- dev_x_sorted,
- dev_radius_sorted,
- dev_contacts,
- dev_distmod) == 1) {
+ topology<<<dimGrid, dimBlock>>>(dev_cellStart,
+ dev_cellEnd,
+ dev_gridParticleIndex,
+ dev_x_sorted,
+ dev_radius_sorted,
+ dev_contacts,
+ dev_distmod);
std::cerr << "Warning! One or more particles have more contacts than a…
- << "Raise NC in datatypes.h to accomodate.\n"
- }
+ << "Raise NC in datatypes.h to accomodate.\n";
// Empty cuPrintf() buffer to console
//cudaThreadSynchronize();
You are viewing proxied material from mx1.adamsgaard.dk. The copyright of proxied material belongs to its original authors. Any comments or complaints in relation to proxied material should be directed to the original authors of the content concerned. Please see the disclaimer for more details.