The alignmentData class inherits from the
alignmentClass class and records information about
a set of alignments of high-throughput sequencing data to a
genome. Details include the alignments themselves, the
chromosomes of the genome to which the data are aligned, and
counts of the aligned tags from each of the libraries from which the
data come.
Objects from the class
Objects can be created by calls of the form new("alignmentData",
...), but more usually by using one of readBAM or
readGeneric functions to generate the object from a set of
alignment files.
Slots
alignments:
Object of class "GRanges".
Stores information about the alignments. See Details.
replicates:
Object of class "factor".
Replicate information for each of the libraries. See Details.
data:
Object of class "matrix". For each
alignment described in the alignments slot, contains the
number of times the alignment is seen in each sample.
libnames:
Object of class "character". The
names of the libraries for which alignment data exists.
libsizes:
Object of class "numeric". The
library sizes (see Details) for each of the libraries.
Details
The alignments slot is the key element of this class. This is a
GRanges object that, in addition to the usual elements defining
the location of aligned objects to a reference genome, also describes
the values ‘tag’, giving the sequence of the tag aligning to the
location, ‘matches’, indicating in how many places that tag matches to
the genome, ‘chunk’, an identifier for the sets of tags that align
close enough together to form a potential locus, and ‘chunkDup’,
indicating whether that tag matches to multiple places within the chunk.
The library sizes, defined in the libsizes slot, provide some
scaling factor for the observed number of counts of a tag in different
samples.
The replicates slot is a vector of factors such that the ith
sample is a replicate of the jth sample if and only if @replicates[i] ==
@replicates[j].
Methods
[
signature(x = "alignmentData"): ...
dim
signature(x = "alignmentData"): ...
initialize
signature(.Object = "alignmentData"): ...
show
signature(object = "alignmentData"): ...
Author(s)
Thomas J. Hardcastle
See Also
alignmentClass, the class from which
'alignmentData' inherits.
readGeneric, which will produce a 'alignmentData'
object from appropriately formatted tab-delimited files.
readBAM, which will produce a 'alignmentData'
object from BAM files.
processAD, which will convert an 'alignmentData'
object into a 'segData' object for segmentation.
Examples
# Define the chromosome lengths for the genome of interest.
chrlens <- c(2e6, 1e6)
# Define the files containing sample information.
datadir <- system.file("extdata", package = "segmentSeq")
libfiles <- c("SL9.txt", "SL10.txt", "SL26.txt", "SL32.txt")
# Establish the library names and replicate structure.
libnames <- c("SL9", "SL10", "SL26", "SL32")
replicates <- c(1,1,2,2)
# Process the files to produce an 'alignmentData' object.
alignData <- readGeneric(file = libfiles, dir = datadir, replicates =
replicates, libnames = libnames, chrs = c(">Chr1", ">Chr2"), chrlens =
chrlens)
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(segmentSeq)
Loading required package: baySeq
Loading required package: GenomicRanges
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: abind
Loading required package: perm
Loading required package: ShortRead
Loading required package: BiocParallel
Loading required package: Biostrings
Loading required package: XVector
Loading required package: Rsamtools
Loading required package: GenomicAlignments
Loading required package: SummarizedExperiment
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")'.
> png(filename="/home/ddbj/snapshot/RGM3/R_BC/result/segmentSeq/alignmentData-class.Rd_%03d_medium.png", width=480, height=480)
> ### Name: alignmentData-class
> ### Title: Class "alignmentData"
> ### Aliases: alignmentData-class alignmentData [,alignmentData-method
> ### [,alignmentData,ANY,ANY-method [,alignmentData,ANY-method
> ### dim,alignmentData-method initialize,alignmentData-method
> ### show,alignmentData-method cbind,alignmentData-method
> ### Keywords: classes
>
> ### ** Examples
>
>
> # Define the chromosome lengths for the genome of interest.
>
> chrlens <- c(2e6, 1e6)
>
> # Define the files containing sample information.
>
> datadir <- system.file("extdata", package = "segmentSeq")
> libfiles <- c("SL9.txt", "SL10.txt", "SL26.txt", "SL32.txt")
>
> # Establish the library names and replicate structure.
>
> libnames <- c("SL9", "SL10", "SL26", "SL32")
> replicates <- c(1,1,2,2)
>
> # Process the files to produce an 'alignmentData' object.
>
> alignData <- readGeneric(file = libfiles, dir = datadir, replicates =
+ replicates, libnames = libnames, chrs = c(">Chr1", ">Chr2"), chrlens =
+ chrlens)
Reading files........done!
Analysing tags...........done!
>
>
>
>
>
>
> dev.off()
null device
1
>