plotGrandLinear {ggbio}R Documentation

Manhattan for GWAS

Description

A Manhattan plot is special scatter plot used to visualize data with a large number of data points, with a distribute of some higher-magnitude values. For example, in the GWAS(genome-wide association studies). Here we mainly focus on GWAS Manhattan plots. X-axis is genomic coordinates and Y-axis is negative logarithm of the associated P-value for each single nucleotide polymorphism. So higher the value, more stronger the association they are.

Usage

plotGrandLinear(obj, ..., facets, space.skip = 0.01, geom = NULL,
                 cutoff = NULL, cutoff.color = "red", cutoff.size = 1,
                 legend = FALSE, xlim, ylim, xlab, ylab, main)

Arguments

obj

GRanges object which contains extra p value, before users pass this object, they need to make sure the pvalue has been changed to -log10(p).

...

extra arguments passed. such as color, size, alpha.

facets

facets formula, such as group ~ .

space.skip

numeric value for skip ratio, between chromosome spaces.default is 0.01.

geom

geometric object, defualt is "point".

cutoff

A numeric vector which used as cutoff for Manhattan plot.

cutoff.color

A character specifying the color used for cutoff. Default is "red".

cutoff.size

A numeric value which used as cutoff line size.

legend

A logical value indicate whether to show legend or not. Default is FALSE which disabled the legend.

xlim

limits for x scale.

ylim

limits for y scale.

xlab

Label for xscale.

ylab

Label for yscale.

main

title.

Details

Please use seqlengths of the object and space.skip arguments to control the layout of the coordiant genome transformation.

aes(y = ...) is requried.

aes(color = ) is used to mapping to data variables, if just pass "color" without aes(), then will recycle the color to represent each chromosomes.please see the example below.

Value

Return a ggplot object.

Author(s)

Tengfei Yin

Examples

##  load
library(ggbio)
data(hg19IdeogramCyto, package = "biovizBase")
data(hg19Ideogram, package = "biovizBase")
library(GenomicRanges)
##  simul_gr
library(biovizBase)
gr <- GRanges(rep(c("chr1", "chr2"), each = 5),
              IRanges(start = rep(seq(1, 100, length = 5), times = 2),
                      width = 50))
autoplot(gr)

plot of chunk unnamed-chunk-1

##  coord:genome
autoplot(gr, coord = "genome")
## using coord:genome to parse x scale
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.

plot of chunk unnamed-chunk-1

gr.t <- transformToGenome(gr)
head(gr.t)
## GRanges with 6 ranges and 2 metadata columns:
##       seqnames     ranges strand |    .start      .end
##          <Rle>  <IRanges>  <Rle> | <numeric> <numeric>
##   [1]     chr1 [  1,  50]      * |         1        50
##   [2]     chr1 [ 25,  74]      * |        25        74
##   [3]     chr1 [ 50,  99]      * |        50        99
##   [4]     chr1 [ 75, 124]      * |        75       124
##   [5]     chr1 [100, 149]      * |       100       149
##   [6]     chr2 [  1,  50]      * |       180       229
##   ---
##   seqlengths:
##    chr1 chr2
##      NA   NA
##  is
is_coord_genome(gr.t)
## [1] TRUE
metadata(gr.t)$coord
## [1] "genome"
##  simul_snp
chrs <- as.character(levels(seqnames(hg19IdeogramCyto)))
seqlths <- seqlengths(hg19Ideogram)[chrs]
set.seed(1)
nchr <- length(chrs)
nsnps <- 100
gr.snp <- GRanges(rep(chrs,each=nsnps),
                  IRanges(start =
                          do.call(c, lapply(chrs, function(chr){
                            N <- seqlths[chr]
                            runif(nsnps,1,N)
                          })), width = 1),
                  SNP=sapply(1:(nchr*nsnps), function(x) paste("rs",x,sep='')),
                  pvalue =  -log10(runif(nchr*nsnps)),
                  group = sample(c("Normal", "Tumor"), size = nchr*nsnps,
                    replace = TRUE)
                  )
