deconvR : Simulation and Deconvolution of Omic Profiles

Introduction

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 or 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.

Installation

The deconvR package can be installed from Bioconductor with:

if (!requireNamespace("BiocManager", quietly = TRUE))
    install.packages("BiocManager")

BiocManager::install("deconvR")

Data

Comprehensive Human Methylome Reference Atlas

The comprehensive human methylome reference atlas created by Moss et al. 1 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.0410

Illumina Infinium MethylationEPIC v1.0 B5 Manifest Probes (hg38)

The 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 seqlengths

Example Workflow For Whole Genome Bisulfate Sequencing Data

Expanding Reference Atlas

As 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.17121393 0.07671448 0.31566169
## 2 cg25913761 0.30988412 0.25511876 0.51702746
## 3 cg26955540 0.18135119 0.12779259 0.38987017
## 4 cg25170017 0.25900278 0.13082197 0.35107494
## 5 cg12827637 0.14201093 0.08595156 0.19779116
## 6 cg19442545 0.08933372 0.03806893 0.09223261

sampleMeta 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
##                  <char>              <char>
## 1:             Sample 1       new cell type
## 2:             Sample 2       new cell type
## 3:             Sample 3       new cell type

Use findSignatures to extend the matrix.

extended_matrix <- findSignatures(samples = samples, 
                                 sampleMeta = sampleMeta, 
                                 atlas = HumanCellTypeMethAtlas,
                                 IDs = "IDs")
## 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_cervix
head(extended_matrix)
##          IDs new_cell_type Monocytes_EPIC B.cells_EPIC CD4T.cells_EPIC
## 1 cg08169020    0.18786337         0.8866       0.2615          0.0149
## 2 cg25913761    0.36067678         0.8363       0.2210          0.2816
## 3 cg26955540    0.23300465         0.7658       0.0222          0.1492
## 4 cg25170017    0.24696657         0.8861       0.5116          0.1021
## 5 cg12827637    0.14191789         0.5212       0.3614          0.0227
## 6 cg19442545    0.07321175         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.0166

WGBS 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 seqlengths

or 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     7

probe_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 seqlengths

Use 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.9166667

This 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.9963
deconvolution$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             0

Constructing tissue specific CpG signature matrix

Alternatively, users can set tissueSpecCpGs as TRUE to construct tissue based methylation signature matrix by using the reference atlas. Here, we used simulated samples to construct tissue specific signature matrix since we don’t have tissue specific samples.

data("HumanCellTypeMethAtlas")
exampleSamples <- simulateCellMix(1,reference = HumanCellTypeMethAtlas)$simulated
exampleMeta <- data.table("Experiment_accession" = "example_sample",
                          "Biosample_term_name" = "example_cell_type")
colnames(exampleSamples) <- c("CpGs", "example_sample")
colnames(HumanCellTypeMethAtlas)[1] <- c("CpGs")

signatures <- findSignatures(
  samples = exampleSamples,
  sampleMeta = exampleMeta,
  atlas = HumanCellTypeMethAtlas,
  IDs = "CpGs", K = 100, tissueSpecCpGs = TRUE)
