Skip to content
Snippets Groups Projects
Commit efe0d0d5 authored by cmaffeo2's avatar cmaffeo2
Browse files

Made tabulated bonded kernels use grid-stride loops

parent ee3383f6
No related branches found
No related tags found
No related merge requests found
...@@ -605,12 +605,10 @@ float ComputeForce::computeTabulated(bool get_energy) { ...@@ -605,12 +605,10 @@ float ComputeForce::computeTabulated(bool get_energy) {
computeTabulatedBonds <<<numBlocks, numThreads>>> ( forceInternal_d, pos_d, num, numParts, sys_d, bondList_d, (numBonds/2), numReplicas, energies_d, get_energy, tableBond_d); computeTabulatedBonds <<<numBlocks, numThreads>>> ( forceInternal_d, pos_d, num, numParts, sys_d, bondList_d, (numBonds/2), numReplicas, energies_d, get_energy, tableBond_d);
} }
if (angleList_d != NULL && tableAngle_d != NULL) if (angleList_d != NULL && tableAngle_d != NULL)
computeTabulatedAngles<<<numBlocks, numThreads>>>(forceInternal_d, pos_d, sys_d, numReplicas, computeTabulatedAngles<<<numBlocks, numThreads>>>(forceInternal_d, pos_d, sys_d, numAngles*numReplicas, angleList_d, tableAngle_d);
numAngles, angleList_d, tableAngle_d);
if (dihedralList_d != NULL && tableDihedral_d != NULL) if (dihedralList_d != NULL && tableDihedral_d != NULL)
computeTabulatedDihedrals<<<numBlocks, numThreads>>>(forceInternal_d, pos_d, sys_d, numReplicas, computeTabulatedDihedrals<<<numBlocks, numThreads>>>(forceInternal_d, pos_d, sys_d, numDihedrals*numReplicas, dihedralList_d, dihedralPotList_d, tableDihedral_d);
numDihedrals, dihedralList_d, dihedralPotList_d, tableDihedral_d);
// Calculate the energy based on the array created by the kernel // Calculate the energy based on the array created by the kernel
......
...@@ -695,9 +695,9 @@ void computeAngles(Vector3 force[], Vector3 pos[], ...@@ -695,9 +695,9 @@ void computeAngles(Vector3 force[], Vector3 pos[],
} }
}*/ }*/
__global__ void computeTabulatedBonds( Vector3* force, Vector3* pos, int num, __global__ void computeTabulatedBonds( Vector3* __restrict__ force, Vector3* __restrict__ pos, int num,
int numParts, BaseGrid* sys, int3* bondList_d, int numParts, BaseGrid* __restrict__ sys, int3* __restrict__ bondList_d,
int numBonds, int numReplicas, float* g_energies, int numBonds, int numReplicas, float* __restrict__ g_energies,
bool get_energy, TabulatedPotential** tableBond) bool get_energy, TabulatedPotential** tableBond)
{ {
// Thread's unique ID. // Thread's unique ID.
...@@ -735,25 +735,21 @@ __global__ void computeTabulatedBonds( Vector3* force, Vector3* pos, int nu ...@@ -735,25 +735,21 @@ __global__ void computeTabulatedBonds( Vector3* force, Vector3* pos, int nu
// TODO: add kernel for energy calculation // TODO: add kernel for energy calculation
__global__ __global__
void computeTabulatedAngles(Vector3* __restrict__ force, Vector3* __restrict__ pos, void computeTabulatedAngles(Vector3* __restrict__ force,
BaseGrid* __restrict__ sys, int numReplicas, Vector3* __restrict__ pos,
int numAngles, int4* angleList_d, TabulatedAnglePotential** tableAngle) { BaseGrid* __restrict__ sys,
int numAngles, int4* __restrict__ angleList_d, TabulatedAnglePotential** tableAngle) {
int currAngle = blockIdx.x * blockDim.x + threadIdx.x; // first particle ID // Loop over ALL angles in ALL replicas
for (int i = threadIdx.x+blockIdx.x*blockDim.x; i<numAngles; i+=blockDim.x*gridDim.x) {
// Loop over ALL angles in ALL replicas int4& ids = angleList_d[i];
if( currAngle < (numAngles * numReplicas) ) { computeAngle(tableAngle[ ids.w ], sys, force, pos, ids.x, ids.y, ids.z);
int4& ids = angleList_d[currAngle];
computeAngle(tableAngle[ ids.w ], sys, force, pos, ids.x, ids.y, ids.z);
// if (get_energy) // if (get_energy)
// { // {
// //TODO: clarification on energy computation needed, consider changing. // //TODO: clarification on energy computation needed, consider changing.
// atomicAdd( &g_energies[i], energy_local); // atomicAdd( &g_energies[i], energy_local);
// //atomicAdd( &g_energies[j], energy_local); // //atomicAdd( &g_energies[j], energy_local);
// } // }
} }
} }
...@@ -847,7 +843,7 @@ void computeDihedrals(Vector3 force[], Vector3 pos[], ...@@ -847,7 +843,7 @@ void computeDihedrals(Vector3 force[], Vector3 pos[],
__global__ __global__
void computeTabulatedDihedrals(Vector3* __restrict__ force, const Vector3* __restrict__ pos, void computeTabulatedDihedrals(Vector3* __restrict__ force, const Vector3* __restrict__ pos,
const BaseGrid* __restrict__ sys, int numReplicas, const BaseGrid* __restrict__ sys,
int numDihedrals, const int4* __restrict__ dihedralList_d, int numDihedrals, const int4* __restrict__ dihedralList_d,
const int* __restrict__ dihedralPotList_d, TabulatedDihedralPotential** tableDihedral) { const int* __restrict__ dihedralPotList_d, TabulatedDihedralPotential** tableDihedral) {
...@@ -855,11 +851,10 @@ void computeTabulatedDihedrals(Vector3* __restrict__ force, const Vector3* __res ...@@ -855,11 +851,10 @@ void computeTabulatedDihedrals(Vector3* __restrict__ force, const Vector3* __res
// Loop over ALL dihedrals in ALL replicas // Loop over ALL dihedrals in ALL replicas
// TODO make grid stride loop // TODO make grid stride loop
if( currDihedral < (numDihedrals * numReplicas) ) { for (int i = threadIdx.x+blockIdx.x*blockDim.x; i < numDihedrals; i+=blockDim.x*gridDim.x) {
const int4& ids = dihedralList_d[currDihedral]; const int4& ids = dihedralList_d[i];
const int& id = dihedralPotList_d[currDihedral]; const int& id = dihedralPotList_d[i];
computeDihedral(tableDihedral[ id ], sys, force, pos, ids.x, ids.y, ids.z, ids.w);
computeDihedral(tableDihedral[ id ], sys, force, pos, ids.x, ids.y, ids.z, ids.w);
// if (get_energy) // if (get_energy)
// { // {
......
...@@ -83,7 +83,7 @@ public: ...@@ -83,7 +83,7 @@ public:
print("INIT\n"); print("INIT\n");
} }
HOST DEVICE inline void print(char* string) HOST DEVICE inline void print(const char* string)
{ {
printf(string); printf(string);
printf("RAND48_SEED = %d, RAND48_MULT = %d, RAND48_ADD = %d, RAND48_MASK = %d\n", RAND48_SEED, RAND48_MULT, RAND48_ADD, RAND48_MASK); printf("RAND48_SEED = %d, RAND48_MULT = %d, RAND48_ADD = %d, RAND48_MASK = %d\n", RAND48_SEED, RAND48_MULT, RAND48_ADD, RAND48_MASK);
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment