deconvR 1.2.0
Recent studies associated the differences of cell-type proportions may be
correlated to certain phenotypes, such as cancer. Therefore, the demand for the
development of computational methods to predict cell type proportions increased.
Hereby, we developed deconvR, a collection of functions designed for analyzing
deconvolution of the bulk sample(s) using an atlas of reference omic signature
profiles and a user-selected model. We wanted to give users an option to extend
their reference atlas. Users can create new reference atlases using
findSignatures extend their atlas by adding more cell types. Additionally, we
included BSMeth2Probe to make mapping whole-genome bisulfite sequencing data
to their probe IDs easier. So users can map WGBS methylation data (as in
methylKit or GRanges object format) to probe IDs, and the results of
this mapping can be used as the bulk samples in the deconvolution. We also
included a comprehensive DNA methylation atlas of 25 different cell types to use
in the main function deconvolute. deconvolute allows the user to specify
their desired deconvolution model (non-negative least squares regression,
support vector regression, quadratic programming, or robust linear regression),
and returns a dataframe which contains predicted cell-type proportions of bulk
methylation profiles, as well as partial R-squared values for each sample.
As an another option, users can generate a simulated table of a desired number
of samples, with either user-specified or random origin proportions using
simulateCellMix. simulateCellMix returns a second data frame called
proportions, which contains the actual cell-type proportions of the simulated
sample. It can be used for testing the accuracy of the deconvolution by
comparing these actual proportions to the proportions predicted by
deconvolute.
deconvolute returns partial R-squares, to check if deconvolution brings
advantages on top of the basic bimodal profiles. The reference matrix usually
follows a bimodal distribution in the case of methylation, and taking the
average of the rows of methylation matrix might give a pretty similar profile
to the bulk methylation profile you are trying to deconvolute. If the
deconvolution is advantageous, partial R-squared is expect to be high.
The deconvR package can be installed from Bioconductor with:
if (!requireNamespace("BiocManager", quietly = TRUE))
    install.packages("BiocManager")
BiocManager::install("deconvR")The comprehensive human methylome reference atlas created by Moss et al. 1 Moss,
J. et al. (2018). Comprehensive human cell-type methylation atlas reveals
origins of circulating cell-free DNA in health and disease. Nature
communications, 9(1), 1-12. https://doi.org/10.1038/s41467-018-07466-6 can be
used as the reference atlas parameter for several functions in this package.
This atlas was modified to remove duplicate CpG loci before being included in
the package as the dataframe. The dataframe is composed of 25 human cell types
and roughly 6000 CpG loci, identified by their Illumina Probe ID. For each cell
type and CpG locus, a methylation value between 0 and 1 is provided. This value
represents the fraction of methylated bases of the CpG locus. The atlas
therefore provides a unique methylation pattern for each cell type and can be
directly used as reference in deconvolute and simulateCellMix, and atlas
in findSignatures. Below is an example dataframe to illustrate the atlas
format.
library(deconvR) 
data("HumanCellTypeMethAtlas")
head(HumanCellTypeMethAtlas[,1:5])##          IDs Monocytes_EPIC B.cells_EPIC CD4T.cells_EPIC NK.cells_EPIC
## 1 cg08169020         0.8866       0.2615          0.0149        0.0777
## 2 cg25913761         0.8363       0.2210          0.2816        0.4705
## 3 cg26955540         0.7658       0.0222          0.1492        0.4005
## 4 cg25170017         0.8861       0.5116          0.1021        0.4363
## 5 cg12827637         0.5212       0.3614          0.0227        0.2120
## 6 cg19442545         0.2013       0.1137          0.0608        0.0410The GRanges object IlluminaMethEpicB5ProbeIDs contains the Illumina probe
IDs of 400000 genomic loci (identified using the “seqnames”, “ranges”, and
“strand” values). This object is based off of the Infinium MethylationEPIC v1.0
B5 Manifest data. Unnecessary columns were removed and rows were truncated to
reduce file size before converting the file to a GRanges object. It can be
used directly as probe_id_locations in BSmeth2Probe.
data("IlluminaMethEpicB5ProbeIDs")
head(IlluminaMethEpicB5ProbeIDs)## GRanges object with 6 ranges and 1 metadata column:
##       seqnames              ranges strand |          ID
##          <Rle>           <IRanges>  <Rle> | <character>
##   [1]    chr19     5236005-5236006      + |  cg07881041
##   [2]    chr20   63216298-63216299      + |  cg18478105
##   [3]     chr1     6781065-6781066      + |  cg23229610
##   [4]     chr2 197438742-197438743      - |  cg03513874
##   [5]     chrX   24054523-24054524      + |  cg09835024
##   [6]    chr14   93114794-93114795      - |  cg05451842
##   -------
##   seqinfo: 25 sequences from an unspecified genome; no seqlengthsAs mentioned in the introduction section, users can extend their reference atlas
to incorporate new data. Or may combine different reference atlases to construct
a more comprehensive one. This can be done using the findSignatures function.
In this example, since we don’t have any additional reference atlas, we will add
simulated data as a new cell type to reference atlas for example purposes.
First, ensure that atlas (the signature matrix to be extended) and samples
(the new data to be added to the signature matrix) are compliant with the
function requirements. Below illustrates the samples format.
samples <- simulateCellMix(3,reference = HumanCellTypeMethAtlas)$simulated
head(samples)##          IDs   Sample 1   Sample 2  Sample 3
## 1 cg08169020 0.27591621 0.14097650 0.1342789
## 2 cg25913761 0.34264959 0.30572697 0.2893796
## 3 cg26955540 0.21397786 0.20976786 0.1922104
## 4 cg25170017 0.33156724 0.18303723 0.1606088
## 5 cg12827637 0.20539436 0.11321829 0.1078107
## 6 cg19442545 0.07518705 0.04936518 0.0463514sampleMeta should include all sample names in samples, and specify the
origins they should be mapped to when added to atlas.
sampleMeta <- data.table("Experiment_accession" = colnames(samples)[-1],
                         "Biosample_term_name" = "new cell type")
