VplotR

Jacques Serizay

2022-11-01

Introduction

Overview

VplotR is an R package streamlining the process of generating V-plots, i.e. two-dimensional paired-end fragment density plots. It contains functions to import paired-end sequencing bam files from any type of DNA accessibility experiments (e.g. ATAC-seq, DNA-seq, MNase-seq) and can produce V-plots and one-dimensional footprint profiles over single or aggregated genomic loci of interest. The R package is well integrated within the Bioconductor environment and easily fits in standard genomic analysis workflows. Integrating V-plots into existing analytical frameworks has already brought additional insights in chromatin organization (Serizay et al., 2020).

The main user-level functions of VplotR are getFragmentsDistribution(), plotVmat(), plotFootprint() and plotProfile().

Installation

VplotR can be installed from Bioconductor:

Importing sequencing datasets

Using importPEBamFiles() function

Paired-end .bam files are imported using the importPEBamFiles() function as follows:

Provided datasets

Several datasets are available for this package:

data(ce11_proms)
ce11_proms
#> GRanges object with 17458 ranges and 3 metadata columns:
#>           seqnames            ranges strand |   TSS.fwd   TSS.rev which.tissues
#>              <Rle>         <IRanges>  <Rle> | <numeric> <numeric>      <factor>
#>       [1]     chrI       11273-11423      + |     11294     11416       Intest.
#>       [2]     chrI       11273-11423      - |     11294     11416       Intest.
#>       [3]     chrI       26903-27053      - |     27038     26901       Ubiq.  
#>       [4]     chrI       36019-36169      - |     36112     36028       Intest.
#>       [5]     chrI       42127-42277      - |     42216     42119       Soma   
#>       ...      ...               ...    ... .       ...       ...           ...
#>   [17454]     chrX 17670496-17670646      + |  17670678  17670478  Muscle      
#>   [17455]     chrX 17684894-17685044      - |  17684919  17684932  Hypod.      
#>   [17456]     chrX 17686030-17686180      - |  17686189  17686064  Unclassified
#>   [17457]     chrX 17694789-17694939      + |  17694962  17694934  Intest.     
#>   [17458]     chrX 17711839-17711989      - |  17711974  17711854  Germline    
#>   -------
#>   seqinfo: 6 sequences from an unspecified genome; no seqlengths
data(ATAC_ce11_Serizay2020)
ATAC_ce11_Serizay2020
#> $Germline
#> GRanges object with 462371 ranges and 0 metadata columns:
#>            seqnames            ranges strand
#>               <Rle>         <IRanges>  <Rle>
#>        [1]     chrI           426-514      +
#>        [2]     chrI         3588-3854      +
#>        [3]     chrI         3640-3798      +
#>        [4]     chrI         3650-3694      +
#>        [5]     chrI         3732-3863      +
#>        ...      ...               ...    ...
#>   [462367]     chrX 17712277-17712469      -
#>   [462368]     chrX 17712279-17712342      -
#>   [462369]     chrX 17712282-17712565      -
#>   [462370]     chrX 17712285-17712384      -
#>   [462371]     chrX 17712287-17712576      -
#>   -------
#>   seqinfo: 7 sequences from an unspecified genome; no seqlengths
#> 
#> $Neurons
#> GRanges object with 367935 ranges and 0 metadata columns:
#>            seqnames            ranges strand
#>               <Rle>         <IRanges>  <Rle>
#>        [1]     chrI         4011-4241      +
#>        [2]     chrI         7397-7614      +
#>        [3]     chrI       11279-11502      +
#>        [4]     chrI       12744-12819      +
#>        [5]     chrI       14381-14433      +
#>        ...      ...               ...    ...
#>   [367931]     chrX 17687948-17687982      -
#>   [367932]     chrX 17699614-17699853      -
#>   [367933]     chrX 17706798-17706923      -
#>   [367934]     chrX 17708264-17708347      -
#>   [367935]     chrX 17709920-17710007      -
#>   -------
#>   seqinfo: 7 sequences from an unspecified genome; no seqlengths

Fragment size distribution

