Contents

1 GenVisR

Intuitively visualizing and interpreting data from high-throughput genomic technologies continues to be challenging. “Genomic Visualizations in R” (GenVisR) attempts to alleviate this burden by providing highly customizable publication-quality graphics supporting multiple species and focused primarily on a cohort level (i.e., multiple samples/patients). GenVisR attempts to maintain a high degree of flexibility while leveraging the abilities of ggplot2 and bioconductor to achieve this goal.

Read the published Bioinformatics paper!

1.1 Install from Bioconductor

For the majority of users we recommend installing GenVisR from the release branch of Bioconductor, Installation instructions using this method can be found on the GenVisR landing page on Bioconductor.

Please note that GenVisR imports a few packages that have “system requirements”, in most cases these requirements will already be installed. If they are not please follow the instructions to install these packages given in the R terminal. Briefly these packages are: “libcurl4-openssl-dev” and “libxml2-dev”

1.2 Development

Development for GenVisR occurs on the griffith lab github repository available here. For users wishing to contribute to development we recommend cloning the GenVisR repo there and submitting a pull request. Please note that development occurs on the R version that will be available at each Bioconductor release cycle. This ensures that GenVisR will be stable for each Bioconductor release but it may necessitate developers download R-devel.

We also encourage users to report bugs and suggest enhancements to GenVisR on the github issue page available here:

1.3 Functions

1.3.1 waterfall (mutation overview graphic)

waterfall provides a method of visualizing the mutational landscape of a cohort. The input to waterfall consists of a data frame derived from either a .maf (version 2.4) file or a file in MGI annotation format (obtained from The Genome Modeling System) specified via the fileType parameter. waterfall will display the mutation occurrence and type in the main panel while showing the mutation burden and the percentage of samples with a mutation in the top and side sub-plots. Conflicts arising from multiple mutations in the same gene/sample cell are resolved by a hierarchical removal of mutations keeping the most deleterious as defined by the order of the “mutation type” legend. Briefly this hierarchy is as follows with the most deleterious defined first:

MAF MGI
Nonsense_Mutation nonsense
Frame_Shift_Ins frame_shift_del
Frame_Shift_Del frame_shift_ins
Translation_Start_Site splice_site_del
Splice_Site splice_site_ins
Nonstop_Mutation splice_site
In_Frame_Ins nonstop
In_Frame_Del in_frame_del
Missense_Mutation in_frame_ins
5’Flank missense
3’Flank splice_region_del
5’UTR splice_region_ins
3’UTR splice_region
RNA 5_prime_flanking_region
Intron 3_prime_flanking_region
IGR 3_prime_untranslated_region
Silent 5_prime_untranslated_region
Targeted_Region rna
intronic
silent

Occasionally a situation may arise in which it may be desireable to run waterfall on an unsupported file type. This can be achieved by setting the fileType parameter to “Custom”. Further the hieararchy of mutations (described above) must be specified with the variant_class_order parameter which expects a character vector describing the mutations observed in order of most to least important. Note that all mutations in the input data must be specified in the variant_class_order parameter. Using this option will require the data frame to contain the following column names: “sample”, “gene”, “variant_class”.

To view the general behavior of waterfall we use the brcaMAF data structure available within GenVisR. This data structure is a truncated MAF file consisting of 50 samples from the TCGA project corresponding to Breast invasive carcinoma (complete data from TCGA public web portal).

# Plot the mutation landscape
waterfall(brcaMAF, fileType="MAF")

This type of view is of limited use without expanding the graphic device given the large number of genes. Often it is beneficial to reduce the number of cells in the plot by limiting the number of genes plotted. There are three ways to accomplish this, the mainRecurCutoff parameter accepts a numeric value between 0 and 1 and will remove genes from the data which do not have at least x proportion of samples mutated. For example if it were desireable to plot those genes with mutations in >= 6% of samples:

# Load GenVisR and set seed
library(GenVisR)
set.seed(383)

# Plot only genes with mutations in 6% or more of samples
waterfall(brcaMAF, mainRecurCutoff = 0.06)

Alternatively one can set a maximum number of genes to plot via the maxGenes parameter which will select the top x recurrently mutated genes. In addition specific genes of interest can be displayed using the plotGenes parameter. This parameter accepts a case insensitive character vector of genes present in the data and will subset the data on those genes. For example, if it was desirable to plot only the following genes “PIK3CA”, “TP53”, “USH2A”, “MLL3”, AND “BRCA1”:

# Plot only the specified genes
waterfall(brcaMAF, plotGenes = c("PIK3CA", "TP53", "USH2A", "MLL3", "BRCA1"))

