"platforms/opencl/vscode:/vscode.git/clone" did not exist on "8eeee16de9bc772321ae2b87672700b076913b56"
pdbxfile.py 21.1 KB
Newer Older
1
"""
2
pdbxfile.py: Used for loading PDBx/mmCIF files.
3
4
5
6
7
8

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.

9
Portions copyright (c) 2015-2023 Stanford University and the Authors.
10
Authors: Peter Eastman
11
Contributors: Jason Swails
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30

Permission is hereby granted, free of charge, to any person obtaining a
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.
"""
31
32
from __future__ import division, absolute_import, print_function

33
__author__ = "Peter Eastman"
34
__version__ = "2.0"
35
36

import sys
37
import math
38
from openmm import Vec3, Platform
39
from datetime import date
40
41
42
from openmm.app.internal.pdbx.reader.PdbxReader import PdbxReader
from openmm.app.internal.unitcell import computePeriodicBoxVectors, computeLengthsAndAngles
from openmm.app import Topology, PDBFile
Peter Eastman's avatar
Peter Eastman committed
43
from openmm.unit import nanometers, angstroms, is_quantity, norm, Quantity
44
from . import element as elem
45
46
47
48
49
50
try:
    import numpy
except:
    pass

class PDBxFile(object):
peastman's avatar
peastman committed
51
    """PDBxFile parses a PDBx/mmCIF file and constructs a Topology and a set of atom positions from it."""
52
53
54
55

    def __init__(self, file):
        """Load a PDBx/mmCIF file.

56
57
        The atom positions and Topology can be retrieved by calling
        getPositions() and getTopology().
58

59
60
61
62
63
        Parameters
        ----------
        file : string
            the name of the file to load.  Alternatively you can pass an open
            file object.
64
65
66
67
68
        """
        top = Topology()
        ## The Topology read from the PDBx/mmCIF file
        self.topology = top
        self._positions = []
69
        PDBFile._loadNameReplacementTables()
70
71
72
73

        # Load the file.

        inputFile = file
74
        ownHandle = False
75
76
        if isinstance(file, str):
            inputFile = open(file)
77
            ownHandle = True
78
79
80
        reader = PdbxReader(inputFile)
        data = []
        reader.read(data)
81
82
        if ownHandle:
            inputFile.close()
83
84
85
86
87
        block = data[0]

        # Build the topology.

        atomData = block.getObj('atom_site')
88
        atomNameCol = atomData.getAttributeIndex('auth_atom_id')
peastman's avatar
peastman committed
89
90
        if atomNameCol == -1:
            atomNameCol = atomData.getAttributeIndex('label_atom_id')
91
        atomIdCol = atomData.getAttributeIndex('id')
92
        resNameCol = atomData.getAttributeIndex('auth_comp_id')
peastman's avatar
peastman committed
93
94
        if resNameCol == -1:
            resNameCol = atomData.getAttributeIndex('label_comp_id')
95
        resNumCol = atomData.getAttributeIndex('auth_seq_id')
peastman's avatar
peastman committed
96
97
        if resNumCol == -1:
            resNumCol = atomData.getAttributeIndex('label_seq_id')
peastman's avatar
peastman committed
98
        resInsertionCol = atomData.getAttributeIndex('pdbx_PDB_ins_code')
99
        chainIdCol = atomData.getAttributeIndex('auth_asym_id')
peastman's avatar
peastman committed
100
101
        if chainIdCol == -1:
            chainIdCol = atomData.getAttributeIndex('label_asym_id')
102
103
104
            altChainIdCol = -1
        else:
            altChainIdCol = atomData.getAttributeIndex('label_asym_id')
105
106
107
108
109
110
111
        if altChainIdCol != -1:
            # Figure out which column is best to use for chain IDs.
            
            idSet = set(row[chainIdCol] for row in atomData.getRowList())
            altIdSet = set(row[altChainIdCol] for row in atomData.getRowList())
            if len(altIdSet) > len(idSet):
                chainIdCol, altChainIdCol = (altChainIdCol, chainIdCol)
112
        elementCol = atomData.getAttributeIndex('type_symbol')
