Skip to content

Commit 1f05314

Browse files
committed
New observable: TopologicalOrderParameter
1 parent 46431f7 commit 1f05314

File tree

5 files changed

+371
-1
lines changed

5 files changed

+371
-1
lines changed

pytim/observables/__init__.py

+2
Original file line numberDiff line numberDiff line change
@@ -11,6 +11,8 @@
1111
from .basic_observables import Position, RelativePosition, Velocity, Force
1212
from .basic_observables import Number, Mass, Charge, NumberOfResidues
1313
from .basic_observables import Distance
14+
from .topological_order_parameter import TopologicalOrderParameter
15+
from .adjacency_matrix import DistanceBasedAdjacencyMatrix, WaterHydrogenBondingAdjacencyMatrix
1416
from .local_frame import LocalReferenceFrame, Curvature
1517

1618
from .intrinsic_distance import IntrinsicDistance

pytim/observables/adjacency_matrix.py

+142
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,142 @@
1+
from __future__ import print_function
2+
from abc import ABCMeta, abstractmethod
3+
import numpy as np
4+
import warnings
5+
from scipy import stats, sparse
6+
from scipy.sparse.csgraph import dijkstra
7+
from scipy.spatial import KDTree
8+
from MDAnalysis.core.groups import Atom, AtomGroup, Residue, ResidueGroup
9+
from MDAnalysis.lib import distances
10+
from ..utilities_geometry import minimum_image
11+
from .observable import Observable
12+
13+
14+
__metaclass__ = ABCMeta
15+
16+
class AdjacencyMatrix(Observable):
17+
"""
18+
Base class for defining different bonding criteria.
19+
"""
20+
@property
21+
@abstractmethod
22+
def base_group(self):
23+
# TODO: comment
24+
pass
25+
26+
class DistanceBasedAdjacencyMatrix(AdjacencyMatrix):
27+
"""
28+
Simple distance-based bonding criterion working on a per-atom basis.
29+
"""
30+
def __init__(self, cutoff_distance=3.5, pbc=True):
31+
self.cutoff_distance = cutoff_distance
32+
self.pbc = pbc
33+
self.group = None
34+
35+
def compute(self, group):
36+
"""
37+
Return an adjacency matrix
38+
39+
:param group: AtomGroup to compute the adjacency matrix of.
40+
41+
:return: a (len(group),len(group)) ndarray
42+
"""
43+
if self.pbc: box = group.dimensions[:3]
44+
else: box = None
45+
self.group = group
46+
self.size = (len(group),len(group))
47+
M = sparse.lil_matrix(self.size)*0
48+
pos = group.pack_into_box(inplace=True)
49+
kdtree = KDTree(pos,boxsize=box)
50+
neighs_list = kdtree.query_ball_tree(kdtree, self.cutoff_distance)
51+
for i,neighs in enumerate(neighs_list): M[i,neighs] = 1
52+
self.M = M
53+
return M
54+
55+
@property
56+
def base_group(self): return self.group
57+
58+
59+
class WaterHydrogenBondingAdjacencyMatrix(AdjacencyMatrix):
60+
"""
61+
Hydrogen bond criterion for water.
62+
:param str method: one of
63+
- 'OO' : Based on O...O distance.
64+
- 'OO_OH' : Based on O...O and O...H distance.
65+
- 'OHO' : Based on O...H...O angle.
66+
67+
:param float OO: the donor-acceptor cutoff distance.
68+
Only donor-acceptor pairs within this cutoff are considered
69+
hydrogen-bonded. Default: 3.5 Angstrom. See also the other
70+
requirement, cut_hydrogen_acceptor. Note that requiring two distance
71+
cutoffs is practically equivalent to requiring a distance and an angle
72+
cutoff, if chosen appropriately.
73+
74+
:param float OH: the maximum distance between a hydrogen atom and its
75+
acceptor to be considered hydrogen bonded. Only used when the chosen
76+
method is 'OO_OH'. Default: 2.5 Angstrom
77+
78+
79+
:param float OHO: the maximum angle for two molecules to be considered
80+
h-bonded. Only used when the chosen method is 'OHO'. Default: 30 deg.
81+
"""
82+
83+
def __init__(self, method="OO_OH", OO=3.5, OH=2.5, OHO=30.0):
84+
self.method = method
85+
self.OO = OO # O...O distance cutoff in Å
86+
self.OH = OH # O...H distance cutoff in Å
87+
self.OHO = OHO # O...H...O angle cutoff in degrees
88+
self.oxygens = None
89+
90+
def compute(self, group):
91+
"""
92+
Return an adjacency matrix.
93+
94+
:param group: AtomGroup to compute the adjacency matrix of.
95+
96+
:return: a (nmol, nmol) ndarray with nmol the number of water
97+
molecules in group.
98+
"""
99+
# we do a soft check using the residue name
100+
assert np.all(group.residues.resnames==group.residues.resnames[0]), 'the group should be homogeneous'
101+
_ = group.pack_into_box(inplace=True)
102+
# number of hydrogens in each donor molecule
103+
self.box = group.dimensions[:3]
104+
self.size = (len(group.residues),len(group.residues))
105+
self.M = sparse.lil_matrix(self.size)*0
106+
self.group = group
107+
self.oxygens = group.atoms[group.atoms.types == 'O']
108+
self.hydrogens = group.atoms[group.atoms.types == 'H']
109+
if self.method == 'OO_OH': return self._adjacency_matrix_OO_OH()
110+
else: raise NotImplemented("methods OO and OHO not implemented")
111+
112+
def _adjacency_matrix_OO_OH(self):
113+
tree = KDTree(self.oxygens.positions,boxsize=self.box)
114+
OOlist = tree.query_ball_tree(tree, self.OO)
115+
# the observable is based on residues, so the max size of the matrix will never be larger than
116+
# the number of residues. We build a lookup table to simplify handling them: to obtain the sequential
117+
# index within the oxygen group, pass its resindex to the table.
118+
# resindices is zero-based, so let's add one. Any residue not in the lookuptable will return -1
119+
lookup = -1 + np.zeros(1+np.max(self.oxygens.resindices),dtype=int)
120+
lookup[self.oxygens.resindices] = np.arange(len(self.oxygens))
121+
for i,neighs in enumerate(OOlist): # for each of the acceptors, the list of neighboring donors
122+
neigh_As = self.oxygens[[a for a in neighs if a != i ]] # group of donors neighbor to the i-th acceptor, exclude self.
123+
sel = neigh_As.residues.atoms # all the atoms in the molecules of the neighboring donors of the i-th acceptor
124+
# displacement between i-th acceptor and all hydrogens in the neighboring donor molecules
125+
disp = minimum_image(sel[sel.types=='H'].positions - self.oxygens[i].position,self.box)
126+
#correpsonding Acceptor-Hydrogen (AH) distances
127+
AH_d = np.linalg.norm(disp,axis=1)
128+
# AH_d: we reshape to extract the closest of the 2 Hs to the acceptor
129+
AH_d = np.min(AH_d.reshape((AH_d.shape[0]//2,2)),axis=1)
130+
# HB criterion on AH distance
131+
cond = np.argwhere(AH_d<self.OH).flatten().tolist()
132+
# Update the list of HB pairs: we use the lookup table
133+
b = lookup[self.oxygens[i].resindex]
134+
assert b==i, "b is not i"+str(i)+' '+str(self.oxygens[i].resindex)+ " b="+str(b)
135+
# remember, neigh_As is the list of donors in the neighborhood of the i-th acceptor
136+
for a in np.unique(lookup[neigh_As[cond].resindices]).tolist(): self.M[b,a], self.M[a,b] = 1,1
137+
return self.M
138+
139+
@property
140+
def base_group(self): return self.oxygens
141+
142+
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,209 @@
1+
# -*- Mode: python; tab-width: 4; indent-tabs-mode:nil; coding: utf-8 -*-
2+
# vim: tabstop=4 expandtab shiftwidth=4 softtabstop=4
3+
""" Module: topological_order_parameter
4+
===================================
5+
6+
Implements the topological order parameter Psi as introduced in
7+
"Correlated Fluctuations of Structural Indicators Close to the
8+
Liquid−Liquid Transition in Supercooled Water",
9+
R. Foffi and F. Sciortino, J. Phys. Chem. B 2023, 127, 378−386
10+
DOI: 10.1021/acs.jpcb.2c07169
11+
12+
"""
13+
from __future__ import print_function
14+
from abc import ABCMeta, abstractmethod
15+
import numpy as np
16+
import warnings
17+
from scipy import stats, sparse
18+
from scipy.sparse.csgraph import dijkstra
19+
from scipy.spatial import KDTree
20+
from MDAnalysis.core.groups import Atom, AtomGroup, Residue, ResidueGroup
21+
from MDAnalysis.lib import distances
22+
from ..utilities_geometry import minimum_image
23+
from .observable import Observable
24+
25+
class TopologicalOrderParameter(Observable):
26+
r"""Bond Topological order parameter. This observable provides
27+
also several helper functions and stores additional information that
28+
is useful for analyzing the topological properties of a set of molecules.
29+
"""
30+
31+
def __init__(self, bonding_criterion, topological_distance=4):
32+
r"""
33+
:param BondingCriterion bonding_criterion: instance of a class to determine whether
34+
atoms in two groups are bonded or not.
35+
:param int topological_distance: which topological distance to consider when
36+
computing the order parameter
37+
38+
"""
39+
40+
Observable.__init__(self, None)
41+
self.bonding_criterion = bonding_criterion
42+
self.topological_distance = topological_distance
43+
self.M = None
44+
45+
def compute(self, inp ,kargs=None):
46+
r"""
47+
Compute the order parameters and Euclidean distances
48+
for atoms in group separated by the chosen topologigal
49+
distance
50+
51+
:param AtomGroup inp: an AtomGroup of including donors,
52+
acceptors and hydrogen
53+
54+
Returns:
55+
56+
psi : an ndarray containing the order parameters
57+
of each atom in group and an array with the indices
58+
of nearest neighbors at the chosen topological
59+
distance.
60+
61+
Note that the observable stores also other information
62+
in the class:
63+
64+
- neighbors: the index of the nearest neighbor at topological
65+
distance N. The distances to these neighbors are the values
66+
of psi.
67+
68+
- dists: the Euclidean distance to all neighbors at the
69+
chosen topological distance N
70+
71+
- predecessors: the matrix of predecessors. This can be
72+
used to obtain the sequence of molecules at topological
73+
distance 1,2, ..., N.
74+
75+
- pairs: the lists of pairs at topological distance N
76+
77+
- dist: the topological distances
78+
79+
Example:
80+
81+
>>> import numpy as np
82+
>>> import MDAnalysis as mda
83+
>>> import pytim
84+
>>> from pytim.observables import TopologicalOrderParameter,WaterHydrogenBondingAdjacencyMatrix
85+
>>> from pytim.datafiles import WATER_GRO
86+
>>>
87+
>>> u = mda.Universe(WATER_GRO)
88+
>>> pos = u.atoms.positions[:]
89+
>>> # First, we need to define the criterion
90+
>>> # used to define the bonding
91+
>>> HB = WaterHydrogenBondingAdjacencyMatrix()
92+
>>> # we create a molecule in the vapor phase,
93+
>>> # its value of the order parameter must be inf
94+
>>> pos[u.atoms.resindices==0]+=75
95+
>>> u.atoms.positions=pos
96+
>>> g = u.select_atoms('name OW')
97+
>>> psi = TopologicalOrderParameter(HB,4)
98+
>>> with np.printoptions(precision=4):
99+
... print(psi.compute(g)[:4])
100+
[ inf 3.2736 3.9129 3.5797]
101+
102+
"""
103+
104+
M = self.bonding_criterion.compute(inp)
105+
# we might use directed=True to speedup because we computed the matrix symmetrically
106+
res = dijkstra(M,
107+
unweighted=True,
108+
limit=self.topological_distance,
109+
directed=True,
110+
return_predecessors=True)
111+
self.dist,self.predecessors = res
112+
113+
pairs = np.argwhere(self.dist==self.topological_distance)
114+
g = self.bonding_criterion.base_group
115+
dists = np.linalg.norm(minimum_image(g[pairs[:,0]].positions - g[pairs[:,1]].positions,
116+
self.bonding_criterion.group.universe.dimensions[:3]),axis=1)
117+
self.dists = dists
118+
# in some case (gas) a molecule might have no 4th neighbor, so we initialize
119+
# psi with inf.
120+
neighbors,psi=-np.ones(len(g),dtype=int),np.inf*np.ones(len(g),dtype=float)
121+
122+
for i in np.arange(len(g)):
123+
# pairs has symmetric entries (eg if [1,2] is present, also [2,1] will)
124+
# because of how M has been built. These are the pairs that are connected
125+
# by a given topological_distance. Now we take the first column, of pairs
126+
# and we are sure, because of symmetry, that these are all the indices of
127+
# molecules that have a linked neighbour. We check for each of the donors
128+
# and acceptors if they are in this list.
129+
cond = np.where(pairs[:,0]==i)[0]
130+
if len(cond) > 0: # cond gives the list of donors or acceptors connected to
131+
# the i-th molecule of the donor+acceptor group.
132+
# select the specific pair with minimum distance
133+
amin = np.argmin(dists[cond])
134+
# the neighbor is the second column in the pair where the distance is
135+
# minimum, and the associated value of psi is the minimum distance
136+
# itself.
137+
neighbors[i],psi[i] = pairs[cond,1][amin], np.min(dists[cond])
138+
# otherwise, there is no 4th neighbor, we keep the default value (inf)
139+
self.psi = psi
140+
# so far neighbors and pairs are indices of the local array, now we store
141+
# them as MDAnalysis indices. To be able to search for paths, we store also
142+
# the indices of the local arrau
143+
self.neighbors = g.indices[neighbors]
144+
self._neighbors_localindex = neighbors
145+
self.pairs = self.bonding_criterion.base_group.indices[pairs]
146+
return self.psi
147+
148+
def path(self, start, end):
149+
r"""
150+
Given the start and end of a path, returns the
151+
molecules within the connected path, if it exists.
152+
Accepts indices (zero-based) of the input group
153+
used to compute the observable, or two Atom objects,
154+
in which case it identifies the molecules to which they
155+
belong and returns an AtomGroup.
156+
157+
:param int|Atom : start, the first atom in the chain
158+
:param int|Atom : end, the last atom in the chain
159+
"""
160+
try: self.psi
161+
except (AttributeError): warnings.warn("The order parameter has not yet been computed",category=UserWarning)
162+
ret_type = list
163+
if isinstance(start,Atom) and isinstance(end,Atom):
164+
ret_type=AtomGroup
165+
if not start in self.group: raise ValueError("the start Atom does not belong to the input group")
166+
if not end in self.group: raise ValueError("the end Atom does not belong to input group")
167+
168+
if ret_type == AtomGroup:
169+
start_id, end_id = start.resindex, end.resindex
170+
else:
171+
start_id, end_id = start, end
172+
p = [end_id]
173+
try:
174+
while p[-1] != start_id:
175+
p.append(self.predecessors[start_id,p[-1]] )
176+
except:
177+
return []
178+
p = np.array(p).tolist()
179+
if ret_type == AtomGroup: return self.bonding_criterion.base_group[p]
180+
else: return p
181+
182+
def path_to_nearest_topological_neighbor(self,start):
183+
r"""
184+
Given a starting atom belonging to the input group
185+
of the observable, return the set of connected neighbors
186+
leading to the nearest neighbor at the chosen
187+
topological distance.
188+
189+
:param int|Atom : start, the first atom in the chain.
190+
If an Atom is passed, this must be part
191+
of either the acceptor or donor groups.
192+
If an index is passed, this is interpreted
193+
as the residue index of one of the acceptor
194+
or donor group members.
195+
196+
"""
197+
try: self.psi
198+
except (AttributeError): warnings.warn("The order parameter has not yet been computed",category=UserWarning)
199+
if isinstance(start, Atom):
200+
if not start in self.bonding_criterion.base_group: raise ValueError("the start Atom does not belong to either of the acceptor or donor groups")
201+
start_id = start.resindex
202+
else:
203+
start_id = start
204+
end_id = self._neighbors_localindex[start_id]
205+
path = self.path(start_id,end_id)
206+
if isinstance(start, Atom): return self.bonding_criterion.base_group[path]
207+
else: return path
208+
209+

pytim/utilities.py

+1-1
Original file line numberDiff line numberDiff line change
@@ -284,4 +284,4 @@ def generate_cube_vertices(box, delta=0.0, jitter=False, dim=3):
284284

285285
return vertices
286286

287-
#
287+

pytim/utilities_geometry.py

+17
Original file line numberDiff line numberDiff line change
@@ -210,3 +210,20 @@ def EulerRotation(phi, theta, psi):
210210
R2 = [-cth * cps * sph - cph * sps, cth * cps * cph - sph * sps, cps * sth]
211211
R3 = [sth * sph, -cph * sth, cth]
212212
return np.array([R1, R2, R3])
213+
214+
def minimum_image(disp, box):
215+
""" wraps displacement vectors so that they satisfy the minimum image
216+
convention
217+
218+
:param ndarray disp: an (N,3) array of displacements
219+
:param ndarray box: (3,) array with the rectangular box edges' length
220+
:return ndarray: the wrapped displacements
221+
222+
"""
223+
cond = np.where(disp>box[:3]/2)
224+
disp[cond] = disp[cond] - box[cond[1]]
225+
cond = np.where(disp<-box[:3]/2)
226+
disp[cond] = disp[cond] + box[cond[1]]
227+
return disp
228+
229+

0 commit comments

Comments
 (0)