EnMCB

Xin Yu

2021-10-26

Introduction

This package is designed to help you to create the methylation correlated blocks using methylation profiles. A stacked ensemble of machine learning models, which combined the Cox regression, support vector regression and elastic-net regression model, can be constructed using this package1. You also can choose one of them to build DNA methylation signatures associated with disease progression.

Note: This package is still under developing. Some of the functions may change.

Followings are brief insturctions for using this package:

You can install and test our package by downloading source package.

Installation

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

Useage

First, you need a methylation data set, currently only most common platform ‘Illumina Infinium Human Methylation 450K’ is supported.

You can use your own datasets,or use our demo data.

You can automatically run following:

suppressPackageStartupMessages(library(EnMCB))

methylation_dataset<-create_demo()

res<-IdentifyMCB(methylation_dataset)
## Loading required package: IlluminaHumanMethylation450kanno.ilmn12.hg19
## 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, 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 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.4    2020-03-24

IdentfyMCB() function will calculated Pearson correlation coefficients between the any two CpGs. A value of Pearson correlation coefficients which under the threshold was used to identify boundaries between any two adjacent markers indicating uncorrelated methylation. Markers not separated by a boundary were combined into the MCB.You can extract the MCB information,

MCB<-res$MCBinformation

and select some of MCBs for further modeling.

MCB<-MCB[MCB[,"CpGs_num"]>=5,]

In order to build models, one may run following:

# sample the dataset into training set and testing set
trainingset<-colnames(methylation_dataset) %in% sample(colnames(methylation_dataset),0.6*length(colnames(methylation_dataset)))

testingset<-!trainingset

#build the models
library(survival)
data(demo_survival_data)

models<-metricMCB(MCB,
                    training_set = methylation_dataset[,trainingset],
                    Surv = demo_survival_data[trainingset],
                    Method = "cox",ci = TRUE)

#select the best
onemodel<-models$best_cox_model$cox_model

Then, you can predict the risk by the model you build:

newcgdata<-data.frame(t(methylation_dataset[,testingset]))
           
prediction_results<-predict(onemodel, newcgdata)

In order to build ensemble model, one may run following:

# You can choose one of MCBs:
select_single_one=1

em<-ensemble_model(t(MCB[select_single_one,]),
                    training_set=methylation_dataset[,trainingset],
                    Surv_training=demo_survival_data[trainingset])

Note that this function only can be used for single MCB only, otherwise the precessing time could be very long.

Then, you can predict the risk by the model you build:

em_prediction_results<-ensemble_prediction(ensemble_model = em,
                    prediction_data = methylation_dataset[,testingset])

This function will return the single vector with risk scores predicted by ensemble model.

For detailed information, you can find at our references.

Session Info

