TestAmberPrmtopFile.py 11 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 *
6
import simtk.openmm.app.element as elem
7

8
9
prmtop1 = AmberPrmtopFile('systems/alanine-dipeptide-explicit.prmtop')
prmtop2 = AmberPrmtopFile('systems/alanine-dipeptide-implicit.prmtop')
10
11
prmtop3 = AmberPrmtopFile('systems/ff14ipq.parm7')
inpcrd3 = AmberInpcrdFile('systems/ff14ipq.rst7')
12

13
14
15
16
17
18
19
20
21
22
23
24
class TestAmberPrmtopFile(unittest.TestCase):

    """Test the AmberPrmtopFile.createSystem() method."""

    def test_NonbondedMethod(self):
        """Test all five options for the nonbondedMethod parameter."""

        methodMap = {NoCutoff:NonbondedForce.NoCutoff, 
                     CutoffNonPeriodic:NonbondedForce.CutoffNonPeriodic, 
                     CutoffPeriodic:NonbondedForce.CutoffPeriodic, 
                     Ewald:NonbondedForce.Ewald, PME: NonbondedForce.PME}
        for method in methodMap:
25
            system = prmtop1.createSystem(nonbondedMethod=method)
26
27
28
29
30
31
32
33
34
            forces = system.getForces()
            self.assertTrue(any(isinstance(f, NonbondedForce) and 
                                f.getNonbondedMethod()==methodMap[method] 
                                for f in forces))

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

        for method in [CutoffNonPeriodic, CutoffPeriodic, Ewald, PME]:
35
36
37
            system = prmtop1.createSystem(nonbondedMethod=method, 
                                          nonbondedCutoff=2*nanometer, 
                                          constraints=HBonds)
