mergeCAGEsets | R Documentation |
Merges two CAGEr
objects into one by combining the CTSS genomic
coordinates and raw tag counts. The resulting object will contain a union
of TSS positions present in the two input objects and raw tag counts for
those TSSs in all samples from both input objects.
mergeCAGEsets(cs1, cs2)
## S4 method for signature 'CAGEexp,CAGEexp'
mergeCAGEsets(cs1, cs2)
cs1 |
A |
cs2 |
A |
Note that merging discards all other information present in the
two CAGEr
objects, that is, the merged object will not contain any
normalised tag counts, CTSS clusters, quantile positions, etc., so these
have to be calculated again by calling the appropriate functions on the
merged object. Also, it is only possible to merge two objects that contain
TSS information for the same reference genome and do not share any sample
names.
Returns a CAGEexp
object, which contains a union of
TSS positions present in the two input objects and raw tag counts for those
TSSs in all samples from both input objects.
Vanja Haberle
Charles Plessy
CAGEexp
library(BSgenome.Drerio.UCSC.danRer7)
pathsToInputFiles <- system.file("extdata", c("Zf.unfertilized.egg.chr17.ctss",
"Zf.30p.dome.chr17.ctss", "Zf.prim6.rep1.chr17.ctss"), package="CAGEr")
ce1 <- CAGEexp(genomeName = "BSgenome.Drerio.UCSC.danRer7",
inputFiles = pathsToInputFiles[1:2], inputFilesType = "ctss", sampleLabels =
c("sample1", "sample2"))
ce1 <- getCTSS(ce1)
ce2 <- CAGEexp(genomeName = "BSgenome.Drerio.UCSC.danRer7",
inputFiles = pathsToInputFiles[3], inputFilesType = "ctss", sampleLabels =
"sample3")
ce2 <- getCTSS(ce2)
ce <- mergeCAGEsets(ce1, ce2)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.