Package 'DEsingle'

Title: DEsingle for detecting three types of differential expression in single-cell RNA-seq data
Description: DEsingle is an R package for differential expression (DE) analysis of single-cell RNA-seq (scRNA-seq) data. It defines and detects 3 types of differentially expressed genes between two groups of single cells, with regard to different expression status (DEs), differential expression abundance (DEa), and general differential expression (DEg). DEsingle employs Zero-Inflated Negative Binomial model to estimate the proportion of real and dropout zeros and to define and detect the 3 types of DE genes. Results showed that DEsingle outperforms existing methods for scRNA-seq DE analysis, and can reveal different types of DE genes that are enriched in different biological functions.
Authors: Zhun Miao <[email protected]>
Maintainer: Zhun Miao <[email protected]>
License: GPL-2
Version: 1.25.0
Built: 2024-09-28 03:42:36 UTC
Source: https://github.com/bioc/DEsingle

Help Index


DEsingle: Detecting differentially expressed genes from scRNA-seq data

Description

This function is used to detect differentially expressed genes between two specified groups of cells in a raw read counts matrix of single-cell RNA-seq (scRNA-seq) data. It takes a non-negative integer matrix of scRNA-seq raw read counts or a SingleCellExperiment object as input. So users should map the reads (obtained from sequencing libraries of the samples) to the corresponding genome and count the reads mapped to each gene according to the gene annotation to get the raw read counts matrix in advance.

Usage

DEsingle(counts, group, parallel = FALSE, BPPARAM = bpparam())

Arguments

counts

A non-negative integer matrix of scRNA-seq raw read counts or a SingleCellExperiment object which contains the read counts matrix. The rows of the matrix are genes and columns are samples/cells.

group

A vector of factor which specifies the two groups to be compared, corresponding to the columns in the counts matrix.

parallel

If FALSE (default), no parallel computation is used; if TRUE, parallel computation using BiocParallel, with argument BPPARAM.

BPPARAM

An optional parameter object passed internally to bplapply when parallel=TRUE. If not specified, bpparam() (default) will be used.

Value

A data frame containing the differential expression (DE) analysis results, rows are genes and columns contain the following items:

  • theta_1, theta_2, mu_1, mu_2, size_1, size_2, prob_1, prob_2: MLE of the zero-inflated negative binomial distribution's parameters of group 1 and group 2.

  • total_mean_1, total_mean_2: Mean of read counts of group 1 and group 2.

  • foldChange: total_mean_1/total_mean_2.

  • norm_total_mean_1, norm_total_mean_2: Mean of normalized read counts of group 1 and group 2.

  • norm_foldChange: norm_total_mean_1/norm_total_mean_2.

  • chi2LR1: Chi-square statistic for hypothesis testing of H0.

  • pvalue_LR2: P value of hypothesis testing of H20 (Used to determine the type of a DE gene).

  • pvalue_LR3: P value of hypothesis testing of H30 (Used to determine the type of a DE gene).

  • FDR_LR2: Adjusted P value of pvalue_LR2 using Benjamini & Hochberg's method (Used to determine the type of a DE gene).

  • FDR_LR3: Adjusted P value of pvalue_LR3 using Benjamini & Hochberg's method (Used to determine the type of a DE gene).

  • pvalue: P value of hypothesis testing of H0 (Used to determine whether a gene is a DE gene).

  • pvalue.adj.FDR: Adjusted P value of H0's pvalue using Benjamini & Hochberg's method (Used to determine whether a gene is a DE gene).

  • Remark: Record of abnormal program information.

Author(s)

Zhun Miao.

See Also

DEtype, for the classification of differentially expressed genes found by DEsingle.

TestData, a test dataset for DEsingle.

Examples

# Load test data for DEsingle
data(TestData)

# Specifying the two groups to be compared
# The sample number in group 1 and group 2 is 50 and 100 respectively
group <- factor(c(rep(1,50), rep(2,100)))

# Detecting the differentially expressed genes
results <- DEsingle(counts = counts, group = group)

