The BSgenomeViews class is a container for storing a set of genomic
positions on a BSgenome object, called the "subject" in this
context.
Usage
## 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).
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><83><c2><a8>s
See Also
The BSgenome class.
The GRanges class in the GenomicRanges
package.
The DNAStringSet class in the Biostrings
package.
The seqinfo and related getters in the
GenomeInfoDb package for getting the sequence information
stored in an object.
TxDb objects in the GenomicFeatures
package.
Examples
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)
Results
R version 3.3.1 (2016-06-21) -- "Bug in Your Hair"
Copyright (C) 2016 The R Foundation for Statistical Computing
Platform: x86_64-pc-linux-gnu (64-bit)
R is free software and comes with ABSOLUTELY NO WARRANTY.
You are welcome to redistribute it under certain conditions.
Type 'license()' or 'licence()' for distribution details.
R is a collaborative project with many contributors.
Type 'contributors()' for more information and
'citation()' on how to cite R or R packages in publications.
Type 'demo()' for some demos, 'help()' for on-line help, or
'help.start()' for an HTML browser interface to help.
Type 'q()' to quit R.
> library(BSgenome)
Loading required package: BiocGenerics
Loading required package: parallel
Attaching package: 'BiocGenerics'
The following objects are masked from 'package:parallel':
clusterApply, clusterApplyLB, clusterCall, clusterEvalQ,
clusterExport, clusterMap, parApply, parCapply, parLapply,
parLapplyLB, parRapply, parSapply, parSapplyLB
The following objects are masked from 'package:stats':
IQR, mad, xtabs
The following objects are masked from 'package:base':
Filter, Find, Map, Position, Reduce, anyDuplicated, append,
as.data.frame, cbind, colnames, do.call, duplicated, eval, evalq,
get, grep, grepl, intersect, is.unsorted, lapply, lengths, mapply,
match, mget, order, paste, pmax, pmax.int, pmin, pmin.int, rank,
rbind, rownames, sapply, setdiff, sort, table, tapply, union,
unique, unsplit
Loading required package: S4Vectors
Loading required package: stats4
Attaching package: 'S4Vectors'
The following objects are masked from 'package:base':
colMeans, colSums, expand.grid, rowMeans, rowSums
Loading required package: IRanges
Loading required package: GenomeInfoDb
Loading required package: GenomicRanges
Loading required package: Biostrings
Loading required package: XVector
Loading required package: rtracklayer
> png(filename="/home/ddbj/snapshot/RGM3/R_BC/result/BSgenome/BSgenomeViews-class.Rd_%03d_medium.png", width=480, height=480)
> ### Name: BSgenomeViews-class
> ### Title: BSgenomeViews objects
> ### Aliases: class:BSgenomeViews BSgenomeViews-class BSgenomeViews
> ### Views,BSgenome-method subject,BSgenomeViews-method
> ### granges,BSgenomeViews-method length,BSgenomeViews-method
> ### names,BSgenomeViews-method seqnames,BSgenomeViews-method
> ### start,BSgenomeViews-method end,BSgenomeViews-method
> ### width,BSgenomeViews-method strand,BSgenomeViews-method
> ### ranges,BSgenomeViews-method elementNROWS,BSgenomeViews-method
> ### seqinfo,BSgenomeViews-method coerce,BSgenomeViews,DNAStringSet-method
> ### coerce,BSgenomeViews,XStringSet-method
> ### as.character,BSgenomeViews-method as.data.frame,BSgenomeViews-method
> ### extractROWS,BSgenomeViews-method getListElement,BSgenomeViews-method
> ### seqtype,BSgenomeViews-method nchar,BSgenomeViews-method
> ### unlist,BSgenomeViews-method alphabetFrequency,BSgenomeViews-method
> ### hasOnlyBaseLetters,BSgenomeViews-method
> ### uniqueLetters,BSgenomeViews-method
> ### letterFrequency,BSgenomeViews-method
> ### oligonucleotideFrequency,BSgenomeViews-method
> ### nucleotideFrequencyAt,BSgenomeViews-method
> ### consensusMatrix,BSgenomeViews-method
> ### consensusString,BSgenomeViews-method show,BSgenomeViews-method
> ### Keywords: methods classes
>
> ### ** Examples
>
> library(BSgenome.Mmusculus.UCSC.mm10)
> genome <- BSgenome.Mmusculus.UCSC.mm10
> library(TxDb.Mmusculus.UCSC.mm10.knownGene)
Loading required package: GenomicFeatures
Loading required package: AnnotationDbi
Loading required package: Biobase
Welcome to Bioconductor
Vignettes contain introductory material; view with
'browseVignettes()'. To cite Bioconductor, see
'citation("Biobase")', and for packages 'citation("pkgname")'.
> txdb <- TxDb.Mmusculus.UCSC.mm10.knownGene
> ex <- exons(txdb, columns=c("exon_id", "tx_name", "gene_id"))
> v <- Views(genome, ex)
> v
BSgenomeViews object with 256721 views and 3 metadata columns:
seqnames ranges strand dna |
<Rle> <IRanges> <Rle> <DNAStringSet> |
[1] chr1 [4807893, 4807982] + [GCACTGTCCG...CACCGCCGCG] |
[2] chr1 [4808455, 4808486] + [GTTATTTTCC...GAGATACAGG] |
[3] chr1 [4828584, 4828649] + [GCATGGATGG...GTCCACATGC] |
[4] chr1 [4830268, 4830315] + [CCCTGTGATG...TGCCTTCTTG] |
[5] chr1 [4832311, 4832381] + [GTTTGATATC...GCAGAAACCG] |
... ... ... ... ... .
[256717] chrUn_JH584304 [55512, 55701] - [GTGTCCCTGT...TTCTTCTGGC] |
[256718] chrUn_JH584304 [55740, 57151] - [GTTGTACTTT...AGTACTCAAA] |
[256719] chrUn_JH584304 [56986, 57151] - [GTTGTACTTT...CCTGAGCAGG] |
[256720] chrUn_JH584304 [58564, 58835] - [CTGTGGTCCT...CAGAGAAATG] |
[256721] chrUn_JH584304 [59592, 59689] - [TCTCTGCTGC...GCCTTCTCAG] |
exon_id tx_name gene_id
<integer> <CharacterList> <CharacterList>
[1] 1 uc007afg.1,uc007afh.1 18777
[2] 2 uc007afg.1,uc007afh.1 18777
[3] 3 uc007afg.1,uc007afh.1 18777
[4] 4 uc007afg.1,uc007afh.1 18777
[5] 5 uc007afg.1,uc007afh.1 18777
... ... ... ...
[256717] 256717 uc029xhi.1 66776
[256718] 256718 uc029xho.1 66776
[256719] 256719 uc029xhj.1 66776
[256720] 256720 uc029xhj.1,uc029xho.1 66776
[256721] 256721 uc029xhi.1,uc029xho.1 66776
-------
seqinfo: 66 sequences (1 circular) from mm10 genome
>
> subject(v)
Mouse genome:
# organism: Mus musculus (Mouse)
# provider: UCSC
# provider version: mm10
# release date: Dec. 2011
# release name: Genome Reference Consortium GRCm38
# 66 sequences:
# chr1 chr2 chr3
# chr4 chr5 chr6
# chr7 chr8 chr9
# chr10 chr11 chr12
# chr13 chr14 chr15
# ... ... ...
# chrUn_GL456372 chrUn_GL456378 chrUn_GL456379
# chrUn_GL456381 chrUn_GL456382 chrUn_GL456383
# chrUn_GL456385 chrUn_GL456387 chrUn_GL456389
# chrUn_GL456390 chrUn_GL456392 chrUn_GL456393
# chrUn_GL456394 chrUn_GL456396 chrUn_JH584304
# (use 'seqnames()' to see all the sequence names, use the '$' or '[[' operator
# to access a given sequence)
> granges(v)
GRanges object with 256721 ranges and 0 metadata columns:
seqnames ranges strand
<Rle> <IRanges> <Rle>
[1] chr1 [4807893, 4807982] +
[2] chr1 [4808455, 4808486] +
[3] chr1 [4828584, 4828649] +
[4] chr1 [4830268, 4830315] +
[5] chr1 [4832311, 4832381] +
... ... ... ...
[256717] chrUn_JH584304 [55512, 55701] -
[256718] chrUn_JH584304 [55740, 57151] -
[256719] chrUn_JH584304 [56986, 57151] -
[256720] chrUn_JH584304 [58564, 58835] -
[256721] chrUn_JH584304 [59592, 59689] -
-------
seqinfo: 66 sequences (1 circular) from mm10 genome
> seqinfo(v)
Seqinfo object with 66 sequences (1 circular) from mm10 genome:
seqnames seqlengths isCircular genome
chr1 195471971 FALSE mm10
chr2 182113224 FALSE mm10
chr3 160039680 FALSE mm10
chr4 156508116 FALSE mm10
chr5 151834684 FALSE mm10
... ... ... ...
chrUn_GL456392 23629 FALSE mm10
chrUn_GL456393 55711 FALSE mm10
chrUn_GL456394 24323 FALSE mm10
chrUn_GL456396 21240 FALSE mm10
chrUn_JH584304 114452 FALSE mm10
> as(v, "DNAStringSet")
A DNAStringSet instance of length 256721
width seq
[1] 90 GCACTGTCCGCCAGCCGGTGGATGTGCGGCA...TGTGCCGGCCGCCCGGAAGGCCACCGCCGCG
[2] 32 GTTATTTTCCTTCACGGATTGGGAGATACAGG
[3] 66 GCATGGATGGGCAGAAGCCTTTGCAGGTATC...GTCCCCACATCAAATACATCTGTCCACATGC
[4] 48 CCCTGTGATGCCAGTCACATTAAATATGAATATGGCTATGCCTTCTTG
[5] 71 GTTTGATATCGTTGGACTTTCACCAGATTCC...GAATCTGGAATTAAACAGGCAGCAGAAACCG
... ... ...
[256717] 190 GTGTCCCTGTATAAAACTGTGCCAACTCGTT...CACTACCACAACCTCAGCGAGTTCTTCTGGC
[256718] 1412 GTTGTACTTTCCCCAGCTGGCCCTATGGCGG...TGCGGGGTTGGAGTATATGAAAGTACTCAAA
[256719] 166 GTTGTACTTTCCCCAGCTGGCCCTATGGCGG...CAGCCAGGTGGGATGGAGGCCCCTGAGCAGG
[256720] 272 CTGTGGTCCTCTAACCTGTGCCTGTGAGCCA...CAGGACCGGAGCTCCAAGCAGCAGAGAAATG
[256721] 98 TCTCTGCTGCCGGAGCAAGCTCAGCTGTCCC...TGGACAGGAAGACTCTCTCCGGCCTTCTCAG
>
> v10 <- v[1:10] # select the first 10 views
> subject(v10) # same as subject(v)
Mouse genome:
# organism: Mus musculus (Mouse)
# provider: UCSC
# provider version: mm10
# release date: Dec. 2011
# release name: Genome Reference Consortium GRCm38
# 66 sequences:
# chr1 chr2 chr3
# chr4 chr5 chr6
# chr7 chr8 chr9
# chr10 chr11 chr12
# chr13 chr14 chr15
# ... ... ...
# chrUn_GL456372 chrUn_GL456378 chrUn_GL456379
# chrUn_GL456381 chrUn_GL456382 chrUn_GL456383
# chrUn_GL456385 chrUn_GL456387 chrUn_GL456389
# chrUn_GL456390 chrUn_GL456392 chrUn_GL456393
# chrUn_GL456394 chrUn_GL456396 chrUn_JH584304
# (use 'seqnames()' to see all the sequence names, use the '$' or '[[' operator
# to access a given sequence)
> granges(v10)
GRanges object with 10 ranges and 0 metadata columns:
seqnames ranges strand
<Rle> <IRanges> <Rle>
[1] chr1 [4807893, 4807982] +
[2] chr1 [4808455, 4808486] +
[3] chr1 [4828584, 4828649] +
[4] chr1 [4830268, 4830315] +
[5] chr1 [4832311, 4832381] +
[6] chr1 [4837001, 4837074] +
[7] chr1 [4839387, 4839488] +
[8] chr1 [4840956, 4841132] +
[9] chr1 [4840956, 4842827] +
[10] chr1 [4844963, 4846735] +
-------
seqinfo: 66 sequences (1 circular) from mm10 genome
> seqinfo(v10) # same as seqinfo(v)
Seqinfo object with 66 sequences (1 circular) from mm10 genome:
seqnames seqlengths isCircular genome
chr1 195471971 FALSE mm10
chr2 182113224 FALSE mm10
chr3 160039680 FALSE mm10
chr4 156508116 FALSE mm10
chr5 151834684 FALSE mm10
... ... ... ...
chrUn_GL456392 23629 FALSE mm10
chrUn_GL456393 55711 FALSE mm10
chrUn_GL456394 24323 FALSE mm10
chrUn_GL456396 21240 FALSE mm10
chrUn_JH584304 114452 FALSE mm10
> as(v10, "DNAStringSet")
A DNAStringSet instance of length 10
width seq
[1] 90 GCACTGTCCGCCAGCCGGTGGATGTGCGGCAAC...GTTGTGCCGGCCGCCCGGAAGGCCACCGCCGCG
[2] 32 GTTATTTTCCTTCACGGATTGGGAGATACAGG
[3] 66 GCATGGATGGGCAGAAGCCTTTGCAGGTATCAAAAGTCCCCACATCAAATACATCTGTCCACATGC
[4] 48 CCCTGTGATGCCAGTCACATTAAATATGAATATGGCTATGCCTTCTTG
[5] 71 GTTTGATATCGTTGGACTTTCACCAGATTCCCA...ATGAATCTGGAATTAAACAGGCAGCAGAAACCG
[6] 74 TAAAAGCCTTGATAGATCAAGAAGTGAAGAATG...TCTAACAGGATTATTTTGGGAGGATTTTCTCAG
[7] 102 GGAGGCGCCTTGTCTTTATACACTGCTCTCACC...TGCTGGCTTCCACTTCGGGCTTCGTTTTCACAG
[8] 177 GGGCCGATCAACAGTGCTAATCGAGATATTTCC...TATGAAGGCATGATGCACAGCTCATGTCAGCAG
[9] 1872 GGGCCGATCAACAGTGCTAATCGAGATATTTCC...GTTAACTGTAACGAAAAGGTGAGGTTGAAGGGT
[10] 1773 GAAATGATGGATGTCAAGCACTTCATTGATAAG...ATCAATGTGAAATAAAATTTAATTTTGGCTTTA
> alphabetFrequency(v10)
A C G T M R W S Y K V H D B N - + .
[1,] 12 36 30 12 0 0 0 0 0 0 0 0 0 0 0 0 0 0
[2,] 7 5 9 11 0 0 0 0 0 0 0 0 0 0 0 0 0 0
[3,] 20 18 14 14 0 0 0 0 0 0 0 0 0 0 0 0 0 0
[4,] 12 11 9 16 0 0 0 0 0 0 0 0 0 0 0 0 0 0
[5,] 23 14 17 17 0 0 0 0 0 0 0 0 0 0 0 0 0 0
[6,] 24 10 17 23 0 0 0 0 0 0 0 0 0 0 0 0 0 0
[7,] 18 31 24 29 0 0 0 0 0 0 0 0 0 0 0 0 0 0
[8,] 50 40 36 51 0 0 0 0 0 0 0 0 0 0 0 0 0 0
[9,] 594 332 423 523 0 0 0 0 0 0 0 0 0 0 0 0 0 0
[10,] 582 285 289 617 0 0 0 0 0 0 0 0 0 0 0 0 0 0
> alphabetFrequency(v10, collapse=TRUE)
A C G T M R W S Y K V H D B N -
1342 782 868 1313 0 0 0 0 0 0 0 0 0 0 0 0
+ .
0 0
>
> v12 <- v[width(v) <= 12] # select the views of 12 nucleotides or less
> head(as.data.frame(v12))
seqnames start end width strand exon_id tx_name gene_id
1 chr1 40008331 40008339 9 + 1052 uc007atr.... 26921
2 chr1 44119184 44119195 12 + 1244 uc007awe.2 246229
3 chr1 58393097 58393103 7 + 1840 uc007bbq.1 66882
4 chr1 87844112 87844121 10 + 3621 uc007bxs.1 20215
5 chr1 135803282 135803285 4 + 5318 uc033flq.1 21952
6 chr1 135841901 135841911 11 + 5339 uc007ctu.... 21956
dna
1 GGAGAAGTG
2 CGCAAATGGCAG
3 GGCGGGG
4 CAAAAGAAAG
5 TGAG
6 GGAACAGGAAG
> trinucleotideFrequency(v12, simplify.as="collapsed")
AAA AAC AAG AAT ACA ACC ACG ACT AGA AGC AGG AGT ATA ATC ATG ATT CAA CAC CAG CAT
139 91 176 74 90 55 29 144 139 79 177 79 96 57 168 103 98 50 129 59
CCA CCC CCG CCT CGA CGC CGG CGT CTA CTC CTG CTT GAA GAC GAG GAT GCA GCC GCG GCT
89 88 20 104 12 9 35 7 127 86 140 92 139 86 177 138 70 71 13 101
GGA GGC GGG GGT GTA GTC GTG GTT TAA TAC TAG TAT TCA TCC TCG TCT TGA TGC TGG TGT
211 74 306 107 98 48 62 67 75 125 81 125 76 98 15 115 112 67 201 48
TTA TTC TTG TTT
92 88 80 91
>
> ## 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
[1] TRUE
> v12[[2]]
12-letter "DNAString" instance
seq: CGCAAATGGCAG
> head(elementNROWS(v12)) # elementNROWS(v) is the same as width(v)
[1] 9 12 7 10 4 11
> unlist(v12)
8065-letter "DNAString" instance
seq: GGAGAAGTGCGCAAATGGCAGGGCGGGGCAAAAGAA...AAACAACAAACAATTGTCCCTGGAGAGATGGTGCCA
> elementType(v12)
[1] "DNAString"
>
>
>
>
>
> dev.off()
null device
1
>