Hi Bob,
On Tue, Dec 30, 2014 at 8:42 PM, Bob Funchess <bfunch...@kelaroo.com> wrote:
>
>
> We’re using the C# wrappers to RDKit for an enumeration engine and we have
> run into an issue with stereochemistry in RDKit.
>
>
>
> The basic problem is that when we use a reaction which adds a stereocenter
> to a molecule, chirality at carbons that were present in the reactants is
> kept, but chirality at carbons introduced by the reaction itself is lost.
>
>
> In the course of trying to figure this out, I created a Python script that
> goes through the steps performed by the enumeration engine for a sample
> illustrating the issue so that I could try to debug the problem more
> easily, only to find that in Python the chirality is retained throughout
> the molecule. I then translated that Python script back into a standalone
> C# program … and the problem reappears.
>
>
>
> The output molecule should have two chiral centers (the carbon attached to
> the fluorine and the carbon attached to the chlorine). The Python version
> shows both of these as chiral; the C# version has only the fluorine carbon
> as chiral.
>
>
Hmm, when I run stereo.py with the 2014.09.2 release of the RDKit I also
only get one stereo-center:
(py34)~/Downloads > python stereo.py
[07:03:18] product atom-mapping number 7 not found in reactants.
[07:03:18] product atom-mapping number 8 not found in reactants.
[07:03:18] product atom-mapping number 9 not found in reactants.
[07:03:18] product atom-mapping number 10 not found in reactants.
[07:03:18] mapped atoms in the reactants were not mapped in the products.
unmapped numbers are: 6
RDKit 2D
9 8 0 0 0 0 0 0 0 0999 V2000
-3.8971 -0.7500 0.0000 C 0 0 0 0 0 0 0 0 0 0 0 0
-2.5981 0.0000 0.0000 C 0 0 0 0 0 0 0 0 0 0 0 0
-2.5981 1.5000 0.0000 F 0 0 0 0 0 0 0 0 0 0 0 0
-1.2990 -0.7500 0.0000 C 0 0 0 0 0 0 0 0 0 0 0 0
0.0000 0.0000 0.0000 C 0 0 0 0 0 0 0 0 0 0 0 0
1.2990 -0.7500 0.0000 C 0 0 0 0 0 0 0 0 0 0 0 0
2.5981 -0.0000 0.0000 C 0 0 0 0 0 0 0 0 0 0 0 0
3.8971 -0.7500 0.0000 C 0 0 0 0 0 0 0 0 0 0 0 0
2.5981 1.5000 0.0000 Cl 0 0 0 0 0 0 0 0 0 0 0 0
2 1 1 1
2 3 1 0
2 4 1 0
4 5 1 0
5 6 1 0
7 6 1 0
7 8 1 0
7 9 1 0
M END
CC(Cl)CCC[C@@H](C)F
This is, I think, more or less what I expect to happen based on the current
state of the code. You need to tell the reaction engine that you want
stereochemistry the stereochemistry you provide in the products to be
maintained. There's an atom property for this, molInversionFlag. Here's an
example of how to set it:
template_p.GetAtomWithIdx(6).SetProp('molInversionFlag','4')
The attached stereo2.py does this and produces the expected results (I
think).
The current behavior seems like a bug, but I'm going to need to think for a
bit about how to fix it. But at least this may be a workaround for now?
>
>
> I’ve attached the Python script and the C# source file. The C# version has
> the Python commands included as comments for easy comparison.
>
>
>
> Does anyone know what’s going on here? The only significant difference
> between the two that I can see is that Python has rxn.Initialize() while C#
> has rxn.initReactantMatchers(). Could this be the cause, and if so what C#
> method should I be using instead?
>
ChemicalReaction.Initialize() from Python is just calling
initReactantMatchers(). I guess you were using different versions of the
underlying C++ code.
-greg
from rdkit import Chem
from rdkit.Chem import AllChem
molblock_r = """
Mrv0541 12291411382D
6 5 0 0 1 0 999 V2000
-2.8875 0.7145 0.0000 Br 0 0 0 0 0 0 0 0 0 6 0 0
-2.0625 0.7145 0.0000 C 0 0 0 0 0 0 0 0 0 5 0 0
-1.6500 0.0000 0.0000 C 0 0 0 0 0 0 0 0 0 4 0 0
-2.0625 -0.7145 0.0000 C 0 0 1 0 0 0 0 0 0 1 0 0
-1.6500 -1.4289 0.0000 C 0 0 0 0 0 0 0 0 0 3 0 0
-2.8875 -0.7145 0.0000 F 0 0 0 0 0 0 0 0 0 2 0 0
1 2 1 0 0 0 0
2 3 1 0 0 0 0
3 4 1 0 0 0 0
4 5 1 0 0 0 0
4 6 1 1 0 0 0
M END"""
molblock_p = """"
Mrv0541 12291411382D
9 8 0 0 1 0 999 V2000
3.1821 0.6188 0.0000 C 0 0 0 0 0 0 0 0 0 3 0 0
3.8966 0.2063 0.0000 C 0 0 2 0 0 0 0 0 0 1 0 0
3.8966 -0.6187 0.0000 F 0 0 0 0 0 0 0 0 0 2 0 0
4.6111 0.6188 0.0000 C 0 0 0 0 0 0 0 0 0 4 0 0
5.3256 0.2062 0.0000 C 0 0 0 0 0 0 0 0 0 5 0 0
6.0400 0.6187 0.0000 C 0 0 0 0 0 0 0 0 0 7 0 0
6.7545 0.2062 0.0000 C 0 0 1 0 0 0 0 0 0 8 0 0
7.4690 0.6187 0.0000 C 0 0 0 0 0 0 0 0 0 9 0 0
6.7545 -0.6188 0.0000 Cl 0 0 0 0 0 0 0 0 0 10 0 0
1 2 1 0 0 0 0
2 3 1 1 0 0 0
2 4 1 0 0 0 0
4 5 1 0 0 0 0
5 6 1 0 0 0 0
7 6 1 0 0 0 0
7 8 1 0 0 0 0
7 9 1 1 0 0 0
M END
"""
rxn = Chem.rdChemReactions.ChemicalReaction()
template_r = Chem.MolFromMolBlock(molblock_r,True,False)
reactant_H = Chem.AddHs(template_r,True)
rxn.AddReactantTemplate(reactant_H)
template_p = Chem.MolFromMolBlock(molblock_p,True,False)
template_p.GetAtomWithIdx(6).SetProp('molInversionFlag','4')
product_H = Chem.AddHs(template_p,True)
rxn.AddProductTemplate(product_H)
rxn.Initialize()
molblock_r0 = """"
Mrv0541 12291411372D
6 5 0 0 1 0 999 V2000
-2.8875 0.7145 0.0000 Br 0 0 0 0 0 0 0 0 0 6 0 0
-2.0625 0.7145 0.0000 C 0 0 0 0 0 0 0 0 0 5 0 0
-1.6500 0.0000 0.0000 C 0 0 0 0 0 0 0 0 0 4 0 0
-2.0625 -0.7145 0.0000 C 0 0 1 0 0 0 0 0 0 1 0 0
-1.6500 -1.4289 0.0000 C 0 0 0 0 0 0 0 0 0 3 0 0
-2.8875 -0.7145 0.0000 F 0 0 0 0 0 0 0 0 0 2 0 0
1 2 1 0 0 0 0
2 3 1 0 0 0 0
3 4 1 0 0 0 0
4 5 1 0 0 0 0
4 6 1 1 0 0 0
M END
"""
r0 = Chem.MolFromMolBlock(molblock_r0)
Chem.SanitizeMol(r0)
r0_H = Chem.AddHs(r0,False)
ps = rxn.RunReactants((r0_H,))
prod = Chem.RemoveHs(ps[0][0],False)
Chem.rdDepictor.Compute2DCoords(prod)
print(Chem.MolToMolBlock(prod,True))
print(Chem.MolToSmiles(prod,True))
------------------------------------------------------------------------------
Dive into the World of Parallel Programming! The Go Parallel Website,
sponsored by Intel and developed in partnership with Slashdot Media, is your
hub for all things parallel software development, from weekly thought
leadership blogs to news, videos, case studies, tutorials and more. Take a
look and join the conversation now. http://goparallel.sourceforge.net
_______________________________________________
Rdkit-discuss mailing list
Rdkit-discuss@lists.sourceforge.net
https://lists.sourceforge.net/lists/listinfo/rdkit-discuss