Skip to content

Commit 495518e

Browse files
author
Charlles Abreu
authored
Improved shortest distance CV (#101)
* Changed shortest distance CV definition * Updated unit test * Improved accuracy
1 parent 0ee9a75 commit 495518e

File tree

2 files changed

+58
-47
lines changed

2 files changed

+58
-47
lines changed

Diff for: cvpack/shortest_distance.py

+47-43
Original file line numberDiff line numberDiff line change
@@ -9,6 +9,7 @@
99

1010
import typing as t
1111

12+
import numpy as np
1213
import openmm
1314
from openmm import unit as mmunit
1415

@@ -18,26 +19,40 @@
1819

1920
class ShortestDistance(CollectiveVariable, openmm.CustomCVForce):
2021
r"""
21-
A smooth approximation of the shortest distance between two atom groups:
22+
A smooth approximation for the shortest distance between two atom groups:
2223
2324
.. math::
24-
r_{\rm min}({\bf r}) = \frac{
25-
\sum_{i \in {\bf g}_1} \sum_{j \in {\bf g}_2}
26-
r_{ij} e^{-\frac{r_{ij}^2}{2 \sigma^2}}
27-
}{
28-
\sum_{i \in {\bf g}_1} \sum_{j \in {\bf g}_2}
29-
e^{-\frac{r_{ij}^2}{2 \sigma^2}}
30-
}
25+
r_{\rm min}({\bf r}) = \frac{S_{rw}({\bf r})}{S_w({\bf r})}
26+
27+
with
28+
29+
.. math::
30+
31+
S_{rw}({\bf r}) = r_c e^{-\beta} +
32+
\sum_{i \in {\bf g}_1} \sum_{j \in {\bf g}_2 \atop r_{ij} < r_c}
33+
r_{ij} e^{-\beta \frac{r_{ij}}{r_c}}
34+
35+
and
36+
37+
.. math::
38+
39+
S_w({\bf r}) = e^{-\beta} +
40+
\sum_{i \in {\bf g}_1} \sum_{j \in {\bf g}_2 \atop r_{ij} < r_c}
41+
r_{ij} e^{-\beta \frac{r_{ij}}{r_c}}
3142
3243
where :math:`r_{ij} = \|{\bf r}_j - {\bf r}_i\|` is the distance between atoms
33-
:math:`i` and :math:`j` and :math:`\sigma` is a parameter that controls the
34-
degree of approximation. The smaller the value of :math:`\sigma`, the closer the
35-
approximation to the true shortest distance.
44+
:math:`i` and :math:`j`, :math:`{\bf g}_1` and :math:`{\bf g}_2` are the sets of
45+
indices of the atoms in the first and second groups, respectively, :math:`r_c` is
46+
the cutoff distance, and :math:`\beta` is a parameter that controls the degree of
47+
approximation.
48+
49+
The larger the value of :math:`\beta`, the closer the approximation to the true
50+
shortest distance. However, an excessively large value may lead to numerical
51+
instability.
3652
37-
In practice, a cutoff distance :math:`r_c` is applied to the atom pairs so that
38-
the collective variable is computed only for pairs of atoms separated by a distance
39-
less than :math:`r_c`. A switching function is also applied to smoothly turn off
40-
the collective variable starting from a distance :math:`r_s < r_c`.
53+
The terms outside the summations guarantee that shortest distance is well-defined
54+
even when there are no atom pairs within the cutoff distance. In this case, the
55+
collective variable evaluates to the cutoff distance.
4156
4257
.. note::
4358
@@ -52,17 +67,12 @@ class ShortestDistance(CollectiveVariable, openmm.CustomCVForce):
5267
The indices of the atoms in the second group.
5368
numAtoms
5469
The total number of atoms in the system, including those that are not in any
55-
of the groups.
56-
sigma
70+
of the groups. This is necessary for the correct initialization of the
71+
collective variable.
72+
beta
5773
The parameter that controls the degree of approximation.
58-
magnitude
59-
The expected order of magnitude of the shortest distance. This parameter does
60-
not affect the value of the collective variable, but a good estimate can
61-
improve numerical stability.
6274
cutoffDistance
6375
The cutoff distance for evaluating atom pairs.
64-
switchDistance
65-
The distance at which the switching function starts to be applied.
6676
pbc
6777
Whether to consider periodic boundary conditions in distance calculations.
6878
name
@@ -84,7 +94,7 @@ class ShortestDistance(CollectiveVariable, openmm.CustomCVForce):
8494
... coords[group1, None, :] - coords[None, group2, :],
8595
... axis=-1,
8696
... ).min()
87-
0.17573...
97+
0.1757...
8898
>>> num_atoms = model.system.getNumParticles()
8999
>>> min_dist = cvpack.ShortestDistance(group1, group2, num_atoms)
90100
>>> min_dist.addToSystem(model.system)
@@ -93,38 +103,34 @@ class ShortestDistance(CollectiveVariable, openmm.CustomCVForce):
93103
>>> context = openmm.Context(model.system, integrator, platform)
94104
>>> context.setPositions(model.positions)
95105
>>> min_dist.getValue(context)
96-
0.17573... nm
106+
0.1785... nm
97107
"""
98108

99109
def __init__(
100110
self,
101111
group1: t.Sequence[int],
102112
group2: t.Sequence[int],
103113
numAtoms: int,
104-
sigma: mmunit.Quantity = Quantity(0.01 * mmunit.nanometers),
105-
magnitude: mmunit.Quantity = Quantity(0.2 * mmunit.nanometers),
106-
cutoffDistance: mmunit.Quantity = Quantity(0.5 * mmunit.nanometers),
107-
switchDistance: mmunit.Quantity = Quantity(0.4 * mmunit.nanometers),
114+
beta: float = 100,
115+
cutoffDistance: mmunit.Quantity = Quantity(1 * mmunit.nanometers),
108116
pbc: bool = True,
109117
name: str = "shortest_distance",
110118
) -> None:
111-
if mmunit.is_quantity(sigma):
112-
sigma = sigma.value_in_unit(mmunit.nanometers)
113-
if mmunit.is_quantity(magnitude):
114-
magnitude = magnitude.value_in_unit(mmunit.nanometers)
115-
weight = f"exp(-0.5*(r^2 - {magnitude**2})/{sigma**2})"
119+
rc = cutoffDistance
120+
if mmunit.is_quantity(rc):
121+
rc = rc.value_in_unit(mmunit.nanometers)
122+
weight = f"exp(-{beta/rc}*r)"
116123
forces = {
117-
"numerator": openmm.CustomNonbondedForce(f"r*{weight}"),
118-
"denominator": openmm.CustomNonbondedForce(weight),
124+
"wrsum": openmm.CustomNonbondedForce(f"{weight}*r"),
125+
"wsum": openmm.CustomNonbondedForce(weight),
119126
}
120-
super().__init__("numerator/denominator")
127+
super().__init__(f"({rc*np.exp(-beta)}+wrsum)/({np.exp(-beta)}+wsum)")
121128
for cv, force in forces.items():
122129
force.setNonbondedMethod(
123130
force.CutoffPeriodic if pbc else force.CutoffNonPeriodic
124131
)
125132
force.setCutoffDistance(cutoffDistance)
126-
force.setUseSwitchingFunction(True)
127-
force.setSwitchingDistance(switchDistance)
133+
force.setUseSwitchingFunction(False)
128134
force.setUseLongRangeCorrection(False)
129135
for _ in range(numAtoms):
130136
force.addParticle([])
@@ -136,10 +142,8 @@ def __init__(
136142
group1=group1,
137143
group2=group2,
138144
numAtoms=numAtoms,
139-
sigma=sigma,
140-
magnitude=magnitude,
141-
cutoffDistance=cutoffDistance,
142-
switchDistance=switchDistance,
145+
beta=beta,
146+
cutoffDistance=rc,
143147
pbc=pbc,
144148
)
145149

Diff for: cvpack/tests/test_cvpack.py

+11-4
Original file line numberDiff line numberDiff line change
@@ -899,7 +899,9 @@ def test_shortest_distance():
899899
group.extend(atom.index for atom in residue.atoms())
900900
coords = model.positions.value_in_unit(unit.nanometers)
901901
num_atoms = model.system.getNumParticles()
902-
min_dist = cvpack.ShortestDistance(group1, group2, num_atoms)
902+
beta = 100
903+
rc = 1
904+
min_dist = cvpack.ShortestDistance(group1, group2, num_atoms, beta, rc)
903905
min_dist.addToSystem(model.system)
904906
platform = openmm.Platform.getPlatformByName("Reference")
905907
integrator = openmm.LangevinMiddleIntegrator(
@@ -914,8 +916,13 @@ def test_shortest_distance():
914916
)
915917
coords = state.getPositions(asNumpy=True).value_in_unit(unit.nanometers)
916918
cv_value = min_dist.getValue(context) / unit.nanometer
917-
compute_value = np.min(
918-
np.linalg.norm(coords[group1][:, None, :] - coords[group2], axis=2)
919+
distances = np.linalg.norm(
920+
coords[group1][:, None, :] - coords[group2], axis=2
921+
).flatten()
922+
distances = np.append(distances, rc)
923+
exponents = beta * (1 - distances / rc)
924+
computed_value = np.exp(
925+
logsumexp(exponents, b=distances) - logsumexp(exponents)
919926
)
920-
assert cv_value == pytest.approx(compute_value, abs=1e-2)
927+
assert cv_value == pytest.approx(computed_value)
921928
perform_common_tests(min_dist, context)

0 commit comments

Comments
 (0)