Re: [Rdkit-discuss] problems with EmbedMol
Open-Source Cheminformatics and Machine Learning
Brought to you by:
glandrum
From: Greg L. <gre...@gm...> - 2018-03-01 08:07:37
|
Hi Brian, On Wed, Feb 28, 2018 at 7:09 PM, Bennion, Brian <ben...@ll...> wrote: > > This is a real newbie question, but what is the use case for this > function? Is it used to generate all possible connections (limited by some > distance) between 3 or more atoms given in a smiles string? > There's a very brief description of how the RDKit uses distance geometry to generate 3D coordinates here: http://www.rdkit.org/docs/GettingStartedInPython.html#working-with-3d-molecules The, perhaps overly short, answer is that the distance bounds matrix describes the minimum (lower triangle) and maximum (upper triangle) possible 3D distances between atoms in the molecule. So you'd expect that all atom-atom distances in every remotely reasonable conformation of a molecule would fall somewhere in the range described by the bounds matrix. Here's a very simple example: In [27]: m = Chem.AddHs(Chem.MolFromSmiles('O=CC=O')) In [28]: print(numpy.array_str(AllChem.GetMoleculeBoundsMatrix(m),precision=3)) [[0. 1.268 2.394 3.482 2.068 3.318] [1.248 0. 1.468 2.394 1.091 2.247] [2.314 1.448 0. 1.268 2.247 1.091] [2.716 2.314 1.248 0. 3.318 2.068] [1.988 1.071 2.167 2.632 0. 3.155] [2.632 2.167 1.071 1.988 2.539 0. ]] (The two Hs that were added are the last two atoms) Notice: - the bounds for atoms that are bonded to each other are pretty tight: so 1.248 <= d(0,1) <= 1.268 - angular constraints are a bit looser: 2.314 <= d(0,2) <= 2.394 - dihedrals are considerably looser: 2.716 <= d(0,3) <= 3.482 Does that help? -greg > Brian > > > ------------------------------ > *From:* Greg Landrum <gre...@gm...> > *Sent:* Wednesday, February 28, 2018 8:53:59 AM > *To:* Jan Halborg Jensen > *Cc:* rdk...@li... > *Subject:* Re: [Rdkit-discuss] problems with EmbedMol > > Hi Jan, > > It took me much longer than it should have to figure this one out... > > The bounds matrix that is returned by GetMoleculeBoundsMatrix() needs to > have triangle bounds smoothing applied to it before it can be embedded. The > bounds smoothing process narrows the possible distance ranges between the > atoms. Here's a quick demo of that. > > We start with your example: > > In [19]: mol = Chem.MolFromSmiles("CCC") > ...: mol = Chem.AddHs(mol) > ...: bounds = AllChem.GetMoleculeBoundsMatrix(mol) > ...: EmbedLib.EmbedMol(mol,bounds) > ...: > ------------------------------------------------------------ > --------------- > ValueError Traceback (most recent call last) > <ipython-input-19-60befaacbe48> in <module>() > 2 mol = Chem.AddHs(mol) > 3 bounds = AllChem.GetMoleculeBoundsMatrix(mol) > ----> 4 EmbedLib.EmbedMol(mol,bounds) > > c:\Users\glandrum\RDKit_git\rdkit\Chem\Pharm3D\EmbedLib.py in > EmbedMol(mol, bm, atomMatch, weight, randomSeed, excludedVolumes) > 183 for i in range(nAts): > 184 weights.append((i, idx, weight)) > --> 185 coords = DG.EmbedBoundsMatrix(bm, weights=weights, > numZeroFail=1, randomSeed=randomSeed) > 186 # for row in coords: > 187 # print(', '.join(['%.2f'%x for x in row])) > > ValueError: could not embed matrix > > > > But if we do the triangle bounds smoothing things embed without problems: > > In [20]: from rdkit import DistanceGeometry > > In [21]: DistanceGeometry.DoTriangleSmoothing(bounds) > Out[21]: True > > In [22]: EmbedLib.EmbedMol(mol,bounds) > > In [23]: > > > There is a good argument to be made for GetMoleculeBoundsMatrix() > returning the smoothed bounds matrix by default. I'll put that on the list > for the next release. > > Best, > -greg > > > > On Wed, Feb 28, 2018 at 10:41 AM, Jan Halborg Jensen <jhj...@ch...> > wrote: > > The following code works fine with ethane (CC) but for propane (CCC) or > anything else I get the following error > ValueError: could not embed matrix > > Any ideas or solutions would be appreciated > > Best regards, Jan > > > from rdkit import Chem > from rdkit.Chem import AllChem > from rdkit.Chem.Pharm3D import EmbedLib > > > mol = Chem.MolFromSmiles("CCC") > mol = Chem.AddHs(mol) > bounds = AllChem.GetMoleculeBoundsMatrix(mol) > > EmbedLib.EmbedMol(mol,bounds) > EmbedLib.OptimizeMol(mol, bounds) > > ------------------------------------------------------------ > ------------------ > Check out the vibrant tech community on one of the world's most > engaging tech sites, Slashdot.org! http://sdm.link/slashdot > _______________________________________________ > Rdkit-discuss mailing list > Rdk...@li... > https://lists.sourceforge.net/lists/listinfo/rdkit-discuss > > > |