"wrappers/python/vscode:/vscode.git/clone" did not exist on "baab6f40376faeeffe7c4e57f9cf57094f3c3388"
TestGromacsTopFile.py 6.6 KB
Newer Older
1
2
3
4
5
import unittest
from validateConstraints import *
from simtk.openmm.app import *
from simtk.openmm import *
from simtk.unit import *
Robert McGibbon's avatar
Robert McGibbon committed
6
from simtk.openmm.app.gromacstopfile import _defaultGromacsIncludeDir
7
import simtk.openmm.app.element as elem
8

Robert McGibbon's avatar
Robert McGibbon committed
9
10
11
GROMACS_INCLUDE = _defaultGromacsIncludeDir()

@unittest.skipIf(not os.path.exists(GROMACS_INCLUDE), 'GROMACS is not installed')
12
13
14
class TestGromacsTopFile(unittest.TestCase):

    """Test the GromacsTopFile.createSystem() method."""
Robert McGibbon's avatar
Robert McGibbon committed
15

16
17
18
19
20
21
22
23
24
    def setUp(self):
        """Set up the tests by loading the input files."""

        # alanine dipeptide with explicit water
        self.top1 = GromacsTopFile('systems/explicit.top', unitCellDimensions=Vec3(6.223, 6.223, 6.223)*nanometers)
        # alanine dipeptide with implicit water
        self.top2 = GromacsTopFile('systems/implicit.top')

    def test_NonbondedMethod(self):
Peter Eastman's avatar
Peter Eastman committed
25
        """Test all six options for the nonbondedMethod parameter."""
26

Robert McGibbon's avatar
Robert McGibbon committed
27
28
29
        methodMap = {NoCutoff:NonbondedForce.NoCutoff,
                     CutoffNonPeriodic:NonbondedForce.CutoffNonPeriodic,
                     CutoffPeriodic:NonbondedForce.CutoffPeriodic,
Peter Eastman's avatar
Peter Eastman committed
30
31
32
                     Ewald:NonbondedForce.Ewald,
                     PME:NonbondedForce.PME,
                     LJPME:NonbondedForce.LJPME}
33
34
35
        for method in methodMap:
            system = self.top1.createSystem(nonbondedMethod=method)
            forces = system.getForces()
Robert McGibbon's avatar
Robert McGibbon committed
36
37
            self.assertTrue(any(isinstance(f, NonbondedForce) and
                                f.getNonbondedMethod()==methodMap[method]
38
39
                                for f in forces))

40
41
42
43
44
45
46
47
48
49
50
51
52
53
    def test_ff99SBILDN(self):
        """ Test Gromacs topology #define replacement as used in ff99SB-ILDN """
        top = GromacsTopFile('systems/aidilnaaaaa.top')
        gro = GromacsGroFile('systems/aidilnaaaaa.gro')
        system = top.createSystem()
        for force in system.getForces():
            if isinstance(force, PeriodicTorsionForce):
                force.setForceGroup(1)
        context = Context(system, VerletIntegrator(1*femtosecond),
                          Platform.getPlatformByName('Reference'))
        context.setPositions(gro.positions)
        ene = context.getState(getEnergy=True, groups=2**1).getPotentialEnergy()
        self.assertAlmostEqual(ene.value_in_unit(kilojoules_per_mole), 341.6905133582857)

54
55
56
    def test_Cutoff(self):
        """Test to make sure the nonbondedCutoff parameter is passed correctly."""

Peter Eastman's avatar
Peter Eastman committed
57
        for method in [CutoffNonPeriodic, CutoffPeriodic, Ewald, PME, LJPME]:
Robert McGibbon's avatar
Robert McGibbon committed
58
59
            system = self.top1.createSystem(nonbondedMethod=method,
                                               nonbondedCutoff=2*nanometer,
60
61
62
63
64
65
66
67
68
69
70
                                               constraints=HBonds)
            cutoff_distance = 0.0*nanometer
            cutoff_check = 2.0*nanometer
            for force in system.getForces():
                if isinstance(force, NonbondedForce):
                    cutoff_distance = force.getCutoffDistance()
            self.assertEqual(cutoff_distance, cutoff_check)

    def test_EwaldErrorTolerance(self):
        """Test to make sure the ewaldErrorTolerance parameter is passed correctly."""

