tagClustersGR | R Documentation |
Extracts tag clusters (TCs) produced by clusterCTSS
function
for a specified CAGE experiment from a CAGEexp
object.
tagClustersGR(
object,
sample = NULL,
returnInterquantileWidth = FALSE,
qLow = NULL,
qUp = NULL
)
## S4 method for signature 'CAGEexp'
tagClustersGR(
object,
sample = NULL,
returnInterquantileWidth = FALSE,
qLow = NULL,
qUp = NULL
)
tagClustersGR(object, sample = NULL) <- value
## S4 replacement method for signature 'CAGEexp,ANY,TagClusters'
tagClustersGR(object, sample = NULL) <- value
## S4 replacement method for signature 'CAGEexp,missing,GRangesList'
tagClustersGR(object, sample = NULL) <- value
object |
A |
sample |
Label of the CAGE dataset (experiment, sample) for which to
extract tag clusters. If |
returnInterquantileWidth |
Return the interquantile width for each tag cluster. |
qLow, qUp |
Position of which quantile should be used as a left (lower)
or right (upper) boundary (for |
value |
A |
Returns a GRangesList
or a GRanges
object with genomic coordinates,
position of dominant TSS, total CAGE signal and additional information for
all TCs from specified CAGE dataset (sample). If
returnInterquantileWidth = TRUE
, interquantile width for each TC is also
calculated using provided quantile positions. The S4Vectors::metadata
slot of the object contains a copy of the CAGEexp
object's column data,
as well as information on the clustering method in a clusteringMethod
element.
Vanja Haberle
Charles Plessy
Other CAGEr accessor methods:
CTSSclusteringMethod()
,
CTSScoordinatesGR()
,
CTSScumulativesTagClusters()
,
CTSSnormalizedTpmDF()
,
CTSStagCountDF()
,
GeneExpDESeq2()
,
GeneExpSE()
,
consensusClustersGR()
,
expressionClasses()
,
genomeName()
,
inputFilesType()
,
inputFiles()
,
librarySizes()
,
sampleLabels()
,
seqNameTotalsSE()
Other CAGEr clusters functions:
CTSSclusteringMethod()
,
CTSScumulativesTagClusters()
,
CustomConsensusClusters()
,
aggregateTagClusters()
,
clusterCTSS()
,
consensusClustersDESeq2()
,
consensusClustersGR()
,
cumulativeCTSSdistribution()
,
plotInterquantileWidth()
,
quantilePositions()
tagClustersGR( exampleCAGEexp, "Zf.high", TRUE, 0.1, 0.9 )
tagClustersGR( exampleCAGEexp, 1
, returnInterquantileWidth = TRUE, qLow = 0.1, qUp = 0.9 )
tagClustersGR( exampleCAGEexp )@metadata$colData
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.