Skip to content

Commit

Permalink
Added switchDistance option to ForceField.createSystem()
Browse files Browse the repository at this point in the history
  • Loading branch information
peastman committed Feb 21, 2018
1 parent d975668 commit fea6d6f
Show file tree
Hide file tree
Showing 3 changed files with 37 additions and 1 deletion.
10 changes: 10 additions & 0 deletions docs-source/usersguide/application.rst
Original file line number Diff line number Diff line change
Expand Up @@ -792,6 +792,16 @@ The error tolerance is roughly equal to the fractional error in the forces due
to truncating the Ewald summation. If you do not specify it, a default value of
0.0005 is used.

Another optional parameter when using a cutoff is :code:`switchDistance`. This
causes Lennard-Jones interactions to smoothly go to zero over some finite range,
rather than being sharply truncated at the cutoff distance. This can improve
energy conservation. To use it, specify a distance at which the interactions
should start being reduced. For example:
::

system = prmtop.createSystem(nonbondedMethod=PME, nonbondedCutoff=1*nanometer,
switchDistance=0.9*nanometer)


Nonbonded Forces for AMOEBA
---------------------------
Expand Down
13 changes: 12 additions & 1 deletion wrappers/python/simtk/openmm/app/forcefield.py
Original file line number Diff line number Diff line change
Expand Up @@ -1011,7 +1011,7 @@ def generateTemplatesForUnmatchedResidues(self, topology):

def createSystem(self, topology, nonbondedMethod=NoCutoff, nonbondedCutoff=1.0*unit.nanometer,
constraints=None, rigidWater=True, removeCMMotion=True, hydrogenMass=None, residueTemplates=dict(),
ignoreExternalBonds=False, **args):
ignoreExternalBonds=False, switchDistance=None, flexibleConstraints=False, **args):
"""Construct an OpenMM System representing a Topology with this force field.
Parameters
Expand Down Expand Up @@ -1047,6 +1047,9 @@ def createSystem(self, topology, nonbondedMethod=NoCutoff, nonbondedCutoff=1.0*u
not terminated properly. This option can create ambiguities where multiple
templates match the same residue. If that happens, use the residueTemplates
argument to specify which one to use.
switchDistance : float=None
The distance at which the potential energy switching function is turned on for
Lennard-Jones interactions. If this is None, no switching function will be used.
flexibleConstraints : boolean=False
If True, parameters for constrained degrees of freedom will be added to the System
args
Expand All @@ -1059,6 +1062,8 @@ def createSystem(self, topology, nonbondedMethod=NoCutoff, nonbondedCutoff=1.0*u
system
the newly created System
"""
args['switchDistance'] = switchDistance
args['flexibleConstraints'] = flexibleConstraints
data = ForceField._SystemData()
data.atoms = list(topology.atoms())
for atom in data.atoms:
Expand Down Expand Up @@ -2322,6 +2327,9 @@ def createForce(self, sys, data, nonbondedMethod, nonbondedCutoff, args):
force.addParticle(values[0], values[1], values[2])
force.setNonbondedMethod(methodMap[nonbondedMethod])
force.setCutoffDistance(nonbondedCutoff)
if args['switchDistance'] is not None:
force.setUseSwitchingFunction(True)
force.setSwitchingDistance(args['switchDistance'])
if 'ewaldErrorTolerance' in args:
force.setEwaldErrorTolerance(args['ewaldErrorTolerance'])
if 'useDispersionCorrection' in args:
Expand Down Expand Up @@ -2439,6 +2447,9 @@ def createForce(self, sys, data, nonbondedMethod, nonbondedCutoff, args):
self.force.setNonbondedMethod(mm.CustomNonbondedForce.CutoffNonPeriodic)
else:
raise AssertionError('Unrecognized nonbonded method [%s]' % nonbondedMethod)
if args['switchDistance'] is not None:
force.setUseSwitchingFunction(True)
force.setSwitchingDistance(args['switchDistance'])

# Add the particles

Expand Down
15 changes: 15 additions & 0 deletions wrappers/python/tests/TestForceField.py
Original file line number Diff line number Diff line change
Expand Up @@ -76,6 +76,21 @@ def test_Cutoff(self):
cutoff_distance = force.getCutoffDistance()
self.assertEqual(cutoff_distance, cutoff_check)

def test_SwitchingDistance(self):
"""Test that the switchDistance parameter is processed correctly."""

for switchDistance in [None, 0.9*nanometers]:
system = self.forcefield1.createSystem(self.pdb1.topology,
nonbondedMethod=PME,
switchDistance=switchDistance)
for force in system.getForces():
if isinstance(force, NonbondedForce):
if switchDistance is None:
self.assertFalse(force.getUseSwitchingFunction())
else:
self.assertTrue(force.getUseSwitchingFunction())
self.assertEqual(switchDistance, force.getSwitchingDistance())

def test_RemoveCMMotion(self):
"""Test both options (True and False) for the removeCMMotion parameter."""
for b in [True, False]:
Expand Down

0 comments on commit fea6d6f

Please sign in to comment.