## CELL TYPES IN EXTENDED ATLAS: 
## example_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_cervix
print(head(signatures[[2]]))
## $hyper
## cg08169020 cg26955540 cg25170017 cg12827637 cg15838173 cg04858631 cg19442545 
## 0.23515522 0.16818944 0.16630292 0.16128419 0.14657616 0.14641148 0.14232089 
## cg10560079 cg00982136 cg22017733 cg26460530 cg13677741 cg09642825 cg24082121 
## 0.13711960 0.13546874 0.12833870 0.12833805 0.12427268 0.12371552 0.12312251 
## cg01424562 cg15376996 cg14789659 cg19300307 cg14785464 cg18856478 cg12474798 
## 0.12310619 0.12004955 0.11722376 0.11601867 0.11592566 0.11578440 0.11574113 
## cg06340704 cg11895835 cg27189395 cg26033520 cg05507654 cg05258935 cg05056497 
## 0.11569486 0.11507444 0.11467647 0.11411044 0.11321372 0.11277035 0.11264467 
## cg06167719 cg20684110 cg00936790 cg08474651 cg24275356 cg27395200 cg08857351 
## 0.11258307 0.11201174 0.11199686 0.11187844 0.11176133 0.11131857 0.11128526 
## cg25913761 cg23220823 cg10509607 cg14855367 cg22897141 cg08425810 cg23247274 
## 0.11119407 0.11115531 0.11084334 0.11046164 0.11022221 0.10964674 0.10959024 
## cg03254916 cg22879098 cg08063160 cg03711944 cg27388962 cg09267773 cg04913246 
## 0.10936571 0.10932215 0.10867351 0.10865236 0.10839221 0.10837209 0.10825230 
## cg02962602 cg10361922 cg08257293 cg26953232 cg08832851 cg26047334 cg16416715 
## 0.10822338 0.10810550 0.10810174 0.10808079 0.10774221 0.10764697 0.10752221 
## cg26757820 cg14057303 cg25158622 cg11975790 cg14845962 cg16127573 cg15014975 
## 0.10717412 0.10640446 0.10621697 0.10608573 0.10599274 0.10573722 0.10533505 
## cg08367326 cg27312961 cg16278496 cg05043794 cg08846870 cg19331221 cg18675610 
## 0.10530563 0.10517904 0.10455454 0.10454068 0.10417915 0.10402661 0.10387438 
## cg22138735 cg09577804 cg04055490 cg07865091 cg11532431 cg16402452 cg24591770 
## 0.10384958 0.10375695 0.10373261 0.10337715 0.10295510 0.10291040 0.10290975 
## cg10661769 cg03063309 cg02297709 cg21181453 cg06398251 cg17847861 cg04462378 
## 0.10261265 0.10254110 0.10235223 0.10232896 0.10222693 0.10206123 0.10195905 
## cg00769843 cg15185986 cg15575375 cg06889535 cg16017089 cg08659179 cg07218880 
## 0.10146514 0.10104857 0.10101545 0.10090343 0.10086865 0.10084209 0.10083496 
## cg01377358 cg11025609 cg01726273 cg12461469 cg16306706 cg07438103 cg13638257 
## 0.10080914 0.10067212 0.10048038 0.10046597 0.10035026 0.10031623 0.09989954 
## cg16829306 cg16988611 
## 0.09975261 0.09968723 
## 
## $hypo
## cg03663120 cg20942286 cg03963853 cg24788483 cg03126713 cg00828556 cg11186858 
##  0.3729023  0.2950421  0.2730039  0.2692113  0.2684764  0.2665202  0.2639320 
## cg05612654 cg13931640 cg15871206 cg22528270 cg15310871 cg13500029 cg10480329 
##  0.2628384  0.2611751  0.2599141  0.2591176  0.2581375  0.2573197  0.2487100 
## cg11231069 cg12655112 cg03549146 cg05923857 cg12458039 cg20610950 cg07737292 
##  0.2475624  0.2410186  0.2374149  0.2373430  0.2339092  0.2326071  0.2313390 
## cg10456459 cg03313271 cg06988336 cg06517984 cg16636767 cg06297318 cg14976569 
##  0.2266382  0.2217562  0.2200453  0.2123761  0.2122725  0.2105193  0.2096647 
## cg27334271 cg04851465 cg11597902 cg18990407 cg15633603 cg22259797 cg11327657 
##  0.2094470  0.2093740  0.2089090  0.2078968  0.2077687  0.2070742  0.2063704 
## cg23952578 cg08425796 cg13403369 cg06373940 cg09322573 cg04354689 cg22185977 
##  0.2062123  0.2032008  0.2028108  0.2026555  0.2016088  0.2007010  0.2003340 
## cg20107506 cg10718056 cg06125903 cg07033722 cg25517015 cg01879591 cg27366072 
##  0.2001573  0.1990423  0.1989825  0.1987527  0.1984383  0.1980188  0.1977045 
## cg12866960 cg26538782 cg03310874 cg13093111 cg19502671 cg02796279 cg08538581 
##  0.1972752  0.1968136  0.1948999  0.1940120  0.1940068  0.1937292  0.1932172 
## cg06978145 cg20429104 cg26923863 cg26783127 cg10967114 cg06721411 cg19380303 
##  0.1923578  0.1920274  0.1915369  0.1911030  0.1910478  0.1910046  0.1905641 
## cg23403750 cg01024962 cg26305504 cg06585734 cg04664897 cg12555086 cg20966357 
##  0.1902843  0.1900853  0.1899795  0.1895845  0.1893594  0.1892563  0.1891488 
## cg23755933 cg07918933 cg11153071 cg05445326 cg24250902 cg04217515 cg18835472 
##  0.1888468  0.1884887  0.1882545  0.1879730  0.1878966  0.1878767  0.1877763 
## cg25596405 cg27340480 cg23334433 cg02244028 cg26298914 cg11802666 cg17117981 
##  0.1877148  0.1868540  0.1866464  0.1865790  0.1860314  0.1854544  0.1853498 
## cg04836151 cg00009088 cg04586126 cg08428868 cg26889953 cg14189391 cg08400494 
##  0.1852400  0.1848756  0.1847763  0.1846981  0.1846446  0.1840610  0.1839995 
## cg15288326 cg13980609 cg17086773 cg22875823 cg07417146 cg00235484 cg22666015 
##  0.1838200  0.1837114  0.1834981  0.1830439  0.1826333  0.1820731  0.1819755 
## cg01622606 cg22662844 
##  0.1817174  0.1817050