It is important to note that the mutation burden sub plot does not change during these subsets, this is calculated directly from the input via the formula: \(mutations\ in\ sample/coverage\ space * 1000000\). The coverage space defaults to the size in base pairs of the “SeqCap EZ Human Exome Library v2.0”. This default can be changed via the parameter coverageSpace. This calculation is only meant to be a rough estimate as actual coverage space can vary from sample to sample, for a more accurate calculation the user has the option to supply an optional argument via the parameter mutBurden supplying the users own calculation of mutation burden for each sample. This should be a data frame with column names ‘sample’, ‘mut_burden’ taking the following form:

sample mut_burden
TCGA-A1-A0SO-01A-22D-A099-09 1.5572403530013
TCGA-A2-A0EU-01A-22W-A071-09 2.19577768355127
TCGA-A2-A0ER-01A-21W-A050-09 1.89335842847617
TCGA-A2-A0EN-01A-13D-A099-09 2.67976843443599
TCGA-A1-A0SI-01A-11D-A142-09 1.64223789887094
TCGA-A2-A0D0-01A-11W-A019-09 2.9426074728573
TCGA-A2-A0D0-01A-11W-A019-09 1.49832578136762
TCGA-A1-A0SI-01A-11D-A142-09 1.55903682620951
TCGA-A2-A0CT-01A-31W-A071-09 2.61283158874499
TCGA-A2-A04U-01A-11D-A10Y-09 1.49772855192887

In addition to specifying the mutation burden the user also has the ability to plot additional clinical data. The clinical data supplied should be a data frame in “long” format with column names “sample”, “variable”, “value”. It is recommended to use the melt function in the package reshape2 to coerce data into this format. Here we add clinical data to be plotted and specify a custom order and colours for these variables putting these values in two columns within the clinical plot legend:

# Create clinical data
subtype <- c("lumA", "lumB", "her2", "basal", "normal")
subtype <- sample(subtype, 50, replace = TRUE)
age <- c("20-30", "31-50", "51-60", "61+")
age <- sample(age, 50, replace = TRUE)
sample <- as.character(unique(brcaMAF$Tumor_Sample_Barcode))
clinical <- as.data.frame(cbind(sample, subtype, age))

# Melt the clinical data into 'long' format.
library(reshape2)
clinical <- melt(clinical, id.vars = c("sample"))

# Run waterfall
waterfall(brcaMAF, clinDat = clinical, clinVarCol = c(lumA = "blue4", lumB = "deepskyblue", 
    her2 = "hotpink2", basal = "firebrick2", normal = "green4", `20-30` = "#ddd1e7", 
    `31-50` = "#bba3d0", `51-60` = "#9975b9", `61+` = "#7647a2"), plotGenes = c("PIK3CA", 
    "TP53", "USH2A", "MLL3", "BRCA1"), clinLegCol = 2, clinVarOrder = c("lumA", 
    "lumB", "her2", "basal", "normal", "20-30", "31-50", "51-60", "61+"))

Occasionally there may be samples not represented within the .maf file (due to a lack of mutations). It may still be desirable to plot these samples. To accomplish this simply add the relevant samples into the appropriate column before loading the data and leave the rest of the columns as NA. Alternatively the user can specify a list of samples to plot via the plotSamples parameter which will accept samples not in the input data.

1.3.2 lolliplot (mutation hotspot graphic)

lolliplot provides a method for visualizing mutation hotspots overlayed on a protein framework. The (basic) input consists of a data frame with required columns “transcript_name”, “gene” and “amino_acid_change” giving the ensembl transcript id, gene name, and the amino acid change in p. notation respectively. The data frame input to lolliplot must contain only one unique transcript name. lolliplot uses the R package biomaRt to obtain sequence and protein domain information and as such needs an active internet connection. lolliplot assumes the species from which to build a mart is hsapiens, this assumption can be changed via the parameter species which which will expect a valid ensembl mart species. Further by default the latest ensembl annotations will be used to build the protein framework. If it is desireable to build the protein framework from an older ensembl annotation the user has the option to change the annotation version by changing the host parameter. For example if the user wanted to use the Dec. 2013 ensembl annotation they would specify host=“dec2013.archive.ensembl.org”, the host parameter in lolliplot will pass it’s value to biomaRt::useMart, see biomaRt doc for using older archived annotations via biomaRt.

It should be noted that to ensure the most accurate graphic representation the ensembl annotation version for the “amino_acid_change” column supplied by the user and that which is used for the biomaRt queries should be identical.

# Create input data
data <- brcaMAF[brcaMAF$Hugo_Symbol == "TP53", c("Hugo_Symbol", "amino_acid_change_WU")]
data <- as.data.frame(cbind(data, "ENST00000269305"))
colnames(data) <- c("gene", "amino_acid_change", "transcript_name")

