AnnotationSketch code examples

Using AnnotationSketch to draw annotations from a file

The following code examples (in C and Lua) illustrate how to produce an image from a given GFF3 file using AnnotationSketch. The result is shown in Fig. 7. In essence, these code examples implement something like a simple version of the gt sketch tool from GenomeTools without most command-line options.

[Example rendering from parsed GFF3 file]

Figure 7: Example rendering of a GFF3 file with default style.

C code

(See src/examples/sketch_parsed.c in the source distribution.)
#include "genometools.h"

static void handle_error(GtError *err)
{
  fprintf(stderr, "error: %s\n", gt_error_get(err));
  exit(EXIT_FAILURE);
}

int main(int argc, char *argv[])
{
  const char *style_file, *png_file, *gff3_file;
  char *seqid;
  GtStyle *style;
  GtFeatureIndex *feature_index;
  GtRange range;
  GtDiagram *diagram;
  GtLayout *layout;
  GtCanvas *canvas;
  GtUword height;
  GtError *err;

  if (argc != 4) {
    fprintf(stderr, "Usage: %s style_file PNG_file GFF3_file\n", argv[0]);
    return EXIT_FAILURE;
  }

  style_file = argv[1];
  png_file = argv[2];
  gff3_file = argv[3];

  /* initialize */
  gt_lib_init();

  /* create error object */
  err = gt_error_new();

  /* create style */
  if (!(style = gt_style_new(err)))
    handle_error(err);

  /* load style file */
  if (gt_style_load_file(style, style_file, err))
    handle_error(err);

  /* create feature index */
  feature_index = gt_feature_index_memory_new();

  /* add GFF3 file to index */
  if (gt_feature_index_add_gff3file(feature_index, gff3_file, err))
    handle_error(err);

  /* create diagram for first sequence ID in feature index */
  if (!(seqid = gt_feature_index_get_first_seqid(feature_index, err))) {
    if (gt_error_is_set(err))
      handle_error(err);
  }
  if (gt_feature_index_get_range_for_seqid(feature_index, &range, seqid, err))
    handle_error(err);
  diagram = gt_diagram_new(feature_index, seqid, &range, style, err);
  gt_free(seqid);
  if (gt_error_is_set(err))
    handle_error(err);

  /* create layout with given width, determine resulting image height */
  layout = gt_layout_new(diagram, 600, style, err);
  if (!layout)
    handle_error(err);
  if (gt_layout_get_height(layout, &height, err))
    handle_error(err);

  /* create PNG canvas */
  canvas = gt_canvas_cairo_file_new(style, GT_GRAPHICS_PNG, 600, height,
                                    NULL, err);
  if (!canvas)
    handle_error(err);

  /* sketch layout on canvas */
  if (gt_layout_sketch(layout, canvas, err))
    handle_error(err);

  /* write canvas to file */
  if (gt_canvas_cairo_file_to_file((GtCanvasCairoFile*) canvas, png_file, err))
    handle_error(err);

  /* free */
  gt_canvas_delete(canvas);
  gt_layout_delete(layout);
  gt_diagram_delete(diagram);
  gt_feature_index_delete(feature_index);
  gt_style_delete(style);
  gt_error_delete(err);
  /* perform static data cleanup */
  gt_lib_clean();
  return EXIT_SUCCESS;
}

Lua code

(See gtscripts/sketch_parsed.lua in the source distribution. This example can be run by the command line gt gtscripts/sketch_parsed.lua <style_file> <PNG_file> <GFF3_file>)
function usage()
  io.stderr:write(string.format("Usage: %s Style_file PNG_file GFF3_file\n", arg[0]))
  io.stderr:write("Create PNG representation of GFF3 annotation file.\n")
  os.exit(1)
end

if #arg == 3 then
  style_file = arg[1]
  png_file   = arg[2]
  gff3_file  = arg[3]
else
  usage()
end

-- load style file
dofile(style_file)

-- create feature index
feature_index = gt.feature_index_memory_new()

-- add GFF3 file to index
feature_index:add_gff3file(gff3_file)

