Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

[HIP] Optimized the spread kernel #9

Open
wants to merge 1 commit into
base: develop-2020
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
20 changes: 13 additions & 7 deletions src/gromacs/ewald/pme_calculate_splines.hip.h
Original file line number Diff line number Diff line change
Expand Up @@ -107,7 +107,7 @@ __device__ __forceinline__ void pme_gpu_stage_atom_data(const PmeGpuHipKernelPar
* \param[out] sm_gridlineIndices Atom gridline indices in the shared memory.
*/

template<const int order, const int atomsPerBlock, const int atomsPerWarp, const bool writeSmDtheta, const bool writeGlobal>
template<const int order, const int atomsPerBlock, const int atomsPerWarp, const bool writeSmDtheta, const bool writeGlobal, const bool useOrderThreads>
__device__ __forceinline__ void calculate_splines(const PmeGpuHipKernelParams kernelParams,
const int atomIndexOffset,
const float3 atomX,
Expand Down Expand Up @@ -137,17 +137,17 @@ __device__ __forceinline__ void calculate_splines(const PmeGpuHipKernelParams ke
/* Atom index w.r.t. global memory */
const int atomIndexGlobal = atomIndexOffset + atomIndexLocal;
/* Spline contribution index in one dimension */
const int threadLocalIdXY = (threadIdx.y * blockDim.x) + threadIdx.x;
const int orderIndex = threadLocalIdXY / DIM;
//const int threadLocalIdXY = (threadIdx.y * blockDim.x) + threadIdx.x;
const int orderIndex = threadIdx.y;//threadLocalIdXY / DIM;
/* Dimension index */
const int dimIndex = threadLocalIdXY % DIM;
const int dimIndex = threadIdx.x;//threadLocalIdXY % DIM;

/* Multi-purpose index of rvec/ivec atom data */
const int sharedMemoryIndex = atomIndexLocal * DIM + dimIndex;

float splineData[order];

const int localCheck = (dimIndex < DIM) && (orderIndex < 1);
const int localCheck = dimIndex < DIM;//(dimIndex < DIM) && (orderIndex < 1);
const int globalCheck = pme_gpu_check_atom_data_index(atomIndexGlobal, kernelParams.atoms.nAtoms);

/* we have 4 threads per atom, but can only use 3 here for the dimensions */
Expand Down Expand Up @@ -253,8 +253,11 @@ __device__ __forceinline__ void calculate_splines(const PmeGpuHipKernelParams ke
if (writeSmDtheta || writeGlobal)
{
/* Differentiation and storing the spline derivatives (dtheta) */
const int ithyMin = useOrderThreads ? 0 : orderIndex;
const int ithyMax = useOrderThreads ? order : orderIndex + 1;
#pragma unroll
for (o = 0; o < order; o++)
for (int o = ithyMin; o < ithyMax; o++)
//for (o = 0; o < order; o++)
{
const int thetaIndex =
getSplineParamIndex<order, atomsPerWarp>(thetaIndexBase, dimIndex, o);
Expand Down Expand Up @@ -286,8 +289,11 @@ __device__ __forceinline__ void calculate_splines(const PmeGpuHipKernelParams ke
splineData[0] = div * (1.0f - dr) * splineData[0];

/* Storing the spline values (theta) */
const int ithyMin = useOrderThreads ? 0 : orderIndex;
const int ithyMax = useOrderThreads ? order : orderIndex + 1;
#pragma unroll
for (o = 0; o < order; o++)
for (int o = ithyMin; o < ithyMax; o++)
//for (o = 0; o < order; o++)
{
const int thetaIndex =
getSplineParamIndex<order, atomsPerWarp>(thetaIndexBase, dimIndex, o);
Expand Down
2 changes: 1 addition & 1 deletion src/gromacs/ewald/pme_gather.hip.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -361,7 +361,7 @@ __launch_bounds__(c_gatherMaxThreadsPerBlock, c_gatherMinBlocksPerMP) __global__
atomX.y = kernelParams.atoms.d_coordinates[atomIndexGlobal * DIM + YY];
atomX.z = kernelParams.atoms.d_coordinates[atomIndexGlobal * DIM + ZZ];
}
calculate_splines<order, atomsPerBlock, atomsPerWarp, true, false>(
calculate_splines<order, atomsPerBlock, atomsPerWarp, true, false, useOrderThreads>(
kernelParams, atomIndexOffset, atomX, atomCharge, sm_theta, sm_dtheta, sm_gridlineIndices);
// __syncwarp();
__all(1);
Expand Down
8 changes: 4 additions & 4 deletions src/gromacs/ewald/pme_spread.hip.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -133,7 +133,7 @@ __device__ __forceinline__ void spread_charges(const PmeGpuHipKernelParams kerne
const int splineIndexY = getSplineParamIndex<order, atomsPerWarp>(splineIndexBase, YY, ithy);
float thetaY = sm_theta[splineIndexY];
const float Val = thetaZ * thetaY * (*atomCharge);
assert(isfinite(Val));
//assert(isfinite(Val));
const int offset = iy * pnz + iz;

#pragma unroll
Expand All @@ -148,8 +148,8 @@ __device__ __forceinline__ void spread_charges(const PmeGpuHipKernelParams kerne
const int splineIndexX =
getSplineParamIndex<order, atomsPerWarp>(splineIndexBase, XX, ithx);
const float thetaX = sm_theta[splineIndexX];
assert(isfinite(thetaX));
assert(isfinite(gm_grid[gridIndexGlobal]));
//assert(isfinite(thetaX));
//assert(isfinite(gm_grid[gridIndexGlobal]));
#if (HIP_VERSION_MAJOR >= 3) && (HIP_VERSION_MINOR > 3)
atomicAddNoRet(gm_grid + gridIndexGlobal, thetaX * Val);
#else
Expand Down Expand Up @@ -254,7 +254,7 @@ __launch_bounds__(c_spreadMaxThreadsPerBlock) CLANG_DISABLE_OPTIMIZATION_ATTRIBU
atomX.y = kernelParams.atoms.d_coordinates[atomIndexGlobal * DIM + YY];
atomX.z = kernelParams.atoms.d_coordinates[atomIndexGlobal * DIM + ZZ];
}
calculate_splines<order, atomsPerBlock, atomsPerWarp, false, writeGlobal>(
calculate_splines<order, atomsPerBlock, atomsPerWarp, false, writeGlobal, useOrderThreads>(
kernelParams, atomIndexOffset, atomX, atomCharge, sm_theta, &dtheta, sm_gridlineIndices);
// __syncwarp();
__all(1);
Expand Down