##  shorter
seqlengths(gr.snp)
##  chr1 chr10 chr11 chr12 chr13 chr14 chr15 chr16 chr17 chr18 chr19  chr2 
##    NA    NA    NA    NA    NA    NA    NA    NA    NA    NA    NA    NA 
## chr20 chr21 chr22  chr3  chr4  chr5  chr6  chr7  chr8  chr9  chrX  chrY 
##    NA    NA    NA    NA    NA    NA    NA    NA    NA    NA    NA    NA
nms <- seqnames(seqinfo(gr.snp))
nms.new <- gsub("chr", "", nms)
names(nms.new) <- nms
gr.snp <- renameSeqlevels(gr.snp, nms.new)
seqlengths(gr.snp)
##  1 10 11 12 13 14 15 16 17 18 19  2 20 21 22  3  4  5  6  7  8  9  X  Y 
## NA NA NA NA NA NA NA NA NA NA NA NA NA NA NA NA NA NA NA NA NA NA NA NA
##  unorder
autoplot(gr.snp, coord = "genome", geom = "point", aes(y = pvalue), space.skip = 0.01)
## using coord:genome to parse x scale
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.

plot of chunk unnamed-chunk-1

##  sort
gr.snp <- keepSeqlevels(gr.snp, c(1:22, "X", "Y"))
autoplot(gr.snp, coord = "genome", geom = "point", aes(y = pvalue), space.skip = 0.01)
## using coord:genome to parse x scale
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.

plot of chunk unnamed-chunk-1

##  with_seql
names(seqlths) <- gsub("chr", "", names(seqlths))
seqlengths(gr.snp) <- seqlths[names(seqlengths(gr.snp))]
autoplot(gr.snp, coord = "genome", geom = "point", aes(y = pvalue), space.skip = 0.01)
## using coord:genome to parse x scale
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.

plot of chunk unnamed-chunk-1

##  line
autoplot(gr.snp, coord = "genome", geom = "line", aes(y = pvalue, group = seqnames,
                                     color = seqnames))
## using coord:genome to parse x scale
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.

plot of chunk unnamed-chunk-1

##  plotGrandLinear
plotGrandLinear(gr.snp, aes(y = pvalue))
## using coord:genome to parse x scale
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.

plot of chunk unnamed-chunk-1

##  morecolor
plotGrandLinear(gr.snp, aes(y = pvalue, color = seqnames))
## using coord:genome to parse x scale
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.

plot of chunk unnamed-chunk-1

plotGrandLinear(gr.snp, aes(y = pvalue), color = c("green", "deepskyblue"))
## using coord:genome to parse x scale
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.

plot of chunk unnamed-chunk-1

plotGrandLinear(gr.snp, aes(y = pvalue), color = c("green", "deepskyblue", "red"))
## using coord:genome to parse x scale
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.

plot of chunk unnamed-chunk-1

plotGrandLinear(gr.snp, aes(y = pvalue), color = "red")
## using coord:genome to parse x scale
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.

plot of chunk unnamed-chunk-1

##  cutoff
plotGrandLinear(gr.snp, aes(y = pvalue), cutoff = 3, cutoff.color = "blue", cutoff.size = 4)
## using coord:genome to parse x scale
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.

plot of chunk unnamed-chunk-1

##  cutoff-low
plotGrandLinear(gr.snp, aes(y = pvalue)) + geom_hline(yintercept = 3, color = "blue", size = 4)
## using coord:genome to parse x scale
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.

plot of chunk unnamed-chunk-1