head(sampleMeta)##    Experiment_accession Biosample_term_name
## 1:             Sample 1       new cell type
## 2:             Sample 2       new cell type
## 3:             Sample 3       new cell typeUse findSignatures to extend the matrix.
extended_matrix <- findSignatures(samples = samples, 
                                 sampleMeta = sampleMeta, 
                                 atlas = HumanCellTypeMethAtlas)## CELL TYPES IN EXTENDED ATLAS: 
## new cell type 
## Monocytes_EPIC 
## B.cells_EPIC 
## CD4T.cells_EPIC 
## NK.cells_EPIC 
## CD8T.cells_EPIC 
## Neutrophils_EPIC 
## Erythrocyte_progenitors 
## Adipocytes 
## Cortical_neurons 
## Hepatocytes 
## Lung_cells 
## Pancreatic_beta_cells 
## Pancreatic_acinar_cells 
## Pancreatic_duct_cells 
## Vascular_endothelial_cells 
## Colon_epithelial_cells 
## Left_atrium 
## Bladder 
## Breast 
## Head_and_neck_larynx 
## Kidney 
## Prostate 
## Thyroid 
## Upper_GI 
## Uterus_cervixhead(extended_matrix)##          IDs new_cell_type Monocytes_EPIC B.cells_EPIC CD4T.cells_EPIC
## 1 cg08169020    0.18372386         0.8866       0.2615          0.0149
## 2 cg25913761    0.31258538         0.8363       0.2210          0.2816
## 3 cg26955540    0.20531871         0.7658       0.0222          0.1492
## 4 cg25170017    0.22507108         0.8861       0.5116          0.1021
## 5 cg12827637    0.14214111         0.5212       0.3614          0.0227
## 6 cg19442545    0.05696788         0.2013       0.1137          0.0608
##   NK.cells_EPIC CD8T.cells_EPIC Neutrophils_EPIC Erythrocyte_progenitors
## 1        0.0777          0.0164           0.8680                  0.9509
## 2        0.4705          0.3961           0.8293                  0.2385
## 3        0.4005          0.3474           0.7915                  0.1374
## 4        0.4363          0.0875           0.7042                  0.9447
## 5        0.2120          0.0225           0.5368                  0.4667
## 6        0.0410          0.0668           0.1952                  0.1601
##   Adipocytes Cortical_neurons Hepatocytes Lung_cells Pancreatic_beta_cells
## 1     0.0336           0.0168      0.0340     0.0416              0.038875
## 2     0.3578           0.3104      0.2389     0.2250              0.132000
## 3     0.1965           0.0978      0.0338     0.0768              0.041725
## 4     0.0842           0.2832      0.2259     0.0544              0.111750
## 5     0.0287           0.1368      0.0307     0.1607              0.065975
## 6     0.0364           0.0222      0.1574     0.0122              0.003825
##   Pancreatic_acinar_cells Pancreatic_duct_cells Vascular_endothelial_cells
## 1                  0.0209                0.0130                     0.0323
## 2                  0.2249                0.1996                     0.3654
## 3                  0.0314                0.0139                     0.2382
## 4                  0.0309                0.0217                     0.0972
## 5                  0.0370                0.0230                     0.0798
## 6                  0.0378                0.0347                     0.0470
##   Colon_epithelial_cells Left_atrium Bladder Breast Head_and_neck_larynx Kidney
## 1                 0.0163      0.0386  0.0462 0.0264               0.0470 0.0269
## 2                 0.2037      0.2446  0.2054 0.1922               0.2045 0.1596
## 3                 0.0193      0.1134  0.1269 0.1651               0.1523 0.1034
## 4                 0.0187      0.0674  0.0769 0.0691               0.0704 0.0604
## 5                 0.0193      0.0432  0.0459 0.0228               0.0687 0.0234
## 6                 0.0193      0.0287  0.0246 0.0081               0.0098 0.0309
##   Prostate Thyroid Upper_GI Uterus_cervix
## 1   0.0353  0.0553   0.0701        0.0344
## 2   0.1557  0.1848   0.1680        0.2026
## 3   0.0686  0.0943   0.1298        0.1075
## 4   0.0369  0.0412   0.0924        0.0697
## 5   0.0508  0.0726   0.0759        0.0196
## 6   0.0055  0.0188   0.0090        0.0166WGBS methylation data first needs to be mapped to probes using BSmeth2Probe
before being deconvoluted. The methylation data WGBS_data in BSmeth2Probe
may be either a GRanges object or a methylKit object.
Format of WGBS_data as GRanges object:
load(system.file("extdata", "WGBS_GRanges.rda",
                                     package = "deconvR"))
