User guide to the dearseq R package

Marine Gauthier, Denis Agniel, Boris Hejblum

2022-07-21

1 Overview of dearseq

dearseq is an R package that performs Differential Expression Analysis of RNA-seq data, while guaranteeing a sound control of false positives (Gauthier et al., 2019). It relies on variance component score test accounting for data heteroscedasticity through precision weights, not unlike voom (Law et al., 2014). It can perform either

In addition, dearseq can deal with various and complex experimental designs such as:

More details are provided Gauthier M, Agniel D, Thiébaut R & Hejblum BP (2020). dearseq: a variance component score test for RNA-Seq differential analysis that effectively controls the false discovery rate, NAR Genomics and Bioinformatics, 2(4):lqaa093.
DOI: 10.1093/nargab/lqaa093

2 Using dearseq for a gene-wise analysis

2.1 Getting started using dear_seq()

Two inputs are required to run dear_seq():

Both can be supplied in the formed of a SummarizedExperiment object for instance. A third optional input is the covariates that will be adjusted on (such as age for instance).

2.2 An example analysis

Tuberculosis (TB) is a disease caused by a bacterium called Mycobacterium tuberculosis. Bacteria typically infect the lungs, but they can also affect other parts of the body. Tuberculosis can remain in a quiescent state in the body. It is called latent tuberculosis infection (LTBI). It is an infection without clinical signs, bacteriological and radiological disease.

Berry et al. first identified a whole blood 393 transcript signature for active TB using microarray analysis. Then, Singhania et al. found a 373-genes signature of active tuberculosis using RNA-Seq, confirming microarray results, that discriminates active tuberculosis from latently infected and healthy individuals.

We sought to investigate how many of the 373 genes Singhania et al. found using edgeR might actually be false positives. As a quick example, we propose a partial re-analysis (see Gauthier et al. for the complete analysis) of the Singhania et al. dataset: we proceed to the Differential Expression Analysis (DEA) of the Active TB group against the Control one, omitting LTBI/Control and ActiveTB/LTBI comparisons. We focused on the Berry London dataset.

It results in 54 patients of whom 21 are active TB patients, 21 are latent TB patients and 12 are healthy control particiants.

For more details, check the article from Berry et al., 2010 here.

2.2.1 Data preparation

Data from Singhania et al. inverstigating active TB are publicly available on GEO website with GEO access number ‘GSE107995’. Thanks to the GEOquery package to get the data files from GEO website (see appendix for more details on GEOquery)

2.2.1.2 Gene expression matrix

This matrix contains the gene expression (in cells) for each gene (in rows) of each sample (in columns) gathered from RNA-seq measurements. The gene expression could already be normalized, or not (e.g. raw counts or pseudo-counts straight from the alignment). In the latter case, the gene expression is normalized then into log(counts) per million (i.e. log-cpm).

The gene expression matrix (called London) is including in the package dearseq. It was extracted from one of the GEO supplementary files (namely the “GSE107991_edgeR_normalized_Berry_London.xlsx” file).

NB: These expression data have already been already normalized using edgeR

We have nrow(London) genes in rows and ncol(London) samples in columns * rownames are Ensembl IDs ENSGxxxxxxxxxxx for each gene (saved in genes) * colnames are the name of each sample Berry_London_Samplex * each data-cell contains the normalized gene expression

The vector genes contains all the gene identifiers.

2.2.2 Identifying differentially expressed genes (DEG) from dearseq variance component score test

The main arguments to the dear_seq() function are:

  • exprmat: a numeric matrix containing the raw RNA-seq counts or preprocessed expressions
  • covariates: a numeric matrix containing the model covariates that needs to be adjusted. Usually, its first column is the intercept (full of 1s).
  • variables2test: a numeric design matrix containing the variables of interest (with whom the expression association is tested)
  • which_test: a character string indicating which method to use to approximate the variance component score test, either “permutation” or “asymptotic”.
  • preprocessed: a logical flag indicating whether the expression data have already been preprocessed (e.g. log2 transformed). Default is FALSE, in which case y is assumed to contain raw counts and is normalized into log(counts) per million.

We use the permutation test because of the relatively small sample size (\(n=33\)).

res_dearseq is a list containing:

  • the input values for which_test, nperm and preprocessed
  • the gene-wise p-values (raw and adjusted - default is correction for multiple testing with the Benjamini-Hochberg procedure)

NB: We excluded the LTBI patients to focus on the comparison between Active TB versus Control patients.

2.2.3 Using gene expression matrix of raw counts

If the gene expression could not already be normalized, the gene expression is normalized then into log(counts) per million (i.e. log-cpm). We set preprocessed equal to FALSE, the default parameter to normalize the raw counts.

The matrix of raw counts (Raw_counts_Berry_London) contains \(58,051\) genes and \(54\) samples. As described in Singhania et al., only genes expressed with counts per million (CPM) \(> 2\) in at least five samples were considered. The filtering is carried out with edgeR. It results in a matrix of raw counts containing \(14,150\) genes and \(54\) samples (London_raw). Of note, many filtering and normalization methods may be considered and the user must be cautious applying them.