##  longer
## let's make a long name
nms <- seqnames(seqinfo(gr.snp))
nms.new <- paste("chr00000", nms, sep = "")
names(nms.new) <- nms
gr.snp <- renameSeqlevels(gr.snp, nms.new)
seqlengths(gr.snp)
##  chr000001  chr000002  chr000003  chr000004  chr000005  chr000006 
##  249250621  243199373  198022430  191154276  180915260  171115067 
##  chr000007  chr000008  chr000009 chr0000010 chr0000011 chr0000012 
##  159138663  146364022  141213431  135534747  135006516  133851895 
## chr0000013 chr0000014 chr0000015 chr0000016 chr0000017 chr0000018 
##  115169878  107349540  102531392   90354753   81195210   78077248 
## chr0000019 chr0000020 chr0000021 chr0000022  chr00000X  chr00000Y 
##   59128983   63025520   48129895   51304566  155270560   59373566
##  rotate
plotGrandLinear(gr.snp, aes(y = pvalue)) + theme(axis.text.x=element_text(angle=-90, hjust=0))
## using coord:genome to parse x scale
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.

plot of chunk unnamed-chunk-1

##  sessionInfo
sessionInfo()
## R Under development (unstable) (2012-07-17 r59871)
## Platform: x86_64-unknown-linux-gnu (64-bit)
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=en_US.UTF-8        LC_COLLATE=en_US.UTF-8    
##  [5] LC_MONETARY=en_US.UTF-8    LC_MESSAGES=en_US.UTF-8   
##  [7] LC_PAPER=C                 LC_NAME=C                 
##  [9] LC_ADDRESS=C               LC_TELEPHONE=C            
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C       
## 
## attached base packages:
## [1] grid      methods   stats     graphics  grDevices utils     datasets 
## [8] base     
## 
## other attached packages:
##  [1] biovizBase_1.5.9                       
##  [2] BSgenome.Hsapiens.UCSC.hg19_1.3.19     
##  [3] BSgenome_1.25.9                        
##  [4] VariantAnnotation_1.3.32               
##  [5] genefilter_1.39.0                      
##  [6] vsn_3.25.0                             
##  [7] rtracklayer_1.17.21                    
##  [8] Rsamtools_1.9.31                       
##  [9] Biostrings_2.25.12                     
## [10] TxDb.Hsapiens.UCSC.hg19.knownGene_2.8.0
## [11] GenomicFeatures_1.9.44                 
## [12] AnnotationDbi_1.19.46                  
## [13] Biobase_2.17.8                         
## [14] GenomicRanges_1.9.66                   
## [15] IRanges_1.15.48                        
## [16] BiocGenerics_0.3.2                     
## [17] ggbio_1.5.20                           
## [18] ggplot2_0.9.2.1                        
## [19] knitr_0.8                              
## 
## loaded via a namespace (and not attached):
##  [1] affy_1.35.1           affyio_1.25.0         annotate_1.35.5      
##  [4] BiocInstaller_1.5.12  biomaRt_2.13.2        bitops_1.0-4.1       
##  [7] cluster_1.14.2        colorspace_1.1-1      DBI_0.2-5            
## [10] dichromat_1.2-4       digest_0.5.1          evaluate_0.4.2       
## [13] formatR_0.6           gridExtra_0.9.1       gtable_0.1.1         
## [16] Hmisc_3.9-3           labeling_0.1          lattice_0.20-6       
## [19] limma_3.13.24         markdown_0.5.2        MASS_7.3-19          
## [22] memoise_0.1           munsell_0.4           parallel_2.16.0      
## [25] plyr_1.7.1            preprocessCore_1.19.0 proto_0.3-9.2        
## [28] RColorBrewer_1.0-5    RCurl_1.95-0          reshape2_1.2.1       
## [31] RSQLite_0.11.2        scales_0.2.2          splines_2.16.0       
## [34] stats4_2.16.0         stringr_0.6.1         survival_2.36-14     
## [37] tools_2.16.0          XML_3.95-0            xtable_1.7-0         
## [40] zlibbioc_1.3.0

[Package ggbio version 1.5.20 ]