Dear Jean-Marc.
Then you only need addChiralHs=False:
Cheers, p. On 27 Jul 2023, at 12:01, Jean-Marc Nuzillard <jm.nuzill...@univ-reims.fr> wrote:
Dear Paolo,
your code works perfectly well.
Running it with the SMILES of trans-decalin,
C1CC[C@H]2CCCC[C@@H]2C1,
I obtained
<F293bnUXx7pCSmb2.png> .
The code with ReapplyMolBlockWedging() inside (that you proposed
me
in a previous post) was intended to have no H atom displayed at
ring junctions.
Sure, the SMILES of trans-decalin I propose, generated by rdkit
from a molblock
devoid of H atoms, has explicit H atoms inside at the ring
junction.
So, starting from an isomeric SMILES, I would like to produce PNG
drawings that retain
the existing configuration information and do not show H atoms at
ring junctions.
Best,
Jean-Marc
Le 27/07/2023 à 10:33, Paolo Tosco a écrit :
Dear Jean-Marc,
You are generating the molecule from SMILES,
therefore it does not have molblock wedging information.
When you call ReapplyMolBlockWedging(), first
existing wedging info will be stripped.
Then, the molblock wedging will be applied, but
there is none.
Hence, you get no wedging.
You may simplify your code and get the expected
wedging as follows:
Cheers.
p.
Hi David,
thank you for your suggestion.
Setting wedgeBonds=True does not change anything,
probably because wedging is performed later by
the call to Chem.ReapplyMolBlockWedging()
Best,
Jean-Marc
Le 26/07/2023 à 22:28, David Cosgrove a écrit :
Hi,
I’m away from my computer at the moment, so
can’t try anything, but I wonder if it’s anything to do
with the ‘wedgeBonds=False’ option you gave when preparing
the drawing.
Dave
Dear all,
I use the following code to produce PNG drawings. I
use RDKit version 2023.03.1 .
The SMILES chain describes a molecule with a single
chiral center of defined configuration.
from rdkit import Chem
from rdkit.Chem import rdCoordGen
from rdkit.Chem.Draw import rdMolDraw2D
from PIL import Image
from io import BytesIO
smi = "C=C[C@@](/C=C/c1ccc(cc1)O)(CCC=C(C)C)C"
filenameOut = "img.png"
mol = Chem.MolFromSmiles(smi)
rdCoordGen.AddCoords(mol)
print(Chem.MolToMolBlock(mol))
d2d = rdMolDraw2D.MolDraw2DCairo(350, 300)
dopts = d2d.drawOptions()
dopts.baseFontSize = 0.6
dopts.prepareMolsBeforeDrawing = False
mol_draw = rdMolDraw2D.PrepareMolForDrawing(mol,
addChiralHs=False, wedgeBonds=False)
Chem.ReapplyMolBlockWedging(mol_draw)
d2d.DrawMolecule(mol_draw, legend='',
highlightAtoms=[])
d2d.FinishDrawing()
bio = BytesIO(d2d.GetDrawingText())
draw_code = Image.open(bio)
draw_code.save(filenameOut)
The resulting image does not show the chirality
wedge:
<wAS5dAnxRREFAFxp.png>
The script prints the MolBlock that comes from the
SMILES and the calculation of the 2D atomic
coordinates:
_________________
RDKit 2D
19 19 0 0 0 0 0 0 0 0999 V2000
2.0515 1.2242 0.0000 C 0 0 0 0 0
0 0 0 0 0 0 0
1.0515 1.2214 0.0000 C 0 0 0 0 0
0 0 0 0 0 0 0
0.5539 0.3540 0.0000 C 0 0 0 0 0
0 0 0 0 0 0 0
-0.4459 0.3512 0.0000 C 0 0 0 0 0
0 0 0 0 0 0 0
-0.9435 -0.5162 0.0000 C 0 0 0 0 0
0 0 0 0 0 0 0
-1.9435 -0.5190 0.0000 C 0 0 0 0 0
0 0 0 0 0 0 0
-2.4411 -1.3866 0.0000 C 0 0 0 0 0
0 0 0 0 0 0 0
-3.4411 -1.3892 0.0000 C 0 0 0 0 0
0 0 0 0 0 0 0
-3.9435 -0.5246 0.0000 C 0 0 0 0 0
0 0 0 0 0 0 0
-3.4459 0.3428 0.0000 C 0 0 0 0 0
0 0 0 0 0 0 0
-2.4459 0.3456 0.0000 C 0 0 0 0 0
0 0 0 0 0 0 0
-4.9435 -0.5274 0.0000 O 0 0 0 0 0
0 0 0 0 0 0 0
1.4215 -0.1436 0.0000 C 0 0 0 0 0
0 0 0 0 0 0 0
2.2861 0.3588 0.0000 C 0 0 0 0 0
0 0 0 0 0 0 0
3.1535 -0.1388 0.0000 C 0 0 0 0 0
0 0 0 0 0 0 0
4.0181 0.3636 0.0000 C 0 0 0 0 0
0 0 0 0 0 0 0
4.0153 1.3636 0.0000 C 0 0 0 0 0
0 0 0 0 0 0 0
4.8855 -0.1338 0.0000 C 0 0 0 0 0
0 0 0 0 0 0 0
0.5567 -0.6460 0.0000 C 0 0 0 0 0
0 0 0 0 0 0 0
1 2 2 0
2 3 1 0
3 4 1 0
4 5 2 0
5 6 1 0
6 7 2 0
7 8 1 0
8 9 2 0
9 10 1 0
10 11 2 0
9 12 1 0
3 13 1 0
13 14 1 0
14 15 1 0
15 16 2 0
16 17 1 0
16 18 1 0
3 19 1 1
11 6 1 0
M END
_____________________
The wedge bond (3-19) is apparently here but is not
drawn as such.
Is there a remedy for that?
Best regards,
Jean-Marc
--
Jean-Marc Nuzillard
Directeur de Recherches au CNRS
Institut de Chimie Moléculaire de Reims
CNRS UMR 7312
Moulin de la Housse
CPCBAI, Bâtiment 18
BP 1039
51687 REIMS Cedex 2
France
ORCID : 0000-0002-5120-2556
Tel : +33 (0)3 26 91 82 10
http://www.univ-reims.fr/icmr
https://nuzillard.github.io/PyLSD
_______________________________________________
Rdkit-discuss mailing list
Rdkit-discuss@lists.sourceforge.net
https://lists.sourceforge.net/lists/listinfo/rdkit-discuss
--
--
Jean-Marc Nuzillard
Directeur de Recherches au CNRS
Institut de Chimie Moléculaire de Reims
CNRS UMR 7312
Moulin de la Housse
CPCBAI, Bâtiment 18
BP 1039
51687 REIMS Cedex 2
France
ORCID : 0000-0002-5120-2556
Tel : +33 (0)3 26 91 82 10
http://www.univ-reims.fr/icmr
https://nuzillard.github.io/PyLSD
_______________________________________________
Rdkit-discuss mailing list
Rdkit-discuss@lists.sourceforge.net
https://lists.sourceforge.net/lists/listinfo/rdkit-discuss
--
Jean-Marc Nuzillard
Directeur de Recherches au CNRS
Institut de Chimie Moléculaire de Reims
CNRS UMR 7312
Moulin de la Housse
CPCBAI, Bâtiment 18
BP 1039
51687 REIMS Cedex 2
France
ORCID : 0000-0002-5120-2556
Tel : +33 (0)3 26 91 82 10
http://www.univ-reims.fr/icmr
https://nuzillard.github.io/PyLSD
|