A preliminary control to check the distribution of fragment sizes (regardless of their location relative to genomic loci) can be performed using the getFragmentsDistribution() function.

df <- getFragmentsDistribution(
    MNase_sacCer3_Henikoff2011, 
    ABF1_sacCer3
)
#> Warning in as.cls(x): NAs introduced by coercion

#> Warning in as.cls(x): NAs introduced by coercion

#> Warning in as.cls(x): NAs introduced by coercion
p <- ggplot(df, aes(x = x, y = y)) + geom_line() + theme_ggplot2()
p
#> Warning: Removed 2 row(s) containing missing values (geom_path).

Vplot(s)

Single Vplot

Once data is imported, a V-plot of paired-end fragments over loci of interest is generated using the plotVmat() function:

Multiple Vplots

The generation of multiple V-plots can be parallelized using a list of parameters:

For instance, ATAC-seq fragment density can be visualized at different classes of ubiquitous and tissue-specific promoters in C. elegans.

Vplots normalization

Different normalization approaches are available using the normFun argument.

Footprints

VplotR also implements a function to profile the footprint from MNase or ATAC-seq over sets of genomic loci. For instance, CTCF is known for its ~40-bp large footprint at its binding loci.

p <- plotFootprint(
    MNase_sacCer3_Henikoff2011,
    ABF1_sacCer3
)
#> - Getting cuts
#> - Getting cut coverage
#> - Getting cut coverage / target
#> - Reformatting data into matrix
#> - Plotting footprint
p

Local fragment distribution

VplotR provides a function to plot the distribution of paired-end fragments over an individual genomic window.

data(MNase_sacCer3_Henikoff2011_subset)
genes_sacCer3 <- GenomicFeatures::genes(TxDb.Scerevisiae.UCSC.sacCer3.sgdGene::
    TxDb.Scerevisiae.UCSC.sacCer3.sgdGene
)
p <- plotProfile(
    fragments = MNase_sacCer3_Henikoff2011_subset,
    window = "chrXV:186,400-187,400", 
    loci = ABF1_sacCer3, 
    annots = genes_sacCer3,
    min = 20, max = 200, alpha = 0.1, size = 1.5
)
#> Filtering and resizing fragments
#> 32276 fragments mapped over 1001 bases
#> Filtering and resizing fragments
#> Generating plot
#> Warning: Removed 49 row(s) containing missing values (geom_path).
#> Warning: Removed 5176 rows containing missing values (geom_point).
#> Warning: Removed 19 row(s) containing missing values (geom_path).
p

Session Info