# Call lolliplot
lolliplot(data)

In an effort to maintain a high degree of flexibility the user has the option of selecting columns on which to fill and label. The parameters fillCol and labelCol allow this behavior by taking column names on which to fill and label respectively. Additionally one can plot the amino acid sidechain information in lieu of protein domains.

# Add additional columns to the data
data$gender <- sample(c("Male", "Female"), 15, replace = TRUE)
data$impact <- sample(c("Low", "Medium", "High"), 15, replace = TRUE)

# Call lolliplot
lolliplot(data, fillCol = "gender", labelCol = "impact", sideChain = TRUE)

The user has the option of plotting an additional track in the area underneath the protein track via the parameter y. Input for this additional layer consists of a data frame with column names “transcript_name” and “amino_acid_change” in p. notation. If input to parameter y is supplied to lolliplot and the fillCol and/or labelCol parameters are specified (see above) lolliplot will look for the columns in both data frames supplied to x and y and act accordingly. Note that input to parameter y must be from the same transcript as specified in the data frame supplied to parameter x.

# Create additional data
data2 <- data.frame(transcript_name = "ENST00000269305", amino_acid_change = "p.Q331*")

# Call lolliplot
lolliplot(data, y = data2, fillCol = "impact", labelCol = "amino_acid_change")

lolliplot uses a force field model from the package FField to repulse and attract data in an attempt to achieve a reasonable degree of separation between points. Suitable defaults have been set for the majority of use cases. On occasion the user may need to manually adjust the force field parameters especially if the number of points to apply the model to is large. This can be done for both upper and lower tracks individually via rep.fact, rep.dist.lmt, attr.fact, adj.max, adj.lmt, iter.max please see documentation for FField::FFieldPtRep for a complete description of these parameters.

1.3.3 genCov (sequence coverage graphic)

genCov provides a methodology for viewing coverage information in relation to a gene track. It takes a named list of data frames with each data frame containing column names “end” and “cov” and rows corresponding to coordinates within the region of interest. Additional required arguments are a GRanges object specifying the region of interest, a BSgenome for gc content calculation, and a TxDb object containing transcription metadata (see the package Granges for more information). genCov will plot a genomic features track and align coverage data in the list to the plot. It is recommended to use bedtools multicov to obtain coverage information for a region of interest. We demonstrate genCov functionality using pseudo-data containing coverage information for the gene PTEN.

# Load transcript meta data
library(TxDb.Hsapiens.UCSC.hg19.knownGene)
txdb <- TxDb.Hsapiens.UCSC.hg19.knownGene

# Load BSgenome
library(BSgenome.Hsapiens.UCSC.hg19)
genome <- BSgenome.Hsapiens.UCSC.hg19

# Define a region of interest
gr <- GRanges(seqnames = c("chr10"), ranges = IRanges(start = c(89622195), end = c(89729532)), 
    strand = strand(c("+")))

# Create Data for input
start <- c(89622194:89729524)
end <- c(89622195:89729525)
chr <- 10
cov <- c(rnorm(1e+05, mean = 40), rnorm(7331, mean = 10))
cov_input_A <- as.data.frame(cbind(chr, start, end, cov))

start <- c(89622194:89729524)
end <- c(89622195:89729525)
chr <- 10
cov <- c(rnorm(50000, mean = 40), rnorm(7331, mean = 10), rnorm(50000, mean = 40))
cov_input_B <- as.data.frame(cbind(chr, start, end, cov))

# Define the data as a list
data <- list(`Sample A` = cov_input_A, `Sample B` = cov_input_B)

# Call genCov
genCov(data, txdb, gr, genome, gene_labelTranscriptSize = 2, transform = NULL, 
    base = NULL)

Often it may be usefull to compress genomic space, genCov will perform such a compression via a log transform for each feature type,‘Intron’,‘CDS’,‘UTR’ specified by the parameter transform. The degree of compression can be set via the parameter base which will perform the appropriate log compression for the features specified in transform. This behavior will occur by default, to turn off compression set the transform and base parameters to NULL. Here we display genCov compression functionality with log-10 compression for intronic space, and log-2 compression for CDS and UTR regions. Further we choose to display a simplified representation of genomic features within the region of interest via the reduce parameter which will merge all genomic features within a region of interest into a single transcript.

# Turn off feature compression and reduce gene transcripts
genCov(data, txdb, gr, genome, transform = c("Intron", "CDS", "UTR"), base = c(10, 
    2, 2), reduce = TRUE)

1.3.4 TvTi (transition/transversion graphic)

