Hello Paolo,
Many thanks for your detailed explanation. I think it is all clear now. I
was not aware that a C with a single @ would always be tagged as CCW in a
Smarts due to the missing implicit H. Pity I cannot add a H  to get the
correct CW tag since it would no longer match the superstructure...
Good that it is less complicated when using Smiles rather than Smarts. For
my purpose, I can use Smiles too so that's what I'll do.
Thanks again.
Best,
Alexis


On Tue, 4 Jan 2022 at 19:51, Paolo Tosco <paolo.tosco.m...@gmail.com> wrote:

> Hi Alexis,
>
> The chiral substructure match does not look at the CIP labels, but rather
> at the atom parities, as you expected.
>
> In the examples below I have reordered your mol_structure as it makes
> things easier to understand.
>
> Let's work with a fragment constructed from SMILES rather than from SMARTS:
>
> from rdkit import Chem
> from rdkit.Chem.Draw import IPythonConsole
> IPythonConsole.drawOptions.addAtomIndices = True
>
> mol_structure = Chem.MolFromSmiles("C[C@](N)1CCCN(C1)C(C)=O")
> mol_structure
> [image: image.png]
> mol_frag = Chem.MolFromSmiles("C[C@H]1CCCN(C1)C(C)=O")
> mol_frag
> [image: image.png]
> In this case, the chiral tags are both CW and we have a match, and the
> fact that the CIP chirality is opposite is not relevant, as you would
> expect:
>
> print("mol_structure", mol_structure.GetAtomWithIdx(1).GetChiralTag(),
> mol_structure.GetAtomWithIdx(1).GetProp("_CIPCode"))
> print("mol_frag", mol_frag.GetAtomWithIdx(1).GetChiralTag(),
> mol_frag.GetAtomWithIdx(1).GetProp("_CIPCode"))
> print(mol_structure.HasSubstructMatch(mol_frag, useChirality=True))
> mol_structure CHI_TETRAHEDRAL_CW R
> mol_frag CHI_TETRAHEDRAL_CW S
> True
>
> If we reverse the parity, the CIP chirality becomes identical, but as
> expected this time the substructure match fails:
>
> mol_frag = Chem.MolFromSmiles("C[C@@H]1CCCN(C1)C(C)=O")
> mol_frag
> [image: image.png]
>
> print("mol_structure", mol_structure.GetAtomWithIdx(1).GetChiralTag(),
> mol_structure.GetAtomWithIdx(1).GetProp("_CIPCode"))
> print("mol_frag", mol_frag.GetAtomWithIdx(1).GetChiralTag(),
> mol_frag.GetAtomWithIdx(1).GetProp("_CIPCode"))
> print(mol_structure.HasSubstructMatch(mol_frag, useChirality=True))
> mol_structure CHI_TETRAHEDRAL_CW R
> mol_frag CHI_TETRAHEDRAL_CCW R
> False
>
> I believe the above matches your expectations.
>
> Let's now see what happens when we switch to a fragment constructed from
> SMARTS:
>
> mol_frag = Chem.MolFromSmarts("C[C@]1CCCN(C1)C(C)=O")
> mol_frag
> [image: image.png]
>
> In this case, the chiral center has only 3 neighbors: in fact, there is no
> implicit H as there was in the fragment constructed from SMILES; the 4th
> substituent is simply not defined by the SMARTS string.
> When we look at the chiral tag assigned by RDKit, we see that this time
> the chiral tag is CCW, while the matching atom in mol_structure is CW;
> this is why the match fails:
>
> print("mol_structure", mol_structure.GetAtomWithIdx(1).GetChiralTag())
> print("mol_frag", mol_frag.GetAtomWithIdx(1).GetChiralTag())
> print(mol_structure.HasSubstructMatch(mol_frag, useChirality=True))
> mol_structure CHI_TETRAHEDRAL_CW
> mol_frag CHI_TETRAHEDRAL_CCW
> False
>
> My guess is that in this case atom 1 is assigned a CCW parity because this
> is what is specified in the SMARTS string with the @ symbol; given that
> there are less than 3 bonds, the parity cannot be determined based on bond
> indices as normally happens, so I guess the parity is simply set based on
> what is specified by the user in the SMARTS string, e.g.:
>
> mol_frag = Chem.MolFromSmarts("[C@@]")
> print("mol_frag", mol_frag.GetAtomWithIdx(0).GetChiralTag())
> mol_frag CHI_TETRAHEDRAL_CW
>
> mol_frag = Chem.MolFromSmarts("[C@]")
> print("mol_frag", mol_frag.GetAtomWithIdx(0).GetChiralTag())
> mol_frag CHI_TETRAHEDRAL_CCW
>
> Maybe Greg can shed more light on this.
>
> Cheers,
> p.
>
>
> On Tue, Jan 4, 2022 at 1:07 PM Alexis Parenty <
> alexis.parenty.h...@gmail.com> wrote:
>
>> Hi Shani,
>> Many thanks for your response. Good observation! I did not realise that
>> the stereoisomer nomenclature had swap on removing the amine. However, is
>> it the intended behaviour of the method?
>> To me, in that case, the S stereoisomer fragment should match the R
>> stereoisomer super structure but it does not. The R stereoisomer fragment
>> matches the R stereoisomer super structure when It should not. What do you
>> think?
>>
>> Best,
>>
>> Alexis
>>
>>
>>
>> On Tue, 4 Jan 2022 at 12:42, Shani Zev <levishan...@gmail.com> wrote:
>>
>>> Hi Alexis,
>>> I think that the reason for this is because in the first case you are
>>> comparing S and R stereoisomer and hence it is false while in the second
>>> you are comparing R and R.
>>> best,
>>> Shani
>>>
>>> [image: q.gif]
>>>
>>> On Tue, Jan 4, 2022 at 1:22 PM Alexis Parenty <
>>> alexis.parenty.h...@gmail.com> wrote:
>>>
>>>> Hi everyone,
>>>>
>>>> Why is it that the following smarts C[C@]1CCCN(C1)C(C)=O does not
>>>> match the following structure C[C@]1(CCCN(C1)C(C)=O)N  when using the
>>>> chirality argument in the HasSubstructMatch method?
>>>>
>>>> mol_frag = Chem.MolFromSmarts("C[C@]1CCCN(C1)C(C)=O")
>>>> mol_structure = Chem.MolFromSmiles("C[C@]1(CCCN(C1)C(C)=O)N")
>>>>
>>>> print(mol_structure.HasSubstructMatch(mol_frag, useChirality=True))
>>>> ==> False
>>>>
>>>>
>>>> worst, the other isomeric fragment matches when it should not:
>>>>
>>>> mol_frag = Chem.MolFromSmarts("C[C@@]1CCCN(C1)C(C)=O")
>>>> mol_structure = Chem.MolFromSmiles("C[C@]1(CCCN(C1)C(C)=O)N")
>>>>
>>>> print(mol_structure.HasSubstructMatch(mol_frag, useChirality=True))
>>>> ==> True
>>>>
>>>> Am I missing something?
>>>>
>>>> Many thanks,
>>>>
>>>> Alexis
>>>> _______________________________________________
>>>> 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
>>
>
_______________________________________________
Rdkit-discuss mailing list
Rdkit-discuss@lists.sourceforge.net
https://lists.sourceforge.net/lists/listinfo/rdkit-discuss

Reply via email to