38
39
40
41
42
43
44
45
46
47
48
            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."""

        for method in [Ewald, PME]:
49
50
51
            system = prmtop1.createSystem(nonbondedMethod=method, 
                                          ewaldErrorTolerance=1e-6, 
                                          constraints=HBonds)
52
53
54
55
56
57
58
59
60
61
62
            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]:
63
            system = prmtop1.createSystem(removeCMMotion=b)
64
65
66
67
68
69
            forces = system.getForces()
            self.assertEqual(any(isinstance(f, CMMotionRemover) for f in forces), b)

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

70
        topology = prmtop1.topology
71
72
        for constraints_value in [None, HBonds, AllBonds, HAngles]:
            for rigidWater_value in [True, False]:
73
74
                system = prmtop1.createSystem(constraints=constraints_value, 
                                              rigidWater=rigidWater_value)
75
76
77
78
79
80
81
82
83
                validateConstraints(self, topology, system, 
                                    constraints_value, rigidWater_value)

    def test_ImplicitSolvent(self):
        """Test the four types of implicit solvents using the implicitSolvent 
        parameter.

        """
        for implicitSolvent_value in [HCT, OBC1, OBC2, GBn]:
84
            system = prmtop2.createSystem(implicitSolvent=implicitSolvent_value)
85
86
87
88
89
90
91
92
93
            forces = system.getForces()
            if implicitSolvent_value in set([HCT, OBC1, GBn]):
                force_type = CustomGBForce
            else:
                force_type = GBSAOBCForce
            
            self.assertTrue(any(isinstance(f, force_type) for f in forces))

    def test_ImplicitSolventParameters(self):
94
95
96
        """Test that parameters are set correctly for the different types of implicit solvent."""
        methodMap = {NoCutoff:NonbondedForce.NoCutoff,
                     CutoffNonPeriodic:NonbondedForce.CutoffNonPeriodic}
97
        for implicitSolvent_value in [HCT, OBC1, OBC2, GBn]:
98
            for method in methodMap:
99
                system = prmtop2.createSystem(implicitSolvent=implicitSolvent_value, 
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
                                    solventDielectric=50.0, soluteDielectric=0.9, nonbondedMethod=method)
                found_matching_solvent_dielectric=False
                found_matching_solute_dielectric=False
                if implicitSolvent_value in set([HCT, OBC1, GBn]):
                    for force in system.getForces():
                        if isinstance(force, CustomGBForce):
                            self.assertEqual(force.getNonbondedMethod(), methodMap[method])
                            for j in range(force.getNumGlobalParameters()):
                                if (force.getGlobalParameterName(j) == 'solventDielectric' and
                                   force.getGlobalParameterDefaultValue(j) == 50.0):
                                    found_matching_solvent_dielectric = True
                                if (force.getGlobalParameterName(j) == 'soluteDielectric' and
                                   force.getGlobalParameterDefaultValue(j) == 0.9):
                                    found_matching_solute_dielectric = True
                        if isinstance(force, NonbondedForce):
                            self.assertEqual(force.getReactionFieldDielectric(), 1.0)
                            self.assertEqual(force.getNonbondedMethod(), methodMap[method])
                    self.assertTrue(found_matching_solvent_dielectric and 
                                    found_matching_solute_dielectric)
                else:
                    for force in system.getForces():
                        if isinstance(force, GBSAOBCForce):
                            self.assertEqual(force.getNonbondedMethod(), methodMap[method])
                            if force.getSolventDielectric() == 50.0:
124
                                found_matching_solvent_dielectric = True
125
                            if force.getSoluteDielectric() == 0.9:
126
                                found_matching_solute_dielectric = True
127
128
129
130
131
                        if isinstance(force, NonbondedForce):
                            self.assertEqual(force.getReactionFieldDielectric(), 1.0)
                            self.assertEqual(force.getNonbondedMethod(), methodMap[method])
                    self.assertTrue(found_matching_solvent_dielectric and 
                                    found_matching_solute_dielectric)
132

133
134
135
    def test_HydrogenMass(self):
        """Test that altering the mass of hydrogens works correctly."""
        
136
        topology = prmtop1.topology
137
        hydrogenMass = 4*amu
138
139
        system1 = prmtop1.createSystem()
        system2 = prmtop1.createSystem(hydrogenMass=hydrogenMass)
140
141
142
143
144
145
146
147
        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)

148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
#   def test_NBFIX_LongRange(self):
#       """Test prmtop files with NBFIX LJ modifications w/ long-range correction"""
#       system = prmtop3.createSystem(nonbondedMethod=PME,
#                                     nonbondedCutoff=8*angstroms)
#       # Check the forces
#       has_nonbond_force = has_custom_nonbond_force = False
#       nonbond_exceptions = custom_nonbond_exclusions = 0
#       for force in system.getForces():
#           if isinstance(force, NonbondedForce):
#               has_nonbond_force = True
#               nonbond_exceptions = force.getNumExceptions()
#           elif isinstance(force, CustomNonbondedForce):
#               has_custom_nonbond_force = True
#               custom_nonbond_exceptions = force.getNumExclusions()
#       self.assertTrue(has_nonbond_force)
#       self.assertTrue(has_custom_nonbond_force)
#       self.assertEqual(nonbond_exceptions, custom_nonbond_exceptions)
#       integrator = VerletIntegrator(1.0*femtoseconds)
#       # Use reference platform, since it should always be present and
#       # 'working', and the system is plenty small so this won't be too slow
#       sim = Simulation(prmtop3.topology, system, integrator, Platform.getPlatformByName('Reference'))
#       # Check that the energy is about what we expect it to be
#       sim.context.setPeriodicBoxVectors(*inpcrd3.boxVectors)
#       sim.context.setPositions(inpcrd3.positions)
#       ene = sim.context.getState(getEnergy=True, enforcePeriodicBox=True).getPotentialEnergy()
#       ene = ene.value_in_unit(kilocalories_per_mole)
#       # Make sure the energy is relatively close to the value we get with
#       # Amber using this force field.
#       self.assertAlmostEqual(-7042.3903307/ene, 1, places=3)

    def test_NBFIX_noLongRange(self):
        """Test prmtop files with NBFIX LJ modifications w/out long-range correction"""
180
181
182
183
184
185
186
187
188
189
190
191
        system = prmtop3.createSystem(nonbondedMethod=PME,
                                      nonbondedCutoff=8*angstroms)
        # Check the forces
        has_nonbond_force = has_custom_nonbond_force = False
        nonbond_exceptions = custom_nonbond_exclusions = 0
        for force in system.getForces():
            if isinstance(force, NonbondedForce):
                has_nonbond_force = True
                nonbond_exceptions = force.getNumExceptions()
            elif isinstance(force, CustomNonbondedForce):
                has_custom_nonbond_force = True
                custom_nonbond_exceptions = force.getNumExclusions()
192
                force.setUseLongRangeCorrection(False)
193
194
195
196
        self.assertTrue(has_nonbond_force)
        self.assertTrue(has_custom_nonbond_force)
        self.assertEqual(nonbond_exceptions, custom_nonbond_exceptions)
        integrator = VerletIntegrator(1.0*femtoseconds)
197
198
199
        # Use reference platform, since it should always be present and
        # 'working', and the system is plenty small so this won't be too slow
        sim = Simulation(prmtop3.topology, system, integrator, Platform.getPlatformByName('Reference'))
200
201
202
203
204
205
206
207
208
        # Check that the energy is about what we expect it to be
        sim.context.setPeriodicBoxVectors(*inpcrd3.boxVectors)
        sim.context.setPositions(inpcrd3.positions)
        ene = sim.context.getState(getEnergy=True, enforcePeriodicBox=True).getPotentialEnergy()
        ene = ene.value_in_unit(kilocalories_per_mole)
        # Make sure the energy is relatively close to the value we get with
        # Amber using this force field.
        self.assertAlmostEqual(-7042.3903307/ene, 1, places=3)

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