scripts/last-dotplot
author Martin C. Frith
Thu Nov 02 20:04:51 2017 +0900 (2017-11-02)
changeset 897 f7bc7d099e2b
parent 896 2238466f36bd
child 898 f6a9c15287ea
permissions -rwxr-xr-x
Refactor last-dotplot
     1 #! /usr/bin/env python
     2 
     3 # Read pair-wise alignments in MAF or LAST tabular format: write an
     4 # "Oxford grid", a.k.a. dotplot.
     5 
     6 # TODO: Currently, pixels with zero aligned nt-pairs are white, and
     7 # pixels with one or more aligned nt-pairs are black.  This can look
     8 # too crowded for large genome alignments.  I tried shading each pixel
     9 # according to the number of aligned nt-pairs within it, but the
    10 # result is too faint.  How can this be done better?
    11 
    12 import gzip
    13 from fnmatch import fnmatchcase
    14 import itertools, optparse, os, re, sys
    15 
    16 # Try to make PIL/PILLOW work:
    17 try: from PIL import Image, ImageDraw, ImageFont, ImageColor
    18 except ImportError: import Image, ImageDraw, ImageFont, ImageColor
    19 
    20 def myOpen(fileName):  # faster than fileinput
    21     if fileName == "-":
    22         return sys.stdin
    23     if fileName.endswith(".gz"):
    24         return gzip.open(fileName)
    25     return open(fileName)
    26 
    27 def warn(message):
    28     if opts.verbose:
    29         prog = os.path.basename(sys.argv[0])
    30         sys.stderr.write(prog + ": " + message + "\n")
    31 
    32 def croppedBlocks(blocks, range1, range2):
    33     cropBeg1, cropEnd1 = range1
    34     cropBeg2, cropEnd2 = range2
    35     if blocks[0][0] < 0: cropBeg1, cropEnd1 = -cropEnd1, -cropBeg1
    36     if blocks[0][1] < 0: cropBeg2, cropEnd2 = -cropEnd2, -cropBeg2
    37     for beg1, beg2, size in blocks:
    38         b1 = max(cropBeg1, beg1)
    39         e1 = min(cropEnd1, beg1 + size)
    40         if b1 >= e1: continue
    41         offset = beg2 - beg1
    42         b2 = max(cropBeg2, b1 + offset)
    43         e2 = min(cropEnd2, e1 + offset)
    44         if b2 >= e2: continue
    45         yield b2 - offset, b2, e2 - b2
    46 
    47 def tabBlocks(beg1, beg2, blocks):
    48     '''Get the gapless blocks of an alignment, from LAST tabular format.'''
    49     for i in blocks.split(","):
    50         if ":" in i:
    51             x, y = i.split(":")
    52             beg1 += int(x)
    53             beg2 += int(y)
    54         else:
    55             size = int(i)
    56             yield beg1, beg2, size
    57             beg1 += size
    58             beg2 += size
    59 
    60 def mafBlocks(beg1, beg2, seq1, seq2):
    61     '''Get the gapless blocks of an alignment, from MAF format.'''
    62     size = 0
    63     for x, y in itertools.izip(seq1, seq2):
    64         if x == "-":
    65             if size:
    66                 yield beg1, beg2, size
    67                 beg1 += size
    68                 beg2 += size
    69                 size = 0
    70             beg2 += 1
    71         elif y == "-":
    72             if size:
    73                 yield beg1, beg2, size
    74                 beg1 += size
    75                 beg2 += size
    76                 size = 0
    77             beg1 += 1
    78         else:
    79             size += 1
    80     if size: yield beg1, beg2, size
    81 
    82 def alignmentInput(lines):
    83     '''Get alignments and sequence lengths, from MAF or tabular format.'''
    84     mafCount = 0
    85     for line in lines:
    86         w = line.split()
    87         if line[0].isdigit():  # tabular format
    88             chr1, beg1, seqlen1 = w[1], int(w[2]), int(w[5])
    89             if w[4] == "-": beg1 -= seqlen1
    90             chr2, beg2, seqlen2 = w[6], int(w[7]), int(w[10])
    91             if w[9] == "-": beg2 -= seqlen2
    92             blocks = tabBlocks(beg1, beg2, w[11])
    93             yield chr1, seqlen1, chr2, seqlen2, blocks
    94         elif line[0] == "s":  # MAF format
    95             if mafCount == 0:
    96                 chr1, beg1, seqlen1, seq1 = w[1], int(w[2]), int(w[5]), w[6]
    97                 if w[4] == "-": beg1 -= seqlen1
    98                 mafCount = 1
    99             else:
   100                 chr2, beg2, seqlen2, seq2 = w[1], int(w[2]), int(w[5]), w[6]
   101                 if w[4] == "-": beg2 -= seqlen2
   102                 blocks = mafBlocks(beg1, beg2, seq1, seq2)
   103                 yield chr1, seqlen1, chr2, seqlen2, blocks
   104                 mafCount = 0
   105 
   106 def seqRequestFromText(text):
   107     if ":" in text:
   108         pattern, interval = text.rsplit(":", 1)
   109         if "-" in interval:
   110             beg, end = interval.rsplit("-", 1)
   111             return pattern, int(beg), int(end)  # beg may be negative
   112     return text, 0, sys.maxsize
   113 
   114 def rangeFromSeqName(seqRequests, name, seqLen):
   115     if not seqRequests: return 0, seqLen
   116     base = name.split(".")[-1]  # allow for names like hg19.chr7
   117     for pat, beg, end in seqRequests:
   118         if fnmatchcase(name, pat) or fnmatchcase(base, pat):
   119             return max(beg, 0), min(end, seqLen)
   120     return None
   121 
   122 def updateSeqs(isTrim, seqNames, seqLimits, seqName, seqRange, blocks, index):
   123     if seqName not in seqLimits:
   124         seqNames.append(seqName)
   125     if isTrim:
   126         beg = blocks[0][index]
   127         end = blocks[-1][index] + blocks[-1][2]
   128         if beg < 0: beg, end = -end, -beg
   129         if seqName in seqLimits:
   130             b, e = seqLimits[seqName]
   131             seqLimits[seqName] = min(b, beg), max(e, end)
   132         else:
   133             seqLimits[seqName] = beg, end
   134     else:
   135         seqLimits[seqName] = seqRange
   136 
   137 def readAlignments(fileName, opts):
   138     '''Get alignments and sequence limits, from MAF or tabular format.'''
   139     seqRequests1 = map(seqRequestFromText, opts.seq1)
   140     seqRequests2 = map(seqRequestFromText, opts.seq2)
   141 
   142     alignments = []
   143     seqNames1 = []
   144     seqNames2 = []
   145     seqLimits1 = {}
   146     seqLimits2 = {}
   147     lines = myOpen(fileName)
   148     for seqName1, seqLen1, seqName2, seqLen2, blocks in alignmentInput(lines):
   149         range1 = rangeFromSeqName(seqRequests1, seqName1, seqLen1)
   150         if not range1: continue
   151         range2 = rangeFromSeqName(seqRequests2, seqName2, seqLen2)
   152         if not range2: continue
   153         b = list(croppedBlocks(list(blocks), range1, range2))
   154         if not b: continue
   155         aln = seqName1, seqName2, b
   156         alignments.append(aln)
   157         updateSeqs(opts.trim1, seqNames1, seqLimits1, seqName1, range1, b, 0)
   158         updateSeqs(opts.trim2, seqNames2, seqLimits2, seqName2, range2, b, 1)
   159     return alignments, seqNames1, seqNames2, seqLimits1, seqLimits2
   160 
   161 def natural_sort_key(my_string):
   162     '''Return a sort key for "natural" ordering, e.g. chr9 < chr10.'''
   163     parts = re.split(r'(\d+)', my_string)
   164     parts[1::2] = map(int, parts[1::2])
   165     return parts
   166 
   167 def textDimensions(imageDraw, font, textRot, text):
   168     x, y = imageDraw.textsize(text, font=font)
   169     return (y, x) if textRot else (x, y)
   170 
   171 def get_text_sizes(my_strings, font, fontsize, image_mode, textRot):
   172     '''Get widths & heights, in pixels, of some strings.'''
   173     if fontsize == 0: return [(0, 0) for i in my_strings]
   174     image_size = 1, 1
   175     im = Image.new(image_mode, image_size)
   176     draw = ImageDraw.Draw(im)
   177     return [textDimensions(draw, font, textRot, i) for i in my_strings]
   178 
   179 def sizeText(size):
   180     suffixes = "bp", "kb", "Mb", "Gb"
   181     for i, x in enumerate(suffixes):
   182         j = 10 ** (i * 3)
   183         if size < j * 10:
   184             return "%.2g" % (1.0 * size / j) + x
   185         if size < j * 1000 or i == len(suffixes) - 1:
   186             return "%.0f" % (1.0 * size / j) + x
   187 
   188 def seqNameAndSizeText(seqName, seqSize):
   189     return seqName + ": " + sizeText(seqSize)
   190 
   191 def getSeqInfo(sortOpt, seqNames, seqLimits,
   192                font, fontsize, image_mode, isShowSize, textRot):
   193     '''Return miscellaneous information about the sequences.'''
   194     if sortOpt == 1:
   195         seqNames.sort(key=natural_sort_key)
   196     seqSizes = [seqLimits[i][1] - seqLimits[i][0] for i in seqNames]
   197     for i in seqNames:
   198         r = seqLimits[i]
   199         out = i, str(r[0]), str(r[1])
   200         warn("\t".join(out))
   201     warn("")
   202     if sortOpt == 2:
   203         seqRecords = sorted(zip(seqSizes, seqNames), reverse=True)
   204         seqSizes = [i[0] for i in seqRecords]
   205         seqNames = [i[1] for i in seqRecords]
   206     if isShowSize:
   207         seqLabels = map(seqNameAndSizeText, seqNames, seqSizes)
   208     else:
   209         seqLabels = seqNames
   210     labelSizes = get_text_sizes(seqLabels, font, fontsize, image_mode, textRot)
   211     margin = max(i[1] for i in labelSizes)
   212     # xxx the margin may be too big, because some labels may get omitted
   213     return seqNames, seqSizes, seqLabels, labelSizes, margin
   214 
   215 def div_ceil(x, y):
   216     '''Return x / y rounded up.'''
   217     q, r = divmod(x, y)
   218     return q + (r != 0)
   219 
   220 def get_bp_per_pix(rangeSizes, pixTweenRanges, maxPixels):
   221     '''Get the minimum bp-per-pixel that fits in the size limit.'''
   222     numOfRanges = len(rangeSizes)
   223     maxPixelsInRanges = maxPixels - pixTweenRanges * (numOfRanges - 1)
   224     if maxPixelsInRanges < numOfRanges:
   225         raise Exception("can't fit the image: too many sequences?")
   226     negLimit = -maxPixelsInRanges
   227     negBpPerPix = sum(rangeSizes) // negLimit
   228     while True:
   229         if sum(i // negBpPerPix for i in rangeSizes) >= negLimit:
   230             return -negBpPerPix
   231         negBpPerPix -= 1
   232 
   233 def get_seq_starts(seq_pix, pixTweenRanges, margin):
   234     '''Get the start pixel for each sequence.'''
   235     seq_starts = []
   236     pix_tot = margin - pixTweenRanges
   237     for i in seq_pix:
   238         pix_tot += pixTweenRanges
   239         seq_starts.append(pix_tot)
   240         pix_tot += i
   241     return seq_starts
   242 
   243 def pixelData(rangeSizes, bp_per_pix, pixTweenRanges, margin):
   244     '''Return pixel information about the sequences.'''
   245     seq_pix = [div_ceil(i, bp_per_pix) for i in rangeSizes]
   246     seq_starts = get_seq_starts(seq_pix, pixTweenRanges, margin)
   247     tot_pix = seq_starts[-1] + seq_pix[-1]
   248     return seq_pix, seq_starts, tot_pix
   249 
   250 def drawLineForward(hits, width, bp_per_pix, beg1, beg2, size):
   251     while True:
   252         q1, r1 = divmod(beg1, bp_per_pix)
   253         q2, r2 = divmod(beg2, bp_per_pix)
   254         hits[q2 * width + q1] |= 1
   255         next_pix = min(bp_per_pix - r1, bp_per_pix - r2)
   256         if next_pix >= size: break
   257         beg1 += next_pix
   258         beg2 += next_pix
   259         size -= next_pix
   260 
   261 def drawLineReverse(hits, width, bp_per_pix, beg1, beg2, size):
   262     beg2 = -1 - beg2
   263     while True:
   264         q1, r1 = divmod(beg1, bp_per_pix)
   265         q2, r2 = divmod(beg2, bp_per_pix)
   266         hits[q2 * width + q1] |= 2
   267         next_pix = min(bp_per_pix - r1, r2 + 1)
   268         if next_pix >= size: break
   269         beg1 += next_pix
   270         beg2 -= next_pix
   271         size -= next_pix
   272 
   273 def alignmentPixels(width, height, alignments, bp_per_pix, origins1, origins2):
   274     hits = [0] * (width * height)  # the image data
   275     for seq1, seq2, blocks in alignments:
   276         ori1 = origins1[seq1]
   277         ori2 = origins2[seq2]
   278         for beg1, beg2, size in blocks:
   279             if beg1 < 0:
   280                 beg1 = -(beg1 + size)
   281                 beg2 = -(beg2 + size)
   282             if beg2 >= 0:
   283                 drawLineForward(hits, width, bp_per_pix,
   284                                 beg1 + ori1, beg2 + ori2, size)
   285             else:
   286                 drawLineReverse(hits, width, bp_per_pix,
   287                                 beg1 + ori1, beg2 - ori2, size)
   288     return hits
   289 
   290 def expandedSeqDict(seqDict):
   291     '''Allow lookup by short sequence names, e.g. chr7 as well as hg19.chr7.'''
   292     newDict = seqDict.copy()
   293     for name, x in seqDict.items():
   294         if "." in name:
   295             base = name.split(".")[-1]
   296             if base in newDict:  # an ambiguous case was found:
   297                 return seqDict   # so give up completely
   298             newDict[base] = x
   299     return newDict
   300 
   301 def readBed(fileName, seqLimits):
   302     if not fileName: return
   303     for line in myOpen(fileName):
   304         w = line.split()
   305         if not w: continue
   306         seqName = w[0]
   307         if seqName not in seqLimits: continue
   308         beg = int(w[1])
   309         end = int(w[2])
   310         layer = 900
   311         color = "#fbf"
   312         if len(w) > 4:
   313             if w[4] != ".":
   314                 layer = float(w[4])
   315             if len(w) > 5:
   316                 if len(w) > 8 and w[8].count(",") == 2:
   317                     color = "rgb(" + w[8] + ")"
   318                 elif w[5] == "+":
   319                     color = "#ffe8e8"
   320                 elif w[5] == "-":
   321                     color = "#e8e8ff"
   322         yield layer, color, seqName, beg, end
   323 
   324 def commaSeparatedInts(text):
   325     return map(int, text.rstrip(",").split(","))
   326 
   327 def readGenePred(opts, fileName, seqLimits):
   328     if not fileName: return
   329     for line in myOpen(fileName):
   330         fields = line.split()
   331         if not fields: continue
   332         if fields[2] not in "+-": fields = fields[1:]
   333         seqName = fields[1]
   334         if seqName not in seqLimits: continue
   335         #strand = fields[2]
   336         cdsBeg = int(fields[5])
   337         cdsEnd = int(fields[6])
   338         exonBegs = commaSeparatedInts(fields[8])
   339         exonEnds = commaSeparatedInts(fields[9])
   340         for beg, end in zip(exonBegs, exonEnds):
   341             yield 300, opts.exon_color, seqName, beg, end
   342             b = max(beg, cdsBeg)
   343             e = min(end, cdsEnd)
   344             if b < e: yield 400, opts.cds_color, seqName, b, e
   345 
   346 def readRmsk(fileName, seqLimits):
   347     if not fileName: return
   348     for line in myOpen(fileName):
   349         fields = line.split()
   350         if len(fields) == 17:  # rmsk.txt
   351             seqName = fields[5]
   352             if seqName not in seqLimits: continue  # do this ASAP for speed
   353             beg = int(fields[6])
   354             end = int(fields[7])
   355             strand = fields[9]
   356             repeatClass = fields[11]
   357         elif len(fields) == 15:  # .out
   358             seqName = fields[4]
   359             if seqName not in seqLimits: continue
   360             beg = int(fields[5]) - 1
   361             end = int(fields[6])
   362             strand = fields[8]
   363             repeatClass = fields[10]
   364         else:
   365             continue
   366         if repeatClass in ("Low_complexity", "Simple_repeat"):
   367             yield 200, "#fbf", seqName, beg, end
   368         elif strand == "+":
   369             yield 100, "#ffe8e8", seqName, beg, end
   370         else:
   371             yield 100, "#e8e8ff", seqName, beg, end
   372 
   373 def isExtraFirstGapField(fields):
   374     return fields[4].isdigit()
   375 
   376 def readGaps(opts, fileName, seqLimits):
   377     '''Read locations of unsequenced gaps, from an agp or gap file.'''
   378     if not fileName: return
   379     for line in myOpen(fileName):
   380         w = line.split()
   381         if not w or w[0][0] == "#": continue
   382         if isExtraFirstGapField(w): w = w[1:]
   383         if w[4] not in "NU": continue
   384         seqName = w[0]
   385         if seqName not in seqLimits: continue
   386         end = int(w[2])
   387         beg = end - int(w[5])  # zero-based coordinate
   388         if w[7] == "yes":
   389             yield 3000, opts.bridged_color, seqName, beg, end
   390         else:
   391             yield 2000, opts.unbridged_color, seqName, beg, end
   392 
   393 def bedBoxes(beds, seqLimits, origins, margin, edge, isTop, bpPerPix):
   394     for layer, color, seqName, bedBeg, bedEnd in beds:
   395         cropBeg, cropEnd = seqLimits[seqName]
   396         beg = max(bedBeg, cropBeg)
   397         end = min(bedEnd, cropEnd)
   398         if beg >= end: continue
   399         ori = origins[seqName]
   400         if layer <= 1000:
   401             # include partly-covered pixels
   402             b = (ori + beg) // bpPerPix
   403             e = div_ceil(ori + end, bpPerPix)
   404         else:
   405             # exclude partly-covered pixels
   406             b = div_ceil(ori + beg, bpPerPix)
   407             e = (ori + end) // bpPerPix
   408             if e <= b: continue
   409         if isTop:
   410             box = b, margin, e, edge
   411         else:
   412             box = margin, b, edge, e
   413         yield layer, color, box
   414 
   415 def drawAnnotations(im, boxes):
   416     # xxx use partial transparency for different-color overlaps?
   417     for layer, color, box in boxes:
   418         im.paste(color, box)
   419 
   420 def make_label(text, text_size, range_start, range_size):
   421     '''Return an axis label with endpoint & sort-order information.'''
   422     text_width, text_height = text_size
   423     label_start = range_start + (range_size - text_width) // 2
   424     label_end   = label_start + text_width
   425     sort_key    = text_width - range_size
   426     return sort_key, label_start, label_end, text, text_height
   427 
   428 def nonoverlappingLabels(labels, minPixTweenLabels):
   429     '''Get a subset of non-overlapping axis labels, greedily.'''
   430     out = []
   431     for i in labels:
   432         beg = i[1] - minPixTweenLabels
   433         end = i[2] + minPixTweenLabels
   434         if all(j[2] <= beg or j[1] >= end for j in out):
   435             out.append(i)
   436     return out
   437 
   438 def axisImage(seqNames, name_sizes, seq_starts, seq_pix, textRot,
   439               textAln, font, image_mode, opts):
   440     '''Make an image of axis labels.'''
   441     beg = seq_starts[0]
   442     end = seq_starts[-1] + seq_pix[-1]
   443     margin = max(i[1] for i in name_sizes)
   444     labels = map(make_label, seqNames, name_sizes, seq_starts, seq_pix)
   445     labels = [i for i in labels if i[1] >= beg and i[2] <= end]
   446     labels.sort()
   447     minPixTweenLabels = 0 if textRot else opts.label_space
   448     labels = nonoverlappingLabels(labels, minPixTweenLabels)
   449     image_size = (margin, end) if textRot else (end, margin)
   450     im = Image.new(image_mode, image_size, opts.margin_color)
   451     draw = ImageDraw.Draw(im)
   452     for i in labels:
   453         base = margin - i[4] if textAln else 0
   454         position = (base, i[1]) if textRot else (i[1], base)
   455         draw.text(position, i[3], font=font, fill=opts.text_color)
   456     return im
   457 
   458 def seqOrigins(seqNames, seq_starts, seqLimits, bp_per_pix):
   459     for i, j in zip(seqNames, seq_starts):
   460         yield i, bp_per_pix * j - seqLimits[i][0]
   461 
   462 def lastDotplot(opts, args):
   463     if opts.fontfile:  font = ImageFont.truetype(opts.fontfile, opts.fontsize)
   464     else:              font = ImageFont.load_default()
   465 
   466     image_mode = 'RGB'
   467     forward_color = ImageColor.getcolor(opts.forwardcolor, image_mode)
   468     reverse_color = ImageColor.getcolor(opts.reversecolor, image_mode)
   469     zipped_colors = zip(forward_color, reverse_color)
   470     overlap_color = tuple([(i + j) // 2 for i, j in zipped_colors])
   471 
   472     warn("reading alignments...")
   473     alignmentInfo = readAlignments(args[0], opts)
   474     alignments, seqNames1, seqNames2, seqLimits1, seqLimits2 = alignmentInfo
   475     warn("done")
   476     if not alignments: raise Exception("there are no alignments")
   477 
   478     textRot1 = "vertical".startswith(opts.rot1)
   479     i1 = getSeqInfo(opts.sort1, seqNames1, seqLimits1,
   480                     font, opts.fontsize, image_mode, opts.lengths1, textRot1)
   481     seqNames1, rangeSizes1, seqLabels1, labelSizes1, tMargin = i1
   482 
   483     textRot2 = "horizontal".startswith(opts.rot2)
   484     i2 = getSeqInfo(opts.sort2, seqNames2, seqLimits2,
   485                     font, opts.fontsize, image_mode, opts.lengths2, textRot2)
   486     seqNames2, rangeSizes2, seqLabels2, labelSizes2, lMargin = i2
   487 
   488     warn("choosing bp per pixel...")
   489     maxPixels1 = opts.width  - lMargin
   490     maxPixels2 = opts.height - tMargin
   491     bpPerPix1 = get_bp_per_pix(rangeSizes1, opts.border_pixels, maxPixels1)
   492     bpPerPix2 = get_bp_per_pix(rangeSizes2, opts.border_pixels, maxPixels2)
   493     bpPerPix = max(bpPerPix1, bpPerPix2)
   494     warn("bp per pixel = " + str(bpPerPix))
   495 
   496     seq_pix1, seq_starts1, width  = pixelData(rangeSizes1, bpPerPix,
   497                                               opts.border_pixels, lMargin)
   498     seq_pix2, seq_starts2, height = pixelData(rangeSizes2, bpPerPix,
   499                                               opts.border_pixels, tMargin)
   500     warn("width:  " + str(width))
   501     warn("height: " + str(height))
   502 
   503     origins1 = dict(seqOrigins(seqNames1, seq_starts1, seqLimits1, bpPerPix))
   504     origins2 = dict(seqOrigins(seqNames2, seq_starts2, seqLimits2, bpPerPix))
   505 
   506     warn("processing alignments...")
   507     hits = alignmentPixels(width, height, alignments, bpPerPix,
   508                            origins1, origins2)
   509     warn("done")
   510 
   511     image_size = width, height
   512     im = Image.new(image_mode, image_size, opts.background_color)
   513 
   514     seqLimits1 = expandedSeqDict(seqLimits1)
   515     seqLimits2 = expandedSeqDict(seqLimits2)
   516     origins1 = expandedSeqDict(origins1)
   517     origins2 = expandedSeqDict(origins2)
   518 
   519     beds1 = itertools.chain(readBed(opts.bed1, seqLimits1),
   520                             readRmsk(opts.rmsk1, seqLimits1),
   521                             readGenePred(opts, opts.genePred1, seqLimits1),
   522                             readGaps(opts, opts.gap1, seqLimits1))
   523     b1 = bedBoxes(beds1, seqLimits1, origins1, tMargin, height, True, bpPerPix)
   524 
   525     beds2 = itertools.chain(readBed(opts.bed2, seqLimits2),
   526                             readRmsk(opts.rmsk2, seqLimits2),
   527                             readGenePred(opts, opts.genePred2, seqLimits2),
   528                             readGaps(opts, opts.gap2, seqLimits2))
   529     b2 = bedBoxes(beds2, seqLimits2, origins2, lMargin, width, False, bpPerPix)
   530 
   531     boxes = sorted(itertools.chain(b1, b2))
   532     drawAnnotations(im, boxes)
   533 
   534     for i in range(height):
   535         for j in range(width):
   536             store_value = hits[i * width + j]
   537             xy = j, i
   538             if   store_value == 1: im.putpixel(xy, forward_color)
   539             elif store_value == 2: im.putpixel(xy, reverse_color)
   540             elif store_value == 3: im.putpixel(xy, overlap_color)
   541 
   542     if opts.fontsize != 0:
   543         axis1 = axisImage(seqLabels1, labelSizes1, seq_starts1, seq_pix1,
   544                           textRot1, False, font, image_mode, opts)
   545         if textRot1:
   546             axis1 = axis1.transpose(Image.ROTATE_90)
   547         axis2 = axisImage(seqLabels2, labelSizes2, seq_starts2, seq_pix2,
   548                           textRot2, textRot2, font, image_mode, opts)
   549         if not textRot2:
   550             axis2 = axis2.transpose(Image.ROTATE_270)
   551         im.paste(axis1, (0, 0))
   552         im.paste(axis2, (0, 0))
   553 
   554     for i in seq_starts1[1:]:
   555         box = i - opts.border_pixels, tMargin, i, height
   556         im.paste(opts.border_color, box)
   557 
   558     for i in seq_starts2[1:]:
   559         box = lMargin, i - opts.border_pixels, width, i
   560         im.paste(opts.border_color, box)
   561 
   562     im.save(args[1])
   563 
   564 if __name__ == "__main__":
   565     usage = """%prog --help
   566    or: %prog [options] maf-or-tab-alignments dotplot.png
   567    or: %prog [options] maf-or-tab-alignments dotplot.gif
   568    or: ..."""
   569     description = "Draw a dotplot of pair-wise sequence alignments in MAF or tabular format."
   570     op = optparse.OptionParser(usage=usage, description=description)
   571     op.add_option("-v", "--verbose", action="count",
   572                   help="show progress messages & data about the plot")
   573     op.add_option("-1", "--seq1", metavar="PATTERN", action="append",
   574                   default=[],
   575                   help="which sequences to show from the 1st genome")
   576     op.add_option("-2", "--seq2", metavar="PATTERN", action="append",
   577                   default=[],
   578                   help="which sequences to show from the 2nd genome")
   579     # Replace "width" & "height" with a single "length" option?
   580     op.add_option("-x", "--width", type="int", default=1000,
   581                   help="maximum width in pixels (default: %default)")
   582     op.add_option("-y", "--height", type="int", default=1000,
   583                   help="maximum height in pixels (default: %default)")
   584     op.add_option("-c", "--forwardcolor", metavar="COLOR", default="red",
   585                   help="color for forward alignments (default: %default)")
   586     op.add_option("-r", "--reversecolor", metavar="COLOR", default="blue",
   587                   help="color for reverse alignments (default: %default)")
   588     op.add_option("--sort1", type="int", default=1, metavar="N",
   589                   help="genome1 sequence order: 0=input order, 1=name order, "
   590                   "2=length order (default=%default)")
   591     op.add_option("--sort2", type="int", default=1, metavar="N",
   592                   help="genome2 sequence order: 0=input order, 1=name order, "
   593                   "2=length order (default=%default)")
   594     op.add_option("--trim1", action="store_true",
   595                   help="trim unaligned sequence flanks from the 1st genome")
   596     op.add_option("--trim2", action="store_true",
   597                   help="trim unaligned sequence flanks from the 2nd genome")
   598     op.add_option("--border-pixels", metavar="INT", type="int", default=1,
   599                   help="number of pixels between sequences (default=%default)")
   600     op.add_option("--border-color", metavar="COLOR", default="black",
   601                   help="color for pixels between sequences (default=%default)")
   602     op.add_option("--margin-color", metavar="COLOR", default="#dcdcdc",
   603                   help="margin color")
   604 
   605     og = optparse.OptionGroup(op, "Text options")
   606     og.add_option("-f", "--fontfile", metavar="FILE",
   607                   help="TrueType or OpenType font file")
   608     og.add_option("-s", "--fontsize", metavar="SIZE", type="int", default=11,
   609                   help="TrueType or OpenType font size (default: %default)")
   610     og.add_option("--rot1", metavar="ROT", default="h",
   611                   help="text rotation for the 1st genome (default=%default)")
   612     og.add_option("--rot2", metavar="ROT", default="v",
   613                   help="text rotation for the 2nd genome (default=%default)")
   614     og.add_option("--lengths1", action="store_true",
   615                   help="show sequence lengths for the 1st (horizontal) genome")
   616     og.add_option("--lengths2", action="store_true",
   617                   help="show sequence lengths for the 2nd (vertical) genome")
   618     op.add_option_group(og)
   619 
   620     og = optparse.OptionGroup(op, "Annotation options")
   621     og.add_option("--bed1", metavar="FILE",
   622                   help="read genome1 annotations from BED file")
   623     og.add_option("--bed2", metavar="FILE",
   624                   help="read genome2 annotations from BED file")
   625     og.add_option("--rmsk1", metavar="FILE", help="read genome1 repeats from "
   626                   "RepeatMasker .out or rmsk.txt file")
   627     og.add_option("--rmsk2", metavar="FILE", help="read genome2 repeats from "
   628                   "RepeatMasker .out or rmsk.txt file")
   629     op.add_option_group(og)
   630 
   631     og = optparse.OptionGroup(op, "Gene options")
   632     og.add_option("--genePred1", metavar="FILE",
   633                   help="read genome1 genes from genePred file")
   634     og.add_option("--genePred2", metavar="FILE",
   635                   help="read genome2 genes from genePred file")
   636     og.add_option("--exon-color", metavar="COLOR", default="PaleGreen",
   637                   help="color for exons (default=%default)")
   638     og.add_option("--cds-color", metavar="COLOR", default="LimeGreen",
   639                   help="color for protein-coding regions (default=%default)")
   640     op.add_option_group(og)
   641 
   642     og = optparse.OptionGroup(op, "Unsequenced gap options")
   643     og.add_option("--gap1", metavar="FILE",
   644                   help="read genome1 unsequenced gaps from agp or gap file")
   645     og.add_option("--gap2", metavar="FILE",
   646                   help="read genome2 unsequenced gaps from agp or gap file")
   647     og.add_option("--bridged-color", metavar="COLOR", default="yellow",
   648                   help="color for bridged gaps (default: %default)")
   649     og.add_option("--unbridged-color", metavar="COLOR", default="orange",
   650                   help="color for unbridged gaps (default: %default)")
   651     op.add_option_group(og)
   652     (opts, args) = op.parse_args()
   653     if len(args) != 2: op.error("2 arguments needed")
   654 
   655     opts.text_color = "black"
   656     opts.background_color = "white"
   657     opts.label_space = 5     # minimum number of pixels between axis labels
   658 
   659     try: lastDotplot(opts, args)
   660     except KeyboardInterrupt: pass  # avoid silly error message
   661     except Exception, e:
   662         prog = os.path.basename(sys.argv[0])
   663         sys.exit(prog + ": error: " + str(e))