File: TestForceField.py

package info (click to toggle)
openmm 8.1.2%2Bdfsg-11
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 119,168 kB
  • sloc: xml: 377,325; cpp: 226,673; ansic: 42,767; python: 32,634; lisp: 2,441; sh: 440; makefile: 254; f90: 233; csh: 19
file content (1484 lines) | stat: -rw-r--r-- 69,424 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
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
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
import unittest
from validateConstraints import *
from openmm.app import *
from openmm import *
from openmm.unit import *
import openmm.app.element as elem
import openmm.app.forcefield as forcefield
import math
import textwrap
try:
    from cStringIO import StringIO
except ImportError:
    from io import StringIO
import os
import warnings

class TestForceField(unittest.TestCase):
    """Test the ForceField.createSystem() method."""

    def setUp(self):
        """Set up the tests by loading the input pdb files and force field
        xml files.

        """
        # alanine dipeptide with explicit water
        self.pdb1 = PDBFile('systems/alanine-dipeptide-explicit.pdb')
        self.forcefield1 = ForceField('amber99sb.xml', 'tip3p.xml')
        self.topology1 = self.pdb1.topology
        self.topology1.setUnitCellDimensions(Vec3(2, 2, 2))

        # alanine dipeptide with implicit water
        self.pdb2 = PDBFile('systems/alanine-dipeptide-implicit.pdb')
        self.forcefield2 = ForceField('amber99sb.xml', 'amber99_obc.xml')


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

        methodMap = {NoCutoff:NonbondedForce.NoCutoff,
                     CutoffNonPeriodic:NonbondedForce.CutoffNonPeriodic,
                     CutoffPeriodic:NonbondedForce.CutoffPeriodic,
                     Ewald:NonbondedForce.Ewald,
                     PME:NonbondedForce.PME,
                     LJPME:NonbondedForce.LJPME}
        for method in methodMap:
            system = self.forcefield1.createSystem(self.pdb1.topology,
                                                  nonbondedMethod=method)
            forces = system.getForces()
            self.assertTrue(any(isinstance(f, NonbondedForce) and
                                f.getNonbondedMethod()==methodMap[method]
                                for f in forces))

    def test_DispersionCorrection(self):
        """Test to make sure that the dispersion/long-range correction is set properly."""
        top = Topology()
        chain = top.addChain()

        for lrc in (True, False):
            xml = textwrap.dedent(
                """
                <ForceField>
                 <LennardJonesForce lj14scale="0.3" useDispersionCorrection="{lrc}">
                  <Atom type="A" sigma="1" epsilon="0.1"/>
                  <Atom type="B" sigma="2" epsilon="0.2"/>
                  <NBFixPair type1="A" type2="B" sigma="2.5" epsilon="1.1"/>
                 </LennardJonesForce>
                 <NonbondedForce coulomb14scale="0.833333" lj14scale="0.5" useDispersionCorrection="{lrc2}">
                  <Atom type="A" sigma="0.315" epsilon="0.635"/>
                 </NonbondedForce>
                </ForceField>
                """
            )
            ff = ForceField(StringIO(xml.format(lrc=lrc, lrc2=lrc)))
            system = ff.createSystem(top)
            checked_nonbonded = False
            checked_custom = False
            for force in system.getForces():
                if isinstance(force, NonbondedForce):
                    self.assertEqual(force.getUseDispersionCorrection(), lrc)
                    checked_nonbonded = True
                elif isinstance(force, CustomNonbondedForce):
                    self.assertEqual(force.getUseLongRangeCorrection(), lrc)
                    checked_custom = True
            self.assertTrue(checked_nonbonded and checked_custom)

            # check that the keyword argument overwrites xml input
            lrc_kwarg = not lrc
            with warnings.catch_warnings(record=True) as w:
                warnings.simplefilter("always")
                system2 = ff.createSystem(top, useDispersionCorrection=lrc_kwarg)
                self.assertTrue(len(w) == 2)
                assert "conflict" in str(w[-1].message).lower()
            checked_nonbonded = False
            checked_custom = False
            for force in system2.getForces():
                if isinstance(force, NonbondedForce):
                    self.assertEqual(force.getUseDispersionCorrection(), lrc_kwarg)
                    checked_nonbonded = True
                elif isinstance(force, CustomNonbondedForce):
                    self.assertEqual(force.getUseLongRangeCorrection(), lrc_kwarg)
                    checked_custom = True
            self.assertTrue(checked_nonbonded and checked_custom)

            # check that no warning is generated when useDispersionCorrection is not in the xml file
            xml = textwrap.dedent(
                """
                <ForceField>
                 <LennardJonesForce lj14scale="0.3">
                  <Atom type="A" sigma="1" epsilon="0.1"/>
                  <Atom type="B" sigma="2" epsilon="0.2"/>
                  <NBFixPair type1="A" type2="B" sigma="2.5" epsilon="1.1"/>
                 </LennardJonesForce>
                 <NonbondedForce coulomb14scale="0.833333" lj14scale="0.5">
                  <Atom type="A" sigma="0.315" epsilon="0.635"/>
                 </NonbondedForce>
                </ForceField>
                """
            )
            ff = ForceField(StringIO(xml))
            system = ff.createSystem(top)
            for lrc_kwarg in [True, False]:
                with warnings.catch_warnings():
                    warnings.simplefilter("error")
                    system2 = ff.createSystem(top, useDispersionCorrection=lrc_kwarg)

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

        for method in [CutoffNonPeriodic, CutoffPeriodic, Ewald, PME, LJPME]:
            system = self.forcefield1.createSystem(self.pdb1.topology,
                                                   nonbondedMethod=method,
                                                   nonbondedCutoff=2*nanometer,
                                                   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_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]:
            system = self.forcefield1.createSystem(self.pdb1.topology,removeCMMotion=b)
            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."""

        topology = self.pdb1.topology
        for constraints_value in [None, HBonds, AllBonds, HAngles]:
            for rigidWater_value in [True, False, None]:
                system = self.forcefield1.createSystem(topology,
                                                       constraints=constraints_value,
                                                       rigidWater=rigidWater_value)
                validateConstraints(self, topology, system,
                                    constraints_value, rigidWater_value != False)

    def test_flexibleConstraints(self):
        """ Test the flexibleConstraints keyword """
        topology = self.pdb1.topology
        system1 = self.forcefield1.createSystem(topology, constraints=HAngles,
                                                rigidWater=True)
        system2 = self.forcefield1.createSystem(topology, constraints=HAngles,
                                                rigidWater=True, flexibleConstraints=True)
        system3 = self.forcefield1.createSystem(topology, constraints=None, rigidWater=False)
        validateConstraints(self, topology, system1, HAngles, True)
        # validateConstraints fails for system2 since by definition atom pairs can be in both bond
        # and constraint lists. So just check that the number of constraints is the same for both
        # system1 and system2
        self.assertEqual(system1.getNumConstraints(), system2.getNumConstraints())
        for force in system1.getForces():
            if isinstance(force, HarmonicBondForce):
                bf1 = force
            elif isinstance(force, HarmonicAngleForce):
                af1 = force
        for force in system2.getForces():
            if isinstance(force, HarmonicBondForce):
                bf2 = force
            elif isinstance(force, HarmonicAngleForce):
                af2 = force
        for force in system3.getForces():
            if isinstance(force, HarmonicAngleForce):
                af3 = force
        # Make sure we picked up extra bond terms with flexibleConstraints
        self.assertGreater(bf2.getNumBonds(), bf1.getNumBonds())
        # Make sure flexibleConstraints yields just as many angles as no constraints
        self.assertEqual(af2.getNumAngles(), af3.getNumAngles())

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

        """

        topology = self.pdb2.topology
        system = self.forcefield2.createSystem(topology)
        forces = system.getForces()
        self.assertTrue(any(isinstance(f, GBSAOBCForce) for f in forces))

    def test_ImplicitSolventParameters(self):
        """Test that solventDielectric and soluteDielectric are passed correctly
        for the different types of implicit solvent.

        """

        topology = self.pdb2.topology
        system = self.forcefield2.createSystem(topology, solventDielectric=50.0,
                                               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)
        self.assertTrue(found_matching_solvent_dielectric and
                        found_matching_solute_dielectric)

    def test_HydrogenMass(self):
        """Test that altering the mass of hydrogens works correctly."""

        topology = self.pdb1.topology
        hydrogenMass = 4*amu
        system1 = self.forcefield1.createSystem(topology)
        system2 = self.forcefield1.createSystem(topology, hydrogenMass=hydrogenMass)
        for atom in topology.atoms():
            if atom.element == elem.hydrogen:
                self.assertNotEqual(hydrogenMass, system1.getParticleMass(atom.index))
                if atom.residue.name == 'HOH':
                    self.assertEqual(system1.getParticleMass(atom.index), system2.getParticleMass(atom.index))
                else:
                    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)

    def test_DrudeMass(self):
        """Test that setting the mass of Drude particles works correctly."""

        forcefield = ForceField('charmm_polar_2013.xml')
        pdb = PDBFile('systems/ala_ala_ala.pdb')
        modeller = Modeller(pdb.topology, pdb.positions)
        modeller.addExtraParticles(forcefield)
        system = forcefield.createSystem(modeller.topology, drudeMass=0)
        trueMass = [system.getParticleMass(i) for i in range(system.getNumParticles())]
        drudeMass = 0.3*amu
        system = forcefield.createSystem(modeller.topology, drudeMass=drudeMass)
        adjustedMass = [system.getParticleMass(i) for i in range(system.getNumParticles())]
        drudeForce = [f for f in system.getForces() if isinstance(f, DrudeForce)][0]
        drudeParticles = set()
        parentParticles = set()
        for i in range(drudeForce.getNumParticles()):
            params = drudeForce.getParticleParameters(i)
            drudeParticles.add(params[0])
            parentParticles.add(params[1])
        for i in range(system.getNumParticles()):
            if i in drudeParticles:
                self.assertEqual(0*amu, trueMass[i])
                self.assertEqual(drudeMass, adjustedMass[i])
            elif i in parentParticles:
                self.assertEqual(trueMass[i]-drudeMass, adjustedMass[i])
            else:
                self.assertEqual(trueMass[i], adjustedMass[i])

    def test_UnusedArgs(self):
        """Test that specifying an argument that is never used throws an exception."""
        topology = self.pdb1.topology
        # Using the default value should not raise an exception.
        self.forcefield1.createSystem(topology, drudeMass=0.4*amu)
        # Specifying a non-default value should.
        with self.assertRaises(ValueError):
            self.forcefield1.createSystem(topology, drudeMass=0.5*amu)
        # Specifying a nonexistant argument should raise an exception.
        with self.assertRaises(ValueError):
            self.forcefield1.createSystem(topology, nonbndedCutoff=1.0*nanometer)

    def test_Forces(self):
        """Compute forces and compare them to ones generated with a previous version of OpenMM to ensure they haven't changed."""

        pdb = PDBFile('systems/lysozyme-implicit.pdb')
        system = self.forcefield2.createSystem(pdb.topology)
        integrator = VerletIntegrator(0.001)
        context = Context(system, integrator)
        context.setPositions(pdb.positions)
        state1 = context.getState(getForces=True)
        with open('systems/lysozyme-implicit-forces.xml') as input:
            state2 = XmlSerializer.deserialize(input.read())
        numDifferences = 0
        for f1, f2, in zip(state1.getForces().value_in_unit(kilojoules_per_mole/nanometer), state2.getForces().value_in_unit(kilojoules_per_mole/nanometer)):
            diff = norm(f1-f2)
            if diff > 0.1 and diff/norm(f1) > 1e-3:
                numDifferences += 1
        self.assertTrue(numDifferences < system.getNumParticles()/20) # Tolerate occasional differences from numerical error

    def test_ImplicitSolventForces(self):
        """Compute forces for different implicit solvent types, and compare them to ones generated with AmberPrmtopFile."""

        solventType = ['hct', 'obc1', 'obc2', 'gbn', 'gbn2']
        nonbondedMethod = [NoCutoff, CutoffNonPeriodic, CutoffNonPeriodic, NoCutoff, NoCutoff]
        kappa = [0.0, 0.0, 1.698295227342757, 1.698295227342757, 0.0]
        file = [None, 'OBC1_NonPeriodic', 'OBC2_NonPeriodic_Salt', None, 'GBn2_NoCutoff']
        for i in range(len(file)):
            forcefield = ForceField('amber96.xml', f'implicit/{solventType[i]}.xml')
            system = forcefield.createSystem(self.pdb2.topology, nonbondedMethod=nonbondedMethod[i], implicitSolventKappa=kappa[i])
            integrator = VerletIntegrator(0.001)
            context = Context(system, integrator, Platform.getPlatformByName("Reference"))
            context.setPositions(self.pdb2.positions)
            state1 = context.getState(getForces=True)
            if file[i] is not None:
                with open('systems/alanine-dipeptide-implicit-forces/'+file[i]+'.xml') as infile:
                    state2 = XmlSerializer.deserialize(infile.read())
                for f1, f2, in zip(state1.getForces().value_in_unit(kilojoules_per_mole/nanometer), state2.getForces().value_in_unit(kilojoules_per_mole/nanometer)):
                    diff = norm(f1-f2)
                    self.assertTrue(diff < 0.1 or diff/norm(f1) < 1e-4)

    def test_ProgrammaticForceField(self):
        """Test building a ForceField programmatically."""

        # Build the ForceField for TIP3P programmatically.
        ff = ForceField()
        ff.registerAtomType({'name':'tip3p-O', 'class':'OW', 'mass':15.99943*daltons, 'element':elem.oxygen})
        ff.registerAtomType({'name':'tip3p-H', 'class':'HW', 'mass':1.007947*daltons, 'element':elem.hydrogen})
        residue = ForceField._TemplateData('HOH')
        residue.atoms.append(ForceField._TemplateAtomData('O', 'tip3p-O', elem.oxygen))
        residue.atoms.append(ForceField._TemplateAtomData('H1', 'tip3p-H', elem.hydrogen))
        residue.atoms.append(ForceField._TemplateAtomData('H2', 'tip3p-H', elem.hydrogen))
        residue.addBond(0, 1)
        residue.addBond(0, 2)
        ff.registerResidueTemplate(residue)
        bonds = forcefield.HarmonicBondGenerator(ff)
        bonds.registerBond({'class1':'OW', 'class2':'HW', 'length':0.09572*nanometers, 'k':462750.4*kilojoules_per_mole/nanometer})
        ff.registerGenerator(bonds)
        angles = forcefield.HarmonicAngleGenerator(ff)
        angles.registerAngle({'class1':'HW', 'class2':'OW', 'class3':'HW', 'angle':1.82421813418*radians, 'k':836.8*kilojoules_per_mole/radian})
        ff.registerGenerator(angles)
        nonbonded = forcefield.NonbondedGenerator(ff, 0.833333, 0.5, True)
        nonbonded.registerAtom({'type':'tip3p-O', 'charge':-0.834, 'sigma':0.31507524065751241*nanometers, 'epsilon':0.635968*kilojoules_per_mole})
        nonbonded.registerAtom({'type':'tip3p-H', 'charge':0.417, 'sigma':1*nanometers, 'epsilon':0*kilojoules_per_mole})
        ff.registerGenerator(nonbonded)

        # Build a water box.
        modeller = Modeller(Topology(), [])
        modeller.addSolvent(ff, boxSize=Vec3(3, 3, 3)*nanometers)

        # Create a system using the programmatic force field as well as one from an XML file.
        system1 = ff.createSystem(modeller.topology)
        ff2 = ForceField('tip3p.xml')
        system2 = ff2.createSystem(modeller.topology)
        self.assertEqual(XmlSerializer.serialize(system1), XmlSerializer.serialize(system2))

    def test_PeriodicBoxVectors(self):
        """Test setting the periodic box vectors."""

        vectors = (Vec3(5, 0, 0), Vec3(-1.5, 4.5, 0), Vec3(0.4, 0.8, 7.5))*nanometers
        self.pdb1.topology.setPeriodicBoxVectors(vectors)
        self.assertEqual(Vec3(5, 4.5, 7.5)*nanometers, self.pdb1.topology.getUnitCellDimensions())
        system = self.forcefield1.createSystem(self.pdb1.topology)
        for i in range(3):
            self.assertEqual(vectors[i], self.pdb1.topology.getPeriodicBoxVectors()[i])
            self.assertEqual(vectors[i], system.getDefaultPeriodicBoxVectors()[i])

    def test_ResidueAttributes(self):
        """Test a ForceField that gets per-particle parameters from residue attributes."""

        xml = """
<ForceField>
 <AtomTypes>
  <Type name="tip3p-O" class="OW" element="O" mass="15.99943"/>
  <Type name="tip3p-H" class="HW" element="H" mass="1.007947"/>
 </AtomTypes>
 <Residues>
  <Residue name="HOH">
   <Atom name="O" type="tip3p-O" charge="-0.834"/>
   <Atom name="H1" type="tip3p-H" charge="0.417"/>
   <Atom name="H2" type="tip3p-H" charge="0.417"/>
   <Bond from="0" to="1"/>
   <Bond from="0" to="2"/>
  </Residue>
 </Residues>
 <NonbondedForce coulomb14scale="0.833333" lj14scale="0.5">
  <UseAttributeFromResidue name="charge"/>
  <Atom type="tip3p-O" sigma="0.315" epsilon="0.635"/>
  <Atom type="tip3p-H" sigma="1" epsilon="0"/>
 </NonbondedForce>
</ForceField>"""
        ff = ForceField(StringIO(xml))

        # Build a water box.
        modeller = Modeller(Topology(), [])
        modeller.addSolvent(ff, boxSize=Vec3(3, 3, 3)*nanometers)

        # Create a system and make sure all nonbonded parameters are correct.
        system = ff.createSystem(modeller.topology)
        nonbonded = [f for f in system.getForces() if isinstance(f, NonbondedForce)][0]
        atoms = list(modeller.topology.atoms())
        for i in range(len(atoms)):
            params = nonbonded.getParticleParameters(i)
            if atoms[i].element == elem.oxygen:
                self.assertEqual(params[0], -0.834*elementary_charge)
                self.assertEqual(params[1], 0.315*nanometers)
                self.assertEqual(params[2], 0.635*kilojoule_per_mole)
            else:
                self.assertEqual(params[0], 0.417*elementary_charge)
                self.assertEqual(params[1], 1.0*nanometers)
                self.assertEqual(params[2], 0.0*kilojoule_per_mole)

    def test_residueMatcher(self):
        """Test using a custom template matcher to select templates."""
        xml = """
<ForceField>
 <AtomTypes>
  <Type name="tip3p-O" class="OW" element="O" mass="15.99943"/>
  <Type name="tip3p-H" class="HW" element="H" mass="1.007947"/>
 </AtomTypes>
 <Residues>
  <Residue name="HOH">
   <Atom name="O" type="tip3p-O" charge="-0.834"/>
   <Atom name="H1" type="tip3p-H" charge="0.417"/>
   <Atom name="H2" type="tip3p-H" charge="0.417"/>
   <Bond from="0" to="1"/>
   <Bond from="0" to="2"/>
  </Residue>
  <Residue name="HOH2">
   <Atom name="O" type="tip3p-O" charge="0.834"/>
   <Atom name="H1" type="tip3p-H" charge="-0.417"/>
   <Atom name="H2" type="tip3p-H" charge="-0.417"/>
   <Bond from="0" to="1"/>
   <Bond from="0" to="2"/>
  </Residue>
 </Residues>
 <NonbondedForce coulomb14scale="0.833333" lj14scale="0.5">
  <UseAttributeFromResidue name="charge"/>
  <Atom type="tip3p-O" sigma="0.315" epsilon="0.635"/>
  <Atom type="tip3p-H" sigma="1" epsilon="0"/>
 </NonbondedForce>
</ForceField>"""
        ff = ForceField(StringIO(xml))

        # Load a water box.
        prmtop = AmberPrmtopFile('systems/water-box-216.prmtop')
        top = prmtop.topology
        
        # Building a System should fail, because two templates match each residue.
        self.assertRaises(Exception, lambda: ff.createSystem(top))
        
        # Register a template matcher that selects a particular one.
        def matcher(ff, res, bondedToAtom, ignoreExternalBonds, ignoreExtraParticles):
            return ff._templates['HOH2']
        ff.registerTemplateMatcher(matcher)
        
        # It should now succeed in building a System.
        system = ff.createSystem(top)
        
        # Make sure it used the correct parameters.
        nb = [f for f in system.getForces() if isinstance(f, NonbondedForce)][0]
        for atom in top.atoms():
            charge, sigma, epsilon = nb.getParticleParameters(atom.index)
            if atom.name == 'O':
                self.assertEqual(0.834*elementary_charge, charge)
            else:
                self.assertEqual(-0.417*elementary_charge, charge)

    def test_residueTemplateGenerator(self):
        """Test the ability to add residue template generators to parameterize unmatched residues."""
        def simpleTemplateGenerator(forcefield, residue):
            """\
            Simple residue template generator.
            This implementation uses the programmatic API to define residue templates.

            NOTE: We presume we have already loaded the force definitions into ForceField.
            """
            # Generate a unique prefix name for generating parameters.
            from uuid import uuid4
            template_name = uuid4()
            # Create residue template.
            from openmm.app.forcefield import _createResidueTemplate
            template = _createResidueTemplate(residue) # use helper function
            template.name = template_name # replace template name
            for (template_atom, residue_atom) in zip(template.atoms, residue.atoms()):
                template_atom.type = 'XXX' # replace atom type
            # Register the template.
            forcefield.registerResidueTemplate(template)

            # Signal that we have successfully parameterized the residue.
            return True

        # Define forcefield parameters used by simpleTemplateGenerator.
        # NOTE: This parameter definition file will currently only work for residues that either have
        # no external bonds or external bonds to other residues parameterized by the simpleTemplateGenerator.
        simple_ffxml_contents = """
<ForceField>
 <AtomTypes>
  <Type name="XXX" class="XXX" element="C" mass="12.0"/>
 </AtomTypes>
 <HarmonicBondForce>
  <Bond type1="XXX" type2="XXX" length="0.1409" k="392459.2"/>
 </HarmonicBondForce>
 <HarmonicAngleForce>
  <Angle type1="XXX" type2="XXX" type3="XXX" angle="2.09439510239" k="527.184"/>
 </HarmonicAngleForce>
 <NonbondedForce coulomb14scale="0.833333" lj14scale="0.5">
  <Atom type="XXX" charge="0.000" sigma="0.315" epsilon="0.635"/>
 </NonbondedForce>
</ForceField>"""

        #
        # Test where we generate parameters for only a ligand.
        #

        # Load the PDB file.
        pdb = PDBFile(os.path.join('systems', 'T4-lysozyme-L99A-p-xylene-implicit.pdb'))
        # Create a ForceField object.
        forcefield = ForceField('amber99sb.xml', 'tip3p.xml', StringIO(simple_ffxml_contents))
        # Add the residue template generator.
        forcefield.registerTemplateGenerator(simpleTemplateGenerator)
        # Parameterize system.
        system = forcefield.createSystem(pdb.topology, nonbondedMethod=NoCutoff)
        # TODO: Test energies are finite?

        #
        # Test for a few systems where we generate all parameters.
        #

        tests = [
            { 'pdb_filename' : 'alanine-dipeptide-implicit.pdb', 'nonbondedMethod' : NoCutoff },
            { 'pdb_filename' : 'lysozyme-implicit.pdb', 'nonbondedMethod' : NoCutoff },
            { 'pdb_filename' : 'alanine-dipeptide-explicit.pdb', 'nonbondedMethod' : CutoffPeriodic },
            ]

        # Test all systems with separate ForceField objects.
        for test in tests:
            # Load the PDB file.
            pdb = PDBFile(os.path.join('systems', test['pdb_filename']))
            # Create a ForceField object.
            forcefield = ForceField(StringIO(simple_ffxml_contents))
            # Add the residue template generator.
            forcefield.registerTemplateGenerator(simpleTemplateGenerator)
            # Parameterize system.
            system = forcefield.createSystem(pdb.topology, nonbondedMethod=test['nonbondedMethod'])
            # TODO: Test energies are finite?

        # Now test all systems with a single ForceField object.
        # Create a ForceField object.
        forcefield = ForceField(StringIO(simple_ffxml_contents))
        # Add the residue template generator.
        forcefield.registerTemplateGenerator(simpleTemplateGenerator)
        for test in tests:
            # Load the PDB file.
            pdb = PDBFile(os.path.join('systems', test['pdb_filename']))
            # Parameterize system.
            system = forcefield.createSystem(pdb.topology, nonbondedMethod=test['nonbondedMethod'])
            # TODO: Test energies are finite?

    def test_getUnmatchedResidues(self):
        """Test retrieval of list of residues for which no templates are available."""

        # Load the PDB file.
        pdb = PDBFile(os.path.join('systems', 'T4-lysozyme-L99A-p-xylene-implicit.pdb'))
        # Create a ForceField object.
        forcefield = ForceField('amber99sb.xml', 'tip3p.xml')
        # Get list of unmatched residues.
        unmatched_residues = forcefield.getUnmatchedResidues(pdb.topology)
        # Check results.
        self.assertEqual(len(unmatched_residues), 1)
        self.assertEqual(unmatched_residues[0].name, 'TMP')
        self.assertEqual(unmatched_residues[0].id, '163')

        # Load the PDB file.
        pdb = PDBFile(os.path.join('systems', 'ala_ala_ala.pdb'))
        # Create a ForceField object.
        forcefield = ForceField('tip3p.xml')
        # Get list of unmatched residues.
        unmatched_residues = forcefield.getUnmatchedResidues(pdb.topology)
        # Check results.
        self.assertEqual(len(unmatched_residues), 3)
        self.assertEqual(unmatched_residues[0].name, 'ALA')
        self.assertEqual(unmatched_residues[0].chain.id, 'X')
        self.assertEqual(unmatched_residues[0].id, '1')

    def test_generateTemplatesForUnmatchedResidues(self):
        """Test generation of blank forcefield residue templates for unmatched residues."""
        #
        # Test where we generate parameters for only a ligand.
        #

        # Load the PDB file.
        pdb = PDBFile(os.path.join('systems', 'nacl-water.pdb'))
        # Create a ForceField object.
        forcefield = ForceField('tip3p.xml')
        # Get list of unmatched residues.
        unmatched_residues = forcefield.getUnmatchedResidues(pdb.topology)
        [templates, residues] = forcefield.generateTemplatesForUnmatchedResidues(pdb.topology)
        # Check results.
        self.assertEqual(len(unmatched_residues), 24)
        self.assertEqual(len(residues), 2)
        self.assertEqual(len(templates), 2)
        unique_names = set([ residue.name for residue in residues ])
        self.assertTrue('HOH' not in unique_names)
        self.assertTrue('NA' in unique_names)
        self.assertTrue('CL' in unique_names)
        template_names = set([ template.name for template in templates ])
        self.assertTrue('HOH' not in template_names)
        self.assertTrue('NA' in template_names)
        self.assertTrue('CL' in template_names)

        # Define forcefield parameters using returned templates.
        # NOTE: This parameter definition file will currently only work for residues that either have
        # no external bonds or external bonds to other residues parameterized by the simpleTemplateGenerator.
        simple_ffxml_contents = """
<ForceField>
 <AtomTypes>
  <Type name="XXX" class="XXX" element="C" mass="12.0"/>
 </AtomTypes>
 <HarmonicBondForce>
  <Bond type1="XXX" type2="XXX" length="0.1409" k="392459.2"/>
 </HarmonicBondForce>
 <HarmonicAngleForce>
  <Angle type1="XXX" type2="XXX" type3="XXX" angle="2.09439510239" k="527.184"/>
 </HarmonicAngleForce>
 <NonbondedForce coulomb14scale="0.833333" lj14scale="0.5">
  <Atom type="XXX" charge="0.000" sigma="0.315" epsilon="0.635"/>
 </NonbondedForce>
</ForceField>"""

        #
        # Test the pre-geenration of missing residue template for a ligand.
        #

        # Load the PDB file.
        pdb = PDBFile(os.path.join('systems', 'T4-lysozyme-L99A-p-xylene-implicit.pdb'))
        # Create a ForceField object.
        forcefield = ForceField('amber99sb.xml', 'tip3p.xml', StringIO(simple_ffxml_contents))
        # Get list of unique unmatched residues.
        [templates, residues] = forcefield.generateTemplatesForUnmatchedResidues(pdb.topology)
        # Add residue templates to forcefield.
        for template in templates:
            # Replace atom types.
            for atom in template.atoms:
                atom.type = 'XXX'
            # Register the template.
            forcefield.registerResidueTemplate(template)
        # Parameterize system.
        system = forcefield.createSystem(pdb.topology, nonbondedMethod=NoCutoff)
        # TODO: Test energies are finite?

    def test_getMatchingTemplates(self):
        """Test retrieval of list of templates that match residues in a topology."""

        # Load the PDB file.
        pdb = PDBFile(os.path.join('systems', 'ala_ala_ala.pdb'))
        # Create a ForceField object.
        forcefield = ForceField('amber99sb.xml')
        # Get list of matching residue templates.
        templates = forcefield.getMatchingTemplates(pdb.topology)
        # Check results.
        residues = [ residue for residue in pdb.topology.residues() ]
        self.assertEqual(len(templates), len(residues))
        self.assertEqual(templates[0].name, 'NALA')
        self.assertEqual(templates[1].name, 'ALA')
        self.assertEqual(templates[2].name, 'CALA')

    def test_Wildcard(self):
        """Test that PeriodicTorsionForces using wildcard ('') for atom types / classes in the ffxml are correctly registered"""

        # Use wildcards in types
        xml = """
<ForceField>
 <AtomTypes>
  <Type name="C" class="C" element="C" mass="12.010000"/>
  <Type name="O" class="O" element="O" mass="16.000000"/>
 </AtomTypes>
 <PeriodicTorsionForce>
  <Proper type1="" type2="C" type3="C" type4="" periodicity1="2" phase1="3.141593" k1="15.167000"/>
  <Improper type1="C" type2="" type3="" type4="O" periodicity1="2" phase1="3.141593" k1="43.932000"/>
 </PeriodicTorsionForce>
</ForceField>"""

        ff = ForceField(StringIO(xml))

        self.assertEqual(len(ff._forces[0].proper), 1)
        self.assertEqual(len(ff._forces[0].improper), 1)

       # Use wildcards in classes
        xml = """
<ForceField>
 <AtomTypes>
  <Type name="C" class="C" element="C" mass="12.010000"/>
  <Type name="O" class="O" element="O" mass="16.000000"/>
 </AtomTypes>
 <PeriodicTorsionForce>
  <Proper class1="" class2="C" class3="C" class4="" periodicity1="2" phase1="3.141593" k1="15.167000"/>
  <Improper class1="C" class2="" class3="" class4="O" periodicity1="2" phase1="3.141593" k1="43.932000"/>
 </PeriodicTorsionForce>
</ForceField>"""

        ff = ForceField(StringIO(xml))

        self.assertEqual(len(ff._forces[0].proper), 1)
        self.assertEqual(len(ff._forces[0].improper), 1)

    def test_ScalingFactorCombining(self):
        """ Tests that FFs can be combined if their scaling factors are very close """
        forcefield = ForceField('amber99sb.xml', os.path.join('systems', 'test_amber_ff.xml'))
        # This would raise an exception if it didn't work

    def test_MultipleFilesandForceTags(self):
        """Test that the order of listing of multiple ffxmls does not matter.
           Tests that one generator per force type is created and that the ffxml
           defining atom types does not have to be listed first"""

        ffxml = """<ForceField>
 <Residues>
  <Residue name="ACE-Test">
   <Atom name="HH31" type="710"/>
   <Atom name="CH3" type="711"/>
   <Atom name="HH32" type="710"/>
   <Atom name="HH33" type="710"/>
   <Atom name="C" type="712"/>
   <Atom name="O" type="713"/>
   <Bond from="0" to="1"/>
   <Bond from="1" to="2"/>
   <Bond from="1" to="3"/>
   <Bond from="1" to="4"/>
   <Bond from="4" to="5"/>
   <ExternalBond from="4"/>
  </Residue>
 </Residues>
 <PeriodicTorsionForce>
  <Proper class1="C" class2="C" class3="C" class4="C" periodicity1="2" phase1="3.14159265359" k1="10.46"/>
  <Improper class1="C" class2="C" class3="C" class4="C" periodicity1="2" phase1="3.14159265359" k1="43.932"/>
 </PeriodicTorsionForce>
</ForceField>"""

        ff1 = ForceField(StringIO(ffxml), 'amber99sbildn.xml')
        ff2 = ForceField('amber99sbildn.xml', StringIO(ffxml))

        self.assertEqual(len(ff1._forces), 4)
        self.assertEqual(len(ff2._forces), 4)

        pertorsion1 = ff1._forces[0]
        pertorsion2 = ff2._forces[2]

        self.assertEqual(len(pertorsion1.proper), 110)
        self.assertEqual(len(pertorsion1.improper), 42)
        self.assertEqual(len(pertorsion2.proper), 110)
        self.assertEqual(len(pertorsion2.improper), 42)

    def test_ResidueTemplateUserChoice(self):
        """Test createSystem does not allow multiple matching templates, unless
           user has specified which template to use via residueTemplates arg"""
        ffxml = """<ForceField>
 <AtomTypes>
  <Type name="Fe2+" class="Fe2+" element="Fe" mass="55.85"/>
  <Type name="Fe3+" class="Fe3+" element="Fe" mass="55.85"/>
 </AtomTypes>
 <Residues>
  <Residue name="FE2">
   <Atom name="FE2" type="Fe2+" charge="2.0"/>
  </Residue>
  <Residue name="FE">
   <Atom name="FE" type="Fe3+" charge="3.0"/>
  </Residue>
 </Residues>
 <NonbondedForce coulomb14scale="0.833333333333" lj14scale="0.5">
  <UseAttributeFromResidue name="charge"/>
  <Atom type="Fe2+" sigma="0.227535532613" epsilon="0.0150312292"/>
  <Atom type="Fe3+" sigma="0.192790482606" epsilon="0.00046095128"/>
 </NonbondedForce>
</ForceField>"""

        pdb_string = "ATOM      1 FE    FE A   1      20.956  27.448 -29.067  1.00  0.00          Fe"
        ff = ForceField(StringIO(ffxml))
        pdb = PDBFile(StringIO(pdb_string))

        self.assertRaises(Exception, lambda: ff.createSystem(pdb.topology))
        sys = ff.createSystem(pdb.topology, residueTemplates={list(pdb.topology.residues())[0] : 'FE2'})
        # confirm charge
        self.assertEqual(sys.getForce(0).getParticleParameters(0)[0]._value, 2.0)
        sys = ff.createSystem(pdb.topology, residueTemplates={list(pdb.topology.residues())[0] : 'FE'})
        # confirm charge
        self.assertEqual(sys.getForce(0).getParticleParameters(0)[0]._value, 3.0)

    def test_ResidueOverriding(self):
        """Test residue overriding via override tag in the XML"""

        ffxml1 = """<ForceField>
 <AtomTypes>
  <Type name="Fe2+_tip3p_HFE" class="Fe2+_tip3p_HFE" element="Fe" mass="55.85"/>
 </AtomTypes>
 <Residues>
  <Residue name="FE2">
   <Atom name="FE2" type="Fe2+_tip3p_HFE" charge="2.0"/>
  </Residue>
 </Residues>
 <NonbondedForce coulomb14scale="0.833333333333" lj14scale="0.5">
  <UseAttributeFromResidue name="charge"/>
  <Atom type="Fe2+_tip3p_HFE" sigma="0.227535532613" epsilon="0.0150312292"/>
 </NonbondedForce>
</ForceField>"""

        ffxml2 = """<ForceField>
 <AtomTypes>
  <Type name="Fe2+_tip3p_standard" class="Fe2+_tip3p_standard" element="Fe" mass="55.85"/>
 </AtomTypes>
 <Residues>
  <Residue name="FE2">
   <Atom name="FE2" type="Fe2+_tip3p_standard" charge="2.0"/>
  </Residue>
 </Residues>
 <NonbondedForce coulomb14scale="0.833333333333" lj14scale="0.5">
  <UseAttributeFromResidue name="charge"/>
  <Atom type="Fe2+_tip3p_standard" sigma="0.241077193129" epsilon="0.03940482832"/>
 </NonbondedForce>
</ForceField>"""

        ffxml3 = """<ForceField>
 <AtomTypes>
  <Type name="Fe2+_tip3p_standard" class="Fe2+_tip3p_standard" element="Fe" mass="55.85"/>
 </AtomTypes>
 <Residues>
  <Residue name="FE2" override="1">
   <Atom name="FE2" type="Fe2+_tip3p_standard" charge="2.0"/>
  </Residue>
 </Residues>
 <NonbondedForce coulomb14scale="0.833333333333" lj14scale="0.5">
  <UseAttributeFromResidue name="charge"/>
  <Atom type="Fe2+_tip3p_standard" sigma="0.241077193129" epsilon="0.03940482832"/>
 </NonbondedForce>
</ForceField>"""

        pdb_string = "ATOM      1 FE    FE A   1      20.956  27.448 -29.067  1.00  0.00          Fe"
        pdb = PDBFile(StringIO(pdb_string))

        self.assertRaises(Exception, lambda: ForceField(StringIO(ffxml1), StringIO(ffxml2)))
        ff = ForceField(StringIO(ffxml1), StringIO(ffxml3))
        self.assertEqual(ff._templates['FE2'].atoms[0].type, 'Fe2+_tip3p_standard')
        ff.createSystem(pdb.topology)

    def test_LennardJonesGenerator(self):
        """ Test the LennardJones generator"""
        warnings.filterwarnings('ignore', category=CharmmPSFWarning)
        psf = CharmmPsfFile('systems/ions.psf')
        pdb = PDBFile('systems/ions.pdb')
        params = CharmmParameterSet('systems/toppar_water_ions.str'
                                    )

        # Box dimensions (found from bounding box)
        psf.setBox(12.009*angstroms,   12.338*angstroms,   11.510*angstroms)

        # Turn off charges so we only test the Lennard-Jones energies
        for a in psf.atom_list:
            a.charge = 0.0

        # Now compute the full energy
        plat = Platform.getPlatformByName('Reference')
        system = psf.createSystem(params, nonbondedMethod=PME,
                                  nonbondedCutoff=5*angstroms)

        con = Context(system, VerletIntegrator(2*femtoseconds), plat)
        con.setPositions(pdb.positions)

        # Now set up system from ffxml.
        xml = """
<ForceField>
 <AtomTypes>
  <Type name="SOD" class="SOD" element="Na" mass="22.98977"/>
  <Type name="CLA" class="CLA" element="Cl" mass="35.45"/>
 </AtomTypes>
 <Residues>
  <Residue name="CLA">
   <Atom name="CLA" type="CLA"/>
  </Residue>
  <Residue name="SOD">
   <Atom name="SOD" type="SOD"/>
  </Residue>
 </Residues>
 <LennardJonesForce lj14scale="1.0" useDispersionCorrection="False">
  <Atom type="CLA" sigma="0.404468018036" epsilon="0.6276"/>
  <Atom type="SOD" sigma="0.251367073323" epsilon="0.1962296"/>
  <NBFixPair type1="CLA" type2="SOD" sigma="0.33239431" epsilon="0.350933"/>
 </LennardJonesForce>
</ForceField> """
        ff = ForceField(StringIO(xml))
        system2 = ff.createSystem(pdb.topology, nonbondedMethod=PME,
                                  nonbondedCutoff=5*angstroms)
        con2 = Context(system2, VerletIntegrator(2*femtoseconds), plat)
        con2.setPositions(pdb.positions)

        state = con.getState(getEnergy=True, enforcePeriodicBox=True)
        ene = state.getPotentialEnergy().value_in_unit(kilocalories_per_mole)
        state2 = con2.getState(getEnergy=True, enforcePeriodicBox=True)
        ene2 = state2.getPotentialEnergy().value_in_unit(kilocalories_per_mole)
        self.assertAlmostEqual(ene, ene2)

    def test_NBFix(self):
        """Test using LennardJonesGenerator to implement NBFix terms."""
        # Create a chain of five atoms.

        top = Topology()
        chain = top.addChain()
        res = top.addResidue('RES', chain)
        top.addAtom('A', elem.oxygen, res)
        top.addAtom('B', elem.carbon, res)
        top.addAtom('C', elem.carbon, res)
        top.addAtom('D', elem.carbon, res)
        top.addAtom('E', elem.nitrogen, res)
        atoms = list(top.atoms())
        top.addBond(atoms[0], atoms[1])
        top.addBond(atoms[1], atoms[2])
        top.addBond(atoms[2], atoms[3])
        top.addBond(atoms[3], atoms[4])

        # Create the force field and system.

        xml = """
<ForceField>
 <AtomTypes>
  <Type name="A" class="A" element="O" mass="1"/>
  <Type name="B" class="B" element="C" mass="1"/>
  <Type name="C" class="C" element="C" mass="1"/>
  <Type name="D" class="D" element="C" mass="1"/>
  <Type name="E" class="E" element="N" mass="1"/>
 </AtomTypes>
 <Residues>
  <Residue name="RES">
   <Atom name="A" type="A"/>
   <Atom name="B" type="B"/>
   <Atom name="C" type="C"/>
   <Atom name="D" type="D"/>
   <Atom name="E" type="E"/>
   <Bond atomName1="A" atomName2="B"/>
   <Bond atomName1="B" atomName2="C"/>
   <Bond atomName1="C" atomName2="D"/>
   <Bond atomName1="D" atomName2="E"/>
  </Residue>
 </Residues>
 <LennardJonesForce lj14scale="0.3">
  <Atom type="A" sigma="1" epsilon="0.1"/>
  <Atom type="B" sigma="2" epsilon="0.2"/>
  <Atom type="C" sigma="3" epsilon="0.3"/>
  <Atom type="D" sigma="4" epsilon="0.4"/>
  <Atom type="E" sigma="4" epsilon="0.4"/>
  <NBFixPair type1="A" type2="D" sigma="2.5" epsilon="1.1"/>
  <NBFixPair type1="A" type2="E" sigma="3.5" epsilon="1.5"/>
 </LennardJonesForce>
</ForceField> """
        ff = ForceField(StringIO(xml))
        system = ff.createSystem(top)

        # Check that it produces the correct energy.

        integrator = VerletIntegrator(0.001)
        context = Context(system, integrator, Platform.getPlatform(0))
        positions = [Vec3(i, 0, 0) for i in range(5)]*nanometers
        context.setPositions(positions)
        def ljEnergy(sigma, epsilon, r):
            return 4*epsilon*((sigma/r)**12-(sigma/r)**6)
        expected = 0.3*ljEnergy(2.5, 1.1, 3) + 0.3*ljEnergy(3.0, sqrt(0.08), 3) + ljEnergy(3.5, 1.5, 4)
        self.assertAlmostEqual(expected, context.getState(getEnergy=True).getPotentialEnergy().value_in_unit(kilojoules_per_mole))

    def test_IgnoreExternalBonds(self):
        """Test the ignoreExternalBonds option"""

        modeller = Modeller(self.pdb2.topology, self.pdb2.positions)
        modeller.delete([next(modeller.topology.residues())])
        self.assertRaises(Exception, lambda: self.forcefield2.createSystem(modeller.topology))
        system = self.forcefield2.createSystem(modeller.topology, ignoreExternalBonds=True)
        templates = self.forcefield2.getMatchingTemplates(modeller.topology, ignoreExternalBonds=True)
        self.assertEqual(2, len(templates))
        self.assertEqual('ALA', templates[0].name)
        self.assertEqual('NME', templates[1].name)

    def test_Includes(self):
        """Test using a ForceField that includes other files."""
        forcefield = ForceField(os.path.join('systems', 'ff_with_includes.xml'))
        self.assertTrue(len(forcefield._atomTypes) > 10)
        self.assertTrue('spce-O' in forcefield._atomTypes)
        self.assertTrue('HOH' in forcefield._templates)

    def test_ImpropersOrdering(self):
        """Test correctness of the ordering of atom indexes in improper torsions
        and the torsion.ordering parameter.
        """

        xml = """
<ForceField>
 <PeriodicTorsionForce ordering="amber">
  <Improper class1="C" class2="" class3="O2" class4="O2" periodicity1="2" phase1="3.14159265359" k1="43.932"/>
 </PeriodicTorsionForce>
</ForceField>
"""
        pdb = PDBFile('systems/impropers_ordering_tetrapeptide.pdb')
        # ff1 uses default ordering of impropers, ff2 uses "amber" for the one
        # problematic improper
        ff1 = ForceField('amber99sbildn.xml')
        ff2 = ForceField(StringIO(xml), 'amber99sbildn.xml')

        system1 = ff1.createSystem(pdb.topology)
        system2 = ff2.createSystem(pdb.topology)

        imp1 = system1.getForce(1).getTorsionParameters(158)
        imp2 = system2.getForce(0).getTorsionParameters(158)

        system1_indexes = [imp1[0], imp1[1], imp1[2], imp1[3]]
        system2_indexes = [imp2[0], imp2[1], imp2[2], imp2[3]]

        self.assertEqual(system1_indexes, [51, 55, 54, 56])
        self.assertEqual(system2_indexes, [51, 55, 54, 56])

    def test_ImpropersOrdering_smirnoff(self):
        """Test correctness of the ordering of atom indexes in improper torsions
        and the torsion.ordering parameter when using the 'smirnoff' mode.
        """

        # SMIRNOFF parameters for formaldehyde
        xml = """
<ForceField>
  <AtomTypes>
    <Type name="[H]C(=O)[H]$C1#0" element="C" mass="12.01078" class="[H]C(=O)[H]$C1#0"/>
    <Type name="[H]C(=O)[H]$O1#1" element="O" mass="15.99943" class="[H]C(=O)[H]$O1#1"/>
    <Type name="[H]C(=O)[H]$H1#2" element="H" mass="1.007947" class="[H]C(=O)[H]$H1#2"/>
    <Type name="[H]C(=O)[H]$H2#3" element="H" mass="1.007947" class="[H]C(=O)[H]$H2#3"/>
  </AtomTypes>
  <PeriodicTorsionForce ordering="smirnoff">
    <Improper class1="[H]C(=O)[H]$C1#0" class2="[H]C(=O)[H]$O1#1" class3="[H]C(=O)[H]$H1#2" class4="[H]C(=O)[H]$H2#3" periodicity1="2" phase1="3.141592653589793" k1="1.5341333333333336"/>
    <Improper class1="[H]C(=O)[H]$C1#0" class2="[H]C(=O)[H]$H1#2" class3="[H]C(=O)[H]$H2#3" class4="[H]C(=O)[H]$O1#1" periodicity1="2" phase1="3.141592653589793" k1="1.5341333333333336"/>
    <Improper class1="[H]C(=O)[H]$C1#0" class2="[H]C(=O)[H]$H2#3" class3="[H]C(=O)[H]$O1#1" class4="[H]C(=O)[H]$H1#2" periodicity1="2" phase1="3.141592653589793" k1="1.5341333333333336"/>
  </PeriodicTorsionForce>
  <Residues>
    <Residue name="[H]C(=O)[H]">
      <Atom name="C1" type="[H]C(=O)[H]$C1#0" charge="0.5632799863815308"/>
      <Atom name="O1" type="[H]C(=O)[H]$O1#1" charge="-0.514739990234375"/>
      <Atom name="H1" type="[H]C(=O)[H]$H1#2" charge="-0.02426999807357788"/>
      <Atom name="H2" type="[H]C(=O)[H]$H2#3" charge="-0.02426999807357788"/>
      <Bond atomName1="C1" atomName2="O1"/>
      <Bond atomName1="C1" atomName2="H1"/>
      <Bond atomName1="C1" atomName2="H2"/>
    </Residue>
  </Residues>
</ForceField>
"""
        pdb = PDBFile('systems/formaldehyde.pdb')
        # ff1 uses default ordering of impropers, ff2 uses "amber" for the one
        # problematic improper
        ff = ForceField(StringIO(xml))

        system = ff.createSystem(pdb.topology)

        # Check that impropers are applied in the correct three-fold trefoil pattern
        forces = { force.__class__.__name__ : force for force in system.getForces() }
        force = forces['PeriodicTorsionForce']
        created_torsions = set()
        for index in range(force.getNumTorsions()):
            i,j,k,l,_,_,_ = force.getTorsionParameters(index)
            created_torsions.add((i,j,k,l))
        expected_torsions = set([(0,3,1,2), (0,1,2,3), (0,2,3,1)])
        self.assertEqual(expected_torsions, created_torsions)

    def test_Disulfides(self):
        """Test that various force fields handle disulfides correctly."""
        pdb = PDBFile('systems/bpti.pdb')
        for ff in ['amber99sb.xml', 'amber14-all.xml', 'charmm36.xml', 'amberfb15.xml', 'amoeba2013.xml']:
            forcefield = ForceField(ff)
            system = forcefield.createSystem(pdb.topology)

    def test_IdenticalTemplates(self):
        """Test a case where patches produce two identical templates."""
        ff = ForceField('charmm36.xml')
        pdb = PDBFile(StringIO("""
ATOM      1  N   HIS     1A   -2.670    -0.476   0.475  1.00  0.00           N
ATOM      2  HT1 HIS     1A   -2.645    -1.336   1.036  1.00  0.00           H
ATOM      3  HT2 HIS     1A   -2.859    -0.751  -0.532  1.00  0.00           H
ATOM      4  HT3 HIS     1A   -3.415     0.201   0.731  1.00  0.00           H
ATOM      5  CA  HIS     1A   -1.347     0.163   0.471  1.00  0.00           C
ATOM      6  HA  HIS     1A   -1.111     0.506   1.479  1.00  0.00           H
ATOM      7  CB  HIS     1A   -0.352    -0.857  -0.040  1.00  0.00           C
ATOM      8  HB1 HIS     1A   -0.360    -1.741   0.636  1.00  0.00           H
ATOM      9  HB2 HIS     1A   -0.640    -1.175  -1.046  1.00  0.00           H
ATOM     10  CG  HIS     1A    1.003    -0.275  -0.063  1.00  0.00           C
ATOM     11  CD2 HIS     1A    2.143    -0.931  -0.476  1.00  0.00           C
ATOM     12  HD2 HIS     1A    2.217    -1.952  -0.840  1.00  0.00           H
ATOM     13  NE2 HIS     1A    3.137    -0.024  -0.328  1.00  0.00           N
ATOM     14  HE2 HIS     1A    4.132    -0.238  -0.565  1.00  0.00           H
ATOM     15  CE1 HIS     1A    2.649     1.130   0.150  1.00  0.00           C
ATOM     16  HE1 HIS     1A    3.233     2.020   0.360  1.00  0.00           H
ATOM     17  ND1 HIS     1A    1.323     0.973   0.314  1.00  0.00           N
ATOM     18  C   HIS     1A   -1.465     1.282  -0.497  1.00  0.00           C
ATOM     19  OT1 HIS     1A   -2.108     2.309  -0.180  1.00  0.00           O
ATOM     20  OT2 HIS     1A   -0.864     1.172  -1.737  1.00  0.00           O
END"""))
        # If the check is not done correctly, this will throw an exception.
        ff.createSystem(pdb.topology)

    @unittest.skipIf(True, 'Skipping test failing due to missing platforms')
    def test_CharmmPolar(self):
        """Test the CHARMM polarizable force field."""
        pdb = PDBFile('systems/ala_ala_ala_drude.pdb')
        pdb.topology.setUnitCellDimensions(Vec3(3, 3, 3))
        ff = ForceField('charmm_polar_2019.xml')
        system = ff.createSystem(pdb.topology, nonbondedMethod=PME, nonbondedCutoff=1.2*nanometers)
        for i,f in enumerate(system.getForces()):
            f.setForceGroup(i)
            if isinstance(f, NonbondedForce):
                f.setPMEParameters(3.4, 64, 64, 64)
        integrator = DrudeLangevinIntegrator(300, 1.0, 1.0, 10.0, 0.001)
        context = Context(system, integrator, Platform.getPlatformByName('Reference'))
        context.setPositions(pdb.positions)

        # Compare the energy to values computed by CHARMM.  Here is what it outputs:

        # ENER ENR:  Eval#     ENERgy      Delta-E         GRMS
        # ENER INTERN:          BONDs       ANGLes       UREY-b    DIHEdrals    IMPRopers
        # ENER CROSS:           CMAPs        PMF1D        PMF2D        PRIMO
        # ENER EXTERN:        VDWaals         ELEC       HBONds          ASP         USER
        # ENER EWALD:          EWKSum       EWSElf       EWEXcl       EWQCor       EWUTil
        #  ----------       ---------    ---------    ---------    ---------    ---------
        # ENER>        0    102.83992      0.00000     13.06415
        # ENER INTERN>       54.72574     40.21459     11.61009     26.10373      0.14113
        # ENER CROSS>        -3.37113      0.00000      0.00000      0.00000
        # ENER EXTERN>       22.74761    -24.21667      0.00000      0.00000      0.00000
        # ENER EWALD>        56.14258  -7279.07968   7197.82192      0.00000      0.00000
        #  ----------       ---------    ---------    ---------    ---------    ---------

        # First check the total energy.
        
        energy = context.getState(getEnergy=True).getPotentialEnergy().value_in_unit(kilocalories_per_mole)
        self.assertAlmostEqual(102.83992, energy, delta=energy*1e-3)

        # Now check individual components.  CHARMM and OpenMM split them up a little differently.  I've tried to
        # match things up, but I think there's still some inconsistency in where forces related to Drude particles
        # are categorized.  That's why the Coulomb and bonds terms match less accurately than the other terms
        # (and less accurately than the total energy, which agrees well).

        coulomb = 0
        vdw = 0
        bonds = 0
        angles = 0
        propers = 0
        impropers = 0
        cmap = 0
        for i,f in enumerate(system.getForces()):
            energy = context.getState(getEnergy=True, groups={i}).getPotentialEnergy().value_in_unit(kilocalories_per_mole)
            if isinstance(f, NonbondedForce):
                coulomb += energy
            elif isinstance(f, CustomNonbondedForce) or isinstance(f, CustomBondForce):
                vdw += energy
            elif isinstance(f, HarmonicBondForce) or isinstance(f, DrudeForce):
                bonds += energy
            elif isinstance(f, HarmonicAngleForce):
                angles += energy
            elif isinstance(f, PeriodicTorsionForce):
                propers += energy
            elif isinstance(f, CustomTorsionForce):
                impropers += energy
            elif isinstance(f, CMAPTorsionForce):
                cmap += energy
        self.assertAlmostEqual(-24.21667+56.14258-7279.07968+7197.82192, coulomb, delta=abs(coulomb)*5e-2) # ELEC+EWKSum+EWSElf+EWEXcl
        self.assertAlmostEqual(22.74761, vdw, delta=vdw*1e-3) # VDWaals
        self.assertAlmostEqual(54.72574+11.61009, bonds, delta=bonds*2e-2) # BONDs+UREY-b
        self.assertAlmostEqual(40.21459, angles, delta=angles*1e-3) # ANGLes
        self.assertAlmostEqual(26.10373, propers, delta=propers*1e-3) # DIHEdrals
        self.assertAlmostEqual(0.14113, impropers, delta=impropers*1e-3) # IMPRopers

    def test_InitializationScript(self):
        """Test that <InitializationScript> tags get executed."""
        xml = """
<ForceField>
  <InitializationScript>
self.scriptExecuted = True
  </InitializationScript>
</ForceField>
"""
        ff = ForceField(StringIO(xml))
        self.assertTrue(ff.scriptExecuted)

    def test_Glycam(self):
        """Test computing energy with GLYCAM."""
        ff = ForceField('amber14/protein.ff14SB.xml', 'amber14/GLYCAM_06j-1.xml')
        pdb = PDBFile('systems/glycopeptide.pdb')
        system = ff.createSystem(pdb.topology)
        for i, f in enumerate(system.getForces()):
            f.setForceGroup(i)
        integrator = VerletIntegrator(0.001)
        context = Context(system, integrator, Platform.getPlatformByName('Reference'))
        context.setPositions(pdb.positions)
        energies = {}
        for i, f in enumerate(system.getForces()):
            energy = context.getState(getEnergy=True, groups={i}).getPotentialEnergy().value_in_unit(kilojoules_per_mole)
            energies[f.getName()] = energy

        # Compare to values computed with ParmEd.

        self.assertAlmostEqual(32.14082401103625, energies['HarmonicBondForce'], 4)
        self.assertAlmostEqual(48.92017455984504, energies['HarmonicAngleForce'], 3)
        self.assertAlmostEqual(291.61241586209286, energies['PeriodicTorsionForce'], 4)
        self.assertAlmostEqual(1547.011267801862, energies['NonbondedForce'], 4)
        self.assertAlmostEqual(1919.6846822348361, sum(list(energies.values())), 3)

    def test_CustomNonbondedGenerator(self):
        """ Test the CustomNonbondedForce generator"""
        pdb = PDBFile('systems/ions.pdb')
        xml = """
<ForceField>
 <AtomTypes>
  <Type name="SOD" class="SOD" element="Na" mass="22.98977"/>
  <Type name="CLA" class="CLA" element="Cl" mass="35.45"/>
 </AtomTypes>
 <Residues>
  <Residue name="CLA">
   <Atom name="CLA" type="CLA"/>
  </Residue>
  <Residue name="SOD">
   <Atom name="SOD" type="SOD"/>
  </Residue>
 </Residues>
 <CustomNonbondedForce energy="scale*epsilon*((sigma/r)^12-(sigma/r)^6); sigma=halfSig1+halfSig2; epsilon=rootEps1*rootEps2" bondCutoff="3">
  <GlobalParameter name="scale" defaultValue="4"/>
  <PerParticleParameter name="sigma"/>
  <PerParticleParameter name="epsilon"/>
  <ComputedValue name="halfSig" expression="0.5*sigma"/>
  <ComputedValue name="rootEps" expression="sqrt(epsilon)"/>
  <Atom type="CLA" sigma="0.404468018036" epsilon="0.6276"/>
  <Atom type="SOD" sigma="0.251367073323" epsilon="0.1962296"/>
 </CustomNonbondedForce>
</ForceField> """
        ff = ForceField(StringIO(xml))
        system = ff.createSystem(pdb.topology)
        context = Context(system, VerletIntegrator(2*femtoseconds), Platform.getPlatformByName('Reference'))
        context.setPositions(pdb.positions)
        energy1 = context.getState(getEnergy=True).getPotentialEnergy().value_in_unit(kilojoules_per_mole)

        # See if it matches an equivalent NonbondedForce.
        
        system = System()
        system.addParticle(1.0)
        system.addParticle(1.0)
        f = NonbondedForce()
        f.addParticle(0, 0.404468018036, 0.6276)
        f.addParticle(0, 0.251367073323, 0.1962296)
        system.addForce(f)
        context = Context(system, VerletIntegrator(2*femtoseconds), Platform.getPlatformByName('Reference'))
        context.setPositions(pdb.positions)
        energy2 = context.getState(getEnergy=True).getPotentialEnergy().value_in_unit(kilojoules_per_mole)
        self.assertAlmostEqual(energy1, energy2)

    def test_OpcEnergy(self):
        pdb = PDBFile('systems/opcbox.pdb')
        topology, positions = pdb.topology, pdb.positions
        self.assertEqual(len(positions), 864)
        forcefield = ForceField('opc.xml')
        system = forcefield.createSystem(
            topology,
            nonbondedMethod=PME,
            nonbondedCutoff=0.7*nanometer,
            constraints=HBonds,
            rigidWater=True,
        )

        integrator = LangevinIntegrator(300*kelvin, 2.0/picoseconds, 2.0*femtoseconds)
        simulation = Simulation(topology, system, integrator)
        context = simulation.context
        context.setPositions(positions)

        # Compare to values computed with Amber (sander).
        energy_amber = -2647.6233 # kcal/mol
        energy_tolerance = 1.0

        state = context.getState(getEnergy=True)
        energy1 = state.getPotentialEnergy().value_in_unit(kilocalorie_per_mole)
        # -2647.2222697324237
        self.assertTrue(abs(energy1 - energy_amber) < energy_tolerance)

        context.applyConstraints(1e-12)
        state = context.getState(getEnergy=True)
        energy2 = state.getPotentialEnergy().value_in_unit(kilocalorie_per_mole)
        # -2647.441600693312
        self.assertTrue(abs(energy1 - energy_amber) < energy_tolerance)
        self.assertTrue(abs(energy1 - energy2) < energy_tolerance)

    def test_Opc3Energy(self):
        pdb = PDBFile('systems/opc3box.pdb')
        topology, positions = pdb.topology, pdb.positions
        self.assertEqual(len(positions), 648)
        forcefield = ForceField('opc3.xml')
        system = forcefield.createSystem(
            topology,
            nonbondedMethod=PME,
            nonbondedCutoff=0.7*nanometer,
            constraints=HBonds,
            rigidWater=True,
        )

        integrator = LangevinIntegrator(300*kelvin, 2.0/picoseconds, 2.0*femtoseconds)
        simulation = Simulation(topology, system, integrator)
        context = simulation.context
        context.setPositions(positions)

        # Compare to values computed with Amber (sander).
        energy_amber = -2532.1414 # kcal/mol
        energy_tolerance = 1.0

        state = context.getState(getEnergy=True)
        energy1 = state.getPotentialEnergy().value_in_unit(kilocalorie_per_mole)
        # -2532.4862082354407
        self.assertTrue(abs(energy1 - energy_amber) < energy_tolerance)

        context.applyConstraints(1e-12)
        state = context.getState(getEnergy=True)
        energy2 = state.getPotentialEnergy().value_in_unit(kilocalorie_per_mole)
        self.assertTrue(abs(energy1 - energy_amber) < energy_tolerance)
        self.assertTrue(abs(energy1 - energy2) < energy_tolerance)