-- create diagram for first sequence ID in feature index
seqid = feature_index:get_first_seqid()
range = feature_index:get_range_for_seqid(seqid)
diagram = gt.diagram_new(feature_index, seqid, range)

-- create layout
layout = gt.layout_new(diagram, 600)
height = layout:get_height()

-- create canvas
canvas = gt.canvas_cairo_file_new_png(600, height, nil)

-- sketch layout on canvas
layout:sketch(canvas)

-- write canvas to file
canvas:to_file(png_file)

Ruby code

(See gtruby/sketch_parsed.rb in the source distribution.)
require 'gtruby'

if ARGV.size != 3 then
  STDERR.puts "Usage: #{$0} style_file PNG_file GFF3_file"
  STDERR.puts "Create PNG representation of GFF3 annotation file."
  exit(1)
end

(stylefile, pngfile, gff3file) = ARGV

# load style file
style = GT::Style.new()
style.load_file(stylefile)

# create feature index
feature_index = GT::FeatureIndexMemory.new()

# add GFF3 file to index
feature_index.add_gff3file(gff3file)

# create diagram for first sequence ID in feature index
seqid = feature_index.get_first_seqid()
range = feature_index.get_range_for_seqid(seqid)
diagram = GT::Diagram.from_index(feature_index, seqid, range, style)

# create layout for given width
layout = GT::Layout.new(diagram, 800, style)

# create canvas with given width and computed height
canvas = GT::CanvasCairoFile.new(style, 800, layout.get_height, nil)

# sketch layout on canvas
layout.sketch(canvas)

# write canvas to file
canvas.to_file(pngfile)

Python code

(See gtpython/sketch_parsed.py in the source distribution.)
#!/usr/bin/python
# -*- coding: utf-8 -*-

from gt.annotationsketch import *
from gt.core.gtrange import Range
import sys

if __name__ == "__main__":
    if len(sys.argv) != 4:
        sys.stderr.write("Usage: " + (sys.argv)[0] +
                         " Style_file PNG_file GFF3_file\n")
        sys.stderr.write("Create PNG representation of GFF3 annotation file.")
        sys.exit(1)

    pngfile = (sys.argv)[2]

  # load style file

    style = Style()
    style.load_file((sys.argv)[1])

  # create feature index

    feature_index = FeatureIndexMemory()

  # add GFF3 file to index

    feature_index.add_gff3file((sys.argv)[3])

  # create diagram for first sequence ID in feature index

    seqid = feature_index.get_first_seqid()
    range = feature_index.get_range_for_seqid(seqid)
    diagram = Diagram.from_index(feature_index, seqid, range, style)

  # create layout

    layout = Layout(diagram, 600, style)
    height = layout.get_height()

  # create canvas

    canvas = CanvasCairoFile(style, 600, height)

  # sketch layout on canvas

    layout.sketch(canvas)

  # write canvas to file

    canvas.to_file(pngfile)

Using AnnotationSketch to draw user-generated annotations

The following C code example illustrates how to produce an image from annotation graphs created by user code. The result is shown in Fig. 8.

[Example rendering from user-generated annotations]

Figure 8: Example rendering of user-generated annotations with default style.

C code

(See src/examples/sketch_constructed.c in the source distribution.)
#include "genometools.h"

static GtArray* create_example_features(void)
{
  GtArray *features;
  GtGenomeNode *gene, *exon, *intron; /* features */
  GtStr *seqid; /* holds the sequence id the features refer to */

  /* construct the example features */
  features = gt_array_new(sizeof (GtGenomeNode*));
  seqid = gt_str_new_cstr("chromosome_21");

  /* construct a gene on the forward strand with two exons */
  gene = gt_feature_node_new(seqid, "gene", 100, 900, GT_STRAND_FORWARD);
  exon = gt_feature_node_new(seqid, "exon", 100, 200, GT_STRAND_FORWARD);
  gt_feature_node_add_child((GtFeatureNode*) gene, (GtFeatureNode*) exon);
  intron = gt_feature_node_new(seqid, "intron", 201, 799, GT_STRAND_FORWARD);
  gt_feature_node_add_child((GtFeatureNode*) gene, (GtFeatureNode*) intron);
  exon = gt_feature_node_new(seqid, "exon", 800, 900, GT_STRAND_FORWARD);
  gt_feature_node_add_child((GtFeatureNode*) gene, (GtFeatureNode*) exon);

  /* store forward gene in feature array */
  gt_array_add(features, gene);

  /* construct a single-exon gene on the reverse strand
     (within the intron of the forward strand gene) */
  gene = gt_feature_node_new(seqid, "gene", 400, 600, GT_STRAND_REVERSE);
  exon = gt_feature_node_new(seqid, "exon", 400, 600, GT_STRAND_REVERSE);
  gt_feature_node_add_child((GtFeatureNode*) gene, (GtFeatureNode*) exon);

  /* store reverse gene in feature array */
  gt_array_add(features, gene);

  /* free */
  gt_str_delete(seqid);

  return features;
}

