plotNhoodCounts | R Documentation |
Plot the number of cells in a neighbourhood per sample and condition
plotNhoodCounts(x, subset.nhoods, design.df, condition, n_col = 3)
x |
A |
subset.nhoods |
A logical, integer or character vector indicating the rows of |
design.df |
A |
condition |
String specifying the condition of interest Has to be a column in the |
n_col |
Number of columns in the output |
A ggplot-class
object
Nick Hirschmüller
require(SingleCellExperiment)
ux.1 <- matrix(rpois(12000, 5), ncol=300)
ux.2 <- matrix(rpois(12000, 4), ncol=300)
ux <- rbind(ux.1, ux.2)
vx <- log2(ux + 1)
pca <- prcomp(t(vx))
sce <- SingleCellExperiment(assays=list(counts=ux, logcounts=vx),
reducedDims=SimpleList(PCA=pca$x))
milo <- Milo(sce)
milo <- buildGraph(milo, k=20, d=10, transposed=TRUE)
milo <- makeNhoods(milo, k=20, d=10, prop=0.3)
milo <- calcNhoodDistance(milo, d=10)
cond <- sample(c("A","B","C"),300,replace=TRUE)
meta.df <- data.frame(Condition=cond, Replicate=c(rep("R1", 100), rep("R2", 100), rep("R3", 100)))
meta.df$SampID <- paste(meta.df$Condition, meta.df$Replicate, sep="_")
milo <- countCells(milo, meta.data=meta.df, samples="SampID")
design.mtx <- data.frame("Condition"=c(rep("A", 3), rep("B", 3), rep("C",3)),
"Replicate"=rep(c("R1", "R2", "R3"), 3))
design.mtx$SampID <- paste(design.mtx$Condition, design.mtx$Replicate, sep="_")
rownames(design.mtx) <- design.mtx$SampID
plotNhoodCounts(x = milo,
subset.nhoods = c(1,2),
design.df = design.mtx,
condition = "Condition")
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.