class AmoebaTestForceField(unittest.TestCase):
    """Test the ForceField.createSystem() method with the AMOEBA forcefield."""

    def setUp(self):
        """Set up the tests by loading the input pdb files and force field
        xml files.

        """

        self.pdb1 = PDBFile('systems/amoeba-ion-in-water.pdb')
        self.forcefield1 = ForceField('amoeba2013.xml')
        self.topology1 = self.pdb1.topology


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

        methodMap = {NoCutoff:AmoebaMultipoleForce.NoCutoff,
                     PME:AmoebaMultipoleForce.PME}

        for method in methodMap:
            system = self.forcefield1.createSystem(self.pdb1.topology,
                                                  nonbondedMethod=method)
            forces = system.getForces()
            self.assertTrue(any(isinstance(f, AmoebaMultipoleForce) and
                                f.getNonbondedMethod()==methodMap[method]
                                for f in forces))
    def test_Cutoff(self):
        """Test to make sure the nonbondedCutoff parameter is passed correctly."""

        cutoff_distance = 0.7*nanometer
        for method in [NoCutoff, PME]:
            system = self.forcefield1.createSystem(self.pdb1.topology,
                                                   nonbondedMethod=method,
                                                   nonbondedCutoff=cutoff_distance,
                                                   constraints=None)

            for force in system.getForces():
                if isinstance(force, AmoebaVdwForce):
                    self.assertEqual(force.getCutoff(), cutoff_distance)
                if isinstance(force, AmoebaMultipoleForce):
                    self.assertEqual(force.getCutoffDistance(), cutoff_distance)

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

        for useDispersionCorrection in [True, False]:
            system = self.forcefield1.createSystem(self.pdb1.topology,
                                                   nonbondedMethod=PME,
                                                   useDispersionCorrection=useDispersionCorrection)

            for force in system.getForces():
                if isinstance(force, AmoebaVdwForce):
                    self.assertEqual(useDispersionCorrection, force.getUseDispersionCorrection())

    def test_RigidWater(self):
        """Test that AMOEBA creates rigid water with the correct geometry."""

        system = self.forcefield1.createSystem(self.pdb1.topology, rigidWater=True)
        constraints = dict()
        for i in range(system.getNumConstraints()):
            p1,p2,dist = system.getConstraintParameters(i)
            if p1 < 3:
                constraints[(min(p1,p2), max(p1,p2))] = dist.value_in_unit(nanometers)
        hoDist = 0.09572
        hohAngle = 108.50*math.pi/180.0
        hohDist = math.sqrt(2*hoDist**2 - 2*hoDist**2*math.cos(hohAngle))
        self.assertAlmostEqual(constraints[(0,1)], hoDist)
        self.assertAlmostEqual(constraints[(0,2)], hoDist)
        self.assertAlmostEqual(constraints[(1,2)], hohDist)
        
        # Check that all values of rigidWater are interpreted correctly.
        
        numWaters = 215
        self.assertEqual(3*numWaters, system.getNumConstraints())
        system = self.forcefield1.createSystem(self.pdb1.topology, rigidWater=False)
        self.assertEqual(0, system.getNumConstraints())
        system = self.forcefield1.createSystem(self.pdb1.topology, rigidWater=None)
        self.assertEqual(0, system.getNumConstraints())

    @unittest.skipIf(True, 'Skipping test failing due to missing platforms')
    def test_Forces(self):
        """Compute forces and compare them to ones generated with a previous version of OpenMM to ensure they haven't changed."""

        pdb = PDBFile('systems/alanine-dipeptide-implicit.pdb')
        forcefield = ForceField('amoeba2013.xml', 'amoeba2013_gk.xml')
        system = forcefield.createSystem(pdb.topology, polarization='direct')
        integrator = VerletIntegrator(0.001)
        context = Context(system, integrator, Platform.getPlatformByName('Reference'))
        context.setPositions(pdb.positions)
        state1 = context.getState(getForces=True)
        with open('systems/alanine-dipeptide-amoeba-forces.xml') as input:
            state2 = XmlSerializer.deserialize(input.read())
        for f1, f2, in zip(state1.getForces().value_in_unit(kilojoules_per_mole/nanometer), state2.getForces().value_in_unit(kilojoules_per_mole/nanometer)):
            diff = norm(f1-f2)
            self.assertTrue(diff < 0.1 or diff/norm(f1) < 1e-3)

    def computeAmoeba18Energies(self, filename):
        pdb = PDBFile(filename)
        forcefield = ForceField('amoeba2018.xml')
        system = forcefield.createSystem(pdb.topology, polarization='mutual', mutualInducedTargetEpsilon=1e-5)
        for i, f in enumerate(system.getForces()):
            f.setForceGroup(i)
        integrator = VerletIntegrator(0.001)
        context = Context(system, integrator, Platform.getPlatformByName('Reference'))
        context.setPositions(pdb.positions)
        energies = {}
        for i, f in enumerate(system.getForces()):
            state = context.getState(getEnergy=True, groups={i})
            energies[f.getName()] = state.getPotentialEnergy().value_in_unit(kilocalories_per_mole)
        return energies

    @unittest.skipIf(True, 'Skipping test failing due to missing platforms')
    def test_Amoeba18BPTI(self):
        """Test that AMOEBA18 computes energies correctly for BPTI."""
        energies = self.computeAmoeba18Energies('systems/bpti.pdb')

        # Compare to values computed with Tinker.

        self.assertAlmostEqual(290.2445, energies['AmoebaBond'], 4)
        self.assertAlmostEqual(496.4300, energies['AmoebaAngle']+energies['AmoebaInPlaneAngle'], 4)
        self.assertAlmostEqual(51.2913, energies['AmoebaOutOfPlaneBend'], 4)
        self.assertAlmostEqual(5.7695, energies['AmoebaStretchBend'], 4)
        self.assertAlmostEqual(75.6890, energies['PeriodicTorsionForce'], 4)
        self.assertAlmostEqual(19.3364, energies['AmoebaPiTorsion'], 4)
        self.assertAlmostEqual(-32.6689, energies['AmoebaTorsionTorsionForce'], 4)
        self.assertAlmostEqual(383.8705, energies['AmoebaVdwForce'], 4)
        self.assertAlmostEqual(-1323.5640-225.3660, energies['AmoebaMultipoleForce'], 2)
        self.assertAlmostEqual(-258.9676, sum(list(energies.values())), 2)

    @unittest.skipIf(True, 'Skipping test failing due to missing platforms')
    def test_Amoeba18Nucleic(self):
        """Test that AMOEBA18 computes energies correctly for DNA and RNA."""
        energies = self.computeAmoeba18Energies('systems/nucleic.pdb')

        # Compare to values computed with Tinker.

        self.assertAlmostEqual(749.6953, energies['AmoebaBond'], 4)
        self.assertAlmostEqual(579.9971, energies['AmoebaAngle']+energies['AmoebaInPlaneAngle'], 4)
        self.assertAlmostEqual(10.6630, energies['AmoebaOutOfPlaneBend'], 4)
        self.assertAlmostEqual(5.2225, energies['AmoebaStretchBend'], 4)
        self.assertAlmostEqual(166.7233, energies['PeriodicTorsionForce'], 4)
        self.assertAlmostEqual(57.2066, energies['AmoebaPiTorsion'], 4)
        self.assertAlmostEqual(-4.2538, energies['AmoebaStretchTorsion'], 4)
        self.assertAlmostEqual(-5.0402, energies['AmoebaAngleTorsion'], 4)
        self.assertAlmostEqual(187.1103, energies['AmoebaVdwForce'], 4)
        self.assertAlmostEqual(1635.1289-236.1484, energies['AmoebaMultipoleForce'], 3)
        self.assertAlmostEqual(3146.3046, sum(list(energies.values())), 3)

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