The hardware and bandwidth for this mirror is donated by dogado GmbH, the Webhosting and Full Service-Cloud Provider. Check out our Wordpress Tutorial.
If you wish to report a bug, or if you are interested in having us mirror your free-software or open-source project, please feel free to contact us at mirror[@]dogado.de.
This vignette assumes that you have a SingleCellExperiment object at the ready. As an example, we create one from the ginhoux
dataset containing 248 dendritic cell progenitors.
library(SCORPIUS)
library(SingleCellExperiment)
data(ginhoux)
<- SingleCellExperiment(
sce assays = list(counts = t(ginhoux$expression)),
colData = ginhoux$sample_info
)
# short hand notation
<- colData(sce)$group_name
group_name
sce
## class: SingleCellExperiment
## dim: 2000 245
## metadata(0):
## assays(1): counts
## rownames(2000): Mpo DQ688647 ... Clcf1 Pip5k1c
## rowData names(0):
## colnames(245): SRR1558744 SRR1558745 ... SRR1558993 SRR1558994
## colData names(1): group_name
## reducedDimNames(0):
## mainExpName: NULL
## altExpNames(0):
<- reduce_dimensionality(t(assay(sce, "counts")), dist = "spearman", ndim = 3) space
The new space is a matrix that can be visualised with or without colouring of the different cell types.
draw_trajectory_plot(
space,progression_group = group_name,
contour = TRUE
)
<- infer_trajectory(space) traj
The result is a list containing the final trajectory path
and the inferred timeline for each sample time
.
draw_trajectory_plot(
space, progression_group = group_name,
path = traj$path,
contour = TRUE
)
<- gene_importances(
gimp t(assay(sce, "counts")),
$time,
trajnum_permutations = 0,
num_threads = 8
)<- gimp[1:50,]
gene_sel <- t(assay(sce, "counts"))[,gene_sel$gene] expr_sel
To visualise the expression of the selected genes, use the draw_trajectory_heatmap
function.
draw_trajectory_heatmap(expr_sel, traj$time, group_name)
Finally, these genes can also be grouped into modules as follows:
<- extract_modules(scale_quantile(expr_sel), traj$time, verbose = FALSE)
modules draw_trajectory_heatmap(expr_sel, traj$time, group_name, modules)
reducedDims(sce) <- SimpleList(MDS = space)
colData(sce)$trajectory_path <- traj$path
colData(sce)$trajectory_pseudotime <- traj$time
rowData(sce)$trajectory_importance <- gimp[match(rownames(sce), gimp$gene),]$importance
sce
## class: SingleCellExperiment
## dim: 2000 245
## metadata(0):
## assays(1): counts
## rownames(2000): Mpo DQ688647 ... Clcf1 Pip5k1c
## rowData names(1): trajectory_importance
## colnames(245): SRR1558744 SRR1558745 ... SRR1558993 SRR1558994
## colData names(3): group_name trajectory_path trajectory_pseudotime
## reducedDimNames(1): MDS
## mainExpName: NULL
## altExpNames(0):
These binaries (installable software) and packages are in development.
They may not be fully stable and should be used with caution. We make no claims about them.
Health stats visible at Monitor.