library(maftools)
#path to TCGA LAML MAF file
laml.maf = system.file('extdata', 'tcga_laml.maf.gz', package = 'maftools')
#clinical information containing survival information and histology. This is optional
laml.clin = system.file('extdata', 'tcga_laml_annot.tsv', package = 'maftools')

laml = read.maf(maf = laml.maf,
                clinicalData = laml.clin,
                verbose = FALSE)

0.1 Including Transition/Transversions into oncoplot

oncoplot(maf = laml, draw_titv = TRUE)

0.2 Changing colors for variant classifications

#One can use any colors, here in this example color palette from RColorBrewer package is used
vc_cols = RColorBrewer::brewer.pal(n = 8, name = 'Paired')
names(vc_cols) = c(
  'Frame_Shift_Del',
  'Missense_Mutation',
  'Nonsense_Mutation',
  'Multi_Hit',
  'Frame_Shift_Ins',
  'In_Frame_Ins',
  'Splice_Site',
  'In_Frame_Del'
)

print(vc_cols)
#>   Frame_Shift_Del Missense_Mutation Nonsense_Mutation         Multi_Hit 
#>         "#A6CEE3"         "#1F78B4"         "#B2DF8A"         "#33A02C" 
#>   Frame_Shift_Ins      In_Frame_Ins       Splice_Site      In_Frame_Del 
#>         "#FB9A99"         "#E31A1C"         "#FDBF6F"         "#FF7F00"

oncoplot(maf = laml, colors = vc_cols, top = 10)

0.3 Including copy number data into oncoplots.

There are two ways one include CN status into MAF. 1. GISTIC results 2. Custom copy number table

0.3.1 GISTIC results

Most widely used tool for copy number analysis from large scale studies is GISTIC and we can simultaneously read gistic results along with MAF. GISTIC generates numerous files but we need mainly four files all_lesions.conf_XX.txt, amp_genes.conf_XX.txt, del_genes.conf_XX.txt, scores.gistic where XX is confidence level. These files contain significantly altered genomic regions along with amplified and deleted genes respectively.

#GISTIC results LAML
all.lesions =
  system.file("extdata", "all_lesions.conf_99.txt", package = "maftools")
amp.genes =
  system.file("extdata", "amp_genes.conf_99.txt", package = "maftools")
del.genes =
  system.file("extdata", "del_genes.conf_99.txt", package = "maftools")
scores.gis =
  system.file("extdata", "scores.gistic", package = "maftools")

#Read GISTIC results along with MAF
laml.plus.gistic = read.maf(
  maf = laml.maf,
  gisticAllLesionsFile = all.lesions,
  gisticAmpGenesFile = amp.genes,
  gisticDelGenesFile = del.genes,
  gisticScoresFile = scores.gis,
  isTCGA = TRUE,
  verbose = FALSE, 
  clinicalData = laml.clin
)
oncoplot(maf = laml.plus.gistic, top = 10)

This plot shows frequent deletions in TP53 gene which is located on one of the significantly deleted locus 17p13.2.

0.3.2 Custom copy-number table

In case there is no GISTIC results available, one can generate a table containing CN status for known genes in known samples. This can be easily created and read along with MAF file.

For example lets create a dummy CN alterations for DNMT3A in random 20 samples.

set.seed(seed = 1024)
barcodes = as.character(getSampleSummary(x = laml)[,Tumor_Sample_Barcode])
#Random 20 samples
dummy.samples = sample(x = barcodes,
                       size = 20,
                       replace = FALSE)

#Genarate random CN status for above samples
cn.status = sample(
  x = c('Amp', 'Del'),
  size = length(dummy.samples),
  replace = TRUE
)

custom.cn.data = data.frame(
  Gene = "DNMT3A",
  Sample_name = dummy.samples,
  CN = cn.status,
  stringsAsFactors = FALSE
)

head(custom.cn.data)
#>     Gene  Sample_name  CN
#> 1 DNMT3A TCGA-AB-2898 Amp
#> 2 DNMT3A TCGA-AB-2879 Amp
#> 3 DNMT3A TCGA-AB-2920 Del
#> 4 DNMT3A TCGA-AB-2866 Amp
#> 5 DNMT3A TCGA-AB-2892 Amp
#> 6 DNMT3A TCGA-AB-2863 Amp

laml.plus.cn = read.maf(maf = laml.maf,
                        cnTable = custom.cn.data,
                        verbose = FALSE)

oncoplot(maf = laml.plus.cn, top = 5)

0.4 Including significance values

This data should be a data.frame or a tsv file with two required columns titled gene and q.

For example, including mutsig q values into oncoplot.

#MutSig results
laml.mutsig = system.file("extdata", "LAML_sig_genes.txt.gz", package = "maftools")

oncoplot(
  maf = laml,
  mutsig = laml.mutsig,
  mutsigQval = 0.01,
)

