Skip to content

Commit

Permalink
Run BoxReciprocalSums on only charged particles
Browse files Browse the repository at this point in the history
  • Loading branch information
LSchwiebert committed Sep 1, 2024
1 parent 417ba24 commit 9f7a9c0
Show file tree
Hide file tree
Showing 2 changed files with 17 additions and 9 deletions.
20 changes: 12 additions & 8 deletions src/Ewald.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -213,10 +213,12 @@ void Ewald::BoxReciprocalSetup(uint box, XYZArray const &molCoords) {
while (thisMol != end) {
MoleculeKind const &thisKind = mols.GetKind(*thisMol);
double lambdaCoef = GetLambdaCoef(*thisMol, box);
for (uint j = 0; j < thisKind.NumAtoms(); j++) {
thisBoxCoords.Set(i, molCoords[mols.MolStart(*thisMol) + j]);
chargeBox.push_back(thisKind.AtomCharge(j) * lambdaCoef);
i++;
for (uint j = 0; j < thisKind.NumAtoms(); ++j) {
if (thisKind.AtomCharge(j) != 0.0) {
thisBoxCoords.Set(i, molCoords[mols.MolStart(*thisMol) + j]);
chargeBox.push_back(thisKind.AtomCharge(j) * lambdaCoef);
i++;
}
}
thisMol++;
}
Expand Down Expand Up @@ -301,10 +303,12 @@ void Ewald::BoxReciprocalSums(uint box, XYZArray const &molCoords) {
while (thisMol != end) {
MoleculeKind const &thisKind = mols.GetKind(*thisMol);
double lambdaCoef = GetLambdaCoef(*thisMol, box);
for (uint j = 0; j < thisKind.NumAtoms(); j++) {
thisBoxCoords.Set(i, molCoords[mols.MolStart(*thisMol) + j]);
chargeBox.push_back(thisKind.AtomCharge(j) * lambdaCoef);
i++;
for (uint j = 0; j < thisKind.NumAtoms(); ++j) {
if (thisKind.AtomCharge(j) != 0.0) {
thisBoxCoords.Set(i, molCoords[mols.MolStart(*thisMol) + j]);
chargeBox.push_back(thisKind.AtomCharge(j) * lambdaCoef);
i++;
}
}
thisMol++;
}
Expand Down
6 changes: 5 additions & 1 deletion src/GPU/CalculateEwaldCUDAKernel.cu
Original file line number Diff line number Diff line change
Expand Up @@ -66,6 +66,8 @@ void CallBoxReciprocalSetupGPU(VariablesCUDA *vars, XYZArray const &coords,
checkLastErrorCUDA(__FILE__, __LINE__);
#endif

// Fewer blocks are needed since we are doing one computation per image
blocksPerGrid = (imageSize + THREADS_PER_BLOCK - 1) / THREADS_PER_BLOCK;
BoxReciprocalGPU<<<blocksPerGrid, threadsPerBlock>>>(
vars->gpu_prefact[box], vars->gpu_sumRnew[box], vars->gpu_sumInew[box],
vars->gpu_recipEnergies, imageSize);
Expand Down Expand Up @@ -111,6 +113,8 @@ void CallBoxReciprocalSumsGPU(VariablesCUDA *vars, XYZArray const &coords,
checkLastErrorCUDA(__FILE__, __LINE__);
#endif

// Fewer blocks are needed since we are doing one computation per image
blocksPerGrid = (imageSize + THREADS_PER_BLOCK - 1) / THREADS_PER_BLOCK;
BoxReciprocalGPU<<<blocksPerGrid, threadsPerBlock>>>(
vars->gpu_prefactRef[box], vars->gpu_sumRnew[box], vars->gpu_sumInew[box],
vars->gpu_recipEnergies, imageSize);
Expand All @@ -133,7 +137,7 @@ __global__ void BoxReciprocalSumsGPU(double *gpu_x, double *gpu_y,
double *gpu_sumRnew, double *gpu_sumInew) {
int image = blockIdx.x;
double sumR = 0.0, sumI = 0.0;
#pragma unroll 8
#pragma unroll 16
for (int particleID = threadIdx.x; particleID < atomNumber; particleID += THREADS_PER_BLOCK) {
double dot = DotProductGPU(gpu_kx[image], gpu_ky[image],
gpu_kz[image], gpu_x[particleID],
Expand Down

0 comments on commit 9f7a9c0

Please sign in to comment.