Statistiques
| Branche: | Tag: | Révision :

dockonsurf / modules / isolated.py @ 0dfa15ff

Historique | Voir | Annoter | Télécharger (6,22 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 5c70d6c6 Carles
remove_C_linked_Hs: Removes hydrogens bonded to a carbon atom from a molecule.
5 92924d24 Carles
gen_confs: Generate a number of conformers in random orientations.
6 92924d24 Carles
get_rmsd: Gets the rmsd matrix of the conformers in a rdkit mol object.
7 92924d24 Carles
get_moments_of_inertia: Computes moments of inertia of the given conformers.
8 92924d24 Carles
mmff_opt_confs: Optimizes the geometry of the given conformers and returns the
9 92924d24 Carles
    new mol object and the energies of its conformers.
10 4614bb6a Carles
run_isolated: directs the execution of functions to achieve the goal
11 4614bb6a Carles
"""
12 fd33bfdb Carles
import os
13 4614bb6a Carles
import logging
14 4614bb6a Carles
15 92924d24 Carles
import numpy as np
16 fff3aab1 Carles
import rdkit.Chem.AllChem as Chem
17 92924d24 Carles
18 4614bb6a Carles
logger = logging.getLogger('DockOnSurf')
19 4614bb6a Carles
20 4614bb6a Carles
21 92924d24 Carles
def remove_C_linked_Hs(mol: Chem.rdchem.Mol):
22 b1f6e69d Carles
    """Removes hydrogen atoms bonded to a carbon atom from a rdkit mol object.
23 5c70d6c6 Carles

24 92924d24 Carles
    @param mol: rdkit mol object of the molecule with hydrogen atoms.
25 92924d24 Carles
    @return: rdkit mol object of the molecule without hydrogen atoms linked to
26 92924d24 Carles
    a carbon atom.
27 92924d24 Carles

28 92924d24 Carles
    The functions removes the hydrogen atoms bonded to carbon atoms while
29 92924d24 Carles
    keeping the ones bonded to other atoms or non-bonded at all.
30 5c70d6c6 Carles
    """
31 5c70d6c6 Carles
32 5c70d6c6 Carles
    mol = Chem.RWMol(mol)
33 5c70d6c6 Carles
    rev_atm_idxs = [atom.GetIdx() for atom in reversed(mol.GetAtoms())]
34 5c70d6c6 Carles
35 5c70d6c6 Carles
    for atm_idx in rev_atm_idxs:
36 5c70d6c6 Carles
        atom = mol.GetAtomWithIdx(atm_idx)
37 5c70d6c6 Carles
        if atom.GetAtomicNum() != 1:
38 5c70d6c6 Carles
            continue
39 5c70d6c6 Carles
        for neigh in atom.GetNeighbors():
40 5c70d6c6 Carles
            if neigh.GetAtomicNum() == 6:
41 5c70d6c6 Carles
                mol.RemoveAtom(atom.GetIdx())
42 5c70d6c6 Carles
    return mol
43 5c70d6c6 Carles
44 5c70d6c6 Carles
45 92924d24 Carles
def gen_confs(mol: Chem.rdchem.Mol, num_confs: int, local_min=True):
46 92924d24 Carles
    """Generate conformers in random orientations.
47 4614bb6a Carles

48 4614bb6a Carles
    @param mol: rdkit mol object of the molecule to be adsorbed.
49 92924d24 Carles
    @param num_confs: number of conformers to randomly generate.
50 b1f6e69d Carles
    @param local_min: bool: if generated conformers should be a local minimum.
51 4614bb6a Carles
    @return: mol: rdkit mol object containing the different conformers.
52 bd573212 Carles
             rmsd_mtx: Matrix with the rmsd values of conformers.
53 b1f6e69d Carles

54 92924d24 Carles
    Using the rdkit library, conformers are randomly generated. If structures 
55 92924d24 Carles
    are required to be local minima, ie. setting the 'local_min' value to 
56 92924d24 Carles
    True, a geometry optimisation using UFF is performed.
57 4614bb6a Carles
    """
58 92924d24 Carles
    logger.debug('Generating Conformers')
59 9e83b5a1 Carles
60 9324a38a Carles
    mol = Chem.AddHs(mol)
61 fab7b517 Carles
    Chem.EmbedMultipleConfs(mol, numConfs=num_confs, numThreads=0)
62 5c70d6c6 Carles
    Chem.AlignMolConformers(mol)
63 92924d24 Carles
    logger.info(f'Generated {len(mol.GetConformers())} conformers')
64 92924d24 Carles
    return mol
65 92924d24 Carles
66 62576f52 Carles
67 92924d24 Carles
def get_moments_of_inertia(mol: Chem.rdchem.Mol):
68 92924d24 Carles
    """Computes the moments of inertia of the given conformers
69 92924d24 Carles

70 92924d24 Carles
    @param mol: rdkit mol object of the relevant molecule.
71 92924d24 Carles
    @return numpy array 2D: The inner array contains the moments of inertia for
72 92924d24 Carles
    the three principal axis of a given conformer. They are ordered by its value
73 92924d24 Carles
    in ascending order. The outer tuple loops over the conformers.
74 92924d24 Carles
    """
75 92924d24 Carles
    from rdkit.Chem.Descriptors3D import PMI1, PMI2, PMI3
76 92924d24 Carles
77 92924d24 Carles
    return np.array([[PMI(mol, confId=conf) for PMI in (PMI1, PMI2, PMI3)]
78 92924d24 Carles
                     for conf in range(mol.GetNumConformers())])
79 92924d24 Carles
80 92924d24 Carles
81 62576f52 Carles
def mmff_opt_confs(mol: Chem.rdchem.Mol, max_iters=2000):
82 92924d24 Carles
    """Optimizes the geometry of the given conformers and returns the new mol
83 92924d24 Carles
    object and the energies of its conformers.
84 92924d24 Carles

85 92924d24 Carles
    @param mol: rdkit mol object of the relevant molecule.
86 62576f52 Carles
    @param max_iters: Maximum number of geometry optimization iterations. With 0
87 62576f52 Carles
    a single point energy calculation is performed and only the conformer
88 62576f52 Carles
    energies are returned.
89 92924d24 Carles
    @return mol: rdkit mol object of the optimized molecule.
90 92924d24 Carles
    @return numpy.ndarray: Array with the energies of the optimized conformers.
91 92924d24 Carles

92 92924d24 Carles
    The MMFF forcefield is used for the geometry optimization in its rdkit
93 62576f52 Carles
    implementation. With max_iters value set to 0, a single point energy
94 62576f52 Carles
    calculation is performed and only the energies are returned. For values
95 62576f52 Carles
    larger than 0, if the geometry does not converge for a certain conformer,
96 92924d24 Carles
    the latter is removed from the list of conformers and its energy is not
97 92924d24 Carles
    included in the returned list.
98 92924d24 Carles
    """
99 92924d24 Carles
    from rdkit.Chem.rdForceFieldHelpers import MMFFOptimizeMoleculeConfs
100 92924d24 Carles
101 92924d24 Carles
    init_num_confs = mol.GetNumConformers()
102 92924d24 Carles
    results = np.array(MMFFOptimizeMoleculeConfs(mol, numThreads=0,
103 62576f52 Carles
                                                 maxIters=max_iters,
104 62576f52 Carles
                                                 nonBondedThresh=10))
105 62576f52 Carles
106 62576f52 Carles
    # Remove non-converged conformers if optimization is on, ie. maxIters > 0
107 62576f52 Carles
    # return all conformers if optimization is switched off, ie. maxIters = 0
108 62576f52 Carles
    if max_iters > 0:
109 62576f52 Carles
        for i, conv in enumerate(results[:, 0]):
110 62576f52 Carles
            if conv != 0:
111 62576f52 Carles
                mol.RemoveConformer(i)
112 62576f52 Carles
        for i, conf in enumerate(mol.GetConformers()):
113 62576f52 Carles
            conf.SetId(i)
114 62576f52 Carles
        if mol.GetNumConformers() < init_num_confs:
115 62576f52 Carles
            logger.warning(f'MMFF Geometry optimization did not comverge for at'
116 62576f52 Carles
                           f'least one conformer. Continuing with '
117 62576f52 Carles
                           f'{mol.GetNumConformers()} converged conformers')
118 826f0462 Carles Martí
        logger.info(f'Pre-optimized conformers with MMFF.')
119 62576f52 Carles
        return mol, np.array([res[1] for res in results if res[0] == 0])
120 62576f52 Carles
    else:
121 62576f52 Carles
        logger.info(f'Computed conformers energy with MMFF.')
122 62576f52 Carles
        return np.array([res[1] for res in results])
123 4614bb6a Carles
124 4614bb6a Carles
125 4614bb6a Carles
def run_isolated(inp_vars):
126 4614bb6a Carles
    """Directs the execution of functions to obtain the conformers to adsorb
127 4614bb6a Carles

128 a9b058a5 Carles
    @param inp_vars: Calculation parameters from input file.
129 4614bb6a Carles
    @return:
130 4614bb6a Carles
    """
131 f3004731 Carles
    from modules.formats import adapt_format, confs_to_mol_list, \
132 f3004731 Carles
        rdkit_mol_to_ase_atoms
133 fff3aab1 Carles
    from modules.clustering import clustering, get_rmsd
134 3d6a9d3c Carles
    from modules.calculation import run_calc
135 54be2a9e Carles
136 4614bb6a Carles
    logger.info('Carrying out procedures for the isolated molecule')
137 4614bb6a Carles
    rd_mol = adapt_format('rdkit', inp_vars['molec_file'])
138 92924d24 Carles
    confs = gen_confs(rd_mol, inp_vars['num_conformers'])
139 fd33bfdb Carles
    if inp_vars['min_confs']:
140 fd33bfdb Carles
        confs, confs_ener = mmff_opt_confs(confs)
141 fab7b517 Carles
    else:
142 fab7b517 Carles
        confs_ener = mmff_opt_confs(confs, max_iters=0)
143 fff3aab1 Carles
    conf_list = confs_to_mol_list(confs)
144 fff3aab1 Carles
    rmsd_mtx = get_rmsd(conf_list)
145 fab7b517 Carles
    confs_moi = get_moments_of_inertia(confs)
146 fd33bfdb Carles
    exemplars = clustering(rmsd_mtx)
147 f3004731 Carles
    mol_list = confs_to_mol_list(confs, exemplars)
148 f3004731 Carles
    ase_atms_list = [rdkit_mol_to_ase_atoms(mol) for mol in mol_list]
149 f3004731 Carles
    run_calc('isolated', inp_vars, ase_atms_list)