pdbfile.py 14 KB
Newer Older
1
2
"""
pdbfile.py: Used for loading PDB files.
3
4
5
6
7
8
9
10
11
12

This is part of the OpenMM molecular simulation toolkit originating from
Simbios, the NIH National Center for Physics-Based Simulation of
Biological Structures at Stanford, funded under the NIH Roadmap for
Medical Research, grant U54 GM072970. See https://simtk.org.

Portions copyright (c) 2012 Stanford University and the Authors.
Authors: Peter Eastman
Contributors:

Justin MacCallum's avatar
Justin MacCallum committed
13
Permission is hereby granted, free of charge, to any person obtaining a
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
copy of this software and associated documentation files (the "Software"),
to deal in the Software without restriction, including without limitation
the rights to use, copy, modify, merge, publish, distribute, sublicense,
and/or sell copies of the Software, and to permit persons to whom the
Software is furnished to do so, subject to the following conditions:

The above copyright notice and this permission notice shall be included in
all copies or substantial portions of the Software.

THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
THE AUTHORS, CONTRIBUTORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM,
DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR
OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE
USE OR OTHER DEALINGS IN THE SOFTWARE.
30
31
32
33
34
35
"""
__author__ = "Peter Eastman"
__version__ = "1.0"

import os
import sys
36
import math
37
38
import xml.etree.ElementTree as etree
from copy import copy
39
40
from datetime import date
from simtk.openmm import Vec3, Platform
41
42
from simtk.openmm.app.internal.pdbstructure import PdbStructure
from simtk.openmm.app import Topology
Peter Eastman's avatar
Peter Eastman committed
43
from simtk.unit import nanometers, angstroms, is_quantity, norm, Quantity
44
import element as elem
45
46
47
48
49
try:
    import numpy
except:
    pass

50
51
class PDBFile(object):
    """PDBFile parses a Protein Data Bank (PDB) file and constructs a Topology and a set of atom positions from it.
Justin MacCallum's avatar
Justin MacCallum committed
52

53
54
55
    This class also provides methods for creating PDB files.  To write a file containing a single model, call
    writeFile().  You also can create files that contain multiple models.  To do this, first call writeHeader(),
    then writeModel() once for each model in the file, and finally writeFooter() to complete the file."""
Justin MacCallum's avatar
Justin MacCallum committed
56

57
58
    _residueNameReplacements = {}
    _atomNameReplacements = {}
Justin MacCallum's avatar
Justin MacCallum committed
59

60
61
    def __init__(self, file):
        """Load a PDB file.
Justin MacCallum's avatar
Justin MacCallum committed
62

63
        The atom positions and Topology can be retrieved by calling getPositions() and getTopology().
Justin MacCallum's avatar
Justin MacCallum committed
64

65
66
67
68
        Parameters:
         - file (string) the name of the file to load
        """
        top = Topology()
69
        ## The Topology read from the PDB file
70
        self.topology = top
Justin MacCallum's avatar
Justin MacCallum committed
71

72
        # Load the PDB file
Justin MacCallum's avatar
Justin MacCallum committed
73

peastman's avatar
peastman committed
74
75
76
77
78
79
80
        if isinstance(file, PdbStructure):
            pdb = file
        else:
            inputfile = file
            if isinstance(file, str):
                inputfile = open(file)
            pdb = PdbStructure(inputfile, load_all_models=True)
81
82
83
84
        PDBFile._loadNameReplacementTables()

        # Build the topology

85
        atomByNumber = {}
86
87
88
89
90
91
92
93
94
95
96
        for chain in pdb.iter_chains():
            c = top.addChain()
            for residue in chain.iter_residues():
                resName = residue.get_name()
                if resName in PDBFile._residueNameReplacements:
                    resName = PDBFile._residueNameReplacements[resName]
                r = top.addResidue(resName, c)
                if resName in PDBFile._atomNameReplacements:
                    atomReplacements = PDBFile._atomNameReplacements[resName]
                else:
                    atomReplacements = {}
97
                for atom in residue.atoms:
98
99
100
101
                    atomName = atom.get_name()
                    if atomName in atomReplacements:
                        atomName = atomReplacements[atomName]
                    atomName = atomName.strip()
102
103
104
                    element = atom.element
                    if element is None:
                        # Try to guess the element.
Justin MacCallum's avatar
Justin MacCallum committed
105

106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
                        upper = atomName.upper()
                        if upper.startswith('CL'):
                            element = elem.chlorine
                        elif upper.startswith('NA'):
                            element = elem.sodium
                        elif upper.startswith('MG'):
                            element = elem.magnesium
                        elif upper.startswith('BE'):
                            element = elem.beryllium
                        elif upper.startswith('LI'):
                            element = elem.lithium
                        elif upper.startswith('K'):
                            element = elem.potassium
                        elif( len( residue ) == 1 and upper.startswith('CA') ):
                            element = elem.calcium
                        else:
                            try:
                                element = elem.get_by_symbol(atomName[0])
                            except KeyError:
                                pass
