On Tue, Jan 25, 2011 at 10:10 AM, [email protected] <[email protected]> wrote:
> On Tue, Jan 25, 2011 at 5:51 AM, Greg Landrum <[email protected]> wrote:
>> The current release of the code doesn't really support this but it's
>> one of those easy-to-add features that make sense to me. I just
>> checked in some changes to the drawing code so that it supports the
>> ignoreHs argument:
>>
>> [5]>>> Draw.ShowMol(mh,ignoreHs=True)
>>
>> or:
>>
>> [6]>>> Draw.MolToFile(mh,'blah.png',ignoreHs=True)
>
> I downloaded the latest revision and tested it briefly.
>
> While this change is doing a good job at removing H atoms from the
> final image, they are still used to compute 2D coordinates so the
> resulting depiction is often distorted; do you think adding the same
> "IgnoreHs" parameter to AllChem.Compute2DCoords() would be much more
> involved?

Yeah, that's definitely more complicated.

This little bit of code should do what you want though:
#--------------------
from rdkit import Chem
from rdkit.Chem import AllChem

def Compute2DCoordsNoHs(mol):
    mh = Chem.RemoveHs(mol)
    AllChem.Compute2DCoords(mh)
    match = mol.GetSubstructMatch(mh)
    conf = Chem.Conformer(mol.GetNumAtoms())
    refconf = mh.GetConformer()
    for i,idx in enumerate(match):
        conf.SetAtomPosition(idx,refconf.GetAtomPosition(i))
    cid=mol.AddConformer(conf,assignId=True)
    return cid

if __name__=='__main__':
    testmol = Chem.AddHs(Chem.MolFromSmiles('c1ccccn1'))
    AllChem.Compute2DCoords(testmol)
    print Chem.MolToMolBlock(testmol)
    cid=Compute2DCoordsNoHs(testmol)
    print Chem.MolToMolBlock(testmol,confId=cid)
#--------------------


>
> I've implemented a workaround in my code (store a map of {indexH:
> indexNoH} pairs and use it to highlight the requested atom in the
> hydrogen-free depiction) but a more proper solution would be probably
> appreciated (not only by me :-) )

take a look at the above and see if it does what you want.

-greg

------------------------------------------------------------------------------
Special Offer-- Download ArcSight Logger for FREE (a $49 USD value)!
Finally, a world-class log management solution at an even better price-free!
Download using promo code Free_Logger_4_Dev2Dev. Offer expires 
February 28th, so secure your free ArcSight Logger TODAY! 
http://p.sf.net/sfu/arcsight-sfd2d
_______________________________________________
Rdkit-discuss mailing list
[email protected]
https://lists.sourceforge.net/lists/listinfo/rdkit-discuss

Reply via email to