getAdjustedMat | R Documentation |
Get Adjusted Matrix with scMerge2 parameter estimated
getAdjustedMat(
exprsMat,
fullalpha,
ctl = rownames(exprsMat),
adjusted_means = NULL,
ruvK = 20,
return_subset_genes = NULL
)
exprsMat |
A gene (row) by cell (column) matrix to be adjusted. |
fullalpha |
A matrix indicates the estimated alpha returned by |
ctl |
A character vector of negative control. It should have a non-empty intersection with the rows of exprsMat. |
adjusted_means |
A rowwise mean of the gene by cell matrix |
ruvK |
An integer indicates the number of unwanted variation factors that are removed, default is 20. |
return_subset_genes |
An optional character vector of indicates the subset of genes will be adjusted. |
Returns the adjusted matrix will be return.
Yingxin Lin
## Loading example data
data('example_sce', package = 'scMerge')
## Previously computed stably expressed genes
data('segList_ensemblGeneID', package = 'scMerge')
## Running an example data with minimal inputs
library(SingleCellExperiment)
scMerge2_res <- scMerge2(exprsMat = logcounts(example_sce),
batch = example_sce$batch,
ctl = segList_ensemblGeneID$mouse$mouse_scSEG,
return_matrix = FALSE)
cosineNorm_mat <- batchelor::cosineNorm(logcounts(example_sce))
adjusted_means <- DelayedMatrixStats::rowMeans2(cosineNorm_mat)
newY <- getAdjustedMat(cosineNorm_mat, scMerge2_res$fullalpha,
ctl = segList_ensemblGeneID$mouse$mouse_scSEG,
ruvK = 20,
adjusted_means = adjusted_means)
assay(example_sce, "scMerge2") <- newY
example_sce = scater::runPCA(example_sce, exprs_values = 'scMerge2')
scater::plotPCA(example_sce, colour_by = 'cellTypes', shape = 'batch')
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.