0.5 Including expression (or any) values

Similar to significance values included as right bar plot, it’s also possible to include expression (or any sort of continuous) data to left side of the plot

#Dummy expression values for top 20 genes 
set.seed(seed = 1024)
exprs_tbl = data.frame(genes = getGeneSummary(x = laml)[1:20, Hugo_Symbol],
                       exprn = rnorm(n = 10, mean = 12, sd = 5))
head(exprs_tbl)
#>    genes     exprn
#> 1   FLT3  8.106686
#> 2 DNMT3A 10.052618
#> 3   NPM1  1.831008
#> 4   IDH2  7.088134
#> 5   IDH1 13.239450
#> 6   TET2  1.480677

oncoplot(maf = laml, exprsTbl = exprs_tbl)

0.6 Including annotations

Annotations are usually stored in clinical.data slot of MAF.

getClinicalData(x = laml)
#>      Tumor_Sample_Barcode FAB_classification days_to_last_followup
#>   1:         TCGA-AB-2802                 M4                   365
#>   2:         TCGA-AB-2803                 M3                   792
#>   3:         TCGA-AB-2804                 M3                  2557
#>   4:         TCGA-AB-2805                 M0                   577
#>   5:         TCGA-AB-2806                 M1                   945
#>  ---                                                              
#> 189:         TCGA-AB-3007                 M3                  1581
#> 190:         TCGA-AB-3008                 M1                   822
#> 191:         TCGA-AB-3009                 M4                   577
#> 192:         TCGA-AB-3011                 M1                  1885
#> 193:         TCGA-AB-3012                 M3                  1887
#>      Overall_Survival_Status
#>   1:                       1
#>   2:                       1
#>   3:                       0
#>   4:                       1
#>   5:                       1
#>  ---                        
#> 189:                       0
#> 190:                       1
#> 191:                       1
#> 192:                       0
#> 193:                       0

Include FAB_classification from clinical data as one of the sample annotations.

oncoplot(maf = laml, clinicalFeatures = 'FAB_classification')

More than one annotations can be included by passing them to the argument clinicalFeatures. Above plot can be further enhanced by sorting according to annotations. Custom colors can be specified as a list of named vectors for each levels.

#Color coding for FAB classification
fabcolors = RColorBrewer::brewer.pal(n = 8,name = 'Spectral')
names(fabcolors) = c("M0", "M1", "M2", "M3", "M4", "M5", "M6", "M7")
fabcolors = list(FAB_classification = fabcolors)

print(fabcolors)
#> $FAB_classification
#>        M0        M1        M2        M3        M4        M5        M6        M7 
#> "#D53E4F" "#F46D43" "#FDAE61" "#FEE08B" "#E6F598" "#ABDDA4" "#66C2A5" "#3288BD"

oncoplot(
  maf = laml,
  clinicalFeatures = 'FAB_classification',
  sortByAnnotation = TRUE,
  annotationColor = fabcolors
)

0.7 Highlighting samples

If you prefer to highlight mutations by a specific attribute, you can use additionalFeature argument.

Example: Highlight all mutations where alt allele is C.

oncoplot(maf = laml,
         additionalFeature = c("Tumor_Seq_Allele2", "C"))

Note that first argument (Tumor_Seq_Allele2) must a be column in MAF file, and second argument (C) is a value in that column. If you want to know what columns are present in the MAF file, use getFields.

getFields(x = laml)
#>  [1] "Hugo_Symbol"            "Entrez_Gene_Id"         "Center"                
#>  [4] "NCBI_Build"             "Chromosome"             "Start_Position"        
#>  [7] "End_Position"           "Strand"                 "Variant_Classification"
#> [10] "Variant_Type"           "Reference_Allele"       "Tumor_Seq_Allele1"     
#> [13] "Tumor_Seq_Allele2"      "Tumor_Sample_Barcode"   "Protein_Change"        
#> [16] "i_TumorVAF_WU"          "i_transcript_name"

0.8 Group by Pathways

0.8.1 Auto grouping

Genes can be auto grouped based on their Biological processess by setting pathways = 'auto'

oncoplot(maf = laml, top = 20, pathways = 'auto', gene_mar = 8, fontSize = 0.6)

0.8.2 Custom pathways

A custom pathway list can also be provided in the form of a two column tsv file or a data.frame containing gene names and their corresponding pathway.

my_path = data.table::data.table(Gene = c("TP53", "WT1", "DNMT3A", "IDH1", "IDH2", "RUNX1", "CEBPA"), Pathway = c("Tumor Suppressors", "DNA methylation", 'DNA methylation', 'DNA methylation', 'DNA methylation', 'Myeloid TF', 'Myeloid TF'))

