{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Example 7: Rotational barrier prediction in KRAS inhibitors\n", "\n", "## Objective\n", "The objective of this task is to use the `molli` workflow to transition from `.CDXML` molecular drawings to full-fledged 3D structures with subsequent CREST//GFN2-XTB conformer generation, ORCA equilibrium structure minimization, as well as ORCA relaxed surface scan.\n", "This workflow is useful to assess semi-quantitatively the possibility of hindered rotation in a number of molecules (see e. g. [*J. Med. Chem.* **2020**, *63*, 1, 52](https://doi.org/10.1021/acs.jmedchem.9b01180))\n", "\n", "## Prerequisites\n", "- `molli` installation\n", "- `pandas` installation (it is NOT a dependency of molli and must be installed separately)\n", "- `OpenBabel` installation (it is NOT a dependency of molli and must be installed separately)\n", " Note: if you experience problems with openbabel refusing to optimize the structures, do not forget to set the `BABEL_DATADIR` environment variable to point to a valid location with openbabel data files. With conda installations, it seems to be consistently pointing to `$CONDA_PREFIX/share/openbabel`\n", "- `barriers.csv` file that has the experimental rotational barriers listed\n", "- `cladosporin.cdxml` file that contains the structural drawings\n", "\n", "## Methods\n", "\n", "### Conformer generation\n", "Rotational barrier profiles are generated using b97-3c method, as reported by Grimme et al." ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [], "source": [ "import molli as ml\n", "from molli.external import openbabel\n", "from math import degrees\n", "\n", "import pandas as pd\n", "\n", "import matplotlib.pyplot as plt\n", "import numpy as np\n", "from scipy.signal import find_peaks\n", "\n", "# This configures on-screen preview of molecules\n", "ml.visual.configure(bgcolor=\"black\", theme=\"dark\", width=500)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Step 1. `CDXML` Parsing\n", "\n", "In this step, all of the molecules are read from the CDXML file and converted into their 3D versions.\n", "This can be done from both the command line interface, but for the sake of exercise, this workflow will show the script that utilizes the key function of `ml.Molecule.join(...)`" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "To make sure that the parsing step was successful, we will perform a basic inspection of the resulting `.mlib` file" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Dihedral 117.02 for Molecule(name='23', formula='C28 H29 Cl1 F1 N5 O2')\n", "Dihedral 113.62 for Molecule(name='25', formula='C27 H27 Cl1 F1 N5 O2')\n", "Dihedral 134.54 for Molecule(name='28', formula='C24 H24 Cl1 F1 N6 O2 S1')\n", "Dihedral 144.75 for Molecule(name='27', formula='C24 H25 Cl1 F1 N7 O2')\n", "Dihedral 131.87 for Molecule(name='18', formula='C27 H27 Cl1 F1 N5 O2')\n", "Dihedral 128.08 for Molecule(name='29', formula='C25 H27 Cl1 F1 N7 O2')\n", "Dihedral 111.87 for Molecule(name='24', formula='C28 H29 Cl1 F1 N5 O2')\n", "Dihedral 139.62 for Molecule(name='26', formula='C26 H26 Cl1 F1 N6 O2')\n", "Dihedral 126.61 for Molecule(name='22', formula='C27 H25 Cl1 F1 N5 O2')\n" ] } ], "source": [ "cdxf = ml.CDXMLFile(\"rot_barrier.cdxml\")\n", "core = cdxf[\"core\"]\n", "core_ap = core.get_attachment_points()[0]\n", "names = set(cdxf.keys())\n", "names.remove(\"core\")\n", "core.add_implicit_hydrogens()\n", "\n", "# This is the destination \n", "lib = ml.MoleculeLibrary(\"kras_inhibitors.mlib\", overwrite=True, readonly=False)\n", "\n", "with lib.writing():\n", " for n in names:\n", " sub = cdxf[n]\n", " sub_ap = sub.get_attachment_points()[0]\n", " \n", " mol = ml.Molecule.join(\n", " core,\n", " sub,\n", " core_ap,\n", " sub_ap,\n", " name=n,\n", " )\n", "\n", " mol.add_implicit_hydrogens()\n", "\n", " openbabel.obabel_optimize(\n", " mol,\n", " ff=\"MMFF94\",\n", " coord_displace=0.01,\n", " inplace=True,\n", " )\n", "\n", " # just a little extra: calculate the preliminary dihedral angle)\n", " d = mol.dihedral(\"1\", \"2\", \"3\", \"4\")\n", " print(f\"Dihedral {degrees(d):0.2f} for {mol}\")\n", "\n", " lib[n] = mol" ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [ { "data": { "application/3dmoljs_load.v0": "
You appear to be running in JupyterLab (or JavaScript failed to load for some other reason). You need to install the 3dmol extension:
\n jupyter labextension install jupyterlab_3dmol
You appear to be running in JupyterLab (or JavaScript failed to load for some other reason). You need to install the 3dmol extension:
\n",
" jupyter labextension install jupyterlab_3dmol
You appear to be running in JupyterLab (or JavaScript failed to load for some other reason). You need to install the 3dmol extension:
\n jupyter labextension install jupyterlab_3dmol
You appear to be running in JupyterLab (or JavaScript failed to load for some other reason). You need to install the 3dmol extension:
\n",
" jupyter labextension install jupyterlab_3dmol
You appear to be running in JupyterLab (or JavaScript failed to load for some other reason). You need to install the 3dmol extension:
\n jupyter labextension install jupyterlab_3dmol
You appear to be running in JupyterLab (or JavaScript failed to load for some other reason). You need to install the 3dmol extension:
\n",
" jupyter labextension install jupyterlab_3dmol
You appear to be running in JupyterLab (or JavaScript failed to load for some other reason). You need to install the 3dmol extension:
\n jupyter labextension install jupyterlab_3dmol
You appear to be running in JupyterLab (or JavaScript failed to load for some other reason). You need to install the 3dmol extension:
\n",
" jupyter labextension install jupyterlab_3dmol
You appear to be running in JupyterLab (or JavaScript failed to load for some other reason). You need to install the 3dmol extension:
\n jupyter labextension install jupyterlab_3dmol
You appear to be running in JupyterLab (or JavaScript failed to load for some other reason). You need to install the 3dmol extension:
\n",
" jupyter labextension install jupyterlab_3dmol
| \n", " | barrier (R->TS) | \n", "barrier (S->TS) | \n", "TS low freq | \n", "TS N imag freq | \n", "EQ1 N imag freq | \n", "EQ2 N imag freq | \n", "
|---|---|---|---|---|---|---|
| 18 | \n", "108.3 | \n", "106.7 | \n", "-32.2 | \n", "1 | \n", "0 | \n", "0 | \n", "
| 22 | \n", "103.2 | \n", "102.6 | \n", "-36.5 | \n", "1 | \n", "0 | \n", "0 | \n", "
| 23 | \n", "132.0 | \n", "133.3 | \n", "-48.3 | \n", "1 | \n", "1 | \n", "0 | \n", "
| 24 | \n", "150.0 | \n", "151.9 | \n", "-30.0 | \n", "1 | \n", "0 | \n", "0 | \n", "
| 25 | \n", "146.1 | \n", "148.4 | \n", "-35.0 | \n", "1 | \n", "0 | \n", "0 | \n", "
| 26 | \n", "92.6 | \n", "89.4 | \n", "-27.8 | \n", "1 | \n", "0 | \n", "0 | \n", "
| 27 | \n", "73.7 | \n", "76.0 | \n", "-38.6 | \n", "1 | \n", "0 | \n", "0 | \n", "
| 28 | \n", "69.9 | \n", "69.2 | \n", "-35.0 | \n", "1 | \n", "0 | \n", "0 | \n", "
| 29 | \n", "101.3 | \n", "103.8 | \n", "-25.5 | \n", "1 | \n", "0 | \n", "0 | \n", "
You appear to be running in JupyterLab (or JavaScript failed to load for some other reason). You need to install the 3dmol extension:
\n jupyter labextension install jupyterlab_3dmol
You appear to be running in JupyterLab (or JavaScript failed to load for some other reason). You need to install the 3dmol extension:
\n",
" jupyter labextension install jupyterlab_3dmol
You appear to be running in JupyterLab (or JavaScript failed to load for some other reason). You need to install the 3dmol extension:
\n jupyter labextension install jupyterlab_3dmol
You appear to be running in JupyterLab (or JavaScript failed to load for some other reason). You need to install the 3dmol extension:
\n",
" jupyter labextension install jupyterlab_3dmol
| \n", " | barrier (R->TS) | \n", "barrier (S->TS) | \n", "TS low freq | \n", "TS N imag freq | \n", "EQ1 N imag freq | \n", "EQ2 N imag freq | \n", "
|---|---|---|---|---|---|---|
| 18 | \n", "108.3 | \n", "106.7 | \n", "-32.2 | \n", "1 | \n", "0 | \n", "0 | \n", "
| 22 | \n", "103.2 | \n", "102.6 | \n", "-36.5 | \n", "1 | \n", "0 | \n", "0 | \n", "
| 23 | \n", "141.0 | \n", "133.3 | \n", "-48.3 | \n", "1 | \n", "0 | \n", "0 | \n", "
| 24 | \n", "150.0 | \n", "151.9 | \n", "-30.0 | \n", "1 | \n", "0 | \n", "0 | \n", "
| 25 | \n", "146.1 | \n", "148.4 | \n", "-35.0 | \n", "1 | \n", "0 | \n", "0 | \n", "
| 26 | \n", "92.6 | \n", "89.4 | \n", "-27.8 | \n", "1 | \n", "0 | \n", "0 | \n", "
| 27 | \n", "73.7 | \n", "76.0 | \n", "-38.6 | \n", "1 | \n", "0 | \n", "0 | \n", "
| 28 | \n", "69.9 | \n", "69.2 | \n", "-35.0 | \n", "1 | \n", "0 | \n", "0 | \n", "
| 29 | \n", "101.3 | \n", "103.8 | \n", "-25.5 | \n", "1 | \n", "0 | \n", "0 | \n", "