Constructing tissue specific DMPs

Alternatively, users can set tissueSpecDMPs as TRUE to obtain tissue based DMPs by using the reference atlas. Here, we used simulated samples since we don’t have tissue specific samples. Note that both tissueSpecCpGs and tissueSpecDMPs can’t be TRUE at the same time.

data("HumanCellTypeMethAtlas")
exampleSamples <- simulateCellMix(1,reference = HumanCellTypeMethAtlas)$simulated
exampleMeta <- data.table("Experiment_accession" = "example_sample",
                          "Biosample_term_name" = "example_cell_type")
colnames(exampleSamples) <- c("CpGs", "example_sample")
colnames(HumanCellTypeMethAtlas)[1] <- c("CpGs")

signatures <- findSignatures(
  samples = exampleSamples,
  sampleMeta = exampleMeta,
  atlas = HumanCellTypeMethAtlas,
  IDs = "CpGs", tissueSpecDMPs = TRUE)
## CELL TYPES IN EXTENDED ATLAS: 
## example_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_cervix
print(head(signatures[[2]]))
##            intercept        f         pval         qval
## cg10480329 -3.287511 211.8281 2.089936e-13 1.259813e-09
## cg06297318 -3.496336 180.6364 1.155970e-12 3.484093e-09
## cg18835472 -3.350244 155.6390 5.567206e-12 1.118637e-08
## cg05923857 -3.515739 147.4746 9.761387e-12 1.471041e-08
## cg27366072 -2.858145 143.0843 1.335151e-11 1.494108e-08
## cg15633603 -1.496293 141.5991 1.487168e-11 1.494108e-08

Example Workflow For RNA Sequencing Data

It 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

