Palantir-Methods.Rd
This function suite uses the Palantir algorithm to first calculate the diffusion map based on pre-calculated dimension reductions in Seurat (e.g., PCA, harmony), adding the diffusion map (dm) and multiscale space (ms) embeddings back to the original Seurat object. Subsequently, it calculates pseudotime to determine the developmental trajectory and fate decisions of cells based on specified starting points.
Palantir.RunDM(
seu,
reduction = "pca",
n_components = 50,
conda_env = "seuratextend"
)
Palantir.Pseudotime(
seu,
start_cell,
terminal_states = NULL,
title = "Pseudotime",
conda_env = "seuratextend"
)
A Seurat object.
The name the dimension reduction to use. Default: 'pca'.
An integer specifying how many dimensions to use as input. Default: 50.
The name of the Conda environment that contains the necessary Python libraries to execute Palantir. Default: 'seuratextend'
A vector of cell identifiers from the Seurat object that marks the starting cells for the trajectory analysis. These cells are typically identified as progenitor or early developmental states in the dataset.
Optional; a named vector of cell identifiers that represent terminal states in the developmental trajectory, with names reflecting the fate labels, such as c("fate1" = "sample1_AAACCCAAGGCCCAAA-1", "fate2" = "sample1_CCTCTAGGTGAACGGT-1"). These terminal states help Palantir more accurately model the trajectory towards cellular differentiation. Default: NULL
A label used to store and retrieve pseudotime information within the `SeuratObj@misc$Palantir` slot of the Seurat object. Default: 'Pseudotime'
The Seurat object is updated to include the diffusion map and multiscale space embeddings in the respective slots. Additionally, it adds the calculated pseudotime information to the `SeuratObj@misc$Palantir` slot of the Seurat object.
These functions integrate the powerful trajectory analysis capabilities of Palantir with the data handling and visualization strengths of Seurat. Initially, `Palantir.RunDM` calculates diffusion maps that represent cellular states and transitions in a high-dimensional space, which are then applied to visually explore cell trajectories in reduced dimensional space. `Palantir.Pseudotime` further leverages these embeddings to model developmental trajectories, allowing for a detailed understanding of cell fate decisions. This comprehensive analysis is crucial for identifying key transitions and stages in cellular development, making it an invaluable tool for developmental biology studies.
library(Seurat)
library(SeuratExtend)
# Download the example Seurat Object
mye_small <- readRDS(url("https://zenodo.org/records/10944066/files/pbmc10k_mye_small_velocyto.rds", "rb"))
# Running Diffusion Map with Palantir
mye_small <- Palantir.RunDM(mye_small)
# View the new dimensional reductions "dm" and "ms" in the Seurat Object:
print(mye_small@reductions)
# View the first 2 ms dimensions:
DimPlot2(mye_small, reduction = "ms", group.by = "cluster", label = TRUE)
# Optional: Running UMAP based on "ms" if having more than 2 dimensions
# Assuming 'ms' has 10 dimensions:
# mye_small <- RunUMAP(mye_small, reduction = "ms", dims = 1:10)
# DimPlot2(mye_small, group.by = "cluster", label = TRUE, reduction = "umap")
# Determining Cell Fates and Calculating Pseudotime. You can manually select "start cells" using the `CellSelector` function.
# p <- DimPlot(mye_small, reduction = "ms", group.by = "cluster")
# cells <- CellSelector(p)
# print(head(cells, 3))
# Calculating pseudotime and updating the Seurat object
mye_small <- Palantir.Pseudotime(mye_small, start_cell = "sample1_GAGAGGTAGCAGTACG-1")
ps <- mye_small@misc$Palantir$Pseudotime
print(head(ps))
# Visualize cell fate on UMAP
colnames(ps)[3:4] <- c("fate1", "fate2")
mye_small@meta.data[,colnames(ps)] <- ps
DimPlot2(mye_small, features = colnames(ps), reduction = "ms",
cols = list(Entropy = "D"))