113
        altIdCol = atomData.getAttributeIndex('label_alt_id')
114
115
116
117
118
        modelCol = atomData.getAttributeIndex('pdbx_PDB_model_num')
        xCol = atomData.getAttributeIndex('Cartn_x')
        yCol = atomData.getAttributeIndex('Cartn_y')
        zCol = atomData.getAttributeIndex('Cartn_z')
        lastChainId = None
119
        lastAltChainId = None
120
        lastResId = None
121
        lastInsertionCode = ''
122
        atomTable = {}
123
        atomsInResidue = set()
124
125
        models = []
        for row in atomData.getRowList():
126
            atomKey = ((row[resNumCol], row[chainIdCol], row[atomNameCol]))
127
128
129
130
131
            model = ('1' if modelCol == -1 else row[modelCol])
            if model not in models:
                models.append(model)
                self._positions.append([])
            modelIndex = models.index(model)
132
133
134
135
            if row[altIdCol] != '.' and atomKey in atomTable and len(self._positions[modelIndex]) > atomTable[atomKey].index:
                # This row is an alternate position for an existing atom, so ignore it.

                continue
136
137
138
            if modelIndex == 0:
                # This row defines a new atom.

139
140
141
142
143
144
                if resInsertionCol == -1:
                    insertionCode = ''
                else:
                    insertionCode = row[resInsertionCol]
                if insertionCode in ('.', '?'):
                    insertionCode = ''
145
                if lastChainId != row[chainIdCol] or (altChainIdCol != -1 and lastAltChainId != row[altChainIdCol]):
146
                    # The start of a new chain.
147
                    chain = top.addChain(row[chainIdCol])
148
149
                    lastChainId = row[chainIdCol]
                    lastResId = None
150
151
                    if altChainIdCol != -1:
                        lastAltChainId = row[altChainIdCol]
152
                if lastResId != row[resNumCol] or lastChainId != row[chainIdCol] or lastInsertionCode != insertionCode or (lastResId == '.' and row[atomNameCol] in atomsInResidue):
153
                    # The start of a new residue.
peastman's avatar
peastman committed
154
                    resId = (None if resNumCol == -1 else row[resNumCol])
155
                    resIC = insertionCode
156
157
158
159
160
161
162
163
                    resName = row[resNameCol]
                    if resName in PDBFile._residueNameReplacements:
                        resName = PDBFile._residueNameReplacements[resName]
                    res = top.addResidue(resName, chain, resId, resIC)
                    if resName in PDBFile._atomNameReplacements:
                        atomReplacements = PDBFile._atomNameReplacements[resName]
                    else:
                        atomReplacements = {}
164
                    lastResId = row[resNumCol]
165
                    lastInsertionCode = insertionCode
166
                    atomsInResidue.clear()
167
168
169
170
171
                element = None
                try:
                    element = elem.get_by_symbol(row[elementCol])
                except KeyError:
                    pass
172
173
174
175
                atomName = row[atomNameCol]
                if atomName in atomReplacements:
                    atomName = atomReplacements[atomName]
                atom = top.addAtom(atomName, element, res, row[atomIdCol])
176
                atomTable[atomKey] = atom
177
                atomsInResidue.add(atomName)
178
179
180
181
182
183
            else:
                # This row defines coordinates for an existing atom in one of the later models.

                try:
                    atom = atomTable[atomKey]
                except KeyError:
184
                    raise ValueError('Unknown atom %s in residue %s %s for model %s' % (row[atomNameCol], row[resNameCol], row[resNumCol], model))
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
                if atom.index != len(self._positions[modelIndex]):
                    raise ValueError('Atom %s for model %s does not match the order of atoms for model %s' % (row[atomIdCol], model, models[0]))
            self._positions[modelIndex].append(Vec3(float(row[xCol]), float(row[yCol]), float(row[zCol]))*0.1)
        for i in range(len(self._positions)):
            self._positions[i] = self._positions[i]*nanometers
        ## The atom positions read from the PDBx/mmCIF file.  If the file contains multiple frames, these are the positions in the first frame.
        self.positions = self._positions[0]
        self.topology.createStandardBonds()
        self._numpyPositions = None

        # Record unit cell information, if present.

        cell = block.getObj('cell')
        if cell is not None and cell.getRowCount() > 0:
            row = cell.getRow(0)