sessionInfo()
## R version 4.4.2 (2024-10-31)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.1 LTS
## 
## Matrix products: default
## BLAS:   /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3 
## LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so;  LAPACK version 3.12.0
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=en_US.UTF-8        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       
## 
## time zone: Etc/UTC
## tzcode source: system (glibc)
## 
## attached base packages:
## [1] parallel  stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
## [1] dplyr_1.1.4       doParallel_1.0.17 iterators_1.0.14  foreach_1.5.2    
## [5] deconvR_1.13.0    data.table_1.16.2 knitr_1.49        BiocStyle_2.35.0 
## 
## loaded via a namespace (and not attached):
##   [1] splines_4.4.2               BiocIO_1.17.1              
##   [3] bitops_1.0-9                tibble_3.2.1               
##   [5] R.oo_1.27.0                 preprocessCore_1.69.0      
##   [7] XML_3.99-0.17               lifecycle_1.0.4            
##   [9] lattice_0.22-6              MASS_7.3-61                
##  [11] base64_2.0.2                scrime_1.3.5               
##  [13] magrittr_2.0.3              minfi_1.53.1               
##  [15] limma_3.63.2                sass_0.4.9                 
##  [17] rmarkdown_2.29              jquerylib_0.1.4            
##  [19] yaml_2.3.10                 robslopes_1.1.3            
##  [21] doRNG_1.8.6                 askpass_1.2.1              
##  [23] minqa_1.2.8                 DBI_1.2.3                  
##  [25] buildtools_1.0.0            RColorBrewer_1.1-3         
##  [27] abind_1.4-8                 zlibbioc_1.52.0            
##  [29] quadprog_1.5-8              GenomicRanges_1.59.1       
##  [31] purrr_1.0.2                 R.utils_2.12.3             
##  [33] BiocGenerics_0.53.3         RCurl_1.98-1.16            
##  [35] GenomeInfoDbData_1.2.13     IRanges_2.41.1             
##  [37] S4Vectors_0.45.2            maketools_1.3.1            
##  [39] rentrez_1.2.3               genefilter_1.89.0          
##  [41] annotate_1.85.0             DelayedMatrixStats_1.29.0  
##  [43] codetools_0.2-20            DelayedArray_0.33.2        
##  [45] xml2_1.3.6                  tidyselect_1.2.1           
##  [47] UCSC.utils_1.3.0            lme4_1.1-35.5              
##  [49] beanplot_1.3.1              matrixStats_1.4.1          
##  [51] stats4_4.4.2                illuminaio_0.49.0          
##  [53] GenomicAlignments_1.43.0    jsonlite_1.8.9             
##  [55] multtest_2.63.0             e1071_1.7-16               
##  [57] survival_3.7-0              bbmle_1.0.25.1             
##  [59] tools_4.4.2                 rsq_2.7                    
##  [61] Rcpp_1.0.13-1               glue_1.8.0                 
##  [63] SparseArray_1.7.2           xfun_0.49                  
##  [65] qvalue_2.39.0               MatrixGenerics_1.19.0      
##  [67] GenomeInfoDb_1.43.2         HDF5Array_1.35.1           
##  [69] numDeriv_2016.8-1.1         BiocManager_1.30.25        
##  [71] fastmap_1.2.0               boot_1.3-31                
##  [73] rhdf5filters_1.19.0         fansi_1.0.6                
##  [75] openssl_2.2.2               digest_0.6.37              
##  [77] R6_2.5.1                    colorspace_2.1-1           
##  [79] gtools_3.9.5                RSQLite_2.3.8              
##  [81] R.methodsS3_1.8.2           utf8_1.2.4                 
##  [83] tidyr_1.3.1                 generics_0.1.3             
##  [85] rtracklayer_1.67.0          class_7.3-22               
##  [87] httr_1.4.7                  S4Arrays_1.7.1             
##  [89] pkgconfig_2.0.3             gtable_0.3.6               
##  [91] blob_1.2.4                  siggenes_1.81.0            
##  [93] XVector_0.47.0              sys_3.4.3                  
##  [95] htmltools_0.5.8.1           scales_1.3.0               
##  [97] Biobase_2.67.0              png_0.1-8                  
##  [99] deming_1.4-1                tzdb_0.4.0                 
## [101] reshape2_1.4.4              rjson_0.2.23               
## [103] nloptr_2.1.1                coda_0.19-4.1              
## [105] nlme_3.1-166                curl_6.0.1                 
## [107] bdsmatrix_1.3-7             bumphunter_1.49.0          
## [109] proxy_0.4-27                cachem_1.1.0               
## [111] rhdf5_2.51.0                stringr_1.5.1              
## [113] AnnotationDbi_1.69.0        restfulr_0.0.15            
## [115] GEOquery_2.75.0             pillar_1.9.0               
## [117] grid_4.4.2                  reshape_0.8.9              
## [119] vctrs_0.6.5                 Deriv_4.1.6                
## [121] xtable_1.8-4                evaluate_1.0.1             
## [123] readr_2.1.5                 GenomicFeatures_1.59.1     
## [125] mvtnorm_1.3-2               cli_3.6.3                  
## [127] locfit_1.5-9.10             compiler_4.4.2             
## [129] Rsamtools_2.23.1            rlang_1.1.4                
## [131] crayon_1.5.3                rngtools_1.5.2             
## [133] nor1mix_1.3-3               mclust_6.1.1               
## [135] emdbook_1.3.13              plyr_1.8.9                 
## [137] stringi_1.8.4               mcr_1.3.3.1                
## [139] BiocParallel_1.41.0         nnls_1.6                   
## [141] assertthat_0.2.1            munsell_0.5.1              
## [143] Biostrings_2.75.1           Matrix_1.7-1               
## [145] hms_1.1.3                   sparseMatrixStats_1.19.0   
## [147] bit64_4.5.2                 ggplot2_3.5.1              
## [149] Rhdf5lib_1.29.0             KEGGREST_1.47.0            
## [151] methylKit_1.33.1            statmod_1.5.0              
## [153] SummarizedExperiment_1.37.0 fastseg_1.53.2             
## [155] memoise_2.0.1               bslib_0.8.0                
## [157] bit_4.5.0
stopCluster(cl)

  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↩︎