Skip to content

Commit

Permalink
Changes for 8.1.1 (#4376)
Browse files Browse the repository at this point in the history
  • Loading branch information
peastman committed Jan 2, 2024
1 parent 43f571d commit ec797ac
Show file tree
Hide file tree
Showing 5 changed files with 19 additions and 2 deletions.
2 changes: 1 addition & 1 deletion CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -165,7 +165,7 @@ ENDIF (NOT CMAKE_CXX_FLAGS_RELEASE)
SET(OPENMM_LIBRARY_NAME OpenMM)
SET(OPENMM_MAJOR_VERSION 8)
SET(OPENMM_MINOR_VERSION 1)
SET(OPENMM_BUILD_VERSION 0)
SET(OPENMM_BUILD_VERSION 1)

ADD_DEFINITIONS(-DOPENMM_LIBRARY_NAME=${OPENMM_LIBRARY_NAME}
-DOPENMM_MAJOR_VERSION=${OPENMM_MAJOR_VERSION}
Expand Down
2 changes: 1 addition & 1 deletion platforms/common/include/openmm/common/CommonKernels.h
Original file line number Diff line number Diff line change
Expand Up @@ -1578,7 +1578,7 @@ class CommonApplyMonteCarloBarostatKernel : public ApplyMonteCarloBarostatKernel
ComputeContext& cc;
bool hasInitializedKernels, rigidMolecules, atomsWereReordered;
int numMolecules;
ComputeArray savedPositions, savedFloatForces, savedLongForces;
ComputeArray savedPositions, savedFloatForces, savedLongForces, savedVelocities;
ComputeArray moleculeAtoms;
ComputeArray moleculeStartIndex;
ComputeKernel kernel;
Expand Down
3 changes: 3 additions & 0 deletions platforms/common/src/CommonKernels.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -7699,6 +7699,7 @@ void CommonApplyMonteCarloBarostatKernel::initialize(const System& system, const
this->rigidMolecules = rigidMolecules;
ContextSelector selector(cc);
savedPositions.initialize(cc, cc.getPaddedNumAtoms(), cc.getUseDoublePrecision() ? sizeof(mm_double4) : sizeof(mm_float4), "savedPositions");
savedVelocities.initialize(cc, cc.getPaddedNumAtoms(), cc.getUseDoublePrecision() || cc.getUseMixedPrecision() ? sizeof(mm_double4) : sizeof(mm_float4), "savedVelocities");
savedLongForces.initialize<long long>(cc, cc.getPaddedNumAtoms()*3, "savedLongForces");
try {
cc.getFloatForceBuffer(); // This will throw an exception on the CUDA platform.
Expand All @@ -7714,6 +7715,7 @@ void CommonApplyMonteCarloBarostatKernel::initialize(const System& system, const
void CommonApplyMonteCarloBarostatKernel::saveCoordinates(ContextImpl& context) {
ContextSelector selector(cc);
cc.getPosq().copyTo(savedPositions);
cc.getVelm().copyTo(savedVelocities);
cc.getLongForceBuffer().copyTo(savedLongForces);
if (savedFloatForces.isInitialized())
cc.getFloatForceBuffer().copyTo(savedFloatForces);
Expand Down Expand Up @@ -7777,6 +7779,7 @@ void CommonApplyMonteCarloBarostatKernel::scaleCoordinates(ContextImpl& context,
void CommonApplyMonteCarloBarostatKernel::restoreCoordinates(ContextImpl& context) {
ContextSelector selector(cc);
savedPositions.copyTo(cc.getPosq());
savedVelocities.copyTo(cc.getVelm());
savedLongForces.copyTo(cc.getLongForceBuffer());
cc.setPosCellOffsets(lastPosCellOffsets);
if (savedFloatForces.isInitialized())
Expand Down
7 changes: 7 additions & 0 deletions platforms/cuda/src/kernels/findInteractingBlocks.cu
Original file line number Diff line number Diff line change
Expand Up @@ -327,6 +327,13 @@ extern "C" __global__ __launch_bounds__(GROUP_SIZE,3) void findBlocksWithInterac
blockDelta.y = max(0.0f, fabs(blockDelta.y)-blockSizeX.y-largeSize.y);
blockDelta.z = max(0.0f, fabs(blockDelta.z)-blockSizeX.z-largeSize.z);
includeLargeBlock = (blockDelta.x*blockDelta.x+blockDelta.y*blockDelta.y+blockDelta.z*blockDelta.z < PADDED_CUTOFF_SQUARED);
#ifdef TRICLINIC
// The calculation to find the nearest periodic copy is only guaranteed to work if the nearest copy is less than half a box width away.
// If there's any possibility we might have missed it, do a detailed check.

if (periodicBoxSize.z/2-blockSizeX.z-largeSize.z < PADDED_CUTOFF || periodicBoxSize.y/2-blockSizeX.y-largeSize.y < PADDED_CUTOFF)
includeLargeBlock = true;
#endif
}
largeBlockFlags = BALLOT(includeLargeBlock);
loadedLargeBlocks = 32;
Expand Down
7 changes: 7 additions & 0 deletions platforms/opencl/src/kernels/findInteractingBlocks.cl
Original file line number Diff line number Diff line change
Expand Up @@ -196,6 +196,13 @@ __kernel void findBlocksWithInteractions(real4 periodicBoxSize, real4 invPeriodi
blockDelta.y = max((real) 0, fabs(blockDelta.y)-blockSizeX.y-largeSize.y);
blockDelta.z = max((real) 0, fabs(blockDelta.z)-blockSizeX.z-largeSize.z);
includeLargeBlock = (blockDelta.x*blockDelta.x+blockDelta.y*blockDelta.y+blockDelta.z*blockDelta.z < PADDED_CUTOFF_SQUARED);
#ifdef TRICLINIC
// The calculation to find the nearest periodic copy is only guaranteed to work if the nearest copy is less than half a box width away.
// If there's any possibility we might have missed it, do a detailed check.

if (periodicBoxSize.z/2-blockSizeX.z-largeSize.z < PADDED_CUTOFF || periodicBoxSize.y/2-blockSizeX.y-largeSize.y < PADDED_CUTOFF)
includeLargeBlock = true;
#endif
}
largeBlockFlags[get_local_id(0)] = includeLargeBlock;
loadedLargeBlocks = 32;
Expand Down

0 comments on commit ec797ac

Please sign in to comment.