126
127
                    newAtom = top.addAtom(atomName, element, r)
                    atomByNumber[atom.serial_number] = newAtom
128
129
130
        self._positions = []
        for model in pdb.iter_models(True):
            coords = []
131
132
133
134
135
            for chain in model.iter_chains():
                for residue in chain.iter_residues():
                    for atom in residue.atoms:
                        pos = atom.get_position().value_in_unit(nanometers)
                        coords.append(Vec3(pos[0], pos[1], pos[2]))
136
137
138
            self._positions.append(coords*nanometers)
        ## The atom positions read from the PDB file.  If the file contains multiple frames, these are the positions in the first frame.
        self.positions = self._positions[0]
139
140
141
        self.topology.setUnitCellDimensions(pdb.get_unit_cell_dimensions())
        self.topology.createStandardBonds()
        self.topology.createDisulfideBonds(self.positions)
142
        self._numpyPositions = None
Justin MacCallum's avatar
Justin MacCallum committed
143

144
        # Add bonds based on CONECT records.
Justin MacCallum's avatar
Justin MacCallum committed
145

146
        connectBonds = []
Peter Eastman's avatar
Bug fix  
Peter Eastman committed
147
        for connect in pdb.models[0].connects:
148
149
            i = connect[0]
            for j in connect[1:]:
peastman's avatar
peastman committed
150
151
                if i in atomByNumber and j in atomByNumber:
                    connectBonds.append((atomByNumber[i], atomByNumber[j]))
152
153
154
155
156
157
158
        if len(connectBonds) > 0:
            # Only add bonds that don't already exist.
            existingBonds = set(top.bonds())
            for bond in connectBonds:
                if bond not in existingBonds and (bond[1], bond[0]) not in existingBonds:
                    top.addBond(bond[0], bond[1])
                    existingBonds.add(bond)
159
160
161
162

    def getTopology(self):
        """Get the Topology of the model."""
        return self.topology
Justin MacCallum's avatar
Justin MacCallum committed
163

164
165
166
    def getNumFrames(self):
        """Get the number of frames stored in the file."""
        return len(self._positions)
Justin MacCallum's avatar
Justin MacCallum committed
167

168
    def getPositions(self, asNumpy=False, frame=0):
169
        """Get the atomic positions.
Justin MacCallum's avatar
Justin MacCallum committed
170

171
172
        Parameters:
         - asNumpy (boolean=False) if true, the values are returned as a numpy array instead of a list of Vec3s
173
         - frame (int=0) the index of the frame for which to get positions
174
175
176
         """
        if asNumpy:
            if self._numpyPositions is None:
177
178
179
180
181
                self._numpyPositions = [None]*len(self._positions)
            if self._numpyPositions[frame] is None:
                self._numpyPositions[frame] = Quantity(numpy.array(self._positions[frame].value_in_unit(nanometers)), nanometers)
            return self._numpyPositions[frame]
        return self._positions[frame]
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

    @staticmethod
    def _loadNameReplacementTables():
        """Load the list of atom and residue name replacements."""
        if len(PDBFile._residueNameReplacements) == 0:
            tree = etree.parse(os.path.join(os.path.dirname(__file__), 'data', 'pdbNames.xml'))
            allResidues = {}
            proteinResidues = {}
            nucleicAcidResidues = {}
            for residue in tree.getroot().findall('Residue'):
                name = residue.attrib['name']
                if name == 'All':
                    PDBFile._parseResidueAtoms(residue, allResidues)
                elif name == 'Protein':
                    PDBFile._parseResidueAtoms(residue, proteinResidues)
                elif name == 'Nucleic':
                    PDBFile._parseResidueAtoms(residue, nucleicAcidResidues)
            for atom in allResidues:
                proteinResidues[atom] = allResidues[atom]
                nucleicAcidResidues[atom] = allResidues[atom]
            for residue in tree.getroot().findall('Residue'):
                name = residue.attrib['name']
                for id in residue.attrib:
                    if id == 'name' or id.startswith('alt'):
                        PDBFile._residueNameReplacements[residue.attrib[id]] = name
                if 'type' not in residue.attrib:
                    atoms = copy(allResidues)
                elif residue.attrib['type'] == 'Protein':
                    atoms = copy(proteinResidues)
                elif residue.attrib['type'] == 'Nucleic':
                    atoms = copy(nucleicAcidResidues)
                else:
                    atoms = copy(allResidues)
                PDBFile._parseResidueAtoms(residue, atoms)
                PDBFile._atomNameReplacements[name] = atoms

    @staticmethod
    def _parseResidueAtoms(residue, map):
        for atom in residue.findall('Atom'):
            name = atom.attrib['name']
            for id in atom.attrib:
                map[atom.attrib[id]] = name
Justin MacCallum's avatar
Justin MacCallum committed
224

