Description Objects from the Class Slots Methods Utilities Coercion Assays Author(s) See Also Examples
A class for representing whole-genome or capture bisulfite sequencing data.
An object from the class links together several pieces of information.
(1) genomic locations stored as a GRanges
object, a location by
samples matrix of M values, a location by samples matrix of Cov
(coverage) values and phenodata information. In addition, there are
slots for representing smoothed data. This class is an extension of
RangedSummarizedExperiment.
trans
:Object of class function
. This
function transforms the coef
slot from the scale the
smoothing was done to the 0-1 methylation scale.
parameters
:Object of class list
. A list of
parameters representing for example how the data was smoothed.
signature(x = "BSseq")
: Subsetting by location
(using integer indices) or sample (using integers or sample
names).
Unlike for RangedSummarizedExperiment
,
length()
is the number of methylation loci (equal to
length(granges(x))
).
Sample names and its replacement
function for the object. This is an alias for colnames
.
Obtain and replace the pData
slot of the
phenoData
slot. This is an alias for colData
.
The show method.
This function combines two BSSeq
objects.
The genomic locations of the new object is the union of the
genomic locations of the individual objects. In addition, the
methylation data matrices are placed next to each other (as
appropriate wrt. the new genomic locations) and zeros are entered
into the matrices as needed.
This class extends RangedSummarizedExperiment and therefore
inherits a number of useful GRanges
methods that operate on the
rowRanges
slot, used for accessing and setting the genomic locations
and also do subsetByOverlaps
.
There are a number of almost methods-like functions for operating on
objects of class BSseq
, including getBSseq
,
collapseBSseq
, and orderBSseq
. They are detailed below.
collapseBSseq(BSseq, columns)
is used to collapse an object of class BSseq
. By
collapsing we simply mean that certain columns (samples) are merge
together by summing up the methylation evidence and coverage.
This is a useful function if you start by reading in a dataset
based on say flowcells and you (after QC) want to simply add a
number of flowcells into one sample. The argument columns
specify which samples are to be merged, in the following way: it
is a character vector of new sample names, and the names of the
column vector indicates which samples in the BSseq
object
are to be collapsed. If columns
have the same length as
the number of rows of BSseq
(and has no names) it is
assumed that the ordering corresponds to the sample ordering in
BSseq
.
orderBSseq(BSseq, seqOrder = NULL)
simply orders an object of class BSseq
according to
(increasing) genomic locations. The seqOrder
vector is a
character vector of seqnames(BSseq)
describing the order of
the chromosomes. This is useful for ordering chr1
before
chr10
.
chrSelectBSseq(BSseq, seqnames = NULL, order = FALSE)
subsets and optionally reorders an object of class BSseq
.
The seqnames
vector is a character vector of
seqnames(BSseq)
describing which chromosomes should be
retained. If order
is TRUE
, the chromosomes are
also re-ordered using orderBSseq
.
getBSseq(BSseq, type = c("Cov", "M", "gr", "coef",
"se.coef", "trans", "parameters"))
is a general accessor: is used to obtain a specific slot of an
object of class BSseq
. It is primarily intended for
internal use in the package, for users we recommend granges
to get the genomic locations, getCoverage
to get the
coverage slots and getMeth
to get the smoothed values (if
they exist).
hasBeenSmoothed(BSseq)
This function returns a logical depending on whether or not the
BSseq
object has been smoothed using BSmooth
.
combineList(list, BACKEND = NULL)
This function function is a faster way of using combine
on
multiple BSseq objects. The input is a list, with each
component an object of class BSseq. The (slower)
alternative is to use Reduce(combine, list)
.
The BACKEND
argument determines which backend should be used for
the 'M' and 'Cov' matrices and, if present, the 'coef' and 'se.coef'
matrices (the latter two can only be combined if all objects have the
same rowRanges). The default, BACKEND = NULL
, corresponds to using
matrix objects. See
?DelayedArray::setAutoRealizationBackend
for alternative
backends.
strandCollapse(BSseq, shift = TRUE)
This function operates on a BSseq
objects which has
stranded loci (i.e. loci where the strand is one of ‘+’ or
‘-’). It will collapse the methylation and coverage
information across the two strands, unstranding the loci in the process
and potentially re-ordering them.
The argument shift
indicates whether the positions for the loci
on the reverse strand should be shifted one (i.e. the positions for
these loci are the positions of the ‘G’ in the
‘CpG’; this is the case for Bismark output for example).
Package versions 1.5.2 and 1.11.1 introduced a new version of representing
‘BSseq’ objects. You can update old serialized (saved)
objects by invoking x <- updateObject(x)
.
This class overrides the default implementation of assays
to
make it faster. Per default, no names are added to the returned data
matrices.
Assay names can conveniently be obtained by the function
assayNames(x)
Kasper Daniel Hansen khansen@jhsph.edu
The package vignette. BSseq
for the constructor
function. RangedSummarizedExperiment
for the underlying class. getBSseq
,
getCoverage
, and getMeth
for accessing the
data stored in the object and finally BSmooth
for
smoothing the bisulfite sequence data.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 | M <- matrix(1:9, 3,3)
colnames(M) <- c("A1", "A2", "A3")
BStest <- BSseq(pos = 1:3, chr = c("chr1", "chr2", "chr1"), M = M, Cov = M + 2)
chrSelectBSseq(BStest, seqnames = "chr1", order = TRUE)
collapseBSseq(BStest, group = c("A", "A", "B"))
#-------------------------------------------------------------------------------
# An example using a HDF5-backed BSseq object
#
hdf5_BStest <- realize(BStest, "HDF5Array")
chrSelectBSseq(hdf5_BStest, seqnames = "chr1", order = TRUE)
collapseBSseq(
BSseq = hdf5_BStest,
group = c("A", "A", "B"),
BACKEND = "HDF5Array",
type = "integer")
|
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.