Skip to content

Commit

Permalink
Merge pull request #47 from raimis/test
Browse files Browse the repository at this point in the history
Add a unit test for the Python wrapper
  • Loading branch information
Raimondas Galvelis authored Feb 3, 2021
2 parents edff9ea + 83e5dda commit db390e2
Showing 1 changed file with 43 additions and 0 deletions.
43 changes: 43 additions & 0 deletions python/tests/TestPlumedForce.py
Original file line number Diff line number Diff line change
@@ -0,0 +1,43 @@
import simtk.openmm as mm
import simtk.unit as unit
from openmmplumed import PlumedForce
import numpy as np
import unittest


class TestTorchForce(unittest.TestCase):

def testForce(self):
# Create a System that applies a force based on the distance between two atoms.

numParticles = 4
system = mm.System()
positions = np.empty((numParticles, 3))
for i in range(numParticles):
system.addParticle(1.0)
positions[i] = [i, 0.1*i, -0.3*i]
script = '''
d: DISTANCE ATOMS=1,3
BIASVALUE ARG=d
'''
force = PlumedForce(script)
system.addForce(force)
integ = mm.LangevinIntegrator(300.0, 1.0, 1.0)
context = mm.Context(system, integ, mm.Platform.getPlatformByName('Reference'))
context.setPositions(positions)

# Compute the forces and energy.

state = context.getState(getEnergy=True, getForces=True)
delta = positions[0] - positions[2]
dist = np.sqrt(np.sum(delta**2))
zero = np.zeros(3)
self.assertAlmostEqual(dist, state.getPotentialEnergy().value_in_unit(unit.kilojoules_per_mole))
self.assertTrue(np.allclose(-delta/dist, state.getForces(asNumpy=True)[0]))
self.assertTrue(np.allclose(zero, state.getForces(asNumpy=True)[1]))
self.assertTrue(np.allclose(delta/dist, state.getForces(asNumpy=True)[2]))
self.assertTrue(np.allclose(zero, state.getForces(asNumpy=True)[3]))


if __name__ == '__main__':
unittest.main()

0 comments on commit db390e2

Please sign in to comment.