R/runTSCAN.R
plotTSCANPseudotimeHeatmap.RdA wrapper function which visualizes outputs from the
runTSCANDEG function. Plots the top genes that change in
expression with increasing pseudotime along the path in the MST.
runTSCANDEG has to be run in advance with using the same
pathIndex of interest.
plotTSCANPseudotimeHeatmap(
inSCE,
pathIndex,
direction = c("both", "increasing", "decreasing"),
topN = 50,
log2fcThreshold = NULL,
useAssay = NULL,
featureDisplay = metadata(inSCE)$featureDisplay
)Input SingleCellExperiment object.
Path index for which the pseudotime values should be used.
Should have being used in runTSCANDEG.
Should we show features with expression increasing or
decreeasing along the increase in TSCAN pseudotime? Choices are
"both", "increasing" or "decreasing".
An integer. Only to plot this number of top genes along the path
in the MST, in terms of FDR value. Use NULL to cancel the top N
subscription. Default 30.
Only output DEGs with the absolute values of log2FC
larger than this value. Default NULL.
A single character to specify a feature expression matrix in
assays slot. The expression of top features from here will be
visualized. Default NULL use the one used for
runTSCANDEG.
Whether to display feature ID and what ID type to
display. Users can set default ID type by setSCTKDisplayRow.
NULL will display when number of features to display is less than 60.
FALSE for no display. Variable name in rowData to indicate ID
type. "rownames" or TRUE for using rownames(inSCE).
A ComplexHeatmap in .ggplot class
data("mouseBrainSubsetSCE", package = "singleCellTK")
mouseBrainSubsetSCE <- runTSCAN(inSCE = mouseBrainSubsetSCE,
useReducedDim = "PCA_logcounts")
#> Tue Jun 28 22:06:02 2022 ... Running 'scran SNN clustering' with 'louvain' algorithm
#> Tue Jun 28 22:06:03 2022 ... Identified 2 clusters
#> Tue Jun 28 22:06:03 2022 ... Running TSCAN to estimate pseudotime
#> Tue Jun 28 22:06:03 2022 ... Clusters involved in path index 2 are: 1, 2
#> Tue Jun 28 22:06:03 2022 ... Number of estimated paths is 1
terminalNodes <- listTSCANTerminalNodes(mouseBrainSubsetSCE)
mouseBrainSubsetSCE <- runTSCANDEG(inSCE = mouseBrainSubsetSCE,
pathIndex = terminalNodes[1])
plotTSCANPseudotimeHeatmap(mouseBrainSubsetSCE,
pathIndex = terminalNodes[1])