200
            (a, b, c) = [float(row[cell.getAttributeIndex(attribute)])*0.1 for attribute in ('length_a', 'length_b', 'length_c')]
201
202
            (alpha, beta, gamma) = [float(row[cell.getAttributeIndex(attribute)])*math.pi/180.0 for attribute in ('angle_alpha', 'angle_beta', 'angle_gamma')]
            self.topology.setPeriodicBoxVectors(computePeriodicBoxVectors(a, b, c, alpha, beta, gamma))
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241

        # Add bonds based on struct_conn records.

        connectData = block.getObj('struct_conn')
        if connectData is not None:
            res1Col = connectData.getAttributeIndex('ptnr1_label_seq_id')
            res2Col = connectData.getAttributeIndex('ptnr2_label_seq_id')
            atom1Col = connectData.getAttributeIndex('ptnr1_label_atom_id')
            atom2Col = connectData.getAttributeIndex('ptnr2_label_atom_id')
            asym1Col = connectData.getAttributeIndex('ptnr1_label_asym_id')
            asym2Col = connectData.getAttributeIndex('ptnr2_label_asym_id')
            typeCol = connectData.getAttributeIndex('conn_type_id')
            connectBonds = []
            for row in connectData.getRowList():
                type = row[typeCol][:6]
                if type in ('covale', 'disulf', 'modres'):
                    key1 = (row[res1Col], row[asym1Col], row[atom1Col])
                    key2 = (row[res2Col], row[asym2Col], row[atom2Col])
                    if key1 in atomTable and key2 in atomTable:
                        connectBonds.append((atomTable[key1], atomTable[key2]))
            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)

    def getTopology(self):
        """Get the Topology of the model."""
        return self.topology

    def getNumFrames(self):
        """Get the number of frames stored in the file."""
        return len(self._positions)

    def getPositions(self, asNumpy=False, frame=0):
        """Get the atomic positions.

Robert McGibbon's avatar
Robert McGibbon committed
242
243
244
245
246
247
248
249
        Parameters
        ----------
        asNumpy : bool=False
            if true, the values are returned as a numpy array instead of a list
            of Vec3s
        frame : int=0
            the index of the frame for which to get positions
        """
250
251
252
253
254
255
256
        if asNumpy:
            if self._numpyPositions is None:
                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]
257
258
259
260

    @staticmethod
    def writeFile(topology, positions, file=sys.stdout, keepIds=False,
                  entry=None):
261
        """Write a PDBx/mmCIF file containing a single model.
262

Robert McGibbon's avatar
Robert McGibbon committed
263
264
265
266
267
268
        Parameters
        ----------
        topology : Topology
            The Topology defining the model to write
        positions : list
            The list of atomic positions to write
269
270
        file : string or file
            the name of the file to write.  Alternatively you can pass an open file object.
Robert McGibbon's avatar
Robert McGibbon committed
271
272
273
274
275
276
277
        keepIds : bool=False
            If True, keep the residue and chain IDs specified in the Topology
            rather than generating new ones.  Warning: It is up to the caller to
            make sure these are valid IDs that satisfy the requirements of the
            PDBx/mmCIF format.  Otherwise, the output file will be invalid.
        entry : str=None
            The entry ID to assign to the CIF file
278
        """
279
280
281
282
283
284
        if isinstance(file, str):
            with open(file, 'w') as output:
                PDBxFile.writeFile(topology, positions, output, keepIds, entry)
        else:
            PDBxFile.writeHeader(topology, file, entry, keepIds)
            PDBxFile.writeModel(topology, positions, file, keepIds=keepIds)
285
286

    @staticmethod
peastman's avatar
peastman committed
287
    def writeHeader(topology, file=sys.stdout, entry=None, keepIds=False):
