Coverage for cvpack/meta_collective_variable.py: 97%
39 statements
« prev ^ index » next coverage.py v7.5.1, created at 2024-05-09 16:14 +0000
« prev ^ index » next coverage.py v7.5.1, created at 2024-05-09 16:14 +0000
1"""
2.. class:: MetaCollectiveVariable
3 :platform: Linux, MacOS, Windows
4 :synopsis: A function of other collective variables
6.. classauthor:: Charlles Abreu <craabreu@gmail.com>
8"""
10from __future__ import annotations
12import typing as t
13from copy import copy
15import openmm
16from openmm import unit as mmunit
18from .collective_variable import CollectiveVariable
19from .units import Quantity, ScalarQuantity, VectorQuantity, in_md_units
20from .utils import compute_effective_mass, get_single_force_state
23class MetaCollectiveVariable(CollectiveVariable, openmm.CustomCVForce):
24 r"""
25 A collective variable that is a function of other collective variables:
27 .. math::
29 f({\bf r}) = F\left(c_1({\bf r}), c_2({\bf r}), \ldots, c_n({\bf r})\right)
31 where :math:`F(c_1,c_2,\ldots,c_n)` is a user-defined function and
32 :math:`c_i({\bf r})` is the value of the :math:`i`-th collective variable at a
33 given configuration :math:`{\bf r}`.
35 The function :math:`F` is defined as a string and can be any expression supported
36 by the :OpenMM:`CustomCVForce` class. If the expression contains named parameters,
37 the value of each parameter can be passed in one of two ways:
39 #. By a semicolon-separated definition in the function string, such as described
40 in the :OpenMM:`CustomCompoundBondForce` documentation. In this case, the
41 parameter value will be the same for all groups of atoms.
43 #. By a scalar passed as a keyword argument to the :class:`AtomicFunction`
44 constructor. In this case, the parameter will apply to all atom groups and will
45 become available for on-the-fly modification during a simulation via the
46 ``setParameter`` method of an :OpenMM:`Context` object. **Warning**: other
47 collective variables or :OpenMM:`Force` objects in the same system will share
48 the same values of equal-named parameters.
50 Parameters
51 ----------
52 function
53 The function to be evaluated. It must be a valid :OpenMM:`CustomCVForce`
54 expression.
55 variables
56 A sequence of :class:`CollectiveVariable` instances that represent the
57 collective variables on which the meta-collective variable depends. The name
58 of each collective variable must be unique and match a symbol used in the
59 function.
60 unit
61 The unit of measurement of the collective variable. It must be compatible
62 with the MD unit system (mass in `daltons`, distance in `nanometers`, time
63 in `picoseconds`, temperature in `kelvin`, energy in `kilojoules_per_mol`,
64 angle in `radians`). If the collective variables does not have a unit, use
65 `unit.dimensionless`
66 periodicBounds
67 The periodic bounds of the collective variable if it is periodic, or `None` if
68 it is not periodic.
70 Keyword Args
71 ------------
72 **parameters
73 The named parameters of the function, if any. They will become settable context
74 parameters when this meta-collective variable is added to an :OpenMM:`System`.
75 The passed objects must be scalar quantities. Their values will be converted to
76 OpenMM's MD unit system to serve as default values for the context parameters.
78 Example
79 -------
80 >>> import cvpack
81 >>> import math
82 >>> import openmm
83 >>> from openmm import unit
84 >>> from openmmtools import testsystems
85 >>> model = testsystems.AlanineDipeptideVacuum()
86 >>> phi = cvpack.Torsion(6, 8, 14, 16, name="phi")
87 >>> driving_force = cvpack.MetaCollectiveVariable(
88 ... f"0.5*kappa*min(delta,{2*math.pi}-delta)^2"
89 ... "; delta=abs(phi-phi0)",
90 ... [phi],
91 ... unit.kilojoules_per_mole,
92 ... kappa = 1e3 * unit.kilojoules_per_mole/unit.radian**2,
93 ... phi0 = 120 * unit.degrees
94 ... )
95 >>> driving_force.addToSystem(model.system)
96 >>> integrator = openmm.VerletIntegrator(0)
97 >>> platform = openmm.Platform.getPlatformByName('Reference')
98 >>> context = openmm.Context(model.system, integrator, platform)
99 >>> context.setPositions(model.positions)
100 >>> driving_force.getParameterDefaultValues()
101 {'kappa': 1000.0 kJ/(mol rad**2), 'phi0': 2.094... rad}
102 >>> driving_force.getParameterValues(context)
103 {'kappa': 1000.0 kJ/(mol rad**2), 'phi0': 2.094... rad}
104 >>> driving_force.getValue(context)
105 548.3... kJ/mol
106 >>> driving_force.getInnerValues(context)
107 {'phi': 3.14... rad}
108 >>> driving_force.getInnerEffectiveMasses(context)
109 {'phi': 0.05119... nm**2 Da/(rad**2)}
110 >>> driving_force.getParameterDerivatives(context)
111 {'kappa': 0.548... rad**2, 'phi0': -1047.19... kJ/(mol rad)}
112 """
114 def __init__(
115 self,
116 function: str,
117 variables: t.Iterable[str],
118 unit: mmunit.Unit,
119 periodicBounds: t.Optional[VectorQuantity] = None,
120 name: str = "meta_collective_variable",
121 **parameters: ScalarQuantity,
122 ) -> None:
123 super().__init__(function)
124 self._cvs = tuple(map(copy, variables))
125 self._parameters = {k: in_md_units(v) for k, v in parameters.items()}
126 for cv in self._cvs:
127 self.addCollectiveVariable(cv.getName(), cv)
128 for parameter, value in self._parameters.items():
129 self.addGlobalParameter(parameter, value / value.unit)
130 self.addEnergyParameterDerivative(parameter)
131 self._registerCV(
132 name,
133 unit,
134 function=function,
135 variables=variables,
136 unit=unit,
137 periodicBounds=periodicBounds,
138 **self._parameters,
139 )
140 if periodicBounds is not None:
141 self._registerPeriodicBounds(*periodicBounds)
143 def getInnerVariables(self) -> t.Tuple[CollectiveVariable]:
144 """
145 Get the collective variables on which the meta-collective variable depends.
147 Returns
148 -------
149 Tuple[CollectiveVariable]
150 A tuple with the collective variables.
151 """
152 return self._cvs
154 def getInnerValues(self, context: openmm.Context) -> t.Dict[str, Quantity]:
155 """
156 Get the values of the collective variables on which the meta-collective variable
157 depends. The values are returned as a dictionary with the names of the
158 collective variables as keys.
160 Parameters
161 ----------
162 context
163 The context in which the collective variables will be evaluated.
165 Returns
166 -------
167 Dict[str, Quantity]
168 A dictionary with the names of the collective variables as keys and their
169 values as values.
170 """
171 values = self.getCollectiveVariableValues(context)
172 return {
173 cv.getName(): Quantity(value, cv.getUnit())
174 for cv, value in zip(self._cvs, values)
175 }
177 def getInnerEffectiveMasses(self, context: openmm.Context) -> t.Dict[str, Quantity]:
178 """
179 Get the effective masses of the collective variables on which the
180 meta-collective variable depends. The effective masses are calculated from the
181 forces acting on the particles that represent the collective variables.
183 Parameters
184 ----------
185 context
186 The context in which the collective variables will be evaluated.
188 Returns
189 -------
190 Dict[str, Quantity]
191 A dictionary with the names of the collective variables as keys and their
192 effective masses as values.
193 """
194 inner_context = self.getInnerContext(context)
195 masses = [
196 compute_effective_mass(force, inner_context)
197 for force in inner_context.getSystem().getForces()
198 ]
199 return {
200 cv.getName(): Quantity(mass, cv.getMassUnit())
201 for cv, mass in zip(self._cvs, masses)
202 }
204 def getParameterDefaultValues(self) -> t.Dict[str, Quantity]:
205 """
206 Get the default values of the named parameters of this meta-collective variable.
208 Returns
209 -------
210 Dict[str, Quantity]
211 A dictionary with the names of the named parameters as keys and their
212 default values as values.
213 """
214 return self._parameters.copy()
216 def getParameterValues(self, context: openmm.Context) -> t.Dict[str, Quantity]:
217 """
218 Get the values of the named parameters of this meta-collective variable. The
219 values are returned as a dictionary with the names of the parameters as keys.
221 Parameters
222 ----------
223 context
224 The context in which the named parameters will be evaluated.
226 Returns
227 -------
228 Dict[str, Quantity]
229 A dictionary with the names of the named parameters as keys and their values
230 as values.
231 """
232 return {
233 name: Quantity(context.getParameter(name), parameter.unit)
234 for name, parameter in self._parameters.items()
235 }
237 def getParameterDerivatives(
238 self,
239 context: openmm.Context,
240 allowReinitialization: bool = False,
241 ) -> t.Dict[str, Quantity]:
242 """
243 Get the derivatives of the named parameters of this meta-collective variable.
245 Returns
246 -------
247 Dict[str, Quantity]
248 A dictionary with the names of the named parameters as keys and their
249 derivatives as values.
250 """
251 state = get_single_force_state(
252 self, context, allowReinitialization, getParameterDerivatives=True
253 )
254 derivatives = state.getEnergyParameterDerivatives()
255 return {
256 name: Quantity(derivatives[name], self._unit / parameter.unit)
257 for name, parameter in self._parameters.items()
258 }
261MetaCollectiveVariable.registerTag("!cvpack.MetaCollectiveVariable")