from mpi4py import MPI from mpi4py.futures import MPICommExecutor import warnings from Bio.PDB import PDBParser, PPBuilder, CaPPBuilder from Bio.PDB.NeighborSearch import NeighborSearch from Bio.PDB.Selection import unfold_entities import numpy as np import dask.array as da from rdkit import Chem from spyrmsd import molecule from spyrmsd import graph import networkx as nx import os import re import sys # all punctuation punctuation_regex = r"""(\(|\)|\.|=|#|-|\+|\\|\/|:|~|@|\?|>>?|\*|\$|\%[0-9]{2}|[0-9])""" # tokenization regex (Schwaller) molecule_regex = r"""(\[[^\]]+]|Br?|Cl?|N|O|S|P|F|I|b|c|n|o|s|p|\(|\)|\.|=|#|-|\+|\\|\/|:|~|@|\?|>>?|\*|\$|\%[0-9]{2}|[0-9])""" max_seq = 2046 # = 2048 - 2 (accounting for [CLS] and [SEP]) max_smiles = 510 # = 512 - 2 chunk_size = '1G' def rot_from_two_vecs(e0_unnormalized, e1_unnormalized): """Create rotation matrices from unnormalized vectors for the x and y-axes. This creates a rotation matrix from two vectors using Gram-Schmidt orthogonalization. Args: e0_unnormalized: vectors lying along x-axis of resulting rotation e1_unnormalized: vectors lying in xy-plane of resulting rotation Returns: Rotations resulting from Gram-Schmidt procedure. """ # Normalize the unit vector for the x-axis, e0. e0 = e0_unnormalized / np.linalg.norm(e0_unnormalized) # make e1 perpendicular to e0. c = np.dot(e1_unnormalized, e0) e1 = e1_unnormalized - c * e0 e1 = e1 / np.linalg.norm(e1) # Compute e2 as cross product of e0 and e1. e2 = np.cross(e0, e1) # local to space frame return np.stack([e0,e1,e2]).T def get_local_frames(mol): # get the two nearest neighbors of every atom on the molecular graph # ties are broken using canonical ordering g = molecule.Molecule.from_rdkit(mol).to_graph() R = [] for node in g: length = nx.single_source_shortest_path_length(g, node) neighbor_a = [n for n,l in length.items() if l==1][0] try: neighbor_b = [n for n,l in length.items() if l==1][1] except: # get next nearest neighbor neighbor_b = [n for n,l in length.items() if l==2][0] xyz = np.array(mol.GetConformer().GetAtomPosition(node)) xyz_a = np.array(mol.GetConformer().GetAtomPosition(neighbor_a)) xyz_b = np.array(mol.GetConformer().GetAtomPosition(neighbor_b)) R.append(rot_from_two_vecs(xyz_a-xyz, xyz_b-xyz)) return R def parse_complex(fn): try: name = os.path.basename(fn) # parse protein sequence and coordinates parser = PDBParser() with warnings.catch_warnings(): warnings.simplefilter("ignore") structure = parser.get_structure('protein',fn+'/'+name+'_protein.pdb') res_frames = [] # extract sequence, Calpha positions and local coordinate frames using the AF2 convention ppb = CaPPBuilder() seq = [] xyz_receptor = [] R_receptor = [] for pp in ppb.build_peptides(structure): seq.append(str(pp.get_sequence())) xyz_receptor += [tuple(a.get_vector()) for a in pp.get_ca_list()] for res in pp: N = np.array(tuple(res['N'].get_vector())) C = np.array(tuple(res['C'].get_vector())) CA = np.array(tuple(res['CA'].get_vector())) R_receptor.append(rot_from_two_vecs(N-CA,C-CA).flatten().tolist()) seq = ''.join(seq) # parse ligand, convert to SMILES and map atoms suppl = Chem.SDMolSupplier(fn+'/'+name+'_ligand.sdf') mol = next(suppl) # bring molecule atoms in canonical order (to determine local frames uniquely) m_neworder = tuple(zip(*sorted([(j, i) for i, j in enumerate(Chem.CanonicalRankAtoms(mol))])))[1] mol = Chem.RenumberAtoms(mol, m_neworder) # position of atoms in SMILES (not counting punctuation) smi = Chem.MolToSmiles(mol) atom_order = [int(s) for s in list(filter(None,re.sub(r'[\[\]]','',mol.GetProp("_smilesAtomOutputOrder")).split(',')))] # tokenize the SMILES tokens = list(filter(None, re.split(molecule_regex, smi))) # remove punctuation masked_tokens = [re.sub(punctuation_regex,'',s) for s in tokens] k = 0 token_pos = [] token_rot = [] frames = get_local_frames(mol) for i,token in enumerate(masked_tokens): if token != '': token_pos.append(tuple(mol.GetConformer().GetAtomPosition(atom_order[k]))) token_rot.append(frames[atom_order[k]].flatten().tolist()) k += 1 else: token_pos.append((np.nan, np.nan, np.nan)) token_rot.append(np.eye(3).flatten().tolist()) return name, seq, smi, xyz_receptor, token_pos, token_rot, R_receptor except Exception as e: print(e) return None if __name__ == '__main__': import glob filenames = glob.glob('data/pdbbind/v2020-other-PL/*') filenames.extend(glob.glob('data/pdbbind/refined-set/*')) filenames = sorted(filenames) comm = MPI.COMM_WORLD with MPICommExecutor(comm, root=0) as executor: if executor is not None: result = executor.map(parse_complex, filenames, chunksize=32) result = list(result) names = [r[0] for r in result if r is not None] seqs = [r[1] for r in result if r is not None] all_smiles = [r[2] for r in result if r is not None] all_xyz_receptor = [r[3] for r in result if r is not None] all_xyz_ligand = [r[4] for r in result if r is not None] all_rot_ligand = [r[5] for r in result if r is not None] all_rot_receptor = [r[6] for r in result if r is not None] import pandas as pd df = pd.DataFrame({'name': names, 'seq': seqs, 'smiles': all_smiles, 'receptor_xyz': all_xyz_receptor, 'ligand_xyz': all_xyz_ligand, 'ligand_rot': all_rot_ligand, 'receptor_rot': all_rot_receptor}) df.to_parquet('data/pdbbind.parquet',index=False)