dockonsurf / modules / isolated.py @ 54be2a9e
Historique | Voir | Annoter | Télécharger (7,9 ko)
1 | 4614bb6a | Carles | """Functions to generate the conformers to be adsorbed and the most stable one.
|
---|---|---|---|
2 | 4614bb6a | Carles |
|
3 | 5c70d6c6 | Carles | functions:
|
4 | 92924d24 | Carles | confs_to_mol_list: Converts the conformers inside a rdkit mol object to a list
|
5 | 92924d24 | Carles | of independent mol objects.
|
6 | 5c70d6c6 | Carles | remove_C_linked_Hs: Removes hydrogens bonded to a carbon atom from a molecule.
|
7 | 92924d24 | Carles | gen_confs: Generate a number of conformers in random orientations.
|
8 | 92924d24 | Carles | get_rmsd: Gets the rmsd matrix of the conformers in a rdkit mol object.
|
9 | 92924d24 | Carles | get_moments_of_inertia: Computes moments of inertia of the given conformers.
|
10 | 92924d24 | Carles | mmff_opt_confs: Optimizes the geometry of the given conformers and returns the
|
11 | 92924d24 | Carles | new mol object and the energies of its conformers.
|
12 | 4614bb6a | Carles | run_isolated: directs the execution of functions to achieve the goal
|
13 | 4614bb6a | Carles | """
|
14 | 4614bb6a | Carles | import logging |
15 | 4614bb6a | Carles | |
16 | 92924d24 | Carles | import numpy as np |
17 | 92924d24 | Carles | from rdkit.Chem import AllChem as Chem |
18 | 92924d24 | Carles | |
19 | 4614bb6a | Carles | logger = logging.getLogger('DockOnSurf')
|
20 | 4614bb6a | Carles | |
21 | 4614bb6a | Carles | |
22 | 92924d24 | Carles | def confs_to_mol_list(mol: Chem.rdchem.Mol): |
23 | 92924d24 | Carles | """Converts the conformers inside a rdkit mol object to a list of
|
24 | 92924d24 | Carles | separate mol objects.
|
25 | 92924d24 | Carles |
|
26 | 92924d24 | Carles | @param mol: rdkit mol object containing at least one conformer.
|
27 | 92924d24 | Carles | @return list: list of separate mol objects.
|
28 | 92924d24 | Carles | """
|
29 | 92924d24 | Carles | return [Chem.MolFromMolBlock(Chem.MolToMolBlock(mol, confId=conf.GetId()))
|
30 | 92924d24 | Carles | for conf in mol.GetConformers()] |
31 | 92924d24 | Carles | |
32 | 92924d24 | Carles | |
33 | 92924d24 | Carles | def remove_C_linked_Hs(mol: Chem.rdchem.Mol): |
34 | b1f6e69d | Carles | """Removes hydrogen atoms bonded to a carbon atom from a rdkit mol object.
|
35 | 5c70d6c6 | Carles |
|
36 | 92924d24 | Carles | @param mol: rdkit mol object of the molecule with hydrogen atoms.
|
37 | 92924d24 | Carles | @return: rdkit mol object of the molecule without hydrogen atoms linked to
|
38 | 92924d24 | Carles | a carbon atom.
|
39 | 92924d24 | Carles |
|
40 | 92924d24 | Carles | The functions removes the hydrogen atoms bonded to carbon atoms while
|
41 | 92924d24 | Carles | keeping the ones bonded to other atoms or non-bonded at all.
|
42 | 5c70d6c6 | Carles | """
|
43 | 5c70d6c6 | Carles | |
44 | 5c70d6c6 | Carles | mol = Chem.RWMol(mol) |
45 | 5c70d6c6 | Carles | rev_atm_idxs = [atom.GetIdx() for atom in reversed(mol.GetAtoms())] |
46 | 5c70d6c6 | Carles | |
47 | 5c70d6c6 | Carles | for atm_idx in rev_atm_idxs: |
48 | 5c70d6c6 | Carles | atom = mol.GetAtomWithIdx(atm_idx) |
49 | 5c70d6c6 | Carles | if atom.GetAtomicNum() != 1: |
50 | 5c70d6c6 | Carles | continue
|
51 | 5c70d6c6 | Carles | for neigh in atom.GetNeighbors(): |
52 | 5c70d6c6 | Carles | if neigh.GetAtomicNum() == 6: |
53 | 5c70d6c6 | Carles | mol.RemoveAtom(atom.GetIdx()) |
54 | 5c70d6c6 | Carles | return mol
|
55 | 5c70d6c6 | Carles | |
56 | 5c70d6c6 | Carles | |
57 | 92924d24 | Carles | def gen_confs(mol: Chem.rdchem.Mol, num_confs: int, local_min=True): |
58 | 92924d24 | Carles | """Generate conformers in random orientations.
|
59 | 4614bb6a | Carles |
|
60 | 4614bb6a | Carles | @param mol: rdkit mol object of the molecule to be adsorbed.
|
61 | 92924d24 | Carles | @param num_confs: number of conformers to randomly generate.
|
62 | b1f6e69d | Carles | @param local_min: bool: if generated conformers should be a local minimum.
|
63 | 4614bb6a | Carles | @return: mol: rdkit mol object containing the different conformers.
|
64 | 92924d24 | Carles | rmsd_mtx: triangular matrix with the rmsd values of conformers.
|
65 | b1f6e69d | Carles |
|
66 | 92924d24 | Carles | Using the rdkit library, conformers are randomly generated. If structures
|
67 | 92924d24 | Carles | are required to be local minima, ie. setting the 'local_min' value to
|
68 | 92924d24 | Carles | True, a geometry optimisation using UFF is performed.
|
69 | 4614bb6a | Carles | """
|
70 | 92924d24 | Carles | logger.debug('Generating Conformers')
|
71 | 9e83b5a1 | Carles | |
72 | 5c70d6c6 | Carles | conf_ids = Chem.EmbedMultipleConfs(mol, numConfs=num_confs, numThreads=0)
|
73 | b1f6e69d | Carles | if local_min:
|
74 | b1f6e69d | Carles | for conf in conf_ids: |
75 | b1f6e69d | Carles | Chem.UFFOptimizeMolecule(mol, confId=conf) |
76 | 5c70d6c6 | Carles | |
77 | 5c70d6c6 | Carles | Chem.AlignMolConformers(mol) |
78 | 92924d24 | Carles | logger.info(f'Generated {len(mol.GetConformers())} conformers')
|
79 | 92924d24 | Carles | return mol
|
80 | 92924d24 | Carles | |
81 | 62576f52 | Carles | |
82 | 92924d24 | Carles | def get_rmsd(mol: Chem.rdchem.Mol, remove_Hs="c"): |
83 | 92924d24 | Carles | """Computes the rmsd matrix of the conformers in a rdkit mol object.
|
84 | 92924d24 | Carles |
|
85 | 92924d24 | Carles | @param mol: rdkit mol object containing at least two conformers.
|
86 | 92924d24 | Carles | @param remove_Hs: bool or str,
|
87 | 92924d24 | Carles | @return rmsd_matrix: Matrix containing the rmsd values of every pair of
|
88 | 92924d24 | Carles | conformers.
|
89 | 92924d24 | Carles |
|
90 | 92924d24 | Carles | The RMSD values of every pair of conformers is computed, stored in matrix
|
91 | 92924d24 | Carles | form and returned back. The calculation of rmsd values can take into
|
92 | 92924d24 | Carles | account all hydrogens, none, or only the ones not linked to carbon atoms.
|
93 | 92924d24 | Carles | """
|
94 | 92924d24 | Carles | if mol.GetNumConformers() < 2: |
95 | 92924d24 | Carles | err = "The provided molecule has less than 2 conformers"
|
96 | 92924d24 | Carles | logger.error(err) |
97 | 92924d24 | Carles | raise ValueError(err) |
98 | 92924d24 | Carles | |
99 | 92924d24 | Carles | if not remove_Hs: |
100 | 92924d24 | Carles | pass
|
101 | 92924d24 | Carles | elif remove_Hs or remove_Hs.lower() == "all": |
102 | 92924d24 | Carles | mol = Chem.RemoveHs(mol) |
103 | 92924d24 | Carles | elif remove_Hs.lower() == "c": |
104 | 92924d24 | Carles | mol = remove_C_linked_Hs(mol) |
105 | 92924d24 | Carles | else:
|
106 | 54be2a9e | Carles | err = "remove_Hs value does not have an acceptable value"
|
107 | 54be2a9e | Carles | logger.error(err) |
108 | 54be2a9e | Carles | raise ValueError(err) |
109 | 92924d24 | Carles | |
110 | 92924d24 | Carles | num_confs = mol.GetNumConformers() |
111 | 92924d24 | Carles | conf_ids = list(range(num_confs)) |
112 | 92924d24 | Carles | rmsd_mtx = np.zeros((num_confs, num_confs)) |
113 | 5c70d6c6 | Carles | for conf1 in conf_ids: |
114 | 5c70d6c6 | Carles | for conf2 in conf_ids[conf1 + 1:]: |
115 | 92924d24 | Carles | rmsd = Chem.GetBestRMS(mol, mol, prbId=conf2, refId=conf1) |
116 | 92924d24 | Carles | rmsd_mtx[conf1][conf2] = rmsd |
117 | 92924d24 | Carles | rmsd_mtx[conf2][conf1] = rmsd |
118 | 92924d24 | Carles | |
119 | 92924d24 | Carles | return rmsd_mtx
|
120 | 4614bb6a | Carles | |
121 | 92924d24 | Carles | |
122 | 92924d24 | Carles | def get_moments_of_inertia(mol: Chem.rdchem.Mol): |
123 | 92924d24 | Carles | """Computes the moments of inertia of the given conformers
|
124 | 92924d24 | Carles |
|
125 | 92924d24 | Carles | @param mol: rdkit mol object of the relevant molecule.
|
126 | 92924d24 | Carles | @return numpy array 2D: The inner array contains the moments of inertia for
|
127 | 92924d24 | Carles | the three principal axis of a given conformer. They are ordered by its value
|
128 | 92924d24 | Carles | in ascending order. The outer tuple loops over the conformers.
|
129 | 92924d24 | Carles | """
|
130 | 92924d24 | Carles | from rdkit.Chem.Descriptors3D import PMI1, PMI2, PMI3 |
131 | 92924d24 | Carles | |
132 | 92924d24 | Carles | return np.array([[PMI(mol, confId=conf) for PMI in (PMI1, PMI2, PMI3)] |
133 | 92924d24 | Carles | for conf in range(mol.GetNumConformers())]) |
134 | 92924d24 | Carles | |
135 | 92924d24 | Carles | |
136 | 62576f52 | Carles | def mmff_opt_confs(mol: Chem.rdchem.Mol, max_iters=2000): |
137 | 92924d24 | Carles | """Optimizes the geometry of the given conformers and returns the new mol
|
138 | 92924d24 | Carles | object and the energies of its conformers.
|
139 | 92924d24 | Carles |
|
140 | 92924d24 | Carles | @param mol: rdkit mol object of the relevant molecule.
|
141 | 62576f52 | Carles | @param max_iters: Maximum number of geometry optimization iterations. With 0
|
142 | 62576f52 | Carles | a single point energy calculation is performed and only the conformer
|
143 | 62576f52 | Carles | energies are returned.
|
144 | 92924d24 | Carles | @return mol: rdkit mol object of the optimized molecule.
|
145 | 92924d24 | Carles | @return numpy.ndarray: Array with the energies of the optimized conformers.
|
146 | 92924d24 | Carles |
|
147 | 92924d24 | Carles | The MMFF forcefield is used for the geometry optimization in its rdkit
|
148 | 62576f52 | Carles | implementation. With max_iters value set to 0, a single point energy
|
149 | 62576f52 | Carles | calculation is performed and only the energies are returned. For values
|
150 | 62576f52 | Carles | larger than 0, if the geometry does not converge for a certain conformer,
|
151 | 92924d24 | Carles | the latter is removed from the list of conformers and its energy is not
|
152 | 92924d24 | Carles | included in the returned list.
|
153 | 92924d24 | Carles | """
|
154 | 92924d24 | Carles | from rdkit.Chem.rdForceFieldHelpers import MMFFOptimizeMoleculeConfs |
155 | 92924d24 | Carles | |
156 | 92924d24 | Carles | init_num_confs = mol.GetNumConformers() |
157 | 92924d24 | Carles | results = np.array(MMFFOptimizeMoleculeConfs(mol, numThreads=0,
|
158 | 62576f52 | Carles | maxIters=max_iters, |
159 | 62576f52 | Carles | nonBondedThresh=10))
|
160 | 62576f52 | Carles | |
161 | 62576f52 | Carles | # Remove non-converged conformers if optimization is on, ie. maxIters > 0
|
162 | 62576f52 | Carles | # return all conformers if optimization is switched off, ie. maxIters = 0
|
163 | 62576f52 | Carles | if max_iters > 0: |
164 | 62576f52 | Carles | for i, conv in enumerate(results[:, 0]): |
165 | 62576f52 | Carles | if conv != 0: |
166 | 62576f52 | Carles | mol.RemoveConformer(i) |
167 | 62576f52 | Carles | for i, conf in enumerate(mol.GetConformers()): |
168 | 62576f52 | Carles | conf.SetId(i) |
169 | 62576f52 | Carles | if mol.GetNumConformers() < init_num_confs:
|
170 | 62576f52 | Carles | logger.warning(f'MMFF Geometry optimization did not comverge for at'
|
171 | 62576f52 | Carles | f'least one conformer. Continuing with '
|
172 | 62576f52 | Carles | f'{mol.GetNumConformers()} converged conformers')
|
173 | 62576f52 | Carles | logger.info(f'Optimized conformers with MMFF.')
|
174 | 62576f52 | Carles | return mol, np.array([res[1] for res in results if res[0] == 0]) |
175 | 62576f52 | Carles | else:
|
176 | 62576f52 | Carles | logger.info(f'Computed conformers energy with MMFF.')
|
177 | 62576f52 | Carles | return np.array([res[1] for res in results]) |
178 | 4614bb6a | Carles | |
179 | 4614bb6a | Carles | |
180 | 4614bb6a | Carles | def run_isolated(inp_vars): |
181 | 4614bb6a | Carles | """Directs the execution of functions to obtain the conformers to adsorb
|
182 | 4614bb6a | Carles |
|
183 | 4614bb6a | Carles | @param inp_vars:
|
184 | 4614bb6a | Carles | @return:
|
185 | 4614bb6a | Carles | """
|
186 | 54be2a9e | Carles | from formats import adapt_format |
187 | 54be2a9e | Carles | |
188 | 34a03ee1 | Carles | # from clustering import *
|
189 | 4614bb6a | Carles | logger.info('Carrying out procedures for the isolated molecule')
|
190 | 4614bb6a | Carles | rd_mol = adapt_format('rdkit', inp_vars['molec_file']) |
191 | 92924d24 | Carles | confs = gen_confs(rd_mol, inp_vars['num_conformers'])
|
192 | 34a03ee1 | Carles | rmsd_mtx = get_rmsd(confs) |
193 | 92924d24 | Carles | |
194 | 34a03ee1 | Carles | if 'moi' in inp_vars['cluster_magns']: |
195 | 34a03ee1 | Carles | confs_moi = get_moments_of_inertia(confs) |
196 | 92924d24 | Carles | |
197 | 34a03ee1 | Carles | if 'energy' in inp_vars['cluster_magns']: |
198 | 34a03ee1 | Carles | if inp_vars['min_confs']: |
199 | 34a03ee1 | Carles | confs, confs_eners = mmff_opt_confs(confs) |
200 | 34a03ee1 | Carles | else:
|
201 | 34a03ee1 | Carles | confs_eners = mmff_opt_confs(confs, max_iters=0)
|
202 | 9e83b5a1 | Carles | |
203 | 34a03ee1 | Carles | # clustering2() |