288
        """Write out the header for a PDBx/mmCIF file.
289

Robert McGibbon's avatar
Robert McGibbon committed
290
291
292
293
294
295
296
297
        Parameters
        ----------
        topology : Topology
            The Topology defining the molecular system being written
        file : file=stdout
            A file to write the file to
        entry : str=None
            The entry ID to assign to the CIF file
peastman's avatar
peastman committed
298
299
300
301
302
        keepIds : bool=False
            If True, keep the residue and chain IDs specified in the Topology
            rather than generating new ones.  Warning: It is up to the caller to
            make sure these are valid IDs that satisfy the requirements of the
            PDBx/mmCIF format.  Otherwise, the output file will be invalid.
303
304
305
306
307
308
309
310
311
        """
        if entry is not None:
            print('data_%s' % entry, file=file)
        else:
            print('data_cell', file=file)
        print("# Created with OpenMM %s, %s" % (Platform.getOpenMMVersion(), str(date.today())), file=file)
        print('#', file=file)
        vectors = topology.getPeriodicBoxVectors()
        if vectors is not None:
312
313
            a, b, c, alpha, beta, gamma = computeLengthsAndAngles(vectors)
            RAD_TO_DEG = 180/math.pi
314
315
316
317
318
319
            print('_cell.length_a     %10.4f' % (a*10), file=file)
            print('_cell.length_b     %10.4f' % (b*10), file=file)
            print('_cell.length_c     %10.4f' % (c*10), file=file)
            print('_cell.angle_alpha  %10.4f' % (alpha*RAD_TO_DEG), file=file)
            print('_cell.angle_beta   %10.4f' % (beta*RAD_TO_DEG), file=file)
            print('_cell.angle_gamma  %10.4f' % (gamma*RAD_TO_DEG), file=file)
peastman's avatar
peastman committed
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
            print('#', file=file)

        # Identify bonds that should be listed in the file.

        bonds = []
        for atom1, atom2 in topology.bonds():
            if atom1.residue.name not in PDBFile._standardResidues or atom2.residue.name not in PDBFile._standardResidues:
                bonds.append((atom1, atom2))
            elif atom1.name == 'SG' and atom2.name == 'SG' and atom1.residue.name == 'CYS' and atom2.residue.name == 'CYS':
                bonds.append((atom1, atom2))
        if len(bonds) > 0:

            # Write the bond information.

            print('loop_', file=file)
            print('_struct_conn.id', file=file)
            print('_struct_conn.conn_type_id', file=file)
            print('_struct_conn.ptnr1_label_asym_id', file=file)
            print('_struct_conn.ptnr1_label_comp_id', file=file)
            print('_struct_conn.ptnr1_label_seq_id', file=file)
            print('_struct_conn.ptnr1_label_atom_id', file=file)
            print('_struct_conn.ptnr2_label_asym_id', file=file)
            print('_struct_conn.ptnr2_label_comp_id', file=file)
            print('_struct_conn.ptnr2_label_seq_id', file=file)
            print('_struct_conn.ptnr2_label_atom_id', file=file)
            chainIds = {}
            resIds = {}
            if keepIds:
                for chain in topology.chains():
                    chainIds[chain] = chain.id
                for res in topology.residues():
                    resIds[res] = res.id
            else:
                for (chainIndex, chain) in enumerate(topology.chains()):
                    chainIds[chain] = chr(ord('A')+chainIndex%26)
                    for (resIndex, res) in enumerate(chain.residues()):
                        resIds[res] = resIndex+1
            for i, (atom1, atom2) in enumerate(bonds):
                if atom1.residue.name == 'CYS' and atom2.residue.name == 'CYS':
                    bondType = 'disulf'
                else:
                    bondType = 'covale'
                line = "bond%d %s %s %-4s %5s %-4s %s %-4s %5s %-4s"
                print(line % (i+1, bondType, chainIds[atom1.residue.chain], atom1.residue.name, resIds[atom1.residue], atom1.name,
                              chainIds[atom2.residue.chain], atom2.residue.name, resIds[atom2.residue], atom2.name), file=file)
            print('#', file=file)

        # Write the header for the atom coordinates.

