scripts/last-dotplot
author Martin C. Frith
Tue Feb 28 16:35:02 2017 +0900 (2017-02-28)
changeset 834 1514199bb31b
parent 756 a0074597cb5d
child 835 cdb1e40c28e4
permissions -rwxr-xr-x
Bugfix: dotplot left border misplaced (for some versions of python/PIL).
     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 fileinput, fnmatch, itertools, optparse, os, re, sys
    13 
    14 # Try to make PIL/PILLOW work:
    15 try: from PIL import Image, ImageDraw, ImageFont, ImageColor
    16 except ImportError: import Image, ImageDraw, ImageFont, ImageColor
    17 
    18 def warn(message):
    19     prog = os.path.basename(sys.argv[0])
    20     sys.stderr.write(prog + ": " + message + "\n")
    21 
    22 def tabBlocks(beg1, beg2, blocks):
    23     '''Get the gapless blocks of an alignment, from LAST tabular format.'''
    24     for i in blocks.split(","):
    25         if ":" in i:
    26             x, y = i.split(":")
    27             beg1 += int(x)
    28             beg2 += int(y)
    29         else:
    30             size = int(i)
    31             yield beg1, beg2, size
    32             beg1 += size
    33             beg2 += size
    34 
    35 def mafBlocks(beg1, beg2, seq1, seq2):
    36     '''Get the gapless blocks of an alignment, from MAF format.'''
    37     size = 0
    38     for x, y in itertools.izip(seq1, seq2):
    39         if x == "-":
    40             if size:
    41                 yield beg1, beg2, size
    42                 beg1 += size
    43                 beg2 += size
    44                 size = 0
    45             beg2 += 1
    46         elif y == "-":
    47             if size:
    48                 yield beg1, beg2, size
    49                 beg1 += size
    50                 beg2 += size
    51                 size = 0
    52             beg1 += 1
    53         else:
    54             size += 1
    55     if size: yield beg1, beg2, size
    56 
    57 def alignmentInput(lines):
    58     '''Get alignments and sequence lengths, from MAF or tabular format.'''
    59     mafCount = 0
    60     for line in lines:
    61         w = line.split()
    62         if line[0].isdigit():  # tabular format
    63             chr1, beg1, seqlen1 = w[1], int(w[2]), int(w[5])
    64             if w[4] == "-": beg1 -= seqlen1
    65             chr2, beg2, seqlen2 = w[6], int(w[7]), int(w[10])
    66             if w[9] == "-": beg2 -= seqlen2
    67             blocks = tabBlocks(beg1, beg2, w[11])
    68             yield chr1, seqlen1, chr2, seqlen2, blocks
    69         elif line[0] == "s":  # MAF format
    70             if mafCount == 0:
    71                 chr1, beg1, seqlen1, seq1 = w[1], int(w[2]), int(w[5]), w[6]
    72                 if w[4] == "-": beg1 -= seqlen1
    73                 mafCount = 1
    74             else:
    75                 chr2, beg2, seqlen2, seq2 = w[1], int(w[2]), int(w[5]), w[6]
    76                 if w[4] == "-": beg2 -= seqlen2
    77                 blocks = mafBlocks(beg1, beg2, seq1, seq2)
    78                 yield chr1, seqlen1, chr2, seqlen2, blocks
    79                 mafCount = 0
    80 
    81 def isWantedSequenceName(name, patterns):
    82     if not patterns: return True
    83     base = name.split(".")[-1]  # allow for names like hg19.chr7
    84     for i in patterns:
    85         if fnmatch.fnmatchcase(name, i): return True
    86         if fnmatch.fnmatchcase(base, i): return True
    87     return False
    88 
    89 def readAlignments(fileName, opts):
    90     '''Get alignments and sequence lengths, from MAF or tabular format.'''
    91     alignments = []
    92     seqLengths1 = {}
    93     seqLengths2 = {}
    94     lines = fileinput.input(fileName)
    95     for chr1, seqlen1, chr2, seqlen2, blocks in alignmentInput(lines):
    96         if not isWantedSequenceName(chr1, opts.seq1): continue
    97         if not isWantedSequenceName(chr2, opts.seq2): continue
    98         aln = chr1, chr2, blocks
    99         alignments.append(aln)
   100         seqLengths1[chr1] = seqlen1
   101         seqLengths2[chr2] = seqlen2
   102     return alignments, seqLengths1, seqLengths2
   103 
   104 def natural_sort_key(my_string):
   105     '''Return a sort key for "natural" ordering, e.g. chr9 < chr10.'''
   106     parts = re.split(r'(\d+)', my_string)
   107     parts[1::2] = map(int, parts[1::2])
   108     return parts
   109 
   110 def get_text_sizes(my_strings, font, fontsize, image_mode):
   111     '''Get widths & heights, in pixels, of some strings.'''
   112     if fontsize == 0: return [(0, 0) for i in my_strings]
   113     image_size = 1, 1
   114     im = Image.new(image_mode, image_size)
   115     draw = ImageDraw.Draw(im)
   116     return [draw.textsize(i, font=font) for i in my_strings]
   117 
   118 def get_seq_info(seq_size_dic, font, fontsize, image_mode):
   119     '''Return miscellaneous information about the sequences.'''
   120     seq_names = seq_size_dic.keys()
   121     seq_names.sort(key=natural_sort_key)
   122     seq_sizes = [seq_size_dic[i] for i in seq_names]
   123     name_sizes = get_text_sizes(seq_names, font, fontsize, image_mode)
   124     margin = max(zip(*name_sizes)[1])  # maximum text height
   125     return seq_names, seq_sizes, name_sizes, margin
   126 
   127 def div_ceil(x, y):
   128     '''Return x / y rounded up.'''
   129     q, r = divmod(x, y)
   130     return q + (r != 0)
   131 
   132 def tot_seq_pix(seq_sizes, bp_per_pix):
   133     '''Return the total pixels needed for sequences of the given sizes.'''
   134     return sum([div_ceil(i, bp_per_pix) for i in seq_sizes])
   135 
   136 def get_bp_per_pix(seq_sizes, pix_tween_seqs, pix_limit):
   137     '''Get the minimum bp-per-pixel that fits in the size limit.'''
   138     seq_num = len(seq_sizes)
   139     seq_pix_limit = pix_limit - pix_tween_seqs * (seq_num - 1)
   140     if seq_pix_limit < seq_num:
   141         raise Exception("can't fit the image: too many sequences?")
   142     lower_bound = div_ceil(sum(seq_sizes), seq_pix_limit)
   143     for bp_per_pix in itertools.count(lower_bound):  # slow linear search
   144         if tot_seq_pix(seq_sizes, bp_per_pix) <= seq_pix_limit: break
   145     return bp_per_pix
   146 
   147 def get_seq_starts(seq_pix, pix_tween_seqs, margin):
   148     '''Get the start pixel for each sequence.'''
   149     seq_starts = []
   150     pix_tot = margin - pix_tween_seqs
   151     for i in seq_pix:
   152         pix_tot += pix_tween_seqs
   153         seq_starts.append(pix_tot)
   154         pix_tot += i
   155     return seq_starts
   156 
   157 def get_pix_info(seq_sizes, bp_per_pix, pix_tween_seqs, margin):
   158     '''Return pixel information about the sequences.'''
   159     seq_pix = [div_ceil(i, bp_per_pix) for i in seq_sizes]
   160     seq_starts = get_seq_starts(seq_pix, pix_tween_seqs, margin)
   161     tot_pix = seq_starts[-1] + seq_pix[-1]
   162     return seq_pix, seq_starts, tot_pix
   163 
   164 def drawLineForward(hits, width, bp_per_pix, origin, beg1, beg2, size):
   165     while True:
   166         q1, r1 = divmod(beg1, bp_per_pix)
   167         q2, r2 = divmod(beg2, bp_per_pix)
   168         hits[origin + q2 * width + q1] |= 1
   169         next_pix = min(bp_per_pix - r1, bp_per_pix - r2)
   170         if next_pix >= size: break
   171         beg1 += next_pix
   172         beg2 += next_pix
   173         size -= next_pix
   174 
   175 def drawLineReverse(hits, width, bp_per_pix, origin, beg1, beg2, size):
   176     beg2 = -1 - beg2
   177     while True:
   178         q1, r1 = divmod(beg1, bp_per_pix)
   179         q2, r2 = divmod(beg2, bp_per_pix)
   180         hits[origin + q2 * width + q1] |= 2
   181         next_pix = min(bp_per_pix - r1, r2 + 1)
   182         if next_pix >= size: break
   183         beg1 += next_pix
   184         beg2 -= next_pix
   185         size -= next_pix
   186 
   187 def alignmentPixels(width, height, alignments, bp_per_pix,
   188                     seq_start_dic1, seq_start_dic2):
   189     hits = [0] * (width * height)  # the image data
   190     for seq1, seq2, blocks in alignments:
   191         seq_start1 = seq_start_dic1[seq1]
   192         seq_start2 = seq_start_dic2[seq2]
   193         origin = seq_start2 * width + seq_start1
   194         for beg1, beg2, size in blocks:
   195             if beg1 < 0:
   196                 beg1 = -(beg1 + size)
   197                 beg2 = -(beg2 + size)
   198             if beg2 >= 0:
   199                 drawLineForward(hits, width, bp_per_pix, origin,
   200                                 beg1, beg2, size)
   201             else:
   202                 drawLineReverse(hits, width, bp_per_pix, origin,
   203                                 beg1, beg2, size)
   204     return hits
   205 
   206 def expandedSeqDict(seqDict):
   207     '''Allow lookup by short sequence names, e.g. chr7 as well as hg19.chr7.'''
   208     newDict = {}
   209     for name, x in seqDict.items():
   210         base = name.split(".")[-1]
   211         newDict[name] = x
   212         newDict[base] = x
   213     return newDict
   214 
   215 def isExtraFirstGapField(fields):
   216     return fields[4].isdigit()
   217 
   218 def readGaps(fileName):
   219     '''Read locations of unsequenced gaps, from an agp or gap file.'''
   220     if not fileName: return
   221     for line in fileinput.input(fileName):
   222         w = line.split()
   223         if not w or w[0][0] == "#": continue
   224         if isExtraFirstGapField(w): w = w[1:]
   225         if w[4] not in "NU": continue
   226         seqName = w[0]
   227         end = int(w[2])
   228         beg = end - int(w[5])  # zero-based coordinate
   229         bridgedText = w[7]
   230         yield seqName, beg, end, bridgedText
   231 
   232 def drawUnsequencedGaps(im, gaps, start_dic, margin, limit, isTop, bridgedText,
   233                         bp_per_pix, color):
   234     '''Draw rectangles representing unsequenced gaps.'''
   235     for seqName, beg, end, b in gaps:
   236         if b != bridgedText: continue
   237         if seqName not in start_dic: continue
   238         origin = start_dic[seqName]
   239         b = div_ceil(beg, bp_per_pix)  # use fully-covered pixels only
   240         e = end // bp_per_pix
   241         if e <= b: continue
   242         if isTop: box = origin + b, margin, origin + e, limit
   243         else:     box = margin, origin + b, limit, origin + e
   244         im.paste(color, box)
   245 
   246 def make_label(text, text_size, range_start, range_size):
   247     '''Return an axis label with endpoint & sort-order information.'''
   248     text_width  = text_size[0]
   249     label_start = range_start + (range_size - text_width) // 2
   250     label_end   = label_start + text_width
   251     sort_key    = text_width - range_size
   252     return sort_key, label_start, label_end, text
   253 
   254 def get_nonoverlapping_labels(labels, label_space):
   255     '''Get a subset of non-overlapping axis labels, greedily.'''
   256     nonoverlapping_labels = []
   257     for i in labels:
   258         if True not in [i[1] < j[2] + label_space and j[1] < i[2] + label_space
   259                         for j in nonoverlapping_labels]:
   260             nonoverlapping_labels.append(i)
   261     return nonoverlapping_labels
   262 
   263 def get_axis_image(seq_names, name_sizes, seq_starts, seq_pix,
   264                    font, image_mode, opts):
   265     '''Make an image of axis labels.'''
   266     min_pos = seq_starts[0]
   267     max_pos = seq_starts[-1] + seq_pix[-1]
   268     height = max(zip(*name_sizes)[1])
   269     labels = [make_label(i, j, k, l) for i, j, k, l in
   270               zip(seq_names, name_sizes, seq_starts, seq_pix)]
   271     labels = [i for i in labels if i[1] >= min_pos and i[2] <= max_pos]
   272     labels.sort()
   273     labels = get_nonoverlapping_labels(labels, opts.label_space)
   274     image_size = max_pos, height
   275     im = Image.new(image_mode, image_size, opts.border_shade)
   276     draw = ImageDraw.Draw(im)
   277     for i in labels:
   278         position = i[1], 0
   279         draw.text(position, i[3], font=font, fill=opts.text_color)
   280     return im
   281 
   282 def lastDotplot(opts, args):
   283     if opts.fontfile:  font = ImageFont.truetype(opts.fontfile, opts.fontsize)
   284     else:              font = ImageFont.load_default()
   285 
   286     image_mode = 'RGB'
   287     forward_color = ImageColor.getcolor(opts.forwardcolor, image_mode)
   288     reverse_color = ImageColor.getcolor(opts.reversecolor, image_mode)
   289     zipped_colors = zip(forward_color, reverse_color)
   290     overlap_color = tuple([(i + j) // 2 for i, j in zipped_colors])
   291 
   292     warn("reading alignments...")
   293     alignments, seq_size_dic1, seq_size_dic2 = readAlignments(args[0], opts)
   294     warn("done")
   295 
   296     if not alignments: raise Exception("there are no alignments")
   297 
   298     seq_info1 = get_seq_info(seq_size_dic1, font, opts.fontsize, image_mode)
   299     seq_info2 = get_seq_info(seq_size_dic2, font, opts.fontsize, image_mode)
   300     seq_names1, seq_sizes1, name_sizes1, margin1 = seq_info1
   301     seq_names2, seq_sizes2, name_sizes2, margin2 = seq_info2
   302 
   303     warn("choosing bp per pixel...")
   304     pix_limit1 = opts.width  - margin1
   305     pix_limit2 = opts.height - margin2
   306     bp_per_pix1 = get_bp_per_pix(seq_sizes1, opts.pix_tween_seqs, pix_limit1)
   307     bp_per_pix2 = get_bp_per_pix(seq_sizes2, opts.pix_tween_seqs, pix_limit2)
   308     bp_per_pix = max(bp_per_pix1, bp_per_pix2)
   309     warn("bp per pixel = " + str(bp_per_pix))
   310 
   311     seq_pix1, seq_starts1, width  = get_pix_info(seq_sizes1, bp_per_pix,
   312                                                  opts.pix_tween_seqs, margin1)
   313     seq_pix2, seq_starts2, height = get_pix_info(seq_sizes2, bp_per_pix,
   314                                                  opts.pix_tween_seqs, margin2)
   315     seq_start_dic1 = dict(zip(seq_names1, seq_starts1))
   316     seq_start_dic2 = dict(zip(seq_names2, seq_starts2))
   317 
   318     warn("processing alignments...")
   319     hits = alignmentPixels(width, height, alignments, bp_per_pix,
   320                            seq_start_dic1, seq_start_dic2)
   321     warn("done")
   322 
   323     image_size = width, height
   324     im = Image.new(image_mode, image_size, opts.background_color)
   325 
   326     start_dic1 = expandedSeqDict(seq_start_dic1)
   327     start_dic2 = expandedSeqDict(seq_start_dic2)
   328     gaps1 = list(readGaps(opts.gap1))
   329     gaps2 = list(readGaps(opts.gap2))
   330     # draw bridged gaps first, then unbridged gaps on top:
   331     drawUnsequencedGaps(im, gaps1, start_dic1, margin2, height, True, "yes",
   332                         bp_per_pix, opts.bridged_color)
   333     drawUnsequencedGaps(im, gaps2, start_dic2, margin1, width, False, "yes",
   334                         bp_per_pix, opts.bridged_color)
   335     drawUnsequencedGaps(im, gaps1, start_dic1, margin2, height, True, "no",
   336                         bp_per_pix, opts.unbridged_color)
   337     drawUnsequencedGaps(im, gaps2, start_dic2, margin1, width, False, "no",
   338                         bp_per_pix, opts.unbridged_color)
   339 
   340     for i in range(height):
   341         for j in range(width):
   342             store_value = hits[i * width + j]
   343             xy = j, i
   344             if   store_value == 1: im.putpixel(xy, forward_color)
   345             elif store_value == 2: im.putpixel(xy, reverse_color)
   346             elif store_value == 3: im.putpixel(xy, overlap_color)
   347 
   348     if opts.fontsize != 0:
   349         axis1 = get_axis_image(seq_names1, name_sizes1, seq_starts1, seq_pix1,
   350                                font, image_mode, opts)
   351         axis2 = get_axis_image(seq_names2, name_sizes2, seq_starts2, seq_pix2,
   352                                font, image_mode, opts)
   353         axis2 = axis2.transpose(Image.ROTATE_270)  # !!! bug hotspot
   354         im.paste(axis1, (0, 0))
   355         im.paste(axis2, (0, 0))
   356 
   357     for i in seq_starts1[1:]:
   358         box = i - opts.pix_tween_seqs, margin2, i, height
   359         im.paste(opts.border_shade, box)
   360 
   361     for i in seq_starts2[1:]:
   362         box = margin1, i - opts.pix_tween_seqs, width, i
   363         im.paste(opts.border_shade, box)
   364 
   365     im.save(args[1])
   366 
   367 if __name__ == "__main__":
   368     usage = """%prog --help
   369    or: %prog [options] maf-or-tab-alignments dotplot.png
   370    or: %prog [options] maf-or-tab-alignments dotplot.gif
   371    or: ..."""
   372     description = "Draw a dotplot of pair-wise sequence alignments in MAF or tabular format."
   373     op = optparse.OptionParser(usage=usage, description=description)
   374     op.add_option("-1", "--seq1", metavar="PATTERN", action="append",
   375                   help="which sequences to show from the 1st genome")
   376     op.add_option("-2", "--seq2", metavar="PATTERN", action="append",
   377                   help="which sequences to show from the 2nd genome")
   378     # Replace "width" & "height" with a single "length" option?
   379     op.add_option("-x", "--width", type="int", default=1000,
   380                   help="maximum width in pixels (default: %default)")
   381     op.add_option("-y", "--height", type="int", default=1000,
   382                   help="maximum height in pixels (default: %default)")
   383     op.add_option("-f", "--fontfile", metavar="FILE",
   384                   help="TrueType or OpenType font file")
   385     op.add_option("-s", "--fontsize", metavar="SIZE", type="int", default=11,
   386                   help="TrueType or OpenType font size (default: %default)")
   387     op.add_option("-c", "--forwardcolor", metavar="COLOR", default="red",
   388                   help="color for forward alignments (default: %default)")
   389     op.add_option("-r", "--reversecolor", metavar="COLOR", default="blue",
   390                   help="color for reverse alignments (default: %default)")
   391     og = optparse.OptionGroup(op, "Unsequenced gap options")
   392     og.add_option("--gap1", metavar="FILE",
   393                   help="read genome1 unsequenced gaps from agp or gap file")
   394     og.add_option("--gap2", metavar="FILE",
   395                   help="read genome2 unsequenced gaps from agp or gap file")
   396     og.add_option("--bridged-color", metavar="COLOR", default="yellow",
   397                   help="color for bridged gaps (default: %default)")
   398     og.add_option("--unbridged-color", metavar="COLOR", default="pink",
   399                   help="color for unbridged gaps (default: %default)")
   400     op.add_option_group(og)
   401     (opts, args) = op.parse_args()
   402     if len(args) != 2: op.error("2 arguments needed")
   403 
   404     opts.text_color = "black"
   405     opts.background_color = "white"
   406     opts.pix_tween_seqs = 2  # number of border pixels between sequences
   407     opts.border_shade = 239, 239, 239  # the shade of grey for border pixels
   408     opts.label_space = 5     # minimum number of pixels between axis labels
   409 
   410     try: lastDotplot(opts, args)
   411     except KeyboardInterrupt: pass  # avoid silly error message
   412     except Exception, e:
   413         prog = os.path.basename(sys.argv[0])
   414         sys.exit(prog + ": error: " + str(e))