TvTi provides a framework for visualizing transversions and transitions for a given cohort. Input consists of a .maf (version 2.4) file containing sample and allele information (see .maf spec). Alternatively the fileType parameter can be set to “MGI” with the input supplied consisting of a data frame with column names “sample”, “reference”, and “variant”. Files for the “MGI” format can be obtained via the Genome Modeling System. TvTi will remove indels and multinucleotide calls from the input and plot the proportion of Transition/Transversion types for each sample specified in the input file.

# Call TvTi
TvTi(brcaMAF, lab_txtAngle=75, fileType="MAF")

TvTi will also plot the observed frequency of each Transition/Transversion type in lieu of proportion if the type parameter is set to “Frequency”. Here we plot the observed frequency from brcaMAF and change the default colors of the plot. When modifying the color palette via the palette parameter specify a character vector of length 6 containing a new color for each Transition/Transversion type.

# Plot the frequency with a different color pallete
TvTi(brcaMAF, type = "Frequency", palette = c("#77C55D", "#A461B4", "#C1524B", 
    "#93B5BB", "#4F433F", "#BFA753"), lab_txtAngle = 75, fileType = "MAF")

If there are prior expectations about the transition/transversion rate the user can specify that information via the parameter y which takes a named vector with names corresponding to each transition/transversion type. The vector must be of length 6 with names “A->C or T->G (TV)”, “A->G or T->C (TI)”, “A->T or T->A (TV)”, “G->A or C->T (TI)”, “G->C or C->G (TV)”, and “G->T or C->A (TV)”. The Resulting plot will contain an additional subplot corresponding to the apriori expectations.

# Create a named vector of apriori expectations
expec <- c(`A->C or T->G (TV)` = 0.066, `A->G or T->C (TI)` = 0.217, `A->T or T->A (TV)` = 0.065, 
    `G->A or C->T (TI)` = 0.4945, `G->C or C->G (TV)` = 0.0645, `G->T or C->A (TV)` = 0.093)

# Call TvTi with the additional data
TvTi(brcaMAF, y = expec, lab_txtAngle = 45, fileType = "MAF")

1.3.5 cnSpec (copy altered cohort graphic)

cnSpec produces a plot displaying copy number segments at a cohort level. Basic input consists of a data frame with column names ‘chromosome’, ‘start’, ‘end’ ‘segmean’ and ‘sample’ with rows denoting segments with copy number alterations. A UCSC genome is also required (defaults to ‘hg19’) to determine chromosomal boundaries. cnSpec will produce a grid faceted on chromosome and sample displaying all CN segment calls in the input. Here we use the attached data set LucCNseg containing copy number segment calls for 4 samples from whole genome sequencing data.

# Call cnSpec with minimum required inputs
cnSpec(LucCNseg, genome = "hg19")
## genome specified is preloaded, retrieving data...

By default a few select genomes are included as part of GenVisR, these are “hg38”, “hg19”, “mm10”, “mm9”, “rn5”. If input into genome is not one of the previously mentioned genomes cnSpec will attempt to query the UCSC sql database to obtain chromosomal boundary information. This has been built in as a convenience, if internet connectivity is an issue, or if copy number segment calls are derived from an assembly not supported by UCSC the user can specify chromosomal boundaries via the argument y. This should take the form of a data frame with column names “chromosome”, “start”, “end” with rows providing positions for each chromosome. An example of this is provided in the included data set hg19chr.

# Call cnSpec with the y parameter
cnSpec(LucCNseg, y = hg19chr)

1.3.6 cnView (copy altered single sample graphic)

cnView provides a method for visualizing raw copy number calls focused on either a single chromosome or all chromosomes. Unlike the majority of plots within GenVisR cnView is intended to be used for a single sample. Input consists of a data frame with column names “chromosome”, “coordinate”, “cn”, and “p_value” (optional) as well as a specification of which chromosome to plot specified via the parameter chr and which genome assembly should be used for chromosome boundaries genome. The algorithm will produce an ideogram on the top track and plot copy number calls beneath. If a “p_value” column is present in the input data cnView will create a transparency value for all calls/observations based on that column with less significant calls having a higher transparency. Eliminating the “p_value” column will terminate this behavior. Here we demonstrate cnView pseudo-data for chromosome 14.

# Create data
chromosome <- "chr14"
coordinate <- sort(sample(0:106455000, size = 2000, replace = FALSE))
cn <- c(rnorm(300, mean = 3, sd = 0.2), rnorm(700, mean = 2, sd = 0.2), rnorm(1000, 
    mean = 3, sd = 0.2))
data <- as.data.frame(cbind(chromosome, coordinate, cn))

# Call cnView with basic input
cnView(data, chr = "chr14", genome = "hg19", ideogram_txtSize = 4)