head(WGBS_GRanges)## GRanges object with 6 ranges and 1 metadata column:
##       seqnames    ranges strand |   sample1
##          <Rle> <IRanges>  <Rle> | <numeric>
##   [1]     chr1     47176      + |  0.833333
##   [2]     chr1     47177      - |  0.818182
##   [3]     chr1     47205      + |  0.681818
##   [4]     chr1     47206      - |  0.583333
##   [5]     chr1     47362      + |  0.416667
##   [6]     chr1     49271      + |  0.733333
##   -------
##   seqinfo: 1 sequence from an unspecified genome; no seqlengthsor as methylKit object:
head(methylKit::methRead(system.file("extdata", "test1.myCpG.txt", 
                                     package = "methylKit"), 
                         sample.id="test", assembly="hg18", 
                         treatment=1, context="CpG", mincov = 0))##     chr   start     end strand coverage numCs numTs
## 1 chr21 9764513 9764513      -       12     0    12
## 2 chr21 9764539 9764539      -       12     3     9
## 3 chr21 9820622 9820622      +       13     0    13
## 4 chr21 9837545 9837545      +       11     0    11
## 5 chr21 9849022 9849022      +      124    90    34
## 6 chr21 9853296 9853296      +       17    10     7probe_id_locations contains information needed to map cellular loci to probe IDs
data("IlluminaMethEpicB5ProbeIDs")
head(IlluminaMethEpicB5ProbeIDs)## GRanges object with 6 ranges and 1 metadata column:
##       seqnames              ranges strand |          ID
##          <Rle>           <IRanges>  <Rle> | <character>
##   [1]    chr19     5236005-5236006      + |  cg07881041
##   [2]    chr20   63216298-63216299      + |  cg18478105
##   [3]     chr1     6781065-6781066      + |  cg23229610
##   [4]     chr2 197438742-197438743      - |  cg03513874
##   [5]     chrX   24054523-24054524      + |  cg09835024
##   [6]    chr14   93114794-93114795      - |  cg05451842
##   -------
##   seqinfo: 25 sequences from an unspecified genome; no seqlengthsUse BSmeth2Probe to map WGBS data to probe IDs.
mapped_WGBS_data <- BSmeth2Probe(probe_id_locations = IlluminaMethEpicB5ProbeIDs, 
                                 WGBS_data = WGBS_GRanges,
                                 multipleMapping = TRUE,
                                 cutoff = 10)
head(mapped_WGBS_data)##          IDs   sample1
## 1 cg00305774 1.0000000
## 2 cg00546078 0.8181818
## 3 cg00546307 0.5454545
## 4 cg00546971 0.5000000
## 5 cg00774867 0.8461538
## 6 cg00830435 0.9166667This mapped data can now be used in deconvolute. Here we will deconvolute it
using the previously extended signature matrix as the reference atlas.
deconvolution <- deconvolute(reference = HumanCellTypeMethAtlas, 
                             bulk = mapped_WGBS_data)##    Min. 1st Qu.  Median    Mean 3rd Qu.    Max. 