# Dividing the differentially expressed genes into 3 categories
results.classified <- DEtype(results = results, threshold = 0.05)

DEtype: Classifying differentially expressed genes from DEsingle

Description

This function is used to classify the differentially expressed genes of single-cell RNA-seq (scRNA-seq) data found by DEsingle. It takes the output data frame from DEsingle as input.

Usage

DEtype(results, threshold)

Arguments

results

A output data frame from DEsingle, which contains the unclassified differential expression analysis results.

threshold

A number of (0,1) to specify the threshold of FDR.

Value

A data frame containing the differential expression (DE) analysis results and DE gene types and states.

  • theta_1, theta_2, mu_1, mu_2, size_1, size_2, prob_1, prob_2: MLE of the zero-inflated negative binomial distribution's parameters of group 1 and group 2.

  • total_mean_1, total_mean_2: Mean of read counts of group 1 and group 2.

  • foldChange: total_mean_1/total_mean_2.

  • norm_total_mean_1, norm_total_mean_2: Mean of normalized read counts of group 1 and group 2.

  • norm_foldChange: norm_total_mean_1/norm_total_mean_2.

  • chi2LR1: Chi-square statistic for hypothesis testing of H0.

  • pvalue_LR2: P value of hypothesis testing of H20 (Used to determine the type of a DE gene).

  • pvalue_LR3: P value of hypothesis testing of H30 (Used to determine the type of a DE gene).

  • FDR_LR2: Adjusted P value of pvalue_LR2 using Benjamini & Hochberg's method (Used to determine the type of a DE gene).

  • FDR_LR3: Adjusted P value of pvalue_LR3 using Benjamini & Hochberg's method (Used to determine the type of a DE gene).

  • pvalue: P value of hypothesis testing of H0 (Used to determine whether a gene is a DE gene).

  • pvalue.adj.FDR: Adjusted P value of H0's pvalue using Benjamini & Hochberg's method (Used to determine whether a gene is a DE gene).

  • Remark: Record of abnormal program information.

  • Type: Types of DE genes. DEs represents different expression status; DEa represents differential expression abundance; DEg represents general differential expression.

  • State: State of DE genes, up represents up-regulated; down represents down-regulated.

Author(s)

Zhun Miao.

See Also

DEsingle, for the detection of differentially expressed genes from scRNA-seq data.

TestData, a test dataset for DEsingle.

Examples

# Load test data for DEsingle
data(TestData)

# Specifying the two groups to be compared
# The sample number in group 1 and group 2 is 50 and 100 respectively
group <- factor(c(rep(1,50), rep(2,100)))

# Detecting the differentially expressed genes
results <- DEsingle(counts = counts, group = group)

# Dividing the differentially expressed genes into 3 categories
results.classified <- DEtype(results = results, threshold = 0.05)

TestData: A test dataset for DEsingle

Description

A toy dataset containing a single-cell RNA-seq (scRNA-seq) read counts matrix and its grouping information.

Usage

data(TestData)

Format

  • counts. A non-negative integer matrix of scRNA-seq raw read counts, rows are genes and columns are cells.

  • group. A vector of factor specifying the two groups to be compared, corresponding to the columns of the counts.

Details

  • counts. A matrix of raw read counts of scRNA-seq data which has 200 genes (rows) and 150 cells (columns).

  • group. A vector of factor specifying the two groups to be compared in counts. Also could be generated by: group <- factor(c(rep(1,50), rep(2,100)))

Source

Petropoulos S, et al. Cell, 2016, 165(4): 1012-1026.

See Also

DEsingle, for the detection of differentially expressed genes from scRNA-seq data.

DEtype, for the classification of differentially expressed genes found by DEsingle.

Examples

# Load test data for DEsingle
data(TestData)

# Specifying the two groups to be compared
# The sample number in group 1 and group 2 is 50 and 100 respectively
group <- factor(c(rep(1,50), rep(2,100)))

# Detecting the differentially expressed genes
results <- DEsingle(counts = counts, group = group)

# Dividing the differentially expressed genes into 3 categories
results.classified <- DEtype(results = results, threshold = 0.05)