3 Using dearseq for gene-set analysis

Here we will use a subsample of the RNA-seq data from Baduel et al. (2016) studying Arabidopsis Arenosa physiology, that is included in the dearseq package.

We investigates 2 gene sets that were derived in a data-driven manner By Baduel et al. (2016).

library(dearseq)
library(SummarizedExperiment)
library(BiocSet)
data("baduel_5gs")
se2 <- SummarizedExperiment(assay = log2(expr_norm_corr + 1), colData = design)
genes_non0var_ind <- which(matrixStats::rowVars(expr_norm_corr) != 0)

KAvsTBG <- dearseq::dgsa_seq(object = se2[genes_non0var_ind, ], covariates = c("Vernalized",
    "AgeWeeks", "Vernalized_Population", "AgeWeeks_Population"), variables2test = "PopulationKA",
    genesets = baduel_gmt$genesets[c(3, 5)], which_test = "permutation", which_weights = "loclin",
    n_perm = 2000, preprocessed = TRUE, verbose = FALSE, parallel_comp = interactive())
summary(KAvsTBG)
#> dearseq identifies 1 significant gene sets out of 2 
#>  - test: permutation 
#>  - FDR significance level: 0.05
KAvsTBG$pvals
#>                rawPval   adjPval
#> VR_TBG     0.911544228 0.9115442
#> DE_KAvsTBG 0.009995002 0.0199900

Cold <- dearseq::dgsa_seq(object = se2[genes_non0var_ind, ], covariates = c("AgeWeeks",
    "PopulationKA", "AgeWeeks_Population"), variables2test = c("Vernalized", "Vernalized_Population"),
    genesets = baduel_gmt$genesets[c(3, 5)], which_test = "permutation", which_weights = "loclin",
    n_perm = 2000, preprocessed = TRUE, verbose = FALSE, parallel_comp = interactive())
summary(Cold)
#> dearseq identifies 2 significant gene sets out of 2 
#>  - test: permutation 
#>  - FDR significance level: 0.05
Cold$pvals
#>               rawPval    adjPval
#> VR_TBG     0.02398801 0.03348326
#> DE_KAvsTBG 0.03348326 0.03348326

We can also visualize the (pseudo-)longitudinal expression within a gene set:

# Remove 'non-vernalized' samples
colData(se2)[["Indiv"]] <- colData(se2)[["Population"]]:colData(se2)[["Replicate"]]
colData(se2)[["Vern"]] <- ifelse(colData(se2)[["Vernalized"]], "Vernalized", "Non-vernalized")

library(ggplot2)
spaghettiPlot1GS(gs_index = 3, gmt = baduel_gmt, expr_mat = assay(se2), design = colData(se2),
    var_time = AgeWeeks, var_indiv = Indiv, sampleIdColname = "sample", var_group = Vern,
    var_subgroup = Population, loess_span = 1.5) + ggplot2::xlab("Age (weeks)")

4 Bibliography

Agniel D & Hejblum BP (2017). Variance component score test for time-course gene set analysis of longitudinal RNA-seq data. Biostatistics 18(4):589-604. DOI: 10.1093/biostatistics/kxx005.

Baduel P, Arnold B, Weisman CM, Hunter B & Bomblies K (2016). Habitat-Associated Life History and Stress-Tolerance Variation in Arabidopsis Arenosa. Plant Physiology, 171(1):437-51. DOI: 10.1104/pp.15.01875.

Berry MP, Graham CM, McNab FW, Xu Z, Bloch SAA, Oni T, Wilkinson KA, Banchereau R, Skinner J, Wilkinson RJ, Quinn C, Blankenship D, Dhawan R, Cush JJ, Mejias A, Ramilo O, Kon OM, Pascual V, Banchereau J, Chaussabel D & O’Garra A (2010). An interferon-inducible neutrophil-driven blood transcriptional signature in human tuberculosis. Nature 466:973–977. DOI: 10.1038/nature09247

Gauthier M, Agniel D, Thiébaut R & Hejblum BP (2020). dearseq: a variance component score test for RNA-Seq differential analysis that effectively controls the false discovery rate. NAR Genomics and Bioinformatics, 2(4):lqaa093. DOI: 10.1093/nargab/lqaa093

Singhania A, Verma R, Graham CM, Lee J, Tran T, Richardson M, Lecine P, Leissner P, Berry MPR, Wilkinson RJ, Kaiser K, Rodrigue M, Woltmann G, Haldar P, O’Garra A, (2018). A modular transcriptional signature identifies phenotypic heterogeneity of human tuberculosis infection. Nature communications 9(1):2308. DOI: 10.1038/s41467-018-04579-w

5 Appendix

5.1 GEOquery package

In case the data you want to analyze is publicly available through Gene Expression Omnibus (GEO), you can access it with the GEOquery package, that can be installed with the following commands:

5.2 readxl

The readxl package allows to easily imports data from Excel format and can be similarly installed with the following commands:

More details can be found on Bioconductor and in Davis S, Meltzer P, (2007) GEOquery: a bridge between the Gene Expression Omnibus (GEO) and Bioconductor Bioinformatics 14:1846-1847.

5.3 Session Info

#> R version 4.2.1 (2022-06-23)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 20.04.4 LTS
#> 
#> Matrix products: default
#> BLAS:   /home/biocbuild/bbs-3.15-bioc/R/lib/libRblas.so
#> LAPACK: /home/biocbuild/bbs-3.15-bioc/R/lib/libRlapack.so
#> 
#> locale:
#>  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
#>  [3] LC_TIME=en_GB              LC_COLLATE=C              
#>  [5] LC_MONETARY=en_US.UTF-8    LC_MESSAGES=en_US.UTF-8   
#>  [7] LC_PAPER=en_US.UTF-8       LC_NAME=C                 
#>  [9] LC_ADDRESS=C               LC_TELEPHONE=C            
#> [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C       
#> 
#> attached base packages:
#> [1] stats4    stats     graphics  grDevices utils     datasets  methods  
#> [8] base     
#> 
#> other attached packages:
#>  [1] ggplot2_3.3.6               BiocSet_1.10.0             
#>  [3] dplyr_1.0.9                 edgeR_3.38.1               
#>  [5] limma_3.52.2                SummarizedExperiment_1.26.1
#>  [7] GenomicRanges_1.48.0        GenomeInfoDb_1.32.2        
#>  [9] IRanges_2.30.0              S4Vectors_0.34.0           
#> [11] MatrixGenerics_1.8.1        matrixStats_0.62.0         
#> [13] dearseq_1.8.4               readxl_1.4.0               
#> [15] GEOquery_2.64.2             Biobase_2.56.0             
#> [17] BiocGenerics_0.42.0        
#> 
#> loaded via a namespace (and not attached):
#>  [1] nlme_3.1-158           bitops_1.0-7           bit64_4.0.5           
#>  [4] httr_1.4.3             tools_4.2.1            bslib_0.4.0           
#>  [7] utf8_1.2.2             R6_2.5.1               KernSmooth_2.23-20    
#> [10] mgcv_1.8-40            DBI_1.1.3              colorspace_2.0-3      
#> [13] withr_2.5.0            tidyselect_1.1.2       ontologyIndex_2.7     
#> [16] bit_4.0.4              curl_4.3.2             compiler_4.2.1        
#> [19] cli_3.3.0              formatR_1.12           xml2_1.3.3            
#> [22] DelayedArray_0.22.0    labeling_0.4.2         sass_0.4.2            
#> [25] scales_1.2.0           readr_2.1.2            pbapply_1.5-0         
#> [28] stringr_1.4.0          digest_0.6.29          rmarkdown_2.14        
#> [31] R.utils_2.12.0         XVector_0.36.0         pkgconfig_2.0.3       
#> [34] htmltools_0.5.3        fastmap_1.1.0          highr_0.9             
#> [37] rlang_1.0.4            RSQLite_2.2.15         BiocIO_1.6.0          
#> [40] jquerylib_0.1.4        generics_0.1.3         farver_2.1.1          
#> [43] jsonlite_1.8.0         R.oo_1.25.0            RCurl_1.98-1.7        
#> [46] magrittr_2.0.3         GenomeInfoDbData_1.2.8 patchwork_1.1.1       
#> [49] Matrix_1.4-1           Rcpp_1.0.9             munsell_0.5.0         
#> [52] fansi_1.0.3            lifecycle_1.0.1        R.methodsS3_1.8.2     
#> [55] stringi_1.7.8          yaml_2.3.5             CompQuadForm_1.4.3    
#> [58] zlibbioc_1.42.0        plyr_1.8.7             blob_1.2.3            
#> [61] grid_4.2.1             parallel_4.2.1         crayon_1.5.1          
#> [64] lattice_0.20-45        Biostrings_2.64.0      splines_4.2.1         
#> [67] KEGGREST_1.36.3        hms_1.1.1              locfit_1.5-9.6        
#> [70] knitr_1.39             pillar_1.8.0           reshape2_1.4.4        
#> [73] glue_1.6.2             evaluate_0.15          mitools_2.4           
#> [76] data.table_1.14.2      png_0.1-7              vctrs_0.4.1           
#> [79] tzdb_0.3.0             cellranger_1.1.0       gtable_0.3.0          
#> [82] purrr_0.3.4            tidyr_1.2.0            assertthat_0.2.1      
#> [85] cachem_1.0.6           xfun_0.31              survey_4.1-1          
#> [88] survival_3.3-1         viridisLite_0.4.0      tibble_3.1.7          
#> [91] memoise_2.0.1          AnnotationDbi_1.58.0   statmod_1.4.36        
#> [94] ellipsis_0.3.2