static void handle_error(GtError *err)
{
  fprintf(stderr, "error writing canvas %s\n", gt_error_get(err));
  exit(EXIT_FAILURE);
}

static void draw_example_features(GtArray *features, const char *style_file,
                                  const char *output_file)
{
  GtRange range = { 1, 1000 }; /* the genomic range to draw */
  GtStyle *style;
  GtDiagram *diagram;
  GtLayout *layout;
  GtCanvas *canvas;
  GtUword height;
  GtError *err = gt_error_new();

  /* create style */
  if (!(style = gt_style_new(err)))
    handle_error(err);

  /* load style file */
  if (gt_style_load_file(style, style_file, err))
    handle_error(err);

  /* create diagram */
  diagram = gt_diagram_new_from_array(features, &range, style);

  /* create layout with given width, determine resulting image height */
  layout = gt_layout_new(diagram, 600, style, err);
  if (!layout)
    handle_error(err);
  if (gt_layout_get_height(layout, &height, err))
    handle_error(err);

  /* create PNG canvas */
  canvas = gt_canvas_cairo_file_new(style, GT_GRAPHICS_PNG, 600, height,
                                    NULL, err);
  if (!canvas)
    handle_error(err);

  /* sketch layout on canvas */
  if (gt_layout_sketch(layout, canvas, err))
    handle_error(err);

  /* write canvas to file */
  if (gt_canvas_cairo_file_to_file((GtCanvasCairoFile*) canvas, output_file,
                                   err)) {
    handle_error(err);
  }

  /* free */
  gt_canvas_delete(canvas);
  gt_layout_delete(layout);
  gt_diagram_delete(diagram);
  gt_style_delete(style);
  gt_error_delete(err);
}

static void delete_example_features(GtArray *features)
{
  GtUword i;
  for (i = 0; i < gt_array_size(features); i++)
    gt_genome_node_delete(*(GtGenomeNode**) gt_array_get(features, i));
  gt_array_delete(features);
}

int main(int argc, char *argv[])
{
  GtArray *features; /* stores the created example features */

  if (argc != 3) {
    fprintf(stderr, "Usage: %s style_file output_file\n", argv[0]);
    return EXIT_FAILURE;
  }

  gt_lib_init();

  features = create_example_features();

  draw_example_features(features, argv[1], argv[2]);

  delete_example_features(features);

  gt_lib_clean();
  return EXIT_SUCCESS;
}

Lua code

(See gtscripts/sketch_constructed.lua in the source distribution. This example can be run by the command line gt gtscripts/sketch_constructed.lua <style_file> <PNG_file>)
function usage()
  io.stderr:write(string.format("Usage: %s Style_file PNG_file\n", arg[0]))
  os.exit(1)
end

if #arg == 2 then
  style_file = arg[1]
  png_file   = arg[2]
else
  usage()
end

-- load style file
dofile(style_file)

-- construct the example features
seqid = "chromosome_21"
nodes = {}

-- construct a gene on the forward strand with two exons
gene   = gt.feature_node_new(seqid, "gene", 100, 900, "+")
exon   = gt.feature_node_new(seqid, "exon", 100, 200, "+")
gene:add_child(exon)
intron = gt.feature_node_new(seqid, "intron", 201, 799, "+")
gene:add_child(intron)
exon   = gt.feature_node_new(seqid, "exon", 800, 900, "+")
gene:add_child(exon)
nodes[1] = gene

