BSgenomeViews-class: BSgenomeViews objects

Description Usage Arguments Constructors Accessors Coercion Subsetting DNAStringSet methods Author(s) See Also Examples

Description

The BSgenomeViews class is a container for storing a set of genomic positions on a BSgenome object, called the "subject" in this context.

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
## Constructor
## ------------

BSgenomeViews(subject, granges)

## Accessors
## ---------

## S4 method for signature 'BSgenomeViews'
subject(x)
## S4 method for signature 'BSgenomeViews'
granges(x, use.mcols=FALSE)

## S4 method for signature 'BSgenomeViews'
length(x)
## S4 method for signature 'BSgenomeViews'
names(x)
## S4 method for signature 'BSgenomeViews'
seqnames(x)
## S4 method for signature 'BSgenomeViews'
start(x)
## S4 method for signature 'BSgenomeViews'
end(x)
## S4 method for signature 'BSgenomeViews'
width(x)
## S4 method for signature 'BSgenomeViews'
strand(x)
## S4 method for signature 'BSgenomeViews'
ranges(x, use.mcols=FALSE)
## S4 method for signature 'BSgenomeViews'
elementNROWS(x)
## S4 method for signature 'BSgenomeViews'
seqinfo(x)

## DNAStringSet methods
## --------------------

## S4 method for signature 'BSgenomeViews'
seqtype(x)

## S4 method for signature 'BSgenomeViews'
nchar(x, type="chars", allowNA=FALSE)

## S4 method for signature 'BSgenomeViews'
unlist(x, recursive=TRUE, use.names=TRUE)

## S4 method for signature 'BSgenomeViews'
alphabetFrequency(x, as.prob=FALSE, collapse=FALSE, baseOnly=FALSE)

## S4 method for signature 'BSgenomeViews'
hasOnlyBaseLetters(x)

## S4 method for signature 'BSgenomeViews'
uniqueLetters(x)

## S4 method for signature 'BSgenomeViews'
letterFrequency(x, letters, OR="|", as.prob=FALSE, collapse=FALSE)

## S4 method for signature 'BSgenomeViews'
oligonucleotideFrequency(x, width, step=1,
                         as.prob=FALSE, as.array=FALSE,
                         fast.moving.side="right", with.labels=TRUE, simplify.as="matrix")

## S4 method for signature 'BSgenomeViews'
nucleotideFrequencyAt(x, at, as.prob=FALSE, as.array=TRUE,
                      fast.moving.side="right", with.labels=TRUE)

## S4 method for signature 'BSgenomeViews'
consensusMatrix(x, as.prob=FALSE, shift=0L, width=NULL, baseOnly=FALSE)

## S4 method for signature 'BSgenomeViews'
consensusString(x, ambiguityMap=IUPAC_CODE_MAP, threshold=0.25,
                shift=0L, width=NULL)

Arguments

subject

A BSgenome object or the name of a reference genome specified in a way that is accepted by the getBSgenome function. In that case the corresponding BSgenome data package needs to be already installed (see ?getBSgenome for the details).

granges

A GRanges object containing ranges relative to the genomic sequences stored in subject.

x

A BSgenomeViews object.

use.mcols

TRUE or FALSE (the default). Whether the metadata columns on x (accessible with mcols(x)) should be propagated to the returned object or not.

type, allowNA, recursive, use.names

Ignored.

as.prob, letters, OR, width

See ?alphabetFrequency and ?oligonucleotideFrequency in the Biostrings package.

collapse, baseOnly

See ?alphabetFrequency in the Biostrings package.

step, as.array, fast.moving.side, with.labels, simplify.as, at

See ?oligonucleotideFrequency in the Biostrings package.

shift, ambiguityMap, threshold

See ?consensusMatrix in the Biostrings package.

Constructors

BSgenomeViews(subject, granges): Make a BSgenomeViews object by putting the views specified by granges on top of the genomic sequences stored in subject. See above for how argument subject and granges should be specified.

Views(subject, granges): Equivalent to BSgenomeViews(subject, granges). Provided for convenience.

Accessors

In the code snippets below, x is a BSgenomeViews object.

subject(x): Return the BSgenome object containing the full genomic sequences on top of which the views in x are defined.

granges(x, use.mcols=FALSE): Return the genomic ranges of the views as a GRanges object. These ranges are relative to the genomic sequences stored in subject(x).

length(x): The number of views in x.

names(x): The names of the views in x.

seqnames(x), start(x), end(x), width(x), strand(x): Equivalent to seqnames(granges(x)), start(granges(x)), end(granges(x)), width(granges(x)), strand(granges(x)), respectively.

ranges(x, use.mcols=FALSE): Equivalent to ranges(granges(x, use.mcols), use.mcols).

elementNROWS(x): Equivalent to width(x).

seqinfo(x): Equivalent to seqinfo(subject(x)) and to seqinfo(granges(x)) (both are guaranteed to be the same). See ?seqinfo in the GenomeInfoDb package for more information.

Coercion

In the code snippets below, x is a BSgenomeViews object.

as(x, "DNAStringSet"): Turn x into a DNAStringSet object by extxracting the DNA sequence corresponding to each view. Alternatively as(x, "XStringSet") can be used for this, and is equivalent to as(x, "DNAStringSet").

as.character(x): Equivalent to as.character(as(x, "DNAStringSet")).

as.data.frame(x): Turn x into a data.frame.

Subsetting

x[i]: Select the views specified by i.

x[[i]]: Extract the one view specified by i.

DNAStringSet methods

For convenience, some methods defined for DNAStringSet objects in the Biostrings package can be used directly on a BSgenomeViews object. In that case, everything happens like if the BSgenomeViews object x was turned into a DNAStringSet object (with as(x, "DNAStringSet")) before it's passed to the method for DNAStringSet objects.

At the moment, the list of such methods is: seqtype, nchar,XStringSet-method, unlist,XStringSet-method, alphabetFrequency, hasOnlyBaseLetters, uniqueLetters, letterFrequency, oligonucleotideFrequency, nucleotideFrequencyAt, consensusMatrix, and consensusString.

See the corresponding man page in the Biostrings package for a description of these methods.

Author(s)

H. Pag<c3><a8>s

See Also

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
library(BSgenome.Mmusculus.UCSC.mm10)
genome <- BSgenome.Mmusculus.UCSC.mm10
library(TxDb.Mmusculus.UCSC.mm10.knownGene)
txdb <- TxDb.Mmusculus.UCSC.mm10.knownGene
ex <- exons(txdb, columns=c("exon_id", "tx_name", "gene_id"))
v <- Views(genome, ex)
v

subject(v)
granges(v)
seqinfo(v)
as(v, "DNAStringSet")

v10 <- v[1:10]  # select the first 10 views
subject(v10)    # same as subject(v)
granges(v10)
seqinfo(v10)    # same as seqinfo(v)
as(v10, "DNAStringSet")
alphabetFrequency(v10)
alphabetFrequency(v10, collapse=TRUE)

v12 <- v[width(v) <= 12]  # select the views of 12 nucleotides or less
head(as.data.frame(v12))
trinucleotideFrequency(v12, simplify.as="collapsed")

## BSgenomeViews objects are list-like objects. That is, the
## BSgenomeViews class derives from List and typical list/List
## operations (e.g. [[, elementNROWS(), unlist(), elementType(),
## etc...) work on these objects:
is(v12, "List")  # TRUE
v12[[2]]
head(elementNROWS(v12))  # elementNROWS(v) is the same as width(v)
unlist(v12) 
elementType(v12)

BSgenome documentation built on Nov. 8, 2020, 7:48 p.m.