print(my_path)
#>      Gene           Pathway
#> 1:   TP53 Tumor Suppressors
#> 2:    WT1   DNA methylation
#> 3: DNMT3A   DNA methylation
#> 4:   IDH1   DNA methylation
#> 5:   IDH2   DNA methylation
#> 6:  RUNX1        Myeloid TF
#> 7:  CEBPA        Myeloid TF

oncoplot(maf = laml, top = 20, pathways = my_path, gene_mar = 6, fontSize = 0.6)

0.9 Combining everything

oncoplot(
  maf = laml.plus.gistic,
  draw_titv = TRUE,
  clinicalFeatures = c('FAB_classification', 'Overall_Survival_Status'),
  additionalFeature = c("Tumor_Seq_Allele2", "C"),
  sortByAnnotation = TRUE,
  mutsig = laml.mutsig,
  exprsTbl = exprs_tbl,
  logColBar = TRUE
)

0.10 SessionInfo

sessionInfo()
#> R version 4.0.0 alpha (2020-03-31 r78116)
#> Platform: x86_64-w64-mingw32/x64 (64-bit)
#> Running under: Windows Server 2012 R2 x64 (build 9600)
#> 
#> Matrix products: default
#> 
#> locale:
#> [1] LC_COLLATE=C                          
#> [2] LC_CTYPE=English_United States.1252   
#> [3] LC_MONETARY=English_United States.1252
#> [4] LC_NUMERIC=C                          
#> [5] LC_TIME=English_United States.1252    
#> 
#> attached base packages:
#> [1] stats4    parallel  stats     graphics  grDevices utils     datasets 
#> [8] methods   base     
#> 
#> other attached packages:
#>  [1] pheatmap_1.0.12                   doParallel_1.0.15                
#>  [3] iterators_1.0.12                  foreach_1.5.0                    
#>  [5] NMF_0.22.0                        Biobase_2.47.3                   
#>  [7] cluster_2.1.0                     rngtools_1.5                     
#>  [9] pkgmaker_0.31.1                   registry_0.5-1                   
#> [11] BSgenome.Hsapiens.UCSC.hg19_1.4.3 BSgenome_1.55.4                  
#> [13] rtracklayer_1.47.0                Biostrings_2.55.7                
#> [15] XVector_0.27.2                    GenomicRanges_1.39.3             
#> [17] GenomeInfoDb_1.23.16              IRanges_2.21.8                   
#> [19] S4Vectors_0.25.15                 BiocGenerics_0.33.3              
#> [21] maftools_2.3.30                  
#> 
#> loaded via a namespace (and not attached):
#>  [1] splines_4.0.0               R.utils_2.9.2              
#>  [3] assertthat_0.2.1            GenomeInfoDbData_1.2.2     
#>  [5] Rsamtools_2.3.7             yaml_2.2.1                 
#>  [7] pillar_1.4.3                lattice_0.20-41            
#>  [9] glue_1.4.0                  digest_0.6.25              
#> [11] RColorBrewer_1.1-2          colorspace_1.4-1           
#> [13] plyr_1.8.6                  htmltools_0.4.0            
#> [15] Matrix_1.2-18               R.oo_1.23.0                
#> [17] pkgconfig_2.0.3             XML_3.99-0.3               
#> [19] bibtex_0.4.2.2              zlibbioc_1.33.1            
#> [21] purrr_0.3.3                 xtable_1.8-4               
#> [23] scales_1.1.0                BiocParallel_1.21.2        
#> [25] tibble_3.0.0                ggplot2_3.3.0              
#> [27] ellipsis_0.3.0              withr_2.1.2                
#> [29] SummarizedExperiment_1.17.5 cli_2.0.2                  
#> [31] survival_3.1-12             magrittr_1.5               
#> [33] crayon_1.3.4                mclust_5.4.6               
#> [35] evaluate_0.14               R.methodsS3_1.8.0          
#> [37] fansi_0.4.1                 tools_4.0.0                
#> [39] data.table_1.12.8           gridBase_0.4-7             
#> [41] lifecycle_0.2.0             matrixStats_0.56.0         
#> [43] stringr_1.4.0               munsell_0.5.0              
#> [45] DelayedArray_0.13.12        compiler_4.0.0             
#> [47] berryFunctions_1.18.2       rlang_0.4.5                
#> [49] grid_4.0.0                  RCurl_1.98-1.1             
#> [51] bitops_1.0-6                rmarkdown_2.1              
#> [53] gtable_0.3.0                codetools_0.2-16           
#> [55] abind_1.4-5                 reshape2_1.4.4             
#> [57] R6_2.4.1                    GenomicAlignments_1.23.2   
#> [59] dplyr_0.8.5                 knitr_1.28                 
#> [61] stringi_1.4.6               Rcpp_1.0.4.6               
#> [63] vctrs_0.2.4                 tidyselect_1.0.0           
#> [65] wordcloud_2.6               xfun_0.12