CalcDEcombn: Performs DE testing between pairs of clusters in sCVdata

CalcDEcombnR Documentation

Performs DE testing between pairs of clusters in sCVdata

Description

Performs differential gene expression tests between each pairwise combination of cluster in an sCVdata object using the gene expression matrix of input data object. Alternatively, this function can be skipped, and existing DE test results can be assigned directly to the sCVdata object.

Usage

CalcDEcombn(sCVd, inD)

## S4 method for signature 'sCVdata'
CalcDEcombn(sCVd, inD)

Arguments

sCVd

An sCVdata object.

inD

The input dataset. An object of class seurat or SingleCellExperiment. Other data classes are not currently supported. Please submit requests for other data objects here!

Details

This function performs Wilcoxon rank sum tests comparing gene expression between the cells of all pairwise combinations of cluster clusters in the input data. Gene expression ratio in log space (logGER) and differences in detection rates (dDR) are reported for all genes in the comparison. Genes are tested if they are detected in at least one of the cluster at a higher proportion than Param(sCVd,"DRthresh"), and both unadjusted p-values and false discovery rates are reported for all genes tested. To help track its progress, this function uses progress bars from pbapply. To disable these, set pboptions(type="none"). To re-enable, set pboptions(type="timer").

If using existing DE test results, assign results of differential gene expression tests for all pairwise combinations of clusters in sCVdata to the DEcombn slot of the sCVdata object. See example and slot documentation.

Value

A named list of data frames, one entry for each pairwise combination of levels in Clusters(sCVd) (with corresponding name where levels are separated by '-'). Each entry is data frame containing gene differential expression stats when comparing the cells of that cluster to all other cells in the input data. Rows represent genes, and variables include logGER (an effect size measure: gene expression ratio in log space, often referred to as logFC), dDR (an effect size measure: difference in detection rate), and FDR (significance measure: false discovery rate). Also included are Wstat and pVal, the test statistic and the p-value of the Wilcoxon rank sum test.

Methods (by class)

  • sCVdata: Calculate DE between cluster pairs

See Also

CalcSCV for wrapper function to calculate all statistics for an sCVdata object, and fx_calcEScombn and fx_calcDEcombn for the internal functions performing the calculations. Wilcox test is now powered by wilcoxauc for super speed.

Examples

## Not run: 
## Example using CalcDEvsRest ##
DEcombn(your_sCV_obj) <- CalcDEcombn(sCVd=your_sCV_obj,
                                     inD=your_scRNAseq_data_object)


## Example using MAST results from Seurat to replace CalcDEcombn ##

MAST_pw <- apply(combn(levels(your_seurat_obj@ident),2),
                 MARGIN=2,
                 function(X) {
                   FindMarkers(your_seurat_obj,
                               ident.1=X[1],
                               ident.2=X[2],
                               logfc.threshold=0,
                               min.pct=0.1,
                               test.use="MAST",
                               latent.vars="nUMI")
                 })
names(MAST_pw) <- apply(combn(levels(your_seurat_obj@ident),2),2,
                        function(X) paste(X,collapse="-"))
# ^ Names must be in "X-Y" format

for (i in names(MAST_pw)) {
  MAST_pw[[i]]$dDR <- MAST_pw[[i]]$pct.1 - MAST_pw[[i]]$pct.2
  # ^ Diff in detect rate (dDR) must be a variable in each dataframe
  names(MAST_pw[[i]])[names(MAST_pw[[i]]) == "avg_logFC"] <- "logGER"
  # ^ Effect size variable must be named 'logGER'
  names(MAST_pw[[i]])[names(MAST_pw[[i]]) == "p_val_adj"] <- "FDR"
  # ^ Significance variable must be named 'FDR'
  # Note: rownames of each dataframe must be gene names, 
  # but FindMarkers should already do this.
}
DEcombn(your_sCV_obj) <- MAST_pw
# ^ Slot MAST results into sCVdata object


## End(Not run)


BaderLab/scClustViz documentation built on Sept. 10, 2023, 11:51 p.m.