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

dockonsurf / modules / isolated.py @ ffa1b366

Historique | Voir | Annoter | Télécharger (7,12 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 4614bb6a Carles
import logging
13 4614bb6a Carles
14 92924d24 Carles
import numpy as np
15 fff3aab1 Carles
import rdkit.Chem.AllChem as Chem
16 92924d24 Carles
17 4614bb6a Carles
logger = logging.getLogger('DockOnSurf')
18 4614bb6a Carles
19 4614bb6a Carles
20 92924d24 Carles
def remove_C_linked_Hs(mol: Chem.rdchem.Mol):
21 b1f6e69d Carles
    """Removes hydrogen atoms bonded to a carbon atom from a rdkit mol object.
22 5c70d6c6 Carles

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

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

47 4614bb6a Carles
    @param mol: rdkit mol object of the molecule to be adsorbed.
48 92924d24 Carles
    @param num_confs: number of conformers to randomly generate.
49 4614bb6a Carles
    @return: mol: rdkit mol object containing the different conformers.
50 bd573212 Carles
             rmsd_mtx: Matrix with the rmsd values of conformers.
51 b1f6e69d Carles

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

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

84 92924d24 Carles
    @param mol: rdkit mol object of the relevant molecule.
85 ad57fd42 Carles Marti
    @param force_field: Force Field to use for the pre-optimization.
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 ad57fd42 Carles Marti
    The MMFF and UFF force fields can be used for the geometry optimization in
93 4670488d Carles Marti
    their rdkit implementation. With max_iters value set to 0, a single point
94 4670488d Carles Marti
    energy calculation is performed and only the energies are returned. For
95 4670488d Carles Marti
    values larger than 0, if the geometry does not converge for a certain
96 4670488d Carles Marti
    conformer, the latter is removed from the list of conformers and its energy
97 4670488d Carles Marti
    is not included in the returned list.
98 92924d24 Carles
    """
99 4670488d Carles Marti
    from rdkit.Chem.rdForceFieldHelpers import MMFFOptimizeMoleculeConfs, \
100 4670488d Carles Marti
        UFFOptimizeMoleculeConfs
101 92924d24 Carles
102 92924d24 Carles
    init_num_confs = mol.GetNumConformers()
103 ad57fd42 Carles Marti
    if force_field == 'mmff':
104 4670488d Carles Marti
        results = np.array(MMFFOptimizeMoleculeConfs(mol, numThreads=0,
105 4670488d Carles Marti
                                                     maxIters=max_iters,
106 4670488d Carles Marti
                                                     nonBondedThresh=10))
107 ad57fd42 Carles Marti
    elif force_field == 'uff':
108 4670488d Carles Marti
        results = np.array(UFFOptimizeMoleculeConfs(mol, numThreads=0,
109 4670488d Carles Marti
                                                    maxIters=max_iters,
110 4670488d Carles Marti
                                                    vdwThresh=10))
111 4670488d Carles Marti
    else:
112 4670488d Carles Marti
        logger.error("force_field parameter must be 'MMFF' or 'UFF'")
113 4670488d Carles Marti
        raise ValueError("force_field parameter must be 'MMFF' or 'UFF'")
114 62576f52 Carles
115 62576f52 Carles
    # Remove non-converged conformers if optimization is on, ie. maxIters > 0
116 62576f52 Carles
    # return all conformers if optimization is switched off, ie. maxIters = 0
117 62576f52 Carles
    if max_iters > 0:
118 62576f52 Carles
        for i, conv in enumerate(results[:, 0]):
119 62576f52 Carles
            if conv != 0:
120 62576f52 Carles
                mol.RemoveConformer(i)
121 62576f52 Carles
        for i, conf in enumerate(mol.GetConformers()):
122 62576f52 Carles
            conf.SetId(i)
123 62576f52 Carles
        if mol.GetNumConformers() < init_num_confs:
124 ad57fd42 Carles Marti
            logger.warning(f'Geometry optimization did not comverge for at'
125 62576f52 Carles
                           f'least one conformer. Continuing with '
126 695dcff8 Carles Marti
                           f'{mol.GetNumConformers()} converged conformers.')
127 ad57fd42 Carles Marti
        logger.info(f'Pre-optimized conformers with {force_field}.')
128 62576f52 Carles
        return mol, np.array([res[1] for res in results if res[0] == 0])
129 62576f52 Carles
    else:
130 ad57fd42 Carles Marti
        logger.info(f'Computed conformers energy with {force_field}.')
131 62576f52 Carles
        return np.array([res[1] for res in results])
132 4614bb6a Carles
133 4614bb6a Carles
134 4614bb6a Carles
def run_isolated(inp_vars):
135 4614bb6a Carles
    """Directs the execution of functions to obtain the conformers to adsorb
136 4614bb6a Carles

137 a9b058a5 Carles
    @param inp_vars: Calculation parameters from input file.
138 4614bb6a Carles
    @return:
139 4614bb6a Carles
    """
140 f3004731 Carles
    from modules.formats import adapt_format, confs_to_mol_list, \
141 f3004731 Carles
        rdkit_mol_to_ase_atoms
142 fff3aab1 Carles
    from modules.clustering import clustering, get_rmsd
143 3d6a9d3c Carles
    from modules.calculation import run_calc
144 54be2a9e Carles
145 695dcff8 Carles Marti
    logger.info('Carrying out procedures for the isolated molecule.')
146 90819cc3 Carles Marti
    rd_mol = adapt_format('rdkit', inp_vars['molec_file'],
147 90819cc3 Carles Marti
                          inp_vars['special_atoms'])
148 92924d24 Carles
    confs = gen_confs(rd_mol, inp_vars['num_conformers'])
149 4670488d Carles Marti
    if inp_vars['pre_opt']:
150 4670488d Carles Marti
        confs, confs_ener = pre_opt_confs(confs, inp_vars['pre_opt'])
151 fab7b517 Carles
    else:
152 4670488d Carles Marti
        confs_ener = pre_opt_confs(confs, max_iters=0)
153 fff3aab1 Carles
    conf_list = confs_to_mol_list(confs)
154 fff3aab1 Carles
    rmsd_mtx = get_rmsd(conf_list)
155 fab7b517 Carles
    confs_moi = get_moments_of_inertia(confs)
156 fd33bfdb Carles
    exemplars = clustering(rmsd_mtx)
157 f3004731 Carles
    mol_list = confs_to_mol_list(confs, exemplars)
158 f3004731 Carles
    ase_atms_list = [rdkit_mol_to_ase_atoms(mol) for mol in mol_list]
159 a44ad3c2 Carles Martí
    if len(ase_atms_list) == 0:
160 a44ad3c2 Carles Martí
        err_msg = "No configurations were generated: Check the parameters in" \
161 a44ad3c2 Carles Martí
                  "dockonsurf.inp"
162 a44ad3c2 Carles Martí
        logger.error(err_msg)
163 a44ad3c2 Carles Martí
        raise ValueError(err_msg)
164 f3004731 Carles
    run_calc('isolated', inp_vars, ase_atms_list)
165 14f39d2a Carles Marti
    logger.info('Finished the procedures for the isolated molecule section.')