Project Pseudotime and Branch Probabilities to Single Cells
projectPseudotimeToCell.Rd
This function projects pseudotime and branch probabilities from pseudobulk
data to single-cell resolution (milo
). The results are stored in the
colData
of the milo
object.
Usage
projectPseudotimeToCell(
milo,
pb_milo,
term_states = NULL,
pseudotime_key = "pseudotime",
suffix = "",
verbose = TRUE
)
Arguments
- milo
A
SingleCellExperiment
orMilo
object. Represents single-cell data where pseudotime and branch probabilities will be projected.- pb_milo
A pseudobulk
Milo
object. Contains aggregated branch probabilities and pseudotime information to be transferred to single cells.- term_states
Named vector of terminal states, with branch probabilities to be transferred. The names should correspond to branches of interest.
- pseudotime_key
Character. The column name in
colData
ofpb_milo
that contains the pseudotime information which was used in themarkovProbability
function. Default is"pseudotime"
.- suffix
Character. A suffix to be added to the new column names in
colData
. Default is an empty string (''
).- verbose
Boolean, whether to print messages/warnings.
Value
subset of milo or SingleCellExperiment object where cell that do not belong to any neighbourhood are removed and projected pseudotime information stored colData
Examples
data(sce_vdj)
# downsample to first 2000 cells
sce_vdj <- sce_vdj[, 1:2000]
sce_vdj <- setupVdjPseudobulk(sce_vdj,
already.productive = FALSE,
allowed_chain_status = c("Single pair", "Extra pair")
)
#> Checking productivity from productive_abT_VDJ, productive_abT_VJ ...
#> 1461 of cells filtered
#> checking allowed chains...
#> 3 of cells filtered
#> VDJ data extraction begin:
#> extract_cols not specified, automatically generate colnames for extraction.
#> Extract main TCR from v_call_abT_VDJ, d_call_abT_VDJ, j_call_abT_VDJ, v_call_abT_VJ, j_call_abT_VJ ...
#> Complete.
#> Filtering cells from v_call_abT_VDJ_main, j_call_abT_VDJ_main, v_call_abT_VJ_main, j_call_abT_VJ_main ...
#> 11 of cells filtered
#> 525 of cells remain.
# Build Milo Object
set.seed(100)
milo_object <- miloR::Milo(sce_vdj)
milo_object <- miloR::buildGraph(milo_object,
k = 50, d = 20,
reduced.dim = "X_scvi"
)
#> Constructing kNN graph with k:50
milo_object <- miloR::makeNhoods(milo_object,
reduced_dims = "X_scvi",
d = 20
)
#> Checking valid object
#> Running refined sampling with reduced_dim
# Construct Pseudobulked VDJ Feature Space
pb.milo <- vdjPseudobulk(milo_object, col_to_take = "anno_lvl_2_final_clean")
pb.milo <- scater::runPCA(pb.milo, assay.type = "Feature_space")
#> Warning: more singular values/vectors requested than available
#> Warning: You're computing too large a percentage of total singular values, use a standard svd instead.
# Define root and branch tips
pca <- t(as.matrix(SingleCellExperiment::reducedDim(pb.milo, type = "PCA")))
branch.tips <- c(which.min(pca[, 2]), which.max(pca[, 2]))
names(branch.tips) <- c("CD8+T", "CD4+T")
root <- which.min(pca[, 1])
# Construct Diffusion Map
dm <- destiny::DiffusionMap(t(pca), n_pcs = 10, n_eigs = 5)
dif.pse <- destiny::DPT(dm, tips = c(root, branch.tips), w_width = 0.1)
# Markov Chain Construction
pb.milo <- markovProbability(
milo = pb.milo,
diffusionmap = dm,
diffusiontime = dif.pse[[paste0("DPT", root)]],
terminal_state = branch.tips,
root_cell = root,
pseudotime_key = "pseudotime"
)
#> Sampling and flocking waypoints...
#> Markov chain construction...
#> Computing fundamental matrix and absorption probabilities...
#> Project probabilites from waypoints to each pseudobulk...
#> Complete.
# Project Pseudobulk Data
projected_milo <- projectPseudotimeToCell(
milo_object,
pb.milo,
branch.tips,
pseudotime_key = "pseudotime"
)
#> 0 number of cells removed due to not belonging to any neighbourhood