225
226
227
    @staticmethod
    def writeFile(topology, positions, file=sys.stdout, modelIndex=None):
        """Write a PDB file containing a single model.
Justin MacCallum's avatar
Justin MacCallum committed
228

229
230
231
232
233
234
235
236
237
238
239
240
        Parameters:
         - topology (Topology) The Topology defining the model to write
         - positions (list) The list of atomic positions to write
         - file (file=stdout) A file to write to
        """
        PDBFile.writeHeader(topology, file)
        PDBFile.writeModel(topology, positions, file)
        PDBFile.writeFooter(topology, file)

    @staticmethod
    def writeHeader(topology, file=sys.stdout):
        """Write out the header for a PDB file.
Justin MacCallum's avatar
Justin MacCallum committed
241

242
243
244
245
        Parameters:
         - topology (Topology) The Topology defining the molecular system being written
         - file (file=stdout) A file to write the file to
        """
246
        print >>file, "REMARK   1 CREATED WITH OPENMM %s, %s" % (Platform.getOpenMMVersion(), str(date.today()))
247
248
249
250
        boxSize = topology.getUnitCellDimensions()
        if boxSize is not None:
            size = boxSize.value_in_unit(angstroms)
            print >>file, "CRYST1%9.3f%9.3f%9.3f  90.00  90.00  90.00 P 1           1 " % size
Justin MacCallum's avatar
Justin MacCallum committed
251

252
253
254
    @staticmethod
    def writeModel(topology, positions, file=sys.stdout, modelIndex=None):
        """Write out a model to a PDB file.
Justin MacCallum's avatar
Justin MacCallum committed
255

256
257
258
259
260
261
262
        Parameters:
         - topology (Topology) The Topology defining the model to write
         - positions (list) The list of atomic positions to write
         - file (file=stdout) A file to write the model to
         - modelIndex (int=None) If not None, the model will be surrounded by MODEL/ENDMDL records with this index
        """
        if len(list(topology.atoms())) != len(positions):
Justin MacCallum's avatar
Justin MacCallum committed
263
            raise ValueError('The number of positions must match the number of atoms')
264
265
        if is_quantity(positions):
            positions = positions.value_in_unit(angstroms)
266
267
268
269
        if any(math.isnan(norm(pos)) for pos in positions):
            raise ValueError('Particle position is NaN')
        if any(math.isinf(norm(pos)) for pos in positions):
            raise ValueError('Particle position is infinite')
270
271
272
273
274
        atomIndex = 1
        posIndex = 0
        if modelIndex is not None:
            print >>file, "MODEL     %4d" % modelIndex
        for (chainIndex, chain) in enumerate(topology.chains()):
Peter Eastman's avatar
Peter Eastman committed
275
            chainName = chr(ord('A')+chainIndex%26)
276
277
278
279
280
281
282
            residues = list(chain.residues())
            for (resIndex, res) in enumerate(residues):
                if len(res.name) > 3:
                    resName = res.name[:3]
                else:
                    resName = res.name
                for atom in res.atoms():
283
                    if len(atom.name) < 4 and atom.name[:1].isalpha() and (atom.element is None or len(atom.element.symbol) < 2):
284
285
286
287
288
289
                        atomName = ' '+atom.name
                    elif len(atom.name) > 4:
                        atomName = atom.name[:4]
                    else:
                        atomName = atom.name
                    coords = positions[posIndex]
Robert McGibbon's avatar
Robert McGibbon committed
290
                    line = "ATOM  %5d %-4s %3s %s%4d    %s%s%s  1.00  0.00          %2s  " % (
291
292
                        atomIndex%100000, atomName, resName, chainName,
                        (resIndex+1)%10000, _format_83(coords[0]),
293
294
                        _format_83(coords[1]), _format_83(coords[2]), atom.element.symbol)
                    assert len(line) == 80, 'Fixed width overflow detected'
295
                    print >>file, line
296
297
298
299
300
301
302
303
304
305
306
                    posIndex += 1
                    atomIndex += 1
                if resIndex == len(residues)-1:
                    print >>file, "TER   %5d      %3s %s%4d" % (atomIndex, resName, chainName, resIndex+1)
                    atomIndex += 1
        if modelIndex is not None:
            print >>file, "ENDMDL"

    @staticmethod
    def writeFooter(topology, file=sys.stdout):
        """Write out the footer for a PDB file.
Justin MacCallum's avatar
Justin MacCallum committed
307

308
309
310
311
312
313
        Parameters:
         - topology (Topology) The Topology defining the molecular system being written
         - file (file=stdout) A file to write the file to
        """
        print >>file, "END"

314
315
316
317
318

def _format_83(f):
    """Format a single float into a string of width 8, with ideally 3 decimal
    places of precision. If the number is a little too large, we can
    gracefully degrade the precision by lopping off some of the decimal
Robert McGibbon's avatar
Robert McGibbon committed
319
    places. If it's much too large, we throw a ValueError"""
320
321
322
323
    if -999.999 < f < 9999.999:
        return '%8.3f' % f
    if -9999999 < f < 99999999:
        return ('%8.3f' % f)[:8]
Robert McGibbon's avatar
Robert McGibbon committed
324
325
    raise ValueError('coordinate "%s" could not be represnted '
                     'in a width-8 field' % f)