plotPCA | R Documentation |
Compute PCA of gene expression for an assay, and plot samples coloring by outlier score
## S4 method for signature 'list'
plotPCA(
object,
assays = names(object),
nPC = 2,
robust = FALSE,
...,
maxOutlierZ = 20,
nrow = 2,
size = 2,
fdr.cutoff = 0.05
)
object |
|
assays |
assays / cell types to analyze |
nPC |
number of PCs to uses for outlier score with |
robust |
use robust covariance method, defaults to |
... |
arguments passed to |
maxOutlierZ |
cap outlier z-scores at this value for plotting to maintain consistent color scale |
nrow |
number of rows in plot |
size |
size passed to |
fdr.cutoff |
FDR cutoff to determine outlier |
outlierByAssay()
library(muscat)
library(SingleCellExperiment)
data(example_sce)
# create pseudobulk for each sample and cell cluster
pb <- aggregateToPseudoBulk(example_sce,
assay = "counts",
cluster_id = "cluster_id",
sample_id = "sample_id",
verbose = FALSE
)
# voom-style normalization
res.proc <- processAssays(pb, ~group_id)
# PCA to identify outliers
# from normalized expression
plotPCA( res.proc, c("B cells", "CD14+ Monocytes"))
# Run on regression residuals
#-----------------------------
# Regression analysis
fit = dreamlet(res.proc, ~ group_id)
# Extract regression residuals
residsObj = residuals(fit)
# PCA on residuals
plotPCA( residsObj, c("B cells", "CD14+ Monocytes"))
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.