On Fri, Mar 11, 2016 at 5:41 AM, Peter S. Shenkin <shen...@gmail.com> wrote:

> Seems that somewhere in the guts of RDKit there might well be code that
> divides atoms into equivalence classes.
>
> In the most common (tetravalent, tetrahedral) chiral situation, if the
> tetrahedral center's four connected atoms fell into 3 equivalence classes,
> the center would be prochiral. Then we'd have something like aC(b)(b)d,
> where a, b, and d represent the equivalence classes of the directly
> attached atoms.
>

Assuming I'm understanding you correctly, that's what the code below is
doing. When you tell CanonicalRankAtoms() not to breakTies, it just finds
the equivalence classes.


> In the past, for a core-hopping application, I've stored a core library as
> the SMILES of the stripped and H-substituted core. I've wished that I could
> use SMILES '@' notation to specify which of the two H's originally held the
> substituent (when only one did).
>
>
You could do that with the current RDKit by using either an isotopic label
or atom map number on the H's that replaced the side chain, *or* by using
dummies to mark the attachment points (this last is what
Chem.ReplaceSidechains does)

The atom map number is used in the canonicalization, as this example
demonstrates:

In [34]: m = Chem.MolFromSmiles('C[C@](F)(Cl)Br')

In [35]: m.GetAtomWithIdx(2).SetAtomicNum(1)

In [36]: m.GetAtomWithIdx(3).SetAtomicNum(1)

In [37]: Chem.MolToSmiles(m,True)
Out[37]: '[H][C@]([H])(C)Br'

In [38]: m.GetAtomWithIdx(2).SetProp('molAtomMapNumber','3')

In [39]: m.GetAtomWithIdx(3).SetProp('molAtomMapNumber','1')

In [40]: Chem.MolToSmiles(m,True)
Out[40]: 'C[C@@](Br)([H:1])[H:3]'

In [41]: m.GetAtomWithIdx(2).SetProp('molAtomMapNumber','1')

In [42]: m.GetAtomWithIdx(3).SetProp('molAtomMapNumber','3')

In [43]: Chem.MolToSmiles(m,True)
Out[43]: 'C[C@](Br)([H:1])[H:3]'


On Thu, Mar 10, 2016 at 10:53 PM, Greg Landrum <greg.land...@gmail.com>
> wrote:
>
>> This isn't an area I've thought much about, so this may be a bit naive.
>>
>> It seems like the interesting atom from the perspective of perception is
>> the carbon that the Hs are attached to, not the Hs themselves; it's the
>> carbon that will become a chiral center.
>>
>> If we neglect dependent stereochemistry for the moment, can't we just
>> find carbon atoms that have two hydrogen substituents and then look to see
>> if they (the carbons) have two differently ranked neighbors? The advantage
>> to this is that it saves the addHs step and will generally be a lot faster
>> to execute.
>>
>> Here's a quick code snippet, it lets you choose between using CIP ranks
>> to distinguish atoms or just using their canonical rankings.
>>
>> def findProchiral(mol,useCIPRanks=False):
>>     if not useCIPRanks:
>>         ranks = Chem.CanonicalRankAtoms(mol,breakTies=False)
>>     else:
>>         Chem.AssignStereochemistry(mol)
>>         ranks = [x.GetProp('_CIPRank') for x in mol.GetAtoms()]
>>     res = []
>>     for atom in mol.GetAtoms():
>>         # only consider atoms with two heavy neighbors and two Hs.
>>         # could probably be further limited to just consider carbons
>>         if atom.GetTotalDegree()!=4 or atom.GetTotalNumHs()!=2:
>>             continue
>>         hvyNbrRanks=[]
>>         for nbr in atom.GetNeighbors():
>>             if nbr.GetAtomicNum()>1:
>>                 hvyNbrRanks.append(ranks[nbr.GetIdx()])
>>                 if len(hvyNbrRanks)==2:
>>                     break
>>         if hvyNbrRanks[0] != hvyNbrRanks[1]:
>>             res.append(atom.GetIdx())
>>     return res
>>
>>
>> Is that headed in the right direction?
>> -greg
>>
>>
>>
>> On Thu, Mar 10, 2016 at 5:30 PM, Peter S. Shenkin <shen...@gmail.com>
>> wrote:
>>
>>> Is the canonical rank of prochiral H's different or the same? (For
>>> example the rank of the H's on C-1 of ethyl chloride.)
>>>
>>> Thanks,
>>> -P.
>>>
>>
>>
>
------------------------------------------------------------------------------
Transform Data into Opportunity.
Accelerate data analysis in your applications with
Intel Data Analytics Acceleration Library.
Click to learn more.
http://pubads.g.doubleclick.net/gampad/clk?id=278785111&iu=/4140
_______________________________________________
Rdkit-discuss mailing list
Rdkit-discuss@lists.sourceforge.net
https://lists.sourceforge.net/lists/listinfo/rdkit-discuss

Reply via email to