deconvR 1.14.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
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.
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.0410
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
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.17109691 0.1769050 0.16327324
## 2 cg25913761 0.31218039 0.3261950 0.30300862
## 3 cg26955540 0.19076966 0.2023801 0.18435168
## 4 cg25170017 0.23086994 0.2358932 0.23022638
## 5 cg12827637 0.13949576 0.1428054 0.14025353
## 6 cg19442545 0.06565946 0.0541015 0.05441349
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.17042505 0.8866 0.2615 0.0149
## 2 cg25913761 0.31379466 0.8363 0.2210 0.2816
## 3 cg26955540 0.19250047 0.7658 0.0222 0.1492
## 4 cg25170017 0.23232985 0.8861 0.5116 0.1021
## 5 cg12827637 0.14085155 0.5212 0.3614 0.0227
## 6 cg19442545 0.05805815 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
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
## cg08425810 cg14845962 cg15014975 cg07865091 cg16402452 cg02297709 cg24591770
## 0.10429018 0.10122011 0.10032584 0.09864293 0.09838995 0.09836362 0.09822001
## cg11532431 cg21181453 cg04462378 cg06398251 cg00769843 cg06889535 cg11025609
## 0.09808618 0.09778135 0.09774968 0.09767473 0.09717395 0.09699778 0.09680839
## cg15185986 cg16017089 cg15575375 cg16306706 cg08659179 cg12461469 cg15059065
## 0.09662739 0.09649041 0.09634534 0.09623870 0.09619830 0.09583274 0.09553131
## cg13852284 cg16389209 cg12946660 cg27016579 cg26401541 cg20227056 cg22109530
## 0.09519153 0.09446071 0.09372705 0.09368668 0.09345740 0.09341975 0.09322179
## cg06493154 cg13506158 cg10576280 cg07275218 cg17084151 cg26860970 cg04962480
## 0.09294297 0.09266027 0.09256645 0.09220819 0.09217500 0.09172374 0.09157816
## cg07520810 cg16401465 cg00928816 cg11402700 cg27269488 cg03803541 cg20746552
## 0.09134166 0.09026783 0.09024196 0.08949115 0.08893629 0.08829809 0.08750974
## cg25843651 cg10906729 cg07748583 cg00270878 cg11630632 cg05337711 cg17179862
## 0.08729710 0.08725741 0.08620712 0.08599091 0.08592595 0.08572027 0.08550487
## cg03484180 cg24019564 cg04173586 cg10672681 cg26264318 cg10807894 cg03419014
## 0.08525131 0.08463745 0.08454423 0.08358678 0.08311039 0.08295728 0.08247893
## cg01615880 cg09942210 cg18842353 cg19147129 cg09367425 cg01005506 cg14419102
## 0.08146106 0.08121097 0.08108209 0.08099492 0.08088796 0.08058658 0.08000364
## cg02670123 cg23231163 cg02486646 cg08979515 cg25607383 cg06653052 cg17936564
## 0.07996189 0.07971349 0.07932359 0.07892877 0.07778166 0.07777632 0.07762934
## cg14515791 cg23977888 cg13512830 cg05103231 cg20510285 cg01409343 cg18684142
## 0.07714342 0.07649495 0.07626389 0.07623266 0.07619434 0.07586378 0.07546813
## cg25124966 cg22706007 cg22365240 cg13909895 cg23051299 cg15935721 cg09226051
## 0.07471482 0.07436235 0.07369521 0.07351274 0.07324142 0.07317815 0.07296740
## cg26313247 cg27292818 cg13442969 cg01403948 cg02676175 cg19828220 cg16439823
## 0.07252414 0.07222293 0.07149571 0.07108186 0.07094632 0.07086620 0.07081143
## cg02016130 cg08061524 cg11006453 cg13298116 cg04778755 cg21869046 cg24441810
## 0.07061565 0.07060385 0.07034352 0.07020023 0.07012499 0.07003364 0.06971749
## cg15926420 cg18353226
## 0.06924026 0.06905023
##
## $hypo
## cg01622606 cg16429499 cg18066690 cg14516183 cg26709988 cg23514375 cg17976229
## 0.16598708 0.16539376 0.16107431 0.15918588 0.15844331 0.15487044 0.15473744
## cg19849557 cg12019801 cg09050670 cg15982099 cg05413628 cg13414750 cg10826999
## 0.15426895 0.15208423 0.15179151 0.15173502 0.15168431 0.15156208 0.15077394
## cg11222173 cg06857116 cg09856274 cg13618516 cg22160073 cg01332882 cg13984928
## 0.15059001 0.15054855 0.14823613 0.14751482 0.14644478 0.14554592 0.14457750
## cg01335658 cg13466002 cg09452568 cg24037651 cg25283432 cg03339817 cg14028598
## 0.14453338 0.14439208 0.14429492 0.14275591 0.14261213 0.14198218 0.13962273
## cg07023317 cg00471371 cg17723958 cg10624395 cg06623899 cg05921138 cg05085169
## 0.13883890 0.13860000 0.13729045 0.13686289 0.13565480 0.13541455 0.13359005
## cg17090611 cg07651316 cg16732654 cg25605731 cg21492378 cg04917472 cg00754552
## 0.13343963 0.13316417 0.13271556 0.13085954 0.12922796 0.12876432 0.12862493
## cg24011073 cg16993108 cg14642045 cg11348106 cg19590591 cg25045219 cg20700740
## 0.12685616 0.12516551 0.12364640 0.12263384 0.12084032 0.11988079 0.11881167
## cg00058449 cg19835973 cg14345154 cg01742370 cg06934654 cg26514623 cg04422741
## 0.11867162 0.11850538 0.11807590 0.11805850 0.11766901 0.11662184 0.11530227
## cg12422199 cg25752703 cg10239163 cg18095675 cg12754495 cg13924996 cg04784315
## 0.11493830 0.11477349 0.11475909 0.11446624 0.11426607 0.11210994 0.11176288
## cg23128584 cg26918117 cg09085220 cg07248377 cg02368812 cg25771026 cg14564351
## 0.11145239 0.11013614 0.10883371 0.10762609 0.10730527 0.10608672 0.10501015
## cg16736889 cg06545464 cg15029631 cg02891728 cg12116137 cg25574765 cg05573654
## 0.10495746 0.10493697 0.10490914 0.10464540 0.10457244 0.10454661 0.10410375
## cg11597277 cg16802439 cg19055828 cg02691035 cg23249922 cg16695253 cg21781157
## 0.10389893 0.10315978 0.10198604 0.10198473 0.10070579 0.10022829 0.09978205
## cg19382157 cg08939850 cg21341821 cg01476222 cg05418105 cg10584478 cg25598890
## 0.09768150 0.09764355 0.09722000 0.09578929 0.09568398 0.09519430 0.09496494
## cg15609237 cg08672140 cg20118717 cg11284147 cg13498610 cg07986378 cg08415156
## 0.09456834 0.09423596 0.09400663 0.09381761 0.09307965 0.09239012 0.09217793
## cg00263619 cg05059607
## 0.09197081 0.09190752
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 171.5821 1.995346e-12 1.202794e-08
## cg06297318 -3.496336 156.0490 5.416128e-12 1.632421e-08
## cg18835472 -3.350244 132.8440 2.866307e-11 5.759366e-08
## cg05923857 -3.515739 127.5169 4.351795e-11 6.558156e-08
## cg27366072 -2.858145 118.4588 9.168791e-11 9.298537e-08
## cg08428868 -4.139304 117.1090 1.028748e-10 9.298537e-08
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()
## R version 4.5.0 RC (2025-04-04 r88126)
## Platform: aarch64-apple-darwin20
## Running under: macOS Ventura 13.7.1
##
## Matrix products: default
## BLAS: /Library/Frameworks/R.framework/Versions/4.5-arm64/Resources/lib/libRblas.0.dylib
## LAPACK: /Library/Frameworks/R.framework/Versions/4.5-arm64/Resources/lib/libRlapack.dylib; LAPACK version 3.12.1
##
## locale:
## [1] C/en_US.UTF-8/en_US.UTF-8/C/en_US.UTF-8/en_US.UTF-8
##
## time zone: America/New_York
## tzcode source: internal
##
## 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.14.0 data.table_1.17.0 knitr_1.50 BiocStyle_2.36.0
##
## loaded via a namespace (and not attached):
## [1] splines_4.5.0 BiocIO_1.18.0
## [3] bitops_1.0-9 tibble_3.2.1
## [5] R.oo_1.27.0 preprocessCore_1.70.0
## [7] XML_3.99-0.18 lifecycle_1.0.4
## [9] Rdpack_2.6.4 lattice_0.22-7
## [11] MASS_7.3-65 base64_2.0.2
## [13] scrime_1.3.5 magrittr_2.0.3
## [15] minfi_1.54.0 limma_3.64.0
## [17] sass_0.4.10 rmarkdown_2.29
## [19] jquerylib_0.1.4 yaml_2.3.10
## [21] robslopes_1.1.3 doRNG_1.8.6.2
## [23] askpass_1.2.1 minqa_1.2.8
## [25] DBI_1.2.3 RColorBrewer_1.1-3
## [27] abind_1.4-8 quadprog_1.5-8
## [29] GenomicRanges_1.60.0 purrr_1.0.4
## [31] R.utils_2.13.0 BiocGenerics_0.54.0
## [33] RCurl_1.98-1.17 GenomeInfoDbData_1.2.14
## [35] IRanges_2.42.0 S4Vectors_0.46.0
## [37] rentrez_1.2.3 genefilter_1.90.0
## [39] annotate_1.86.0 DelayedMatrixStats_1.30.0
## [41] codetools_0.2-20 DelayedArray_0.34.0
## [43] xml2_1.3.8 tidyselect_1.2.1
## [45] UCSC.utils_1.4.0 lme4_1.1-37
## [47] beanplot_1.3.1 matrixStats_1.5.0
## [49] stats4_4.5.0 illuminaio_0.50.0
## [51] GenomicAlignments_1.44.0 jsonlite_2.0.0
## [53] multtest_2.64.0 e1071_1.7-16
## [55] survival_3.8-3 bbmle_1.0.25.1
## [57] tools_4.5.0 rsq_2.7
## [59] Rcpp_1.0.14 glue_1.8.0
## [61] SparseArray_1.8.0 xfun_0.52
## [63] qvalue_2.40.0 MatrixGenerics_1.20.0
## [65] GenomeInfoDb_1.44.0 HDF5Array_1.36.0
## [67] numDeriv_2016.8-1.1 BiocManager_1.30.25
## [69] fastmap_1.2.0 boot_1.3-31
## [71] rhdf5filters_1.20.0 openssl_2.3.2
## [73] digest_0.6.37 R6_2.6.1
## [75] colorspace_2.1-1 gtools_3.9.5
## [77] RSQLite_2.3.9 R.methodsS3_1.8.2
## [79] h5mread_1.0.0 tidyr_1.3.1
## [81] generics_0.1.3 rtracklayer_1.68.0
## [83] class_7.3-23 httr_1.4.7
## [85] S4Arrays_1.8.0 pkgconfig_2.0.3
## [87] gtable_0.3.6 blob_1.2.4
## [89] siggenes_1.82.0 XVector_0.48.0
## [91] htmltools_0.5.8.1 bookdown_0.43
## [93] scales_1.3.0 Biobase_2.68.0
## [95] png_0.1-8 reformulas_0.4.0
## [97] deming_1.4-1 tzdb_0.5.0
## [99] reshape2_1.4.4 rjson_0.2.23
## [101] nloptr_2.2.1 coda_0.19-4.1
## [103] nlme_3.1-168 curl_6.2.2
## [105] bdsmatrix_1.3-7 bumphunter_1.50.0
## [107] proxy_0.4-27 cachem_1.1.0
## [109] rhdf5_2.52.0 stringr_1.5.1
## [111] AnnotationDbi_1.70.0 restfulr_0.0.15
## [113] GEOquery_2.76.0 pillar_1.10.2
## [115] grid_4.5.0 reshape_0.8.9
## [117] vctrs_0.6.5 Deriv_4.1.6
## [119] xtable_1.8-4 evaluate_1.0.3
## [121] readr_2.1.5 GenomicFeatures_1.60.0
## [123] mvtnorm_1.3-3 cli_3.6.4
## [125] locfit_1.5-9.12 compiler_4.5.0
## [127] Rsamtools_2.24.0 rlang_1.1.6
## [129] crayon_1.5.3 rngtools_1.5.2
## [131] nor1mix_1.3-3 mclust_6.1.1
## [133] emdbook_1.3.13 plyr_1.8.9
## [135] stringi_1.8.7 mcr_1.3.3.1
## [137] BiocParallel_1.42.0 nnls_1.6
## [139] assertthat_0.2.1 munsell_0.5.1
## [141] Biostrings_2.76.0 Matrix_1.7-3
## [143] hms_1.1.3 sparseMatrixStats_1.20.0
## [145] bit64_4.6.0-1 ggplot2_3.5.2
## [147] Rhdf5lib_1.30.0 KEGGREST_1.48.0
## [149] methylKit_1.34.0 statmod_1.5.0
## [151] SummarizedExperiment_1.38.0 fastseg_1.54.0
## [153] rbibutils_2.3 memoise_2.0.1
## [155] bslib_0.9.0 bit_4.6.0
stopCluster(cl)