-- construct a single-exon gene on the reverse strand
-- (within the intron of the forward strand gene)
reverse_gene = gt.feature_node_new(seqid, "gene", 400, 600, "-")
reverse_exon = gt.feature_node_new(seqid, "exon", 400, 600, "-")
reverse_gene:add_child(reverse_exon)
nodes[2] = reverse_gene

-- create diagram
diagram = gt.diagram_new_from_array(nodes, 1, 1000)
layout = gt.layout_new(diagram, 600)
height = layout:get_height()

-- create canvas
canvas = gt.canvas_cairo_file_new_png(600, height, nil)

-- sketch layout on canvas
layout:sketch(canvas)

-- write canvas to file
canvas:to_file(png_file)

Ruby code

(See gtruby/sketch_constructed.rb in the source distribution.)
require 'gtruby'

if ARGV.size != 2 then
  STDERR.puts "Usage: #{$0} style_file PNG_file"
  exit(1)
end

seqid = "chromosome_21"

# construct a gene on the forward strand with two exons
gene   = GT::FeatureNode.create(seqid, "gene", 100, 900, "+")
exon   = GT::FeatureNode.create(seqid, "exon", 100, 200, "+")
gene.add_child(exon)
intron = GT::FeatureNode.create(seqid, "intron", 201, 799, "+")
gene.add_child(intron)
exon   = GT::FeatureNode.create(seqid, "exon", 800, 900, "+")
gene.add_child(exon)

# construct a single-exon gene on the reverse strand
# (within the intron of the forward strand gene)
reverse_gene = GT::FeatureNode.create(seqid, "gene", 400, 600, "-")
reverse_exon = GT::FeatureNode.create(seqid, "exon", 400, 600, "-")
reverse_gene.add_child(reverse_exon)

pngfile = ARGV[1]

style = GT::Style.new()
style.load_file(ARGV[0])

rng = GT::Range.new(1, 1000)

diagram = GT::Diagram.from_array([gene, reverse_gene], rng, style)

layout = GT::Layout.new(diagram, 600, style)
canvas = GT::CanvasCairoFile.new(style, 600, layout.get_height, nil)
layout.sketch(canvas)

canvas.to_file(pngfile)

Python code

(See gtpython/sketch_constructed.py in the source distribution.)
#!/usr/bin/python
# -*- coding: utf-8 -*-

from gt.core import *
from gt.extended import *
from gt.annotationsketch import *
from gt.annotationsketch.custom_track import CustomTrack
from gt.core.gtrange import Range
import sys

if __name__ == "__main__":
    if len(sys.argv) != 3:
        sys.stderr.write("Usage: " + (sys.argv)[0] +
                         " style_file PNG_file\n")
        sys.exit(1)

    seqid = "chromosome_21"
    nodes = []

  # construct a gene on the forward strand with two exons

    gene = FeatureNode.create_new(seqid, "gene", 100, 900, "+")
    exon = FeatureNode.create_new(seqid, "exon", 100, 200, "+")
    gene.add_child(exon)
    intron = FeatureNode.create_new(seqid, "intron", 201, 799, "+")
    gene.add_child(intron)
    exon = FeatureNode.create_new(seqid, "exon", 800, 900, "+")
    gene.add_child(exon)

  # construct a single-exon gene on the reverse strand
  # (within the intron of the forward strand gene)

    reverse_gene = FeatureNode.create_new(seqid, "gene", 400, 600, "-")
    reverse_exon = FeatureNode.create_new(seqid, "exon", 400, 600, "-")
    reverse_gene.add_child(reverse_exon)

    pngfile = (sys.argv)[2]

    style = Style()
    style.load_file((sys.argv)[1])

    diagram = Diagram.from_array([gene, reverse_gene], Range(1, 1000),
                                 style)

    layout = Layout(diagram, 600, style)
    height = layout.get_height()
    canvas = CanvasCairoFile(style, 600, height)
    layout.sketch(canvas)

    canvas.to_file(pngfile)