1
2
3
4
5
6
7 import Chem
8 import RDConfig
9 import math
10
11 elemDict={
12 7:(0,0,1),
13 8:(1,0,0),
14 9:(.2,.8,.8),
15 15:(1,.5,0),
16 16:(.8,.8,0),
17 17:(0,.8,0),
18 35:(.5,.3,.1),
19 0:(.5,.5,.5),
20 }
21
32
34 dotsPerAngstrom = 30
35
36 atomLabelFontFace = "sans"
37 atomLabelFontSize = 12
38 atomLabelMinFontSize = 10
39
40 bondLineWidth = 1.2
41 dblBondOffset = .2
42 dblBondLengthFrac = .8
43
44 defaultColor = (1,0,0)
45 selectColor = (1,0,0)
46
47 colorBonds=True
48 noCarbonSymbols=True
49 includeAtomNumbers=False
50 atomNumberOffset=0
51
52 dash = (1,1)
53 atomPs = None
54 canvas = None
55 canvasSize=None
56
57 wedgeDashedBonds=False
58
64
71
72
87
91 if not lenFrac:
92 lenFrac = self.dblBondLengthFrac
93
94 dx = p2[0]-p1[0]
95 dy = p2[1]-p1[1]
96
97
98
99
100 fracP1 = p1[0]+offsetX,p1[1]+offsetY
101
102
103 frac = (1.-lenFrac)/2
104 fracP1 = fracP1[0]+dx*frac,\
105 fracP1[1]+dy*frac
106
107 fracP2 = fracP1[0]+dx*lenFrac,\
108 fracP1[1]+dy*lenFrac
109 return fracP1,fracP2
110
111 - def _offsetDblBond(self,p1,p2,bond,a1,a2,conf,dir=1,
112 lenFrac=None):
148
152 perp,offsetX,offsetY = self._getBondOffset(pos,nbrPos)
153 offsetX *=.75
154 offsetY *=.75
155 poly = ((pos[0],pos[1]),
156 (nbrPos[0]+offsetX,nbrPos[1]+offsetY),
157 (nbrPos[0]-offsetX,nbrPos[1]-offsetY))
158
159 if not dash:
160 addCanvasPolygon(canvas,poly,color=color)
161 elif self.wedgeDashedBonds and addCanvasDashedWedge:
162 addCanvasDashedWedge(canvas,poly[0],poly[1],poly[2],color=color)
163 else:
164 addCanvasLine(canvas,pos,nbrPos,linewidth=width*2,color=color,
165 dashes=dash)
166
167 - def _drawBond(self,canvas,bond,atom,nbr,pos,nbrPos,conf,
168 width=bondLineWidth,color=defaultColor,color2=None):
169 if bond.GetBondType() == Chem.BondType.SINGLE:
170 bDir = bond.GetBondDir()
171 if bDir in (Chem.BondDir.BEGINWEDGE,Chem.BondDir.BEGINDASH):
172
173 if bond.GetBeginAtom().GetChiralTag() in (Chem.ChiralType.CHI_TETRAHEDRAL_CW,
174 Chem.ChiralType.CHI_TETRAHEDRAL_CCW):
175 p1,p2 = pos,nbrPos
176 else:
177 p2,p1 = pos,nbrPos
178 if bDir==Chem.BondDir.BEGINWEDGE:
179 self._drawWedgedBond(canvas,bond,p1,p2,color=(0,0,0),width=width)
180 elif bDir==Chem.BondDir.BEGINDASH:
181 self._drawWedgedBond(canvas,bond,p1,p2,color=(0,0,0),width=width,
182 dash=self.dash)
183 else:
184 addCanvasLine(canvas,pos,nbrPos,linewidth=width,color=color,color2=color2)
185 elif bond.GetBondType() == Chem.BondType.DOUBLE:
186 if bond.IsInRing() or (atom.GetDegree()!=1 and bond.GetOtherAtom(atom).GetDegree()!=1):
187 addCanvasLine(canvas,pos,nbrPos,linewidth=width,color=color,color2=color2)
188 fp1,fp2 = self._offsetDblBond(pos,nbrPos,bond,atom,nbr,conf)
189 addCanvasLine(canvas,fp1,fp2,linewidth=width,color=color,color2=color2)
190 else:
191 fp1,fp2 = self._offsetDblBond(pos,nbrPos,bond,atom,nbr,conf,dir=.5,
192 lenFrac=1.0)
193 addCanvasLine(canvas,fp1,fp2,linewidth=width,color=color,color2=color2)
194 fp1,fp2 = self._offsetDblBond(pos,nbrPos,bond,atom,nbr,conf,dir=-.5,
195 lenFrac=1.0)
196 addCanvasLine(canvas,fp1,fp2,linewidth=width,color=color,color2=color2)
197 elif bond.GetBondType() == Chem.BondType.AROMATIC:
198 addCanvasLine(canvas,pos,nbrPos,linewidth=width,color=color,color2=color2)
199 fp1,fp2 = self._offsetDblBond(pos,nbrPos,bond,atom,nbr,conf)
200 addCanvasLine(canvas,fp1,fp2,linewidth=width,color=color,color2=color2,
201 dash=self.dash)
202
203 elif bond.GetBondType() == Chem.BondType.TRIPLE:
204 addCanvasLine(canvas,pos,nbrPos,linewidth=width,color=color,color2=color2)
205 fp1,fp2 = self._offsetDblBond(pos,nbrPos,bond,atom,nbr,conf)
206 addCanvasLine(canvas,fp1,fp2,linewidth=width,color=color,color2=color2)
207 fp1,fp2 = self._offsetDblBond(pos,nbrPos,bond,atom,nbr,conf,dir=-1)
208 addCanvasLine(canvas,fp1,fp2,linewidth=width,color=color,color2=color2)
209
211 canvasSize = self.canvasSize
212 xAccum = 0
213 yAccum = 0
214 nAts = 0
215 minX = 1e8
216 minY = 1e8
217 maxX = -1e8
218 maxY = -1e8
219
220 for atom in mol.GetAtoms():
221 pos = conf.GetAtomPosition(atom.GetIdx())
222 xAccum += pos[0]
223 yAccum += pos[1]
224 minX = min(minX,pos[0])
225 minY = min(minY,pos[1])
226 maxX = max(maxX,pos[0])
227 maxY = max(maxY,pos[1])
228 nAts += 1
229
230 dx = abs(maxX-minX)
231 dy = abs(maxY-minY)
232 xSize = dx*self.dotsPerAngstrom
233 ySize = dy*self.dotsPerAngstrom
234
235 if coordCenter:
236 molTrans = -xAccum/nAts,-yAccum/nAts
237 else:
238 molTrans = -(minX+(maxX-minX)/2),-(minY+(maxY-minY)/2)
239 self.dotsPerAngstrom=30.0
240 self.molTrans = molTrans
241 self.drawingTrans=(0,0)
242 tMin = self.transformPoint((minX,minY))
243 tMax = self.transformPoint((maxX,maxY))
244
245 if xSize>=.95*canvasSize[0]:
246 scale = .9*canvasSize[0]/xSize
247 xSize*=scale
248 ySize*=scale
249 self.dotsPerAngstrom*=scale
250 self.atomLabelFontSize = max(self.atomLabelFontSize*scale,
251 self.atomLabelMinFontSize)
252 if ySize>=.95*canvasSize[1]:
253 scale = .9*canvasSize[1]/ySize
254 xSize*=scale
255 ySize*=scale
256 self.dotsPerAngstrom*=scale
257 self.atomLabelFontSize = max(self.atomLabelFontSize*scale,
258 self.atomLabelMinFontSize)
259 drawingTrans = canvasSize[0]/2,canvasSize[1]/2
260 self.drawingTrans = drawingTrans
261 tMax = self.transformPoint((maxX,maxY))
262 tMin = self.transformPoint((minX,minY))
263
264 - def _drawLabel(self,canvas,label,pos,font,color=None,
265 highlightIt=False):
274
275 - def AddMol(self,mol,canvas=None,centerIt=True,molTrans=(0,0),drawingTrans=(0,0),
276 highlightAtoms=[],confId=-1):
277 """
278
279 Notes:
280 - specifying centerIt will cause molTrans and drawingTrans to be ignored
281
282 """
283 try:
284 dl = addCanvasLine
285 except NameError:
286 registerCanvas('sping')
287 if canvas is None:
288 canvas = self.canvas
289 else:
290 self.canvas = canvas
291 self.canvasSize=canvas.size
292
293 conf = mol.GetConformer(confId)
294
295 if centerIt:
296 self._scaleAndCenter(mol,conf)
297 else:
298 self.molTrans = molTrans
299 self.drawingTrans = drawingTrans
300 font = Font(face=self.atomLabelFontFace,size=self.atomLabelFontSize)
301
302 if not mol.HasProp('_drawingBondsWedged'):
303 Chem.WedgeMolBonds(mol,conf)
304 mol.SetProp('_drawingBondsWedged','set',True)
305
306 self.atomPs[mol] = {}
307 self.activeMol = mol
308 self.bondRings = mol.GetRingInfo().BondRings()
309 for atom in mol.GetAtoms():
310 idx = atom.GetIdx()
311 pos = self.atomPs[mol].get(idx,None)
312 if pos is None:
313 pos = self.transformPoint(conf.GetAtomPosition(idx))
314 self.atomPs[mol][idx] = pos
315 nbrSum = [0,0]
316 for bond in atom.GetBonds():
317 nbr = bond.GetOtherAtom(atom)
318 nbrIdx = nbr.GetIdx()
319 if nbrIdx > idx:
320 nbrPos = self.atomPs[mol].get(nbrIdx,None)
321 if nbrPos is None:
322 nbrPos = self.transformPoint(conf.GetAtomPosition(nbrIdx))
323 self.atomPs[mol][nbrIdx] = nbrPos
324
325 if highlightAtoms and idx in highlightAtoms and nbrIdx in highlightAtoms:
326 width=2.0*self.bondLineWidth
327 color = self.selectColor
328 color2 = self.selectColor
329 else:
330 width=self.bondLineWidth
331 if self.colorBonds:
332 color = elemDict.get(atom.GetAtomicNum(),(0,0,0))
333 color2 = elemDict.get(nbr.GetAtomicNum(),(0,0,0))
334 else:
335 color = self.defaultColor
336 color2= color
337
338
339
340 if atom.GetIdx()==bond.GetBeginAtomIdx():
341 self._drawBond(canvas,bond,atom,nbr,pos,nbrPos,conf,
342 color=color,width=width,color2=color2)
343 else:
344 self._drawBond(canvas,bond,nbr,atom,nbrPos,pos,conf,
345 color=color2,width=width,color2=color)
346 else:
347 nbrPos = self.atomPs[mol][nbrIdx]
348 nbrSum[0] += nbrPos[0]-pos[0]
349 nbrSum[1] += nbrPos[1]-pos[1]
350
351
352 base = atom.GetSymbol()
353 nHs = atom.GetTotalNumHs()
354 if nHs>0:
355 if nHs>1:
356 hs='H%d'%nHs
357 else:
358 hs ='H'
359 else:
360 hs = ''
361 chg = atom.GetFormalCharge()
362 if chg!=0:
363 if chg==1:
364 chg = '+'
365 elif chg==-1:
366 chg = '-'
367 elif chg>1:
368 chg = '+%d'%chg
369 elif chg<-1:
370 chg = '-%d'%chg
371 else:
372 chg = ''
373 if nbrSum[0]<=0:
374 symbol = '%s%s%s'%(base,hs,chg)
375 else:
376 symbol = '%s%s%s'%(chg,hs,base)
377
378
379 labelIt= not self.noCarbonSymbols or \
380 atom.GetAtomicNum()!=6 or \
381 atom.GetFormalCharge()!=0 or \
382 self.includeAtomNumbers
383 if self.includeAtomNumbers:
384 symbol = str(atom.GetIdx())
385 if labelIt:
386 color = elemDict.get(atom.GetAtomicNum(),(0,0,0))
387 self._drawLabel(canvas,symbol,pos,font,color=color,
388 highlightIt=(highlightAtoms and idx in highlightAtoms))
389
391 g= globals()
392 if canvasNm in ('sping','SPING'):
393 from spingCanvas import addCanvasLine,addCanvasText,addCanvasPolygon,addCanvasDashedWedge
394 elif canvasNm in ('agg','AGG'):
395 from aggCanvas import addCanvasLine,addCanvasText,addCanvasPolygon,addCanvasDashedWedge
396 elif canvasNm in ('mpl','MPL'):
397 from mplCanvas import addCanvasLine,addCanvasText,addCanvasPolygon
398 addCanvasDashedWedge=None
399 else:
400 raise ValueError,'unrecognized canvas type'
401 g['addCanvasLine']=addCanvasLine
402 g['addCanvasText']=addCanvasText
403 g['addCanvasPolygon']=addCanvasPolygon
404 g['addCanvasDashedWedge']=addCanvasDashedWedge
405
406 if __name__=='__main__':
407 import sys
408 if len(sys.argv)<2:
409 mol = Chem.MolFromSmiles('O=C1C([C@@H](F)C=CN[C@H](Cl)Br)C(c2c(O)c(NN)ccc2)=C1C#N')
410 else:
411 mol = Chem.MolFromSmiles(sys.argv[1])
412
413
414 Chem.Kekulize(mol)
415 from Chem import rdDepictor
416 rdDepictor.Compute2DCoords(mol)
417
418 if 1:
419 from aggdraw import Draw
420 registerCanvas('agg')
421 from PIL import Image
422 img = Image.new("RGBA",(300,300),"white")
423 canvas=Draw(img)
424 canvas.setantialias(True)
425 drawer = MolDrawing(canvas)
426 drawer.AddMol(mol)
427 canvas.flush()
428 img.save("foo.png")
429 elif 0:
430 from matplotlib.figure import Figure
431 from matplotlib.backends.backend_agg import FigureCanvasAgg
432 registerCanvas('mpl')
433 fig = Figure(figsize=(3,3))
434 ax = fig.add_axes((0,0,1,1),xticks=[],yticks=[],frame_on=False)
435 xd = ax.get_xlim()
436 xd = xd[1]-xd[0]
437 yd = ax.get_ylim()
438 yd = yd[1]-yd[0]
439 ax.size=(xd,yd)
440 drawer = MolDrawing(ax)
441 drawer.AddMol(mol)
442 canv = FigureCanvasAgg(fig)
443 canv.print_figure("foo.png",dpi=80)
444 else:
445 from sping.PDF.pidPDF import PDFCanvas as Canvas
446 canvas = Canvas(size=(300,300),name='test.pdf')
447 registerCanvas('sping')
448 drawer = MolDrawing(canvas)
449 drawer.AddMol(mol)
450 canvas.save()
451