% Generated by roxygen2: do not edit by hand % Please edit documentation in R/runDEAnalysis.R \name{runDEAnalysis} \alias{runDEAnalysis} \alias{runDESeq2} \alias{runLimmaDE} \alias{runANOVA} \alias{runMAST} \alias{runWilcox} \title{Perform differential expression analysis on SCE object} \usage{ runDEAnalysis(inSCE, method = "wilcox", ...) runDESeq2( inSCE, useAssay = "counts", useReducedDim = NULL, index1 = NULL, index2 = NULL, class = NULL, classGroup1 = NULL, classGroup2 = NULL, analysisName, groupName1, groupName2, covariates = NULL, fullReduced = TRUE, onlyPos = FALSE, log2fcThreshold = NULL, fdrThreshold = NULL, minGroup1MeanExp = NULL, maxGroup2MeanExp = NULL, minGroup1ExprPerc = NULL, maxGroup2ExprPerc = NULL, overwrite = FALSE, verbose = TRUE ) runLimmaDE( inSCE, useAssay = "logcounts", useReducedDim = NULL, index1 = NULL, index2 = NULL, class = NULL, classGroup1 = NULL, classGroup2 = NULL, analysisName, groupName1, groupName2, covariates = NULL, onlyPos = FALSE, log2fcThreshold = NULL, fdrThreshold = NULL, minGroup1MeanExp = NULL, maxGroup2MeanExp = NULL, minGroup1ExprPerc = NULL, maxGroup2ExprPerc = NULL, overwrite = FALSE, verbose = TRUE ) runANOVA( inSCE, useAssay = "logcounts", useReducedDim = NULL, index1 = NULL, index2 = NULL, class = NULL, classGroup1 = NULL, classGroup2 = NULL, analysisName, groupName1, groupName2, covariates = NULL, onlyPos = FALSE, log2fcThreshold = NULL, fdrThreshold = NULL, minGroup1MeanExp = NULL, maxGroup2MeanExp = NULL, minGroup1ExprPerc = NULL, maxGroup2ExprPerc = NULL, overwrite = FALSE, verbose = TRUE ) runMAST( inSCE, useAssay = "logcounts", useReducedDim = NULL, index1 = NULL, index2 = NULL, class = NULL, classGroup1 = NULL, classGroup2 = NULL, analysisName, groupName1, groupName2, covariates = NULL, onlyPos = FALSE, log2fcThreshold = NULL, fdrThreshold = NULL, minGroup1MeanExp = NULL, maxGroup2MeanExp = NULL, minGroup1ExprPerc = NULL, maxGroup2ExprPerc = NULL, overwrite = FALSE, check_sanity = TRUE, verbose = TRUE ) runWilcox( inSCE, useAssay = "logcounts", useReducedDim = NULL, index1 = NULL, index2 = NULL, class = "cluster", classGroup1 = c(1), classGroup2 = c(2), analysisName = "cluster1_VS_2", groupName1 = "cluster1", groupName2 = "cluster2", covariates = NULL, onlyPos = FALSE, log2fcThreshold = NULL, fdrThreshold = NULL, minGroup1MeanExp = NULL, maxGroup2MeanExp = NULL, minGroup1ExprPerc = NULL, maxGroup2ExprPerc = NULL, overwrite = FALSE, verbose = TRUE ) } \arguments{ \item{inSCE}{\linkS4class{SingleCellExperiment} inherited object.} \item{method}{Character. Specify which method to use when using \code{runDEAnalysis()}. Choose from \code{"wilcox"}, \code{"MAST"}, \code{"DESeq2"}, \code{"Limma"}, \code{"ANOVA"}. Default \code{"wilcox"}.} \item{...}{Arguments to pass to specific methods when using the generic \code{runDEAnalysis()}.} \item{useAssay}{character. A string specifying which assay to use for the DE regression. Ignored when \code{useReducedDim} is specified. Default \code{"counts"} for DESeq2, \code{"logcounts"} for other methods.} \item{useReducedDim}{character. A string specifying which reducedDim to use for DE analysis. Will treat the dimensions as features. Default \code{NULL}.} \item{index1}{Any type of indices that can subset a \linkS4class{SingleCellExperiment} inherited object by cells. Specifies which cells are of interests. Default \code{NULL}.} \item{index2}{Any type of indices that can subset a \linkS4class{SingleCellExperiment} inherited object by cells. specifies the control group against those specified by \code{index1}. If \code{NULL} when using index specification, \code{index1} cells will be compared with all other cells. Default \code{NULL}.} \item{class}{A vector/factor with \code{ncol(inSCE)} elements, or a character scalar that specifies a column name of \code{colData(inSCE)}. Default \code{"cluster"}.} \item{classGroup1}{a vector specifying which "levels" given in \code{class} are of interests. Default \code{c(1)}.} \item{classGroup2}{a vector specifying which "levels" given in \code{class} is the control group against those specified by \code{classGroup1}. If \code{NULL} when using annotation specification, \code{classGroup1} cells will be compared with all other cells. Default \code{c(2)}.} \item{analysisName}{A character scalar naming the DEG analysis. Default \code{"cluster1_VS_2"}.} \item{groupName1}{A character scalar naming the group of interests. Default \code{"cluster1"}.} \item{groupName2}{A character scalar naming the control group. Default \code{"cluster2"}.} \item{covariates}{A character vector of additional covariates to use when building the model. All covariates must exist in \code{names(colData(inSCE))}. Default \code{NULL}.} \item{fullReduced}{Logical, DESeq2 only argument. Whether to apply LRT (Likelihood ratio test) with a 'full' model. Default \code{TRUE}.} \item{onlyPos}{Whether to only output DEG with positive log2_FC value. Default \code{FALSE}.} \item{log2fcThreshold}{Only out put DEGs with the absolute values of log2FC greater than this value. Default \code{NULL}.} \item{fdrThreshold}{Only out put DEGs with FDR value less than this value. Default \code{NULL}.} \item{minGroup1MeanExp}{Only out put DEGs with mean expression in group1 greater then this value. Default \code{NULL}.} \item{maxGroup2MeanExp}{Only out put DEGs with mean expression in group2 less then this value. Default \code{NULL}.} \item{minGroup1ExprPerc}{Only out put DEGs expressed in greater then this fraction of cells in group1. Default \code{NULL}.} \item{maxGroup2ExprPerc}{Only out put DEGs expressed in less then this fraction of cells in group2. Default \code{NULL}.} \item{overwrite}{A logical scalar. Whether to overwrite result if exists. Default \code{FALSE}.} \item{verbose}{A logical scalar. Whether to show messages. Default \code{TRUE}.} \item{check_sanity}{Logical, MAST only argument. Whether to perform MAST's sanity check to see if the counts are logged. Default \code{TRUE}.} } \value{ The input \linkS4class{SingleCellExperiment} object, where \code{metadata(inSCE)$diffExp} is updated with a list named by \code{analysisName}, with elements of: \item{$groupNames}{the naming of the two conditions} \item{$useAssay, $useReducedDim}{the matrix name that was used for calculation} \item{$select}{the cell selection indices (logical) for each condition} \item{$result}{a \code{data.frame} of the DEGs table} \item{$method}{the method used} } \description{ Perform differential expression analysis on SCE object } \details{ SCTK provides Limma, MAST, DESeq2, ANOVA and Wilcoxon test for differential expression analysis, where DESeq2 expects non-negtive integer assay input while others expect logcounts. Condition specification allows two methods: 1. Index level selection. Only use arguments \code{index1} and \code{index2}. 2. Annotation level selection. Only use arguments \code{class}, \code{classGroup1} and \code{classGroup2}. } \examples{ data(scExample, package = "singleCellTK") sce <- subsetSCECols(sce, colData = "type != 'EmptyDroplet'") sce <- scaterlogNormCounts(sce, assayName = "logcounts") sce <- runDEAnalysis(method = "Limma", inSCE = sce, groupName1 = "group1", groupName2 = "group2", index1 = seq(20), index2 = seq(21,40), analysisName = "Limma") } \seealso{ See \code{\link{plotDEGHeatmap}}, \code{\link{plotDEGRegression}}, \code{\link{plotDEGViolin}} and \code{\link{plotDEGVolcano}} for visualization method after running DE analysis. }