-
Notifications
You must be signed in to change notification settings - Fork 16
/
Copy pathpolymerutils.py
224 lines (172 loc) · 7.62 KB
/
polymerutils.py
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
"""
Loading and saving individual conformations
===========================================
The module :py:mod:`polychrom.polymerutils` provides tools for saving and loading individual conformations. Note that
saving and loading trajectories should generally be done using :py:mod:`polychrom.hdf5_format` module. This module
provides tools for loading/saving invividual conformations, or for working with projects that have both old-style
and new-style trajectories.
For projects using both old-style and new-style trajectories(e.g. in a project that was switched to polychrom,
and new files were added), a function :py:func:`polychrom.polymerutils.fetch_block` can be helpful as it provides the
same interface for fetching a conformation from both old-style and new-style trajectory. Note however that it is not
the fastest way to iterate over conformations in the new-style trajectory, and the
:py:func:`polychrom.hdf5_format.list_URIs` is faster.
A typical workflow with the new-style trajectories should be:
.. code-block:: python
URIs = polychrom.hdf5_format.list_URIs(folder)
for URI in URIs:
data = polychrom.hdf5_format.load_URI(URI)
xyz = data["pos"]
"""
from __future__ import absolute_import, division, print_function, unicode_literals
import glob
import io
import os
import joblib
import numpy as np
import six
from polychrom.hdf5_format import load_URI
from . import hdf5_format
def load(filename):
"""Universal load function for any type of data file It always returns just XYZ
positions - use fetch_block or hdf5_format.load_URI for loading the whole metadata
Accepted file types
-------------------
New-style URIs (HDF5 based storage)
Text files in openmm-polymer format
joblib files in openmm-polymer format
Parameters
----------
filename: str
filename to load or a URI
"""
if "::" in filename:
return hdf5_format.load_URI(filename)["pos"]
if not os.path.exists(filename):
raise IOError("File not found :( \n %s" % filename)
try: # loading from a joblib file here
return dict(joblib.load(filename)).pop("data")
except Exception: # checking for a text file
data_file = open(filename)
line0 = data_file.readline()
try:
N = int(line0)
except (ValueError, UnicodeDecodeError):
raise TypeError("Could not read the file. Not text or joblib.")
data = [list(map(float, i.split())) for i in data_file.readlines()]
if len(data) != N:
raise ValueError("N does not correspond to the number of lines!")
return np.array(data)
def fetch_block(folder, ind, full_output=False):
"""
A more generic function to fetch block number "ind" from a trajectory in a folder
This function is useful both if you want to load both "old style" trajectories (block1.dat),
and "new style" trajectories ("blocks_1-50.h5")
It will be used in files "show"
Parameters
----------
folder: str, folder with a trajectory
ind: str or int, number of a block to fetch
full_output: bool (default=False)
If set to true, outputs a dict with positions, eP, eK, time etc.
if False, outputs just the conformation
(relevant only for new-style URIs, so default is False)
Returns
-------
data, Nx3 numpy array
if full_output==True, then dict with data and metadata; XYZ is under key "pos"
"""
blocksh5 = glob.glob(os.path.join(folder, "blocks*.h5"))
blocksdat = glob.glob(os.path.join(folder, "block*.dat"))
ind = int(ind)
if (len(blocksh5) > 0) and (len(blocksdat) > 0):
raise ValueError("both .h5 and .dat files found in folder - exiting")
if (len(blocksh5) == 0) and (len(blocksdat) == 0):
raise ValueError("no blocks found")
if len(blocksh5) > 0:
fnames = [os.path.split(i)[-1] for i in blocksh5]
inds = [i.split("_")[-1].split(".")[0].split("-") for i in fnames]
exists = [(int(i[0]) <= ind) and (int(i[1]) >= ind) for i in inds]
if True not in exists:
raise ValueError(f"block {ind} not found in files")
if exists.count(True) > 1:
raise ValueError("Cannot find the file uniquely: names are wrong")
pos = exists.index(True)
block = load_URI(blocksh5[pos] + f"::{ind}")
if not full_output:
block = block["pos"]
if len(blocksdat) > 0:
block = load(os.path.join(folder, f"block{ind}.dat"))
return block
def save(data, filename, mode="txt", pdbGroups=None):
"""
Basically unchanged polymerutils.save function from openmm-polymer
It can save into txt or joblib formats used by old openmm-polymer
It is also very useful for saving files to PDB format to make them compatible
with nglview, pymol_show and others
"""
data = np.asarray(data, dtype=np.float32)
if mode.lower() == "joblib":
joblib.dump({"data": data}, filename=filename, compress=9)
return
if mode.lower() == "txt":
lines = [str(len(data)) + "\n"]
for particle in data:
lines.append("{0:.3f} {1:.3f} {2:.3f}\n".format(*particle))
if filename is None:
return lines
elif isinstance(filename, six.string_types):
with open(filename, "w") as myfile:
myfile.writelines(lines)
elif hasattr(filename, "writelines"):
filename.writelines(lines)
else:
raise ValueError("Not sure what to do with filename {0}".format(filename))
elif mode == "pdb":
data = data - np.minimum(np.min(data, axis=0), np.zeros(3, float) - 100)[None, :]
retret = ""
def add(st, n):
if len(st) > n:
return st[:n]
else:
return st + " " * (n - len(st))
if pdbGroups is None:
pdbGroups = ["A" for i in range(len(data))]
else:
pdbGroups = [str(int(i)) for i in pdbGroups]
for i, line, group in zip(list(range(len(data))), data, pdbGroups):
atomNum = (i + 1) % 9000
segmentNum = (i + 1) // 9000 + 1
line = [float(j) for j in line]
ret = add("ATOM", 6)
ret = add(ret + "{:5d}".format(atomNum), 11)
ret = ret + " "
ret = add(ret + "CA", 17)
ret = add(ret + "ALA", 21)
ret = add(ret + group[0], 22)
ret = add(ret + str(atomNum), 26)
ret = add(ret + " ", 30)
# ret = add(ret + "%i" % (atomNum), 30)
ret = add(ret + ("%8.3f" % line[0]), 38)
ret = add(ret + ("%8.3f" % line[1]), 46)
ret = add(ret + ("%8.3f" % line[2]), 54)
ret = add(ret + (" 1.00"), 61)
ret = add(ret + str(float(i % 8 > 4)), 67)
ret = add(ret, 73)
ret = add(ret + str(segmentNum), 77)
retret += ret + "\n"
with open(filename, "w") as f:
f.write(retret)
f.flush()
elif mode == "pyxyz":
with open(filename, "w") as f:
for i in data:
filename.write("C {0} {1} {2}".format(*i))
else:
raise ValueError("Unknown mode : %s, use h5dict, joblib, txt or pdb" % mode)
def rotation_matrix(rotate):
"""Calculates rotation matrix based on three rotation angles"""
tx, ty, tz = rotate
Rx = np.array([[1, 0, 0], [0, np.cos(tx), -np.sin(tx)], [0, np.sin(tx), np.cos(tx)]])
Ry = np.array([[np.cos(ty), 0, -np.sin(ty)], [0, 1, 0], [np.sin(ty), 0, np.cos(ty)]])
Rz = np.array([[np.cos(tz), -np.sin(tz), 0], [np.sin(tz), np.cos(tz), 0], [0, 0, 1]])
return np.dot(Rx, np.dot(Ry, Rz))