Perfect, thank you!

On Fri, Sep 15, 2017 at 9:47 AM, Greg Landrum <[email protected]> wrote:
>
>
> On Fri, Sep 15, 2017 at 9:25 AM, Michał Nowotka <[email protected]> wrote:
>>
>> Thanks Greg, very helpful!
>> Can you tell me how should I modify my code to provide the list of
>> bonds to be highlighted, so I get the image generated by DrawMolecules
>> looking the same way as produced by DrawMolecule?
>
>
> Sure, here's an example GIST with the code:
> https://gist.github.com/greglandrum/431483ac1f9edb03b09c8577031c10e0
>
>
> -greg
>
>> On Thu, Sep 14, 2017 at 4:36 PM, Greg Landrum <[email protected]>
>> wrote:
>> > Hi Michal,
>> >
>> > There are a couple of things in here.
>> >
>> > On Fri, Aug 25, 2017 at 11:42 AM, Michał Nowotka <[email protected]>
>> > wrote:
>> >>
>> >> Hi,
>> >>
>> >> I finally decided to try the new C++ drawing code and I found some
>> >> issues with it. I'll try to descibe my problems.
>> >>
>> >> First, lets start with a code that works perfectly:
>> >>
>> >> from rdkit import Chem
>> >> from rdkit.Chem.AllChem import Compute2DCoords
>> >> from rdkit.Chem.Draw import rdMolDraw2D
>> >>
>> >> m = Chem.MolFromSmiles('O=C(C)Oc1ccccc1C(=O)O')
>> >> Compute2DCoords(m)
>> >> drawer = rdMolDraw2D.MolDraw2DCairo(500, 500)
>> >>
>> >>
>> >> drawer.DrawMolecule(m,highlightAtoms=m.GetSubstructMatch(Chem.MolFromSmarts('c1ccccc1')))
>> >> drawer.FinishDrawing()
>> >> with open('aspirin_a.png','wb') as f:
>> >>     f.write(drawer.GetDrawingText())
>> >>
>> >> This code produces the attached 'aspirin_a.png' image, that looks
>> >> perfect.
>> >
>> >
>> > Yay! :-)
>> >
>> >>
>> >> Now, apart from the `DrawMolecule` there is also `DrawMolecules`
>> >> function exposed by the `MolDraw2DCairo` module. So what will happen
>> >> when I use it to render the same molecule?
>> >>
>> >> from rdkit import Chem
>> >> from rdkit.Chem.AllChem import Compute2DCoords
>> >> from rdkit.Chem.Draw import rdMolDraw2D
>> >>
>> >> m = Chem.MolFromSmiles('O=C(C)Oc1ccccc1C(=O)O')
>> >> Compute2DCoords(m)
>> >> drawer = rdMolDraw2D.MolDraw2DCairo(500, 500)
>> >>
>> >>
>> >> drawer.DrawMolecules([m],highlightAtoms=[m.GetSubstructMatch(Chem.MolFromSmarts('c1ccccc1'))])
>> >> drawer.FinishDrawing()
>> >> with open('aspirin_b.png','wb') as f:
>> >>     f.write(drawer.GetDrawingText())
>> >>
>> >>
>> >> The image is different - only atoms are highlighted, bonds between
>> >> atoms are not (see the attached 'aspirin_b.png') image.
>> >
>> >
>> > Short answer: DrawMolecules() requires you to provide the list of bonds
>> > to
>> > be highlighted too. This is just a couple of lines of code.
>> >
>> >>
>> >> OK, but who would use `DrawMolecules` to render a single compound
>> >> anyway? I think it's meant to render multiple compounds at once. The
>> >> 'old' drawing code has a function called `MolsToGridImage` to do this.
>> >> And it has a `molsPerRow` parameter which makes it easy to define a
>> >> shape of the grid. I couldn't find a similar function in the 'new'
>> >> drawing code so it looks like using `DrawMolecules` is the only way to
>> >> render multiple mols at once. Let's try that:
>> >>
>> >> from rdkit import Chem
>> >> from rdkit.Chem.AllChem import Compute2DCoords
>> >> from rdkit.Chem.Draw import rdMolDraw2D
>> >>
>> >> m1 = Chem.MolFromSmiles('O=C(C)Oc1ccccc1C(=O)O')
>> >> Compute2DCoords(m1)
>> >> m2 = Chem.MolFromSmiles('c1cccnc1O')
>> >> Compute2DCoords(m2)
>> >> drawer = rdMolDraw2D.MolDraw2DCairo(500, 500)
>> >> drawer.DrawMolecules([m1, m2])
>> >> drawer.FinishDrawing()
>> >> with open('mols.png','wb') as f:
>> >>     f.write(drawer.GetDrawingText())
>> >
>> >
>> > The MolDraw2D constructors can be called with additional arguments that
>> > provide the size of the panes used to render grids of molecules.
>> > For example, here's the code to draw two molecules side by side derived
>> > from
>> > your example:
>> > from rdkit import Chem
>> > from rdkit.Chem.Draw import rdMolDraw2D
>> >
>> > m1 = Chem.MolFromSmiles('O=C(C)Oc1ccccc1C(=O)O')
>> > pm1 = rdMolDraw2D.PrepareMolForDrawing(m1)
>> > m2 = Chem.MolFromSmiles('c1cccnc1O')
>> > pm2 = rdMolDraw2D.PrepareMolForDrawing(m2)
>> > drawer = rdMolDraw2D.MolDraw2DCairo(600,300,300,300)
>> > drawer.DrawMolecules([m1, m2])
>> > drawer.FinishDrawing()
>> > with open('mols.png','wb') as f:
>> >     f.write(drawer.GetDrawingText())
>> >
>> >
>> > I've switched to use PrepareMolForDrawing() (because I think it's nice
>> > to
>> > have kekule structures and wedged bonds) and provided a canvas size of
>> > 600x300 divided into 300x300 panels. The code figures out that this
>> > means 2
>> > mols per row.
>> >
>> > It's somewhat curious (and I don't mean that in a good way) that this
>> > hasn't
>> > made it onto the RDKit master yet (instead of using the hack that's
>> > currently there). I will try to get that done before the next release.
>> >
>> > To summarize:
>> >>
>> >>
>> >> - `DrawMolecules` doesn't highlight bonds in the same way as the
>> >> `DrawMolecule` does
>> >
>> >
>> > correct. That's by design (correct method described above)
>> >
>> >> - There is no equivalent of the `MolsToGridImage` where I can define
>> >> the shape of the grid of molecules at least nothing documented in
>> >>
>> >>
>> >> http://www.rdkit.org/Python_Docs/rdkit.Chem.Draw.rdMolDraw2D.MolDraw2D-class.html
>> >
>> >
>> > Right, it's not documented, but hopefully what I have above helps you
>> > see
>> > how to do it.
>> >>
>> >>
>> >> - `DrawMolecules` does nothing to layout molecules in such a way that
>> >> they don't obscure each other
>> >
>> >
>> > Correct. It does what you tell it to do, but hopefully the rest of the
>> > answer tells you how to tell it the right thing.
>> >
>> > -greg
>> >
>> >
>> >>
>> >>
>> >> Regards,
>> >>
>> >> Michał Nowotka
>> >>
>> >>
>> >>
>> >> ------------------------------------------------------------------------------
>> >> Check out the vibrant tech community on one of the world's most
>> >> engaging tech sites, Slashdot.org! http://sdm.link/slashdot
>> >> _______________________________________________
>> >> Rdkit-discuss mailing list
>> >> [email protected]
>> >> https://lists.sourceforge.net/lists/listinfo/rdkit-discuss
>> >>
>> >
>
>

------------------------------------------------------------------------------
Check out the vibrant tech community on one of the world's most
engaging tech sites, Slashdot.org! http://sdm.link/slashdot
_______________________________________________
Rdkit-discuss mailing list
[email protected]
https://lists.sourceforge.net/lists/listinfo/rdkit-discuss

Reply via email to