Coverage for cvpack/torsion_similarity.py: 100%

17 statements  

« prev     ^ index     » next       coverage.py v7.5.1, created at 2024-05-09 16:14 +0000

1""" 

2.. class:: DihedralSimilarity 

3 :platform: Linux, MacOS, Windows 

4 :synopsis: The degree of similarity between pairs of torsion angles. 

5 

6.. classauthor:: Charlles Abreu <craabreu@gmail.com> 

7 

8""" 

9 

10import numpy as np 

11import openmm 

12from numpy.typing import ArrayLike 

13from openmm import unit as mmunit 

14 

15from .collective_variable import CollectiveVariable 

16 

17 

18class TorsionSimilarity(CollectiveVariable, openmm.CustomCompoundBondForce): 

19 r""" 

20 The degree of similarity between `n` pairs of torsion angles: 

21 

22 .. math:: 

23 

24 s({\bf r}) = \frac{n}{2} + \frac{1}{2} \sum_{i=1}^n \cos\Big( 

25 \phi^{\rm 1st}_i({\bf r}) - \phi^{\rm 2nd}_i({\bf r}) 

26 \Big) 

27 

28 where :math:`\phi^{\rm kth}_i` is the torsion angle at position :math:`i` in the 

29 :math:`k`-th list. 

30 

31 .. note:: 

32 

33 In `PLUMED <https://www.plumed.org/doc-v2.8/user-doc/html/_colvar.html>`_, this 

34 collective variable is called ``DIHCOR``. 

35 

36 Parameters 

37 ---------- 

38 firstList 

39 A list of :math:`n` tuples of four atom indices defining the first torsion 

40 angle in each pair. 

41 secondList 

42 A list of :math:`n` tuples of four atom indices defining the second torsion 

43 angle in each pair. 

44 pbc 

45 Whether to use periodic boundary conditions in distance calculations. 

46 name 

47 The name of the collective variable. 

48 

49 Example 

50 ------- 

51 >>> import cvpack 

52 >>> import mdtraj 

53 >>> import openmm 

54 >>> import openmm.unit as mmunit 

55 >>> from openmmtools import testsystems 

56 >>> model = testsystems.LysozymeImplicit() 

57 >>> traj = mdtraj.Trajectory( 

58 ... model.positions, mdtraj.Topology.from_openmm(model.topology) 

59 ... ) 

60 >>> phi_atoms, _ = mdtraj.compute_phi(traj) 

61 >>> valid_atoms = traj.top.select("resid 59 to 79 and backbone") 

62 >>> phi_atoms = [ 

63 ... phi 

64 ... for phi in phi_atoms 

65 ... if all(atom in valid_atoms for atom in phi) 

66 ... ] 

67 >>> torsion_similarity = cvpack.TorsionSimilarity( 

68 ... phi_atoms[1:], phi_atoms[:-1] 

69 ... ) 

70 >>> torsion_similarity.addToSystem(model.system) 

71 >>> integrator = openmm.VerletIntegrator(0) 

72 >>> platform = openmm.Platform.getPlatformByName("Reference") 

73 >>> context = openmm.Context(model.system, integrator, platform) 

74 >>> context.setPositions(model.positions) 

75 >>> torsion_similarity.getValue(context) 

76 18.659... dimensionless 

77 """ 

78 

79 def __init__( 

80 self, 

81 firstList: ArrayLike, 

82 secondList: ArrayLike, 

83 pbc: bool = False, 

84 name: str = "torsion_similarity", 

85 ) -> None: 

86 firstList = [list(map(int, first)) for first in firstList] 

87 secondList = [list(map(int, second)) for second in secondList] 

88 function = f"0.5*(1 + cos(min(delta, {2*np.pi} - delta)))" 

89 definition = "delta = dihedral(p1, p2, p3, p4) - dihedral(p5, p6, p7, p8)" 

90 super().__init__(8, f"{function}; {definition}") 

91 for first, second in zip(firstList, secondList): 

92 self.addBond([*first, *second], []) 

93 self.setUsesPeriodicBoundaryConditions(pbc) 

94 self._registerCV( 

95 name, 

96 mmunit.dimensionless, 

97 firstList=firstList, 

98 secondList=secondList, 

99 ) 

100 

101 

102TorsionSimilarity.registerTag("!cvpack.TorsionSimilarity")