FlowSorted.Blood.EPIC

Lucas A Salas

2023-10-26

Illumina Human Methylation data from EPIC on immunomagnetic sorted adult blood cell populations. The FlowSorted.Blood.EPIC package contains Illumina HumanMethylationEPIC (EPIC)) DNA methylation microarray data from the immunomethylomics group (manuscript submitted), consisting of 37 magnetic sorted blood cell references and 12 samples, formatted as an RGChannelSet object for integration and normalization using most of the existing Bioconductor packages.

This package contains data similar to the FlowSorted.Blood.450k package consisting of data from peripheral blood samples generated from adult men and women. However, when using the newer EPIC microarray minfi estimates of cell type composition using the FlowSorted.Blood.450k package are less precise compared to actual cell counts. Hence, this package consists of appropriate data for deconvolution of adult blood samples used in for example EWAS relying in the newer EPIC technology.

Researchers may find this package useful as these samples represent different cellular populations ( T lymphocytes (CD4+ and CD8+), B cells (CD19+), monocytes (CD14+), NK cells (CD56+) and Neutrophils of cell sorted blood generated with high purity estimates. As a test of accuracy 12 experimental mixtures were reconstructed using fixed amounts of DNA from purified cells.

Objects included:
1. FlowSorted.Blood.EPIC is the RGChannelSet object containing the reference library

library(FlowSorted.Blood.EPIC)
#> Loading required package: minfi
#> Loading required package: BiocGenerics
#> 
#> Attaching package: 'BiocGenerics'
#> The following objects are masked from 'package:stats':
#> 
#>     IQR, mad, sd, var, xtabs
#> The following objects are masked from 'package:base':
#> 
#>     Filter, Find, Map, Position, Reduce, anyDuplicated, aperm, append,
#>     as.data.frame, basename, cbind, colnames, dirname, do.call,
#>     duplicated, eval, evalq, get, grep, grepl, intersect, is.unsorted,
#>     lapply, mapply, match, mget, order, paste, pmax, pmax.int, pmin,
#>     pmin.int, rank, rbind, rownames, sapply, setdiff, sort, table,
#>     tapply, union, unique, unsplit, which.max, which.min
#> Loading required package: GenomicRanges
#> Loading required package: stats4
#> Loading required package: S4Vectors
#> 
#> Attaching package: 'S4Vectors'
#> The following object is masked from 'package:utils':
#> 
#>     findMatches
#> The following objects are masked from 'package:base':
#> 
#>     I, expand.grid, unname
#> Loading required package: IRanges
#> Loading required package: GenomeInfoDb
#> Loading required package: SummarizedExperiment
#> Loading required package: MatrixGenerics
#> Loading required package: matrixStats
#> 
#> Attaching package: 'MatrixGenerics'
#> The following objects are masked from 'package:matrixStats':
#> 
#>     colAlls, colAnyNAs, colAnys, colAvgsPerRowSet, colCollapse,
#>     colCounts, colCummaxs, colCummins, colCumprods, colCumsums,
#>     colDiffs, colIQRDiffs, colIQRs, colLogSumExps, colMadDiffs,
#>     colMads, colMaxs, colMeans2, colMedians, colMins, colOrderStats,
#>     colProds, colQuantiles, colRanges, colRanks, colSdDiffs, colSds,
#>     colSums2, colTabulates, colVarDiffs, colVars, colWeightedMads,
#>     colWeightedMeans, colWeightedMedians, colWeightedSds,
#>     colWeightedVars, rowAlls, rowAnyNAs, rowAnys, rowAvgsPerColSet,
#>     rowCollapse, rowCounts, rowCummaxs, rowCummins, rowCumprods,
#>     rowCumsums, rowDiffs, rowIQRDiffs, rowIQRs, rowLogSumExps,
#>     rowMadDiffs, rowMads, rowMaxs, rowMeans2, rowMedians, rowMins,
#>     rowOrderStats, rowProds, rowQuantiles, rowRanges, rowRanks,
#>     rowSdDiffs, rowSds, rowSums2, rowTabulates, rowVarDiffs, rowVars,
#>     rowWeightedMads, rowWeightedMeans, rowWeightedMedians,
#>     rowWeightedSds, rowWeightedVars
#> Loading required package: Biobase
#> Welcome to Bioconductor
#> 
#>     Vignettes contain introductory material; view with
#>     'browseVignettes()'. To cite Bioconductor, see
#>     'citation("Biobase")', and for packages 'citation("pkgname")'.
#> 
#> Attaching package: 'Biobase'
#> The following object is masked from 'package:MatrixGenerics':
#> 
#>     rowMedians
#> The following objects are masked from 'package:matrixStats':
#> 
#>     anyMissing, rowMedians
#> Loading required package: Biostrings
#> Loading required package: XVector
#> 
#> Attaching package: 'Biostrings'
#> The following object is masked from 'package:base':
#> 
#>     strsplit
#> Loading required package: bumphunter
#> Loading required package: foreach
#> Loading required package: iterators
#> Loading required package: parallel
#> Loading required package: locfit
#> locfit 1.5-9.8    2023-06-11
#> Setting options('download.file.method.GEOquery'='auto')
#> Setting options('GEOquery.inmemory.gpl'=FALSE)
#> Loading required package: ExperimentHub
#> Loading required package: AnnotationHub
#> Loading required package: BiocFileCache
#> Loading required package: dbplyr
#> 
#> Attaching package: 'AnnotationHub'
#> The following object is masked from 'package:Biobase':
#> 
#>     cache
FlowSorted.Blood.EPIC <- libraryDataGet('FlowSorted.Blood.EPIC')
#> see ?FlowSorted.Blood.EPIC and browseVignettes('FlowSorted.Blood.EPIC') for documentation
#> loading from cache
FlowSorted.Blood.EPIC
#> class: RGChannelSet 
#> dim: 1051815 49 
#> metadata(0):
#> assays(2): Green Red
#> rownames(1051815): 1600101 1600111 ... 99810990 99810992
#> rowData names(0):
#> colnames(49): 201868500150_R01C01 201868500150_R03C01 ...
#>   201870610111_R06C01 201870610111_R07C01
#> colData names(32): Sample_Plate Sample_Well ... filenames normalmix
#> Annotation
#>   array: IlluminaHumanMethylationEPIC
#>   annotation: ilm10b4.hg19

The raw dataset is hosted in both ExperimentHub (EH1136) and GEO GSE110554

  1. IDOLOptimizedCpGs the IDOL L-DMR library for EPIC arrays
data("IDOLOptimizedCpGs") 
head(IDOLOptimizedCpGs)  
  1. IDOLOptimizedCpGs450klegacy the IDOL L-DMR library for legacy 450k arrays
data("IDOLOptimizedCpGs450klegacy") 
head(IDOLOptimizedCpGs450klegacy)  

See the object help files for additional information

estimateCellCounts2 function for cell type deconvolution:

We offer the function estimateCellCounts2 a modification of the popular estimatesCellCounts function in minfi. Our function corrected small glitches when dealing with combining the DataFrame objects with the reference libraries. We allow the use of MethylSet objects such as those from GEO. However, we offer only Quantile normalization for those datasets (assuming that they have not been previously normalized). The estimates are calculated using the CP/QP method described in Houseman et al. 2012. and adapted in minfi. CIBERSORT and RPC are allowed using external packages.
see ?estimateCellCounts2 for details

# Step 1: Load the reference library to extract the artificial mixtures
FlowSorted.Blood.EPIC <- libraryDataGet('FlowSorted.Blood.EPIC')
#> see ?FlowSorted.Blood.EPIC and browseVignettes('FlowSorted.Blood.EPIC') for documentation
#> loading from cache
FlowSorted.Blood.EPIC
#> class: RGChannelSet 
#> dim: 1051815 49 
#> metadata(0):
#> assays(2): Green Red
#> rownames(1051815): 1600101 1600111 ... 99810990 99810992
#> rowData names(0):
#> colnames(49): 201868500150_R01C01 201868500150_R03C01 ...
#>   201870610111_R06C01 201870610111_R07C01
#> colData names(32): Sample_Plate Sample_Well ... filenames normalmix
#> Annotation
#>   array: IlluminaHumanMethylationEPIC
#>   annotation: ilm10b4.hg19

library(minfi)
# Note: If your machine does not allow access to internet you can download
# and save the file. Then load it manually into the R environment


# Step 2 separate the reference from the testing dataset if you want to run 
# examples for estimations for this function example

RGsetTargets <- FlowSorted.Blood.EPIC[,
             FlowSorted.Blood.EPIC$CellType == "MIX"]
             
sampleNames(RGsetTargets) <- paste(RGsetTargets$CellType,
                            seq_len(dim(RGsetTargets)[2]), sep = "_")
RGsetTargets
#> class: RGChannelSet 
#> dim: 1051815 12 
#> metadata(0):
#> assays(2): Green Red
#> rownames(1051815): 1600101 1600111 ... 99810990 99810992
#> rowData names(0):
#> colnames(12): MIX_1 MIX_2 ... MIX_11 MIX_12
#> colData names(32): Sample_Plate Sample_Well ... filenames normalmix
#> Annotation
#>   array: IlluminaHumanMethylationEPIC
#>   annotation: ilm10b4.hg19

# Step 3: use your favorite package for deconvolution.
# Deconvolve a target data set consisting of EPIC DNA methylation 
# data profiled in blood, using your prefered method.

# You can use our IDOL optimized DMR library for the EPIC array.  This object
# contains a vector of length 450 consisting of the IDs of the IDOL optimized
# CpG probes.  These CpGs are used as the backbone for deconvolution and were
# selected because their methylation signature differs across the six normal 
# leukocyte subtypes. Use the option "IDOL"

head (IDOLOptimizedCpGs)
#> [1] "cg08769189" "cg07661835" "cg00219921" "cg13468685" "cg04329870"
#> [6] "cg14085952"
# If you need to deconvolve a 450k legacy dataset use 
# IDOLOptimizedCpGs450klegacy instead

# We recommend using Noob processMethod = "preprocessNoob" in minfi for the 
# target and reference datasets. 
# Cell types included are "CD8T", "CD4T", "NK", "Bcell", "Mono", "Neu"

# To use the IDOL optimized list of CpGs (IDOLOptimizedCpGs) use 
# estimateCellCounts2 an adaptation of the popular estimateCellCounts in 
# minfi. This function also allows including customized reference arrays. 

# Do not run with limited RAM the normalization step requires a big amount 
# of memory resources

 propEPIC<-estimateCellCounts2(RGsetTargets, compositeCellType = "Blood", 
                                processMethod = "preprocessNoob",
                                probeSelect = "IDOL", 
                                cellTypes = c("CD8T", "CD4T", "NK", "Bcell", 
                                "Mono", "Neu"))
#> see ?FlowSorted.Blood.EPIC and browseVignettes('FlowSorted.Blood.EPIC') for documentation
#> loading from cache
#> [estimateCellCounts2] Combining user data with reference (flow sorted) data.
#> [estimateCellCounts2] Processing user and reference data together.
#> Loading required package: IlluminaHumanMethylationEPICmanifest
#> Loading required package: IlluminaHumanMethylationEPICanno.ilm10b4.hg19
#> [estimateCellCounts2] Using IDOL L-DMR probes for composition estimation.
#> [estimateCellCounts2] Estimating proportion composition (prop), if you provide cellcounts those will be provided as counts in the composition estimation.
                                
print(head(propEPIC$prop))
#>         CD8T   CD4T     NK  Bcell   Mono    Neu
#> MIX_1 0.1916 0.0708 0.1515 0.1904 0.1901 0.2118
#> MIX_2 0.0463 0.1763 0.0179 0.0442 0.0582 0.6699
#> MIX_3 0.0675 0.1012 0.0043 0.0223 0.1092 0.7033
#> MIX_4 0.1208 0.1781 0.0200 0.0167 0.0764 0.6032
#> MIX_5 0.2891 0.1603 0.1525 0.0718 0.2227 0.1126
#> MIX_6 0.0963 0.1576 0.0269 0.0243 0.0703 0.6401
percEPIC<-round(propEPIC$prop*100,1)

Advanced user deconvolution CP/QP, CIBERSORT and/or RPC deconvolution


noobset<- preprocessNoob(RGsetTargets)
#or from estimateCellCounts2 returnAll=TRUE

 propEPIC<-projectCellType_CP (
 getBeta(noobset)[IDOLOptimizedCpGs,],
 IDOLOptimizedCpGs.compTable, contrastWBC=NULL, nonnegative=TRUE,
 lessThanOne=FALSE)

print(head(propEPIC))
#>         CD8T   CD4T     NK  Bcell   Mono    Neu
#> MIX_1 0.1916 0.0708 0.1515 0.1904 0.1901 0.2118
#> MIX_2 0.0463 0.1763 0.0179 0.0442 0.0582 0.6699
#> MIX_3 0.0675 0.1012 0.0043 0.0223 0.1092 0.7033
#> MIX_4 0.1208 0.1781 0.0200 0.0167 0.0764 0.6032
#> MIX_5 0.2891 0.1603 0.1525 0.0718 0.2227 0.1126
#> MIX_6 0.0963 0.1576 0.0269 0.0243 0.0703 0.6401
percEPIC<-round(propEPIC*100,1)


# If you prefer CIBERSORT or RPC deconvolution use EpiDISH or similar

# Example not to run

# library(EpiDISH)
# RPC <- epidish(getBeta(noobset)[IDOLOptimizedCpGs,],
# IDOLOptimizedCpGs.compTable, method = "RPC")
# RPC$estF#RPC proportion estimates
# percEPICRPC<-round(RPC$estF*100,1)#percentages
# 
# CBS <- epidish(getBeta(noobset)[IDOLOptimizedCpGs,],
# IDOLOptimizedCpGs.compTable, method = "CBS")
# CBS$estF#CBS proportion estimates
# percEPICCBS<-round(CBS$estF*100,1)#percentages

Umbilical Cord Blood

# # UMBILICAL CORD BLOOD DECONVOLUTION
# 
# library (FlowSorted.CordBloodCombined.450k)
# # Step 1: Load the reference library to extract the umbilical cord samples
# FlowSorted.CordBloodCombined.450k <- 
#     libraryDataGet('FlowSorted.CordBloodCombined.450k') 
# 
# FlowSorted.CordBloodCombined.450k
# 
# # Step 2 separate the reference from the testing dataset if you want to run 
# # examples for estimations for this function example
# 
# RGsetTargets <- FlowSorted.CordBloodCombined.450k[,
# FlowSorted.CordBloodCombined.450k$CellType == "WBC"]
# sampleNames(RGsetTargets) <- paste(RGsetTargets$CellType,
#                               seq_len(dim(RGsetTargets)[2]), sep = "_")
# RGsetTargets
# 
# # Step 3: use your favorite package for deconvolution.
# # Deconvolve a target data set consisting of 450K DNA methylation 
# # data profiled in blood, using your prefered method.
# # You can use our IDOL optimized DMR library for the Cord Blood,  This object
# # contains a vector of length 517 consisting of the IDs of the IDOL optimized
# # CpG probes.  These CpGs are used as the backbone for deconvolution and were
# # selected because their methylation signature differs across the six normal 
# # leukocyte subtypes plus the nucleated red blood cells.
# 
# # We recommend using Noob processMethod = "preprocessNoob" in minfi for the 
# # target and reference datasets. 
# # Cell types included are "CD8T", "CD4T", "NK", "Bcell", "Mono", "Gran", 
# # "nRBC"
# # To use the IDOL optimized list of CpGs (IDOLOptimizedCpGsCordBlood) use 
# # estimateCellCounts2 from FlowSorted.Blood.EPIC. 
# # Do not run with limited RAM the normalization step requires a big amount 
# # of memory resources. Use the parameters as specified below for 
# # reproducibility.
# # 
# 
#     propUCB<-estimateCellCounts2(RGsetTargets,
#                                     compositeCellType =
#                                                "CordBloodCombined",
#                                     processMethod = "preprocessNoob",
#                                     probeSelect = "IDOL",
#                                     cellTypes = c("CD8T", "CD4T", "NK",
#                                     "Bcell", "Mono", "Gran", "nRBC"))
# 
#     head(propUCB$prop)
#     percUCB<-round(propUCB$prop*100,1)

Using cell counts instead of proportions. Note: These are random numbers, not the actual cell counts of the experiment

# library(FlowSorted.Blood.450k)
# RGsetTargets2 <- FlowSorted.Blood.450k[,
#                              FlowSorted.Blood.450k$CellType == "WBC"]
# sampleNames(RGsetTargets2) <- paste(RGsetTargets2$CellType,
#                              seq_len(dim(RGsetTargets2)[2]), sep = "_")
# RGsetTargets2
# propEPIC2<-estimateCellCounts2(RGsetTargets2, compositeCellType = "Blood",
#                              processMethod = "preprocessNoob",
#                              probeSelect = "IDOL",
#                              cellTypes = c("CD8T", "CD4T", "NK", "Bcell",
#                              "Mono", "Neu"), cellcounts = rep(10000,6))
# head(propEPIC2$prop)
# head(propEPIC2$counts)
# percEPIC2<-round(propEPIC2$prop*100,1)

Blood Extended deconvolution

# # Blood Extended deconvolution or any external reference
# #please contact <Technology.Transfer@dartmouth.edu>
# 
# # Do not run
# library (FlowSorted.BloodExtended.EPIC)
# # 
# # Step 1: Extract the mix samples
# 
# FlowSorted.Blood.EPIC <- libraryDataGet('FlowSorted.Blood.EPIC')
# 
# # Step 2 separate the reference from the testing dataset if you want to run 
# # examples for estimations for this function example
# 
# RGsetTargets <- FlowSorted.Blood.EPIC[,
# FlowSorted.Blood.EPIC$CellType == "MIX"]
# sampleNames(RGsetTargets) <- paste(RGsetTargets$CellType,
#                               seq_len(dim(RGsetTargets)[2]), sep = "_")
# RGsetTargets
# 
# # Step 3: use your favorite package for deconvolution.
# # Deconvolve the target data set 450K or EPIC blood DNA methylation. 
# # We recommend ONLY the IDOL method, the automatic method can lead to severe
# # biases.
# 
# # We recommend using Noob processMethod = "preprocessNoob" in minfi for the 
# # target and reference datasets. 
# # Cell types included are "Bas", "Bmem", "Bnv", "CD4mem", "CD4nv", 
# # "CD8mem", "CD8nv", "Eos", "Mono", "Neu", "NK", and "Treg"
# # Use estimateCellCounts2 from FlowSorted.Blood.EPIC. 
# # Do not run with limited RAM the normalization step requires a big amount 
# # of memory resources. Use the parameters as specified below for 
# # reproducibility.
# # 
# 
#     prop_ext<-estimateCellCounts2(RGsetTargets,
#                                     compositeCellType =
#                                                "BloodExtended",
#                                     processMethod = "preprocessNoob",
#                                     probeSelect = "IDOL",
#                                     cellTypes = c("Bas", "Bmem", "Bnv",
#                                                "CD4mem", "CD4nv",
#                                               "CD8mem", "CD8nv", "Eos",
#                                               "Mono", "Neu", "NK", "Treg"),
#     CustomCpGs =if(RGsetTargets@annotation[1]=="IlluminaHumanMethylationEPIC"){
#     IDOLOptimizedCpGsBloodExtended}else{IDOLOptimizedCpGsBloodExtended450k})
# 
#    perc_ext<-round(prop_ext$prop*100,1)
#    head(perc_ext)
sessionInfo()
#> R version 4.3.1 (2023-06-16)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 22.04.3 LTS
#> 
#> Matrix products: default
#> BLAS:   /home/biocbuild/bbs-3.18-bioc/R/lib/libRblas.so 
#> LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.0
#> 
#> 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       
#> 
#> time zone: America/New_York
#> tzcode source: system (glibc)
#> 
#> attached base packages:
#> [1] parallel  stats4    stats     graphics  grDevices utils     datasets 
#> [8] methods   base     
#> 
#> other attached packages:
#>  [1] IlluminaHumanMethylationEPICanno.ilm10b4.hg19_0.6.0
#>  [2] IlluminaHumanMethylationEPICmanifest_0.3.0         
#>  [3] FlowSorted.Blood.EPIC_2.6.0                        
#>  [4] ExperimentHub_2.10.0                               
#>  [5] AnnotationHub_3.10.0                               
#>  [6] BiocFileCache_2.10.0                               
#>  [7] dbplyr_2.3.4                                       
#>  [8] minfi_1.48.0                                       
#>  [9] bumphunter_1.44.0                                  
#> [10] locfit_1.5-9.8                                     
#> [11] iterators_1.0.14                                   
#> [12] foreach_1.5.2                                      
#> [13] Biostrings_2.70.1                                  
#> [14] XVector_0.42.0                                     
#> [15] SummarizedExperiment_1.32.0                        
#> [16] Biobase_2.62.0                                     
#> [17] MatrixGenerics_1.14.0                              
#> [18] matrixStats_1.0.0                                  
#> [19] GenomicRanges_1.54.0                               
#> [20] GenomeInfoDb_1.38.0                                
#> [21] IRanges_2.36.0                                     
#> [22] S4Vectors_0.40.0                                   
#> [23] BiocGenerics_0.48.0                                
#> 
#> loaded via a namespace (and not attached):
#>   [1] RColorBrewer_1.1-3            jsonlite_1.8.7               
#>   [3] magrittr_2.0.3                GenomicFeatures_1.54.0       
#>   [5] rmarkdown_2.25                BiocIO_1.12.0                
#>   [7] zlibbioc_1.48.0               vctrs_0.6.4                  
#>   [9] multtest_2.58.0               memoise_2.0.1                
#>  [11] Rsamtools_2.18.0              DelayedMatrixStats_1.24.0    
#>  [13] RCurl_1.98-1.12               askpass_1.2.0                
#>  [15] htmltools_0.5.6.1             S4Arrays_1.2.0               
#>  [17] progress_1.2.2                curl_5.1.0                   
#>  [19] Rhdf5lib_1.24.0               SparseArray_1.2.0            
#>  [21] rhdf5_2.46.0                  sass_0.4.7                   
#>  [23] nor1mix_1.3-0                 bslib_0.5.1                  
#>  [25] plyr_1.8.9                    cachem_1.0.8                 
#>  [27] GenomicAlignments_1.38.0      mime_0.12                    
#>  [29] lifecycle_1.0.3               pkgconfig_2.0.3              
#>  [31] Matrix_1.6-1.1                R6_2.5.1                     
#>  [33] fastmap_1.1.1                 shiny_1.7.5.1                
#>  [35] GenomeInfoDbData_1.2.11       digest_0.6.33                
#>  [37] siggenes_1.76.0               reshape_0.8.9                
#>  [39] AnnotationDbi_1.64.0          RSQLite_2.3.1                
#>  [41] base64_2.0.1                  filelock_1.0.2               
#>  [43] fansi_1.0.5                   httr_1.4.7                   
#>  [45] abind_1.4-5                   compiler_4.3.1               
#>  [47] beanplot_1.3.1                rngtools_1.5.2               
#>  [49] withr_2.5.1                   bit64_4.0.5                  
#>  [51] BiocParallel_1.36.0           DBI_1.1.3                    
#>  [53] HDF5Array_1.30.0              biomaRt_2.58.0               
#>  [55] MASS_7.3-60                   openssl_2.1.1                
#>  [57] rappdirs_0.3.3                DelayedArray_0.28.0          
#>  [59] rjson_0.2.21                  tools_4.3.1                  
#>  [61] interactiveDisplayBase_1.40.0 httpuv_1.6.12                
#>  [63] glue_1.6.2                    quadprog_1.5-8               
#>  [65] restfulr_0.0.15               promises_1.2.1               
#>  [67] nlme_3.1-163                  rhdf5filters_1.14.0          
#>  [69] grid_4.3.1                    generics_0.1.3               
#>  [71] tzdb_0.4.0                    preprocessCore_1.64.0        
#>  [73] tidyr_1.3.0                   data.table_1.14.8            
#>  [75] hms_1.1.3                     xml2_1.3.5                   
#>  [77] utf8_1.2.4                    BiocVersion_3.18.0           
#>  [79] pillar_1.9.0                  stringr_1.5.0                
#>  [81] limma_3.58.0                  later_1.3.1                  
#>  [83] genefilter_1.84.0             splines_4.3.1                
#>  [85] dplyr_1.1.3                   lattice_0.22-5               
#>  [87] survival_3.5-7                rtracklayer_1.62.0           
#>  [89] bit_4.0.5                     GEOquery_2.70.0              
#>  [91] annotate_1.80.0               tidyselect_1.2.0             
#>  [93] knitr_1.44                    xfun_0.40                    
#>  [95] scrime_1.3.5                  statmod_1.5.0                
#>  [97] stringi_1.7.12                yaml_2.3.7                   
#>  [99] evaluate_0.22                 codetools_0.2-19             
#> [101] tibble_3.2.1                  BiocManager_1.30.22          
#> [103] cli_3.6.1                     xtable_1.8-4                 
#> [105] jquerylib_0.1.4               Rcpp_1.0.11                  
#> [107] png_0.1-8                     XML_3.99-0.14                
#> [109] ellipsis_0.3.2                readr_2.1.4                  
#> [111] blob_1.2.4                    prettyunits_1.2.0            
#> [113] mclust_6.0.0                  doRNG_1.8.6                  
#> [115] sparseMatrixStats_1.14.0      bitops_1.0-7                 
#> [117] illuminaio_0.44.0             purrr_1.0.2                  
#> [119] crayon_1.5.2                  rlang_1.1.1                  
#> [121] KEGGREST_1.42.0

References

LA Salas et al. (2018). An optimized library for reference-based deconvolution of whole-blood biospecimens assayed using the Illumina HumanMethylationEPIC BeadArray. Genome Biology 19, 64. doi: 10.1186/s13059-018-1448-7.

LA Salas et al. (2022). . Nat Comm 13, 761 (2022). doi: 10.1038/s41467-021-27864-7.

DC Koestler et al. (2016). Improving cell mixture deconvolution by identifying optimal DNA methylation libraries (IDOL). BMC bioinformatics. 17, 120. doi: 10.1186/s12859-016-0943-7.

K Gervin, LA Salas et al. (2019) . Clin Epigenetics 11,125. doi: 10.1186/s13148-019-0717-y.

EA Houseman et al. (2012) DNA methylation arrays as surrogate measures of cell mixture distribution. BMC Bioinformatics 13, 86. doi: 10.1186/1471-2105-13-86.

minfi Tools to analyze & visualize Illumina Infinium methylation arrays.