Peter Eastman's avatar
Peter Eastman committed
71
        for method in [Ewald, PME, LJPME]:
Robert McGibbon's avatar
Robert McGibbon committed
72
73
            system = self.top1.createSystem(nonbondedMethod=method,
                                               ewaldErrorTolerance=1e-6,
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
                                               constraints=HBonds)
            tolerance = 0
            tolerance_check = 1e-6
            for force in system.getForces():
                if isinstance(force, NonbondedForce):
                    tolerance = force.getEwaldErrorTolerance()
            self.assertEqual(tolerance, tolerance_check)

    def test_RemoveCMMotion(self):
        """Test both options (True and False) for the removeCMMotion parameter."""

        for b in [True, False]:
            system = self.top1.createSystem(removeCMMotion=b)
            self.assertEqual(any(isinstance(f, CMMotionRemover) for f in system.getForces()), b)

    def test_RigidWaterAndConstraints(self):
        """Test all eight options for the constraints and rigidWater parameters."""

        topology = self.top1.topology
        for constraints_value in [None, HBonds, AllBonds, HAngles]:
            for rigidWater_value in [True, False]:
Robert McGibbon's avatar
Robert McGibbon committed
95
                system = self.top1.createSystem(constraints=constraints_value,
96
                                                   rigidWater=rigidWater_value)
Robert McGibbon's avatar
Robert McGibbon committed
97
                validateConstraints(self, topology, system,
98
99
100
101
102
103
104
105
106
107
108
109
110
                                    constraints_value, rigidWater_value)

    def test_ImplicitSolvent(self):
        """Test implicit solvent using the implicitSolvent parameter.

        """
        system = self.top2.createSystem(implicitSolvent=OBC2)
        self.assertTrue(any(isinstance(f, GBSAOBCForce) for f in system.getForces()))

    def test_ImplicitSolventParameters(self):
        """Test that solventDielectric and soluteDielectric are passed correctly.

        """
Robert McGibbon's avatar
Robert McGibbon committed
111
112
        system = self.top2.createSystem(implicitSolvent=OBC2,
                                           solventDielectric=50.0,
113
114
115
116
117
118
119
120
121
122
123
                                           soluteDielectric = 0.9)
        found_matching_solvent_dielectric=False
        found_matching_solute_dielectric=False
        for force in system.getForces():
            if isinstance(force, GBSAOBCForce):
                if force.getSolventDielectric() == 50.0:
                    found_matching_solvent_dielectric = True
                if force.getSoluteDielectric() == 0.9:
                    found_matching_solute_dielectric = True
            if isinstance(force, NonbondedForce):
                self.assertEqual(force.getReactionFieldDielectric(), 1.0)
Robert McGibbon's avatar
Robert McGibbon committed
124
        self.assertTrue(found_matching_solvent_dielectric and
125
126
                        found_matching_solute_dielectric)

127
128
    def test_HydrogenMass(self):
        """Test that altering the mass of hydrogens works correctly."""
Robert McGibbon's avatar
Robert McGibbon committed
129

130
131
132
133
134
135
136
137
138
139
140
141
        topology = self.top1.topology
        hydrogenMass = 4*amu
        system1 = self.top1.createSystem()
        system2 = self.top1.createSystem(hydrogenMass=hydrogenMass)
        for atom in topology.atoms():
            if atom.element == elem.hydrogen:
                self.assertNotEqual(hydrogenMass, system1.getParticleMass(atom.index))
                self.assertEqual(hydrogenMass, system2.getParticleMass(atom.index))
        totalMass1 = sum([system1.getParticleMass(i) for i in range(system1.getNumParticles())]).value_in_unit(amu)
        totalMass2 = sum([system2.getParticleMass(i) for i in range(system2.getNumParticles())]).value_in_unit(amu)
        self.assertAlmostEqual(totalMass1, totalMass2)

142
143
144
if __name__ == '__main__':
    unittest.main()