Plot single-molecule footprinting data for a single genomic region.
Source:R/plotRegion.R
plotRegion.Rd
This function will visualize read-level or collapsed single-molecule
footprinting data, such as data imported using readModkitExtract
or readBedMethyl
.
Arguments
- se
A
SummarizedExperiment
object with read-level or collapsed single-molecule footprinting data (positions in rows and reads or samples in columns).- region
A
GRanges
object with a single region. Only data fromse
overlapping this region will be plotted. Alternatively, the region can be specified as a character scalar (e.g. "chr1:1200-1300") that can be coerced into aGRanges
object. IfNULL
(the default), all the data on the first sequence inse
will be visualized.- tracks
A list of named lists, representing the tracks to generate. Each element of the outer list defines one track, and has to contain at least list entries named 'trackData' (the name of a suitable assay in
se
) and 'trackType' (the type of plot), plus any additional arguments to the respective plot function. Currently supported plot types are"Point"
: A point plot displaying values in the assay.
"Smooth"
: A smoothed line plot displaying values in the assay.
"PointSmooth"
: A point and smoothed line plot displaying values in the assay.
"Lollipop"
: Lollipop plot (filled circles with the color representing the values in the assay).
"Heatmap"
: Heatmap plot (tiles with the color representing the values in the assay).
- modbaseSpace
A logical scalar. If
TRUE
, the x-axis will be shown in the space of modified bases and contain only the positions at which there are modified bases in the data without any gaps between them. IfFALSE
, the x-axis will show the genomic coordinate on which the modified bases are typically irregularly spaced.- sequenceContext
A character vector with sequence context(s) to plot. Only positions that match one of the provided sequence contexts will be included in the plot. Sequence contexts can be provided using IUPAC redundancy codes. The sequence contexts of modified bases are obtained from
rowData(se)$sequenceContext
and thus requires thatse
contains the appropriate information, for example by setting thesequenceContext
andsequenceReference
arguments ofreadBedMethyl
when it was generated, or by adding it usingaddSeqContext
.
Value
A ggplot
object with tracks selected by
tracks
.
See also
readModBam
, readModkitExtract
and
readBedMethyl
for reading read-level and summarized
footprinting data.
Examples
# summarized data (5mC)
bmfiles <- system.file("extdata",
c("modkit_pileup_1.bed.gz", "modkit_pileup_2.bed.gz"),
package = "footprintR")
reffile <- system.file("extdata", "reference.fa.gz", package = "footprintR")
seA <- readBedMethyl(bmfiles, modbase = "m",
sequenceContextWidth = 3, sequenceReference = reffile)
plotRegion(seA, region = "chr1:6940000-6955000", sequenceContext = "GCH")
plotRegion(seA, region = "chr1:6940000-6955000", sequenceContext = "HCG")
plotRegion(seA, region = "chr1:6940000-6955000",
tracks = list(list(trackData = "Nvalid", trackType = "Smooth")))
# read-level data (6mA)
extractfiles <- system.file("extdata",
c("modkit_extract_rc_6mA_1.tsv.gz",
"modkit_extract_rc_6mA_2.tsv.gz"),
package = "footprintR")
seB <- readModkitExtract(extractfiles, modbase = "a", filter = "modkit")
# Lollipop plot
plotRegion(seB, region = "chr1:6935800-6935900",
tracks = list(list(trackData = "mod_prob", trackType = "Lollipop")))
#> Warning: the standard deviation is zero
# Heatmap plots (observed only or filled)
plotRegion(seB, region = "chr1:6935800-6935900",
tracks = list(list(trackData = "mod_prob", trackType = "Heatmap")))
#> Warning: the standard deviation is zero
plotRegion(seB, region = "chr1:6935800-6935900",
tracks = list(list(trackData = "mod_prob", trackType = "Heatmap",
interpolate = TRUE)))
#> Warning: the standard deviation is zero
# multiple plots
plotRegion(seB, region = "chr1:6935400-6935450",
tracks = list(list(trackData = "mod_prob", trackType = "Lollipop",
size = 4),
list(trackData = "mod_prob", trackType = "Heatmap")),
modbaseSpace = TRUE)