sessionInfo()
#> R version 4.2.1 (2022-06-23)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 20.04.5 LTS
#> 
#> Matrix products: default
#> BLAS:   /home/biocbuild/bbs-3.16-bioc/R/lib/libRblas.so
#> LAPACK: /home/biocbuild/bbs-3.16-bioc/R/lib/libRlapack.so
#> 
#> locale:
#>  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
#>  [3] LC_TIME=en_GB              LC_COLLATE=C              
#>  [5] LC_MONETARY=en_US.UTF-8    LC_MESSAGES=en_US.UTF-8   
#>  [7] LC_PAPER=en_US.UTF-8       LC_NAME=C                 
#>  [9] LC_ADDRESS=C               LC_TELEPHONE=C            
#> [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C       
#> 
#> attached base packages:
#> [1] stats4    stats     graphics  grDevices utils     datasets  methods  
#> [8] base     
#> 
#> other attached packages:
#> [1] VplotR_1.8.0         magrittr_2.0.3       ggplot2_3.3.6       
#> [4] GenomicRanges_1.50.0 GenomeInfoDb_1.34.0  IRanges_2.32.0      
#> [7] S4Vectors_0.36.0     BiocGenerics_0.44.0 
#> 
#> loaded via a namespace (and not attached):
#>  [1] bitops_1.0-7                               
#>  [2] matrixStats_0.62.0                         
#>  [3] bit64_4.0.5                                
#>  [4] filelock_1.0.2                             
#>  [5] RColorBrewer_1.1-3                         
#>  [6] progress_1.2.2                             
#>  [7] httr_1.4.4                                 
#>  [8] tools_4.2.1                                
#>  [9] bslib_0.4.0                                
#> [10] utf8_1.2.2                                 
#> [11] R6_2.5.1                                   
#> [12] DBI_1.1.3                                  
#> [13] colorspace_2.0-3                           
#> [14] withr_2.5.0                                
#> [15] tidyselect_1.2.0                           
#> [16] prettyunits_1.1.1                          
#> [17] bit_4.0.4                                  
#> [18] curl_4.3.3                                 
#> [19] compiler_4.2.1                             
#> [20] cli_3.4.1                                  
#> [21] Biobase_2.58.0                             
#> [22] xml2_1.3.3                                 
#> [23] DelayedArray_0.24.0                        
#> [24] rtracklayer_1.58.0                         
#> [25] labeling_0.4.2                             
#> [26] sass_0.4.2                                 
#> [27] scales_1.2.1                               
#> [28] rappdirs_0.3.3                             
#> [29] stringr_1.4.1                              
#> [30] digest_0.6.30                              
#> [31] Rsamtools_2.14.0                           
#> [32] rmarkdown_2.17                             
#> [33] XVector_0.38.0                             
#> [34] pkgconfig_2.0.3                            
#> [35] htmltools_0.5.3                            
#> [36] MatrixGenerics_1.10.0                      
#> [37] dbplyr_2.2.1                               
#> [38] fastmap_1.1.0                              
#> [39] highr_0.9                                  
#> [40] rlang_1.0.6                                
#> [41] RSQLite_2.2.18                             
#> [42] TxDb.Scerevisiae.UCSC.sacCer3.sgdGene_3.2.2
#> [43] jquerylib_0.1.4                            
#> [44] BiocIO_1.8.0                               
#> [45] farver_2.1.1                               
#> [46] generics_0.1.3                             
#> [47] zoo_1.8-11                                 
#> [48] jsonlite_1.8.3                             
#> [49] BiocParallel_1.32.0                        
#> [50] dplyr_1.0.10                               
#> [51] RCurl_1.98-1.9                             
#> [52] GenomeInfoDbData_1.2.9                     
#> [53] Matrix_1.5-1                               
#> [54] Rcpp_1.0.9                                 
#> [55] munsell_0.5.0                              
#> [56] fansi_1.0.3                                
#> [57] lifecycle_1.0.3                            
#> [58] stringi_1.7.8                              
#> [59] yaml_2.3.6                                 
#> [60] SummarizedExperiment_1.28.0                
#> [61] zlibbioc_1.44.0                            
#> [62] plyr_1.8.7                                 
#> [63] BiocFileCache_2.6.0                        
#> [64] grid_4.2.1                                 
#> [65] blob_1.2.3                                 
#> [66] parallel_4.2.1                             
#> [67] crayon_1.5.2                               
#> [68] lattice_0.20-45                            
#> [69] Biostrings_2.66.0                          
#> [70] cowplot_1.1.1                              
#> [71] GenomicFeatures_1.50.0                     
#> [72] hms_1.1.2                                  
#> [73] KEGGREST_1.38.0                            
#> [74] knitr_1.40                                 
#> [75] pillar_1.8.1                               
#> [76] rjson_0.2.21                               
#> [77] reshape2_1.4.4                             
#> [78] codetools_0.2-18                           
#> [79] biomaRt_2.54.0                             
#> [80] XML_3.99-0.12                              
#> [81] glue_1.6.2                                 
#> [82] evaluate_0.17                              
#> [83] png_0.1-7                                  
#> [84] vctrs_0.5.0                                
#> [85] gtable_0.3.1                               
#> [86] assertthat_0.2.1                           
#> [87] cachem_1.0.6                               
#> [88] xfun_0.34                                  
#> [89] restfulr_0.0.15                            
#> [90] tibble_3.1.8                               
#> [91] GenomicAlignments_1.34.0                   
#> [92] AnnotationDbi_1.60.0                       
#> [93] memoise_2.0.1                              
#> [94] ellipsis_0.3.2