369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
        print('loop_', file=file)
        print('_atom_site.group_PDB', file=file)
        print('_atom_site.id', file=file)
        print('_atom_site.type_symbol', file=file)
        print('_atom_site.label_atom_id', file=file)
        print('_atom_site.label_alt_id', file=file)
        print('_atom_site.label_comp_id', file=file)
        print('_atom_site.label_asym_id', file=file)
        print('_atom_site.label_entity_id', file=file)
        print('_atom_site.label_seq_id', file=file)
        print('_atom_site.pdbx_PDB_ins_code', file=file)
        print('_atom_site.Cartn_x', file=file)
        print('_atom_site.Cartn_y', file=file)
        print('_atom_site.Cartn_z', file=file)
        print('_atom_site.occupancy', file=file)
        print('_atom_site.B_iso_or_equiv', file=file)
        print('_atom_site.Cartn_x_esd', file=file)
        print('_atom_site.Cartn_y_esd', file=file)
        print('_atom_site.Cartn_z_esd', file=file)
        print('_atom_site.occupancy_esd', file=file)
        print('_atom_site.B_iso_or_equiv_esd', file=file)
        print('_atom_site.pdbx_formal_charge', file=file)
        print('_atom_site.auth_seq_id', file=file)
        print('_atom_site.auth_comp_id', file=file)
        print('_atom_site.auth_asym_id', file=file)
        print('_atom_site.auth_atom_id', file=file)
        print('_atom_site.pdbx_PDB_model_num', file=file)

    @staticmethod
    def writeModel(topology, positions, file=sys.stdout, modelIndex=1, keepIds=False):
399
        """Write out a model to a PDBx/mmCIF file.
400

Robert McGibbon's avatar
Robert McGibbon committed
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
        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=1
            The model number of this frame
        keepIds : bool=False
            If True, keep the residue and chain IDs specified in the Topology
            rather than generating new ones.  Warning: It is up to the caller to
            make sure these are valid IDs that satisfy the requirements of the
            PDBx/mmCIF format.  Otherwise, the output file will be invalid.
416
417
418
419
420
        """
        if len(list(topology.atoms())) != len(positions):
            raise ValueError('The number of positions must match the number of atoms')
        if is_quantity(positions):
            positions = positions.value_in_unit(angstroms)
421
        import numpy as np
422
        positions = np.asarray(positions)
423
        if np.isnan(positions).any():
424
            raise ValueError('Particle position is NaN.  For more information, see https://github.com/openmm/openmm/wiki/Frequently-Asked-Questions#nan')
425
        if np.isinf(positions).any():
426
            raise ValueError('Particle position is infinite.  For more information, see https://github.com/openmm/openmm/wiki/Frequently-Asked-Questions#nan')
427
428
        nonHeterogens = PDBFile._standardResidues[:]
        nonHeterogens.remove('HOH')
429
430
431
432
433
434
435
436
437
438
439
        atomIndex = 1
        posIndex = 0
        for (chainIndex, chain) in enumerate(topology.chains()):
            if keepIds:
                chainName = chain.id
            else:
                chainName = chr(ord('A')+chainIndex%26)
            residues = list(chain.residues())
            for (resIndex, res) in enumerate(residues):
                if keepIds:
                    resId = res.id
440
                    resIC = (res.insertionCode if res.insertionCode.strip() else '.')
441
442
                else:
                    resId = resIndex + 1
peastman's avatar
peastman committed
443
                    resIC = '.'
444
445
446
447
                if res.name in nonHeterogens:
                    recordName = "ATOM"
                else:
                    recordName = "HETATM"
448
449
450
451
452
453
                for atom in res.atoms():
                    coords = positions[posIndex]
                    if atom.element is not None:
                        symbol = atom.element.symbol
                    else:
                        symbol = '?'
454
455
                    line = "%s  %5d %-3s %-4s . %-4s %s ? %5s %s %10.4f %10.4f %10.4f  0.0  0.0  ?  ?  ?  ?  ?  .  %5s %4s %s %4s %5d"
                    print(line % (recordName, atomIndex, symbol, atom.name, res.name, chainName, resId, resIC, coords[0], coords[1], coords[2],
456
457
458
                                  resId, res.name, chainName, atom.name, modelIndex), file=file)
                    posIndex += 1
                    atomIndex += 1