GeneTrendHeatmap-Methods.Rd
This function suite provides tools for visualizing gene expression dynamics along pseudotime trajectories using heatmaps, leveraging results from either Palantir or Slingshot pseudotime analyses. GeneTrendHeatmap.Palantir() and GeneTrendHeatmap.Slingshot() each cater to the specific formatting and requirements of the pseudotime data generated by their respective tools.
GeneTrendHeatmap.Palantir(
seu,
features,
pseudotime.data = "Pseudotime",
lineage,
magic = FALSE,
scale = TRUE,
sort = TRUE,
color_scheme = "D",
conda_env = "seuratextend"
)
GeneTrendHeatmap.Slingshot(
seu,
features,
pseudotime.data = "PCA",
lineage,
magic = FALSE,
scale = TRUE,
sort = TRUE,
color_scheme = "D",
conda_env = "seuratextend"
)
A Seurat object that has undergone pseudotime analysis with either Palantir or Slingshot, containing necessary gene expression data and pseudotime calculations.
A vector of gene names for which expression dynamics will be visualized in the heatmap.
The name of the slot within `seu@misc` where pseudotime data is stored or a data.frame containing pseudotime information; for Palantir, this is typically 'Pseudotime', and for Slingshot, 'PCA'.
The specific lineage or trajectory to visualize in the heatmap, which corresponds to the clusters or trajectories defined within the pseudotime analysis.
Logical indicating whether to use denoised gene expression data for plotting, obtained via MAGIC. Requires prior execution of Palantir.Magic() if set to TRUE. Default: FALSE.
Boolean indicating whether to scale gene expression data between 0 and 1, which can help in normalizing expression levels for better visual contrast in the heatmap. Default: TRUE.
Boolean indicating whether to sort genes by expression peak along pseudotime from left to right. Default: TRUE.
The color scheme to use for the heatmap. Default: 'D'.
Name of the Conda environment where necessary Python dependencies are installed for operations that require Python. Default: 'seuratextend'.
A ggplot object, displaying gene expression trends along computed pseudotime trajectories by heatmap.
These visualization tools are designed to highlight gene expression changes over pseudotime, providing insights into cellular behaviors and gene regulation during differentiation processes. Each function is tailored to the specific pseudotime analysis, ensuring that visualizations are closely aligned with the underlying data and pseudotime results.
library(Seurat)
library(SeuratExtend)
# Load an example Seurat Object
mye_small <- readRDS(url("https://zenodo.org/records/10944066/files/pbmc10k_mye_small_velocyto.rds", "rb"))
# Calculate diffusion map and pseudotime using Palantir
mye_small <- Palantir.RunDM(mye_small)
mye_small <- Palantir.Pseudotime(mye_small, start_cell = "sample1_GAGAGGTAGCAGTACG-1")
# Customizing gene trend heatmap with Palantir
ps <- mye_small@misc$Palantir$Pseudotime
colnames(ps)[3:4] <- c("fate1", "fate2")
GeneTrendHeatmap.Palantir(
mye_small,
features = c("CD14", VariableFeatures(mye_small)[1:10]),
pseudotime.data = ps,
lineage = "fate1"
)
# Run Slingshot
mye_small <- RunSlingshot(mye_small, group.by = "cluster", start.clus = "Mono CD14")
# Using Slingshot for similar visualizations
GeneTrendHeatmap.Slingshot(
mye_small,
features = c("CD14", VariableFeatures(mye_small)[1:10]),
lineage = "slingPseudotime_2"
)