VplotR

knitr::opts_chunk$set(
    collapse = TRUE,
    comment = "#>"
)
suppressPackageStartupMessages({
    library(GenomicRanges)
    library(ggplot2)
    library(magrittr)
    library(VplotR)
})

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:

if(!requireNamespace("BiocManager", quietly = TRUE))
install.packages("BiocManager")
BiocManager::install("VplotR")
library("VplotR")

Importing sequencing datasets

Using importPEBamFiles() function

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

library(VplotR)
bamfile <- system.file("extdata", "ex1.bam", package = "Rsamtools")
fragments <- importPEBamFiles(
    bamfile, 
    shift_ATAC_fragments = TRUE
)
fragments

Provided datasets

Several datasets are available for this package:

data(ce11_proms)
ce11_proms
data(ATAC_ce11_Serizay2020)
ATAC_ce11_Serizay2020
data(ABF1_sacCer3)
ABF1_sacCer3
data(MNase_sacCer3_Henikoff2011)
MNase_sacCer3_Henikoff2011

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
)
p <- ggplot(df, aes(x = x, y = y)) + geom_line() + theme_ggplot2()
p

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:

p <- plotVmat(x = MNase_sacCer3_Henikoff2011, granges = ABF1_sacCer3)
p

Multiple Vplots

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

list_params <- list(
    "MNase\n@ ABF1" = list(MNase_sacCer3_Henikoff2011, ABF1_sacCer3), 
    "MNase\n@ random loci" = list(
        MNase_sacCer3_Henikoff2011, sampleGRanges(ABF1_sacCer3)
    )
)
p <- plotVmat(
    list_params, 
    cores = 1
)
p

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

list_params <- list(
    "Germline ATACseq\n@ Ubiq. proms" = list(
        ATAC_ce11_Serizay2020[['Germline']], 
        ce11_proms[ce11_proms$which.tissues == 'Ubiq.']
    ), 
    "Germline ATACseq\n@ Germline proms" = list(
        ATAC_ce11_Serizay2020[['Germline']], 
        ce11_proms[ce11_proms$which.tissues == 'Germline']
    ),
    "Neuron ATACseq\n@ Ubiq. proms" = list(
        ATAC_ce11_Serizay2020[['Neurons']], 
        ce11_proms[ce11_proms$which.tissues == 'Ubiq.']
    ), 
    "Neuron ATACseq\n@ Neuron proms" = list(
        ATAC_ce11_Serizay2020[['Neurons']], 
        ce11_proms[ce11_proms$which.tissues == 'Neurons']
    )
)
p <- plotVmat(
    list_params, 
    cores = 1,
    nrow = 2, ncol = 5
)
p

Vplots normalization

Different normalization approaches are available using the normFun argument.

# No normalization 
p <- plotVmat(
    list_params, 
    cores = 1, 
    nrow = 2, ncol = 5, 
    verbose = FALSE,
    normFun = 'none'
)
p
# Library depth + number of loci of interest (default)
p <- plotVmat(
    list_params, 
    cores = 1, 
    nrow = 2, ncol = 5, 
    verbose = FALSE,
    normFun = 'libdepth+nloci'
)
p
# Zscore
p <- plotVmat(
    list_params, 
    cores = 1, 
    nrow = 2, ncol = 5, 
    verbose = FALSE,
    normFun = 'zscore'
)
p
# Quantile
p <- plotVmat(
    list_params, 
    cores = 1, 
    nrow = 2, ncol = 5, 
    verbose = FALSE,
    normFun = 'quantile', 
    s = 0.99
)
p

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
)
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
)
p

Session Info

sessionInfo()


Try the VplotR package in your browser

Any scripts or data that you put into this service are public.

VplotR documentation built on Nov. 8, 2020, 7:50 p.m.