##  0.9963  0.9963  0.9963  0.9963  0.9963  0.9963deconvolution$proportions##         Monocytes_EPIC B.cells_EPIC CD4T.cells_EPIC NK.cells_EPIC
## sample1              0            0               0             0
##         CD8T.cells_EPIC Neutrophils_EPIC Erythrocyte_progenitors Adipocytes
## sample1               0                0                       0  0.5011789
##         Cortical_neurons Hepatocytes Lung_cells Pancreatic_beta_cells
## sample1        0.2917357           0          0                     0
##         Pancreatic_acinar_cells Pancreatic_duct_cells
## sample1                       0                     0
##         Vascular_endothelial_cells Colon_epithelial_cells Left_atrium Bladder
## sample1                          0            0.001012524           0       0
##         Breast Head_and_neck_larynx Kidney Prostate Thyroid Upper_GI
## sample1      0            0.2060729      0        0       0        0
##         Uterus_cervix
## sample1             0It is possible to use RNA-seq data for deconvolution via deconvR package.
Beware that you have to set IDs column that contains Gene names to run
deconvR functions. Therefore you can simulate bulk RNA-seq data via
simulateCellMix and, extend RNA-seq reference atlas via findSignatures.
sessionInfo()## R version 4.2.0 RC (2022-04-19 r82224)
## 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] parallel  stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
## [1] dplyr_1.0.8       doParallel_1.0.17 iterators_1.0.14  foreach_1.5.2    
## [5] deconvR_1.2.0     data.table_1.14.2 knitr_1.38        BiocStyle_2.24.0 
## 
## loaded via a namespace (and not attached):
##   [1] nlme_3.1-157                rsq_2.2                    
##   [3] bitops_1.0-7                matrixStats_0.62.0         
##   [5] GenomeInfoDb_1.32.0         numDeriv_2016.8-1.1        
##   [7] Deriv_4.1.3                 tools_4.2.0                
##   [9] bslib_0.3.1                 utf8_1.2.2                 
##  [11] R6_2.5.1                    DBI_1.1.2                  
##  [13] BiocGenerics_0.42.0         colorspace_2.0-3           
##  [15] fastseg_1.42.0              tidyselect_1.1.2           
##  [17] compiler_4.2.0              cli_3.3.0                  
##  [19] Biobase_2.56.0              DelayedArray_0.22.0        
##  [21] rtracklayer_1.56.0          bookdown_0.26              
##  [23] sass_0.4.1                  scales_1.2.0               
##  [25] nnls_1.4                    mvtnorm_1.1-3              
##  [27] quadprog_1.5-8              proxy_0.4-26               
##  [29] stringr_1.4.0               digest_0.6.29              
##  [31] Rsamtools_2.12.0            minqa_1.2.4                
##  [33] rmarkdown_2.14              R.utils_2.11.0             
##  [35] XVector_0.36.0              pkgconfig_2.0.3            
##  [37] htmltools_0.5.2             lme4_1.1-29                
##  [39] MatrixGenerics_1.8.0        fastmap_1.1.0              
##  [41] bbmle_1.0.24                limma_3.52.0               
##  [43] rlang_1.0.2                 jquerylib_0.1.4            
##  [45] BiocIO_1.6.0                generics_0.1.2             
##  [47] jsonlite_1.8.0              mclust_5.4.9               
##  [49] BiocParallel_1.30.0         gtools_3.9.2               
##  [51] R.oo_1.24.0                 RCurl_1.98-1.6             
##  [53] magrittr_2.0.3              GenomeInfoDbData_1.2.8     
##  [55] Matrix_1.4-1                Rcpp_1.0.8.3               
##  [57] munsell_0.5.0               S4Vectors_0.34.0           
##  [59] fansi_1.0.3                 lifecycle_1.0.1            
##  [61] R.methodsS3_1.8.1           stringi_1.7.6              
##  [63] yaml_2.3.5                  SummarizedExperiment_1.26.0
##  [65] MASS_7.3-57                 zlibbioc_1.42.0            
##  [67] plyr_1.8.7                  qvalue_2.28.0              
##  [69] grid_4.2.0                  bdsmatrix_1.3-4            
##  [71] crayon_1.5.1                lattice_0.20-45            
##  [73] Biostrings_2.64.0           splines_4.2.0              
##  [75] methylKit_1.22.0            pillar_1.7.0               
##  [77] GenomicRanges_1.48.0        boot_1.3-28                
##  [79] rjson_0.2.21                reshape2_1.4.4             
##  [81] codetools_0.2-18            stats4_4.2.0               
##  [83] XML_3.99-0.9                glue_1.6.2                 
##  [85] evaluate_0.15               BiocManager_1.30.17        
##  [87] nloptr_2.0.0                vctrs_0.4.1                
##  [89] tidyr_1.2.0                 gtable_0.3.0               
##  [91] purrr_0.3.4                 assertthat_0.2.1           
##  [93] ggplot2_3.3.5               emdbook_1.3.12             
##  [95] xfun_0.30                   restfulr_0.0.13            
##  [97] e1071_1.7-9                 coda_0.19-4                
##  [99] class_7.3-20                tibble_3.1.6               
## [101] GenomicAlignments_1.32.0    IRanges_2.30.0             
## [103] ellipsis_0.3.2stopCluster(cl)