sessionInfo()
## R version 4.1.1 (2021-08-10)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.3 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.14-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.14-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  stats4    stats     graphics  grDevices utils     datasets 
## [8] methods   base     
## 
## other attached packages:
##  [1] survival_3.2-13                                   
##  [2] IlluminaHumanMethylation450kanno.ilmn12.hg19_0.6.0
##  [3] minfi_1.40.0                                      
##  [4] bumphunter_1.36.0                                 
##  [5] locfit_1.5-9.4                                    
##  [6] iterators_1.0.13                                  
##  [7] foreach_1.5.1                                     
##  [8] Biostrings_2.62.0                                 
##  [9] XVector_0.34.0                                    
## [10] SummarizedExperiment_1.24.0                       
## [11] Biobase_2.54.0                                    
## [12] MatrixGenerics_1.6.0                              
## [13] matrixStats_0.61.0                                
## [14] GenomicRanges_1.46.0                              
## [15] GenomeInfoDb_1.30.0                               
## [16] IRanges_2.28.0                                    
## [17] S4Vectors_0.32.0                                  
## [18] BiocGenerics_0.40.0                               
## [19] EnMCB_1.6.0                                       
## 
## loaded via a namespace (and not attached):
##   [1] utf8_1.2.2                survivalsvm_0.0.5        
##   [3] rms_6.2-0                 tidyselect_1.1.1         
##   [5] RSQLite_2.2.8             AnnotationDbi_1.56.0     
##   [7] htmlwidgets_1.5.4         grid_4.1.1               
##   [9] BiocParallel_1.28.0       munsell_0.5.0            
##  [11] codetools_0.2-18          preprocessCore_1.56.0    
##  [13] colorspace_2.0-2          filelock_1.0.2           
##  [15] knitr_1.36                rstudioapi_0.13          
##  [17] GenomeInfoDbData_1.2.7    bit64_4.0.5              
##  [19] rhdf5_2.38.0              TH.data_1.1-0            
##  [21] vctrs_0.3.8               generics_0.1.1           
##  [23] xfun_0.27                 BiocFileCache_2.2.0      
##  [25] R6_2.5.1                  doParallel_1.0.16        
##  [27] illuminaio_0.36.0         bitops_1.0-7             
##  [29] rhdf5filters_1.6.0        cachem_1.0.6             
##  [31] reshape_0.8.8             DelayedArray_0.20.0      
##  [33] assertthat_0.2.1          BiocIO_1.4.0             
##  [35] scales_1.1.1              multcomp_1.4-17          
##  [37] nnet_7.3-16               gtable_0.3.0             
##  [39] conquer_1.0.2             sandwich_3.0-1           
##  [41] rlang_0.4.12              MatrixModels_0.5-0       
##  [43] genefilter_1.76.0         splines_4.1.1            
##  [45] rtracklayer_1.54.0        GEOquery_2.62.0          
##  [47] checkmate_2.0.0           yaml_2.2.1               
##  [49] GenomicFeatures_1.46.0    backports_1.2.1          
##  [51] Hmisc_4.6-0               inum_1.0-4               
##  [53] tools_4.1.1               nor1mix_1.3-0            
##  [55] ggplot2_3.3.5             ellipsis_0.3.2           
##  [57] stabs_0.6-4               jquerylib_0.1.4          
##  [59] RColorBrewer_1.1-2        siggenes_1.68.0          
##  [61] Rcpp_1.0.7                plyr_1.8.6               
##  [63] base64enc_0.1-3           sparseMatrixStats_1.6.0  
##  [65] progress_1.2.2            zlibbioc_1.40.0          
##  [67] purrr_0.3.4               RCurl_1.98-1.5           
##  [69] prettyunits_1.1.1         rpart_4.1-15             
##  [71] openssl_1.4.5             zoo_1.8-9                
##  [73] cluster_2.1.2             magrittr_2.0.1           
##  [75] data.table_1.14.2         SparseM_1.81             
##  [77] mvtnorm_1.1-3             mboost_2.9-5             
##  [79] hms_1.1.1                 evaluate_0.14            
##  [81] xtable_1.8-4              XML_3.99-0.8             
##  [83] jpeg_0.1-9                mclust_5.4.7             
##  [85] gridExtra_2.3             shape_1.4.6              
##  [87] compiler_4.1.1            biomaRt_2.50.0           
##  [89] tibble_3.1.5              crayon_1.4.1             
##  [91] htmltools_0.5.2           tzdb_0.1.2               
##  [93] Formula_1.2-4             tidyr_1.1.4              
##  [95] libcoin_1.0-9             DBI_1.1.1                
##  [97] dbplyr_2.1.1              MASS_7.3-54              
##  [99] rappdirs_0.3.3            boot_1.3-28              
## [101] Matrix_1.3-4              readr_2.0.2              
## [103] quadprog_1.5-8            pkgconfig_2.0.3          
## [105] GenomicAlignments_1.30.0  foreign_0.8-81           
## [107] xml2_1.3.2                annotate_1.72.0          
## [109] bslib_0.3.1               rngtools_1.5.2           
## [111] multtest_2.50.0           beanplot_1.2             
## [113] doRNG_1.8.2               scrime_1.3.5             
## [115] stringr_1.4.0             digest_0.6.28            
## [117] pracma_2.3.3              rmarkdown_2.11           
## [119] base64_2.0                htmlTable_2.3.0          
## [121] DelayedMatrixStats_1.16.0 restfulr_0.0.13          
## [123] curl_4.3.2                Rsamtools_2.10.0         
## [125] quantreg_5.86             rjson_0.2.20             
## [127] lifecycle_1.0.1           nlme_3.1-153             
## [129] jsonlite_1.7.2            Rhdf5lib_1.16.0          
## [131] survivalROC_1.0.3         askpass_1.1              
## [133] limma_3.50.0              fansi_0.5.0              
## [135] pillar_1.6.4              lattice_0.20-45          
## [137] KEGGREST_1.34.0           fastmap_1.1.0            
## [139] httr_1.4.2                glue_1.4.2               
## [141] png_0.1-7                 glmnet_4.1-2             
## [143] bit_4.0.4                 stringi_1.7.5            
## [145] sass_0.4.0                HDF5Array_1.22.0         
## [147] nnls_1.4                  blob_1.2.2               
## [149] polspline_1.1.19          partykit_1.2-15          
## [151] latticeExtra_0.6-29       memoise_2.0.0            
## [153] dplyr_1.0.7

References


  1. Xin Yu et al. 2019 Predicting disease progression in lung adenocarcinoma patients based on methylation correlated blocks using ensemble machine learning classifiers (under review)