Hi Alfredo,

There are two things going on here:

1) If stereocenters are not specified, the RDKit will not take them into
account when generating conformers. This means you'll get a "random"
configuration for each stereoisomer in each conformation (I put random in
quotes because that's not strictly true, but that's a longer discussion).
The fact that you happened to get an S,S configurations when you didnt'
specify the stereochemistry is really just luck.

2) The conformation generation includes a random process. If you make two
separate calls to the EmbedMolecule for exactly the same molecule, you'll
typically get different results.

Here's a demonstration of that for your molecule with fully specified
stereochemistry:

In [81]: def run_smiles(smi,seed):
    ...:     m = Chem.MolFromSmiles(smi)
    ...:     ma = Chem.AddHs(m)
    ...:     mb = Chem.AddHs(m)
    ...:     params = AllChem.ETKDGv2()
    ...:     if seed is not None:
    ...:         params.randomSeed = seed
    ...:     AllChem.EmbedMolecule(ma,params)
    ...:     AllChem.EmbedMolecule(mb,params)

    ...:     # use GetBestRMS() to show how different the conformations are:

    ...:     print(AllChem.GetBestRMS(Chem.RemoveHs(ma),Chem.RemoveHs(mb)))
    ...:     Chem.AssignStereochemistryFrom3D(ma)
    ...:     Chem.AssignStereochemistryFrom3D(mb)
    ...:
print(Chem.FindMolChiralCenters(ma),Chem.FindMolChiralCenters(mb))
    ...:


In [86]: run_smiles('CCCC[C@@H](C(=O)COc1c(F)c(F)cc(F)c1F)n1cc([C@
@](C)(NCc2ccn3nccc3n2)C2CCCC2)nn1',None)
2.1825663314354418
[(4, 'S'), (22, 'S')] [(4, 'S'), (22, 'S')]


If you want to get reproducible results from the conformation generation,
you need to provide a random number seed to the embedding call:

In [87]: 
run_smiles('CCCC[C@@H](C(=O)COc1c(F)c(F)cc(F)c1F)n1cc([C@@](C)(NCc2ccn3nccc3n2)C2CCCC2)nn1',0xf00d)

1.4715514179429237e-07
[(4, 'S'), (22, 'S')] [(4, 'S'), (22, 'S')]


Here's a demo of the same effect for the version of your molecule with
partially specified stereo:

In [89]: run_smiles('CCCCC(C(=O)COc1c(F)c(F)cc(F)c1F)n1cc([C@
@](C)(NCc2ccn3nccc3n2)C2CCCC2)nn1',None)
2.889576407915814
[(4, 'R'), (22, 'S')] [(4, 'S'), (22, 'S')]

In [90]: run_smiles('CCCCC(C(=O)COc1c(F)c(F)cc(F)c1F)n1cc([C@
@](C)(NCc2ccn3nccc3n2)C2CCCC2)nn1',0xf00d)
1.0405439864921206e-07
[(4, 'R'), (22, 'S')] [(4, 'R'), (22, 'S')]

In [91]: run_smiles('CCCCC(C(=O)COc1c(F)c(F)cc(F)c1F)n1cc([C@
@](C)(NCc2ccn3nccc3n2)C2CCCC2)nn1',23)
0.0
[(4, 'S'), (22, 'S')] [(4, 'S'), (22, 'S')]



I hope this helps,
-greg

On Thu, May 23, 2019 at 8:34 PM Alfredo Quevedo <maquevedo....@gmail.com>
wrote:

> Good afternoon,
>
> I am generating 3D structures in .pdb format starting from smiles
> notation. I am working with chiral molecules, and noted that the 3D
> structure obtained when generating the molecule explicitly indicating
> the configuration of the stereocenter is sligthly different from that
> obtained when the configuration of the stereocenter is not explicitly
> indicated. A pair of examples are as follows:
>
> #No explicit stereochemistry in the 1st chiral center
> m =
> Chem.MolFromSmiles('CCCCC(C(=O)COc1c(F)c(F)cc(F)c1F)n1cc([C@
> @](C)(NCc2ccn3nccc3n2)C2CCCC2)nn1´)
> m3d = Chem.AddHs(m)
> AllChem.EmbedMolecule(m3d, AllChem.ETKDG())
> print(Chem.MolToPDBBlock(m3d, flavor=2),file = open(' molecule1.pdb',
> 'w+'))
>
>
> #Explicit stereochemistry in the 1st chiral center
> m =
> Chem.MolFromSmiles('CCCC[C@@H](C(=O)COc1c(F)c(F)cc(F)c1F)n1cc([C@
> @](C)(NCc2ccn3nccc3n2)C2CCCC2)nn1´)
> m3d = Chem.AddHs(m)
> AllChem.EmbedMolecule(m3d, AllChem.ETKDG())
> print(Chem.MolToPDBBlock(m3d, flavor=2),file = open(' molecule2.pdb',
> 'w+'))
>
>
> So in both cases I get the S,S enantiomers, but molecule1 is not
> perfectly superimposable on molecule2, which leads to further resp
> fitted charges differences.
>
> Any idea why expicitly indicating the configuration of the first
> stereocenter is giving slighlty different results regarding the 3D
> conformation compared to not indicating it?
>
> Thanks in advance for the help,
>
> best regards
>
> Alfredo
>
>
>
>
>
>
>
> _______________________________________________
> Rdkit-discuss mailing list
> Rdkit-discuss@lists.sourceforge.net
> https://lists.sourceforge.net/lists/listinfo/rdkit-discuss
>
_______________________________________________
Rdkit-discuss mailing list
Rdkit-discuss@lists.sourceforge.net
https://lists.sourceforge.net/lists/listinfo/rdkit-discuss

Reply via email to