Project Pseudotime and Branch Probabilities to Single Cells
projectPseudotimeToCell.Rd
This function projects pseudotime and branch probabilities from pseudobulk (pb.milo
)
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 = ""
)
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
A 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 (''
).
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)
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 ...
#> 7279 of cells filtered
#> checking allowed chain status...
#> 12 of cells filtered
#> VDJ data extraction begin:
#> Parameter extract_cols do not provided, automatically geneterate colnames for extraction.
#> Detect whether colData v_call_abT_VDJ, d_call_abT_VDJ, j_call_abT_VDJ, v_call_abT_VJ, j_call_abT_VJ already exist...
#> 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 ...
#> 63 of cells filtered
#> 2646 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")
# Define root and branch tips
pca <- t(as.matrix(SingleCellExperiment::reducedDim(pb.milo, type = "PCA")))
branch.tips <- c(189, 198) # which.min(pca[, 2]) and which.max(pca[, 2])
names(branch.tips) <- c("CD8+T", "CD4+T")
root <- 177 # which.max(pca[, 1])
# Construct Diffusion Map
dm <- destiny::DiffusionMap(t(pca), n_pcs = 50, n_eigs = 10)
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,
terminal_state = branch.tips,
diffusiontime = dif.pse[[paste0("DPT", root)]],
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...
# Project Pseudobulk Data
projected_milo <- projectPseudotimeToCell(
milo_object,
pb.milo,
branch.tips,
pseudotime_key = "pseudotime"
)
#> 3 number of cells removed due to not belonging to any neighbourhood