Contents

1 Using RaMWAS with other methylation platforms or data types

RaMWAS is primarily designed for studies of methylation measurements from enrichment platforms.

However, RaMWAS can also be useful for the analysis of methylation measurements from other platforms (e.g. Illumina HumanMethylation450K array) or other data types such as gene expression levels or genotype information. RaMWAS can perform several analysis steps on such data including: principal component analysis (PCA), association testing (MWAS, TWAS, GWAS), and multimarker analysis with cross validation using the elastic net.

1.1 Import data from other sources

Without external data source at hand, we show how to create and fill data matrices with artificial data. Importing real data can be done in a similar way, with random data generation replaced with reading data from existing sources.

We create data files in the same format as produced by Step 3 of RaMWAS.

These files include

  • CpG_locations.* – filematrix with the location of the CpGs / SNPs / gene trascription start sites.
    It must have two columns with integer values – chromosome number and location (chr and position).
  • CpG_chromosome_names.txt – file with chromosome names (factor levels) for the integer column chr in the location filematrix.
  • Coverage.* – filematrix with the data for all samples and all locations.
    Each row has data for a single sample. Row names must be sample names.
    Each column has data for a single location (CpG / SNP / gene trascription start site). Columns must match rows of the location filematrix.

First, we load the package and set up a working directory. The project directory dr can be set to a more convenient location when running the code.

library(ramwas)

# work in a temporary directory
dr = paste0(tempdir(), "/simulated_matrix_data")
dir.create(dr, showWarnings = FALSE)
cat(dr,"\n")
## /tmp/Rtmp11V0Mo/simulated_matrix_data

Let the sample data matrix have 200 samples and 100,000 variables.

nsamples = 200
nvariables = 100000

For these 200 samples we generate a data frame with age and sex phenotypes and a batch effect covariate.

covariates = data.frame(
    sample = paste0("Sample_",seq_len(nsamples)),
    sex = seq_len(nsamples) %% 2,
    age = runif(nsamples, min = 20, max = 80),
    batch = paste0("batch",(seq_len(nsamples) %% 3))
)
pander(head(covariates))
sample sex age batch
Sample_1 1 71.5 batch1
Sample_2 0 35.8 batch2
Sample_3 1 60.4 batch0
Sample_4 0 64.5 batch1
Sample_5 1 28.4 batch2
Sample_6 0 26.3 batch0

Next, we create the genomic locations for 100,000 variables.

temp = cumsum(sample(20e7 / nvariables, nvariables, replace = TRUE) + 0)
chr      = as.integer(temp %/% 1e7) + 1L
position = as.integer(temp %% 1e7)

locmat = cbind(chr = chr, position = position)
chrnames = paste0("chr", 1:10)
pander(head(locmat))
chr position
1 958
1 1850
1 2916
1 4390
1 5386
1 6104

Now we save locations in a filematrix and create a text file with chromosome names.

fmloc = fm.create.from.matrix(
            filenamebase = paste0(dr,"/CpG_locations"),
            mat = locmat)
close(fmloc)
writeLines(
        con = paste0(dr,"/CpG_chromosome_names.txt"),
        text = chrnames)

Finally, we create data matrix. We include sex effect in 225 variables and age effect in 16 variables out of each 2000. Each variable is also affected by noise and batch effects.

fm = fm.create(paste0(dr,"/Coverage"), nrow = nsamples, ncol = nvariables)

# Row names of the matrix are set to sample names
rownames(fm) = as.character(covariates$sample)

# The matrix is filled, 2000 variables at a time
byrows = 2000
for( i in seq_len(nvariables/byrows) ){ # i=1
    slice = matrix(runif(nsamples*byrows), nrow = nsamples, ncol = byrows)
    slice[,  1:225] = slice[,  1:225] + covariates$sex / 30 / sd(covariates$sex)
    slice[,101:116] = slice[,101:116] + covariates$age / 10 / sd(covariates$age)
    slice = slice + ((as.integer(factor(covariates$batch))+i) %% 3) / 40
    fm[,(1:byrows) + byrows*(i-1)] = slice
}
close(fm)

1.2 Principal Component Analysis (PCA)

To run PCA with RaMWAS we specify three parameters:

  • dircoveragenorm – directory with the data matrix
  • covariates – data frame with covariates
  • modelcovariates – names of covariates to regress out
param = ramwasParameters(
    dircoveragenorm = dr,
    covariates = covariates,
    modelcovariates = NULL
)

Now we run PCA.

ramwas4PCA(param)

The top several PCs are marginally distinct from the rest.

pfull = parameterPreprocess(param)
eigenvalues = fm.load(paste0(pfull$dirpca, "/eigenvalues"))
eigenvectors = fm.open(
                filenamebase = paste0(pfull$dirpca, "/eigenvectors"),
                readonly = TRUE)
plotPCvalues(eigenvalues)

plotPCvectors(eigenvectors[,1], 1)

plotPCvectors(eigenvectors[,2], 2)

plotPCvectors(eigenvectors[,3], 3)

plotPCvectors(eigenvectors[,4], 4)

close(eigenvectors)

There are strong correlations between top PCs with sex, age, and batch covariates.
Note, for the categorical covariate (batch) the table shows R2 instead of correlations.

# Get the directory with PCA results
pfull = parameterPreprocess(param)
tblcr = read.table(
            file = paste0(pfull$dirpca, "/PC_vs_covs_corr.txt"),
            header = TRUE,
            sep = "\t")
pander(head(tblcr, 5))
name sex age batch_R2
PC1 -0.0278 -0.0991 0.984
PC2 0.0326 0.0372 0.986
PC3 -0.938 -0.163 0.00167
PC4 0.286 -0.942 0.000988
PC5 -0.0126 0.0021 8.94e-05

The p-values for these correlations and R2 show that the top two PCs are correlated with sex and age while a number of other PCs are affected by sample batch effects.

pfull = parameterPreprocess(param)
tblpv = read.table(
            file = paste0(pfull$dirpca, "/PC_vs_covs_pvalue.txt"),
            header = TRUE,
            sep = "\t")
pander(head(tblpv, 5))
name sex age batch_R2
PC1 0.696 0.163 1.11e-178
PC2 0.647 0.601 1.6e-183
PC3 2.55e-93 0.0211 0.848
PC4 4.06e-05 6.37e-96 0.907
PC5 0.86 0.976 0.991

1.3 PCA with batch regressed out

It is common to regress out batch and lab-technical effects from the data in the analysis.

Let’s regress out batch in our example by changing modelcovariates parameter.

param$modelcovariates = "batch"

ramwas4PCA(param)

The p-values for association between PCs and covariates changed slightly:

# Get the directory with PCA results
pfull = parameterPreprocess(param)
tblpv = read.table(
            file = paste0(pfull$dirpca, "/PC_vs_covs_pvalue.txt"),
            header = TRUE,
            sep = "\t")
pander(head(tblpv, 5))
name sex age batch_R2
PC1 4.54e-93 0.0185 NA
PC2 4.41e-05 1.25e-98 NA
PC3 0.86 0.997 NA
PC4 0.852 0.584 NA
PC5 0.883 0.692 NA

Note that the PCs are now orthogonal to the batch effects and thus the corresponding p-values all equal to 1.

1.4 Association testing

Let us test for association between variables in the data matrix and the sex covariate (modeloutcome parameter) correcting for batch effects (modelcovariates parameter). Save top 20 results (toppvthreshold parameter) in a text file.

param$modelcovariates = "batch"
param$modeloutcome = "sex"
param$toppvthreshold = 20

ramwas5MWAS(param)

The QQ-plot shows mild enrichment among a large number of variables, which is consistent with how the data was generated – 22% of variables are affected by sex.

mwas = getMWAS(param)
qqPlotFast(mwas$`p-value`)
title(pfull$qqplottitle)

The top finding saved in the text file are:

# Get the directory with testing results
pfull = parameterPreprocess(param)
toptbl = read.table(
            file = paste0(pfull$dirmwas,"/Top_tests.txt"),
            header = TRUE,
            sep = "\t")
pander(head(toptbl, 5))
chr start end cor t.test p.value q.value beta
chr9 3943528 3943529 0.373 5.63 6.12e-08 0.00612 0.631
chr2 7934644 7934645 0.347 5.18 5.46e-07 0.0162 0.582
chr9 5974532 5974533 0.345 5.15 6.45e-07 0.0162 0.607
chr1 207416 207417 0.345 5.15 6.47e-07 0.0162 0.629
chr2 5946217 5946218 0.342 5.09 8.34e-07 0.0163 0.589

1.5 Further steps of RaMWAS pipeline

Steps 6 and 7 of RaMWAS pipeline can also be applied to the data matrix exactly as described in the overview vignette.

1.6 Cleanup

Here we remove all the files created by the code above.

unlink(paste0(dr,"/*"), recursive=TRUE)

2 Version information

sessionInfo()
## R version 4.3.0 RC (2023-04-13 r84269)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 22.04.2 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.17-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] stats4    stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
##  [1] BSgenome.Ecoli.NCBI.20080805_1.3.1000 BSgenome_1.68.0                      
##  [3] rtracklayer_1.60.0                    Biostrings_2.68.0                    
##  [5] XVector_0.40.0                        GenomicRanges_1.52.0                 
##  [7] GenomeInfoDb_1.36.0                   IRanges_2.34.0                       
##  [9] S4Vectors_0.38.0                      BiocGenerics_0.46.0                  
## [11] ramwas_1.24.0                         filematrix_1.3                       
## [13] pander_0.6.5                          knitr_1.42                           
## [15] BiocStyle_2.28.0                     
## 
## loaded via a namespace (and not attached):
##  [1] tidyselect_1.2.0            dplyr_1.1.2                
##  [3] blob_1.2.4                  filelock_1.0.2             
##  [5] bitops_1.0-7                fastmap_1.1.1              
##  [7] RCurl_1.98-1.12             BiocFileCache_2.8.0        
##  [9] GenomicAlignments_1.36.0    XML_3.99-0.14              
## [11] digest_0.6.31               lifecycle_1.0.3            
## [13] survival_3.5-5              KEGGREST_1.40.0            
## [15] RSQLite_2.3.1               magrittr_2.0.3             
## [17] compiler_4.3.0              rlang_1.1.0                
## [19] sass_0.4.5                  progress_1.2.2             
## [21] tools_4.3.0                 utf8_1.2.3                 
## [23] yaml_2.3.7                  prettyunits_1.1.1          
## [25] curl_5.0.0                  bit_4.0.5                  
## [27] DelayedArray_0.26.0         xml2_1.3.3                 
## [29] BiocParallel_1.34.0         KernSmooth_2.23-20         
## [31] grid_4.3.0                  fansi_1.0.4                
## [33] iterators_1.0.14            biomaRt_2.56.0             
## [35] SummarizedExperiment_1.30.0 cli_3.6.1                  
## [37] rmarkdown_2.21              crayon_1.5.2               
## [39] generics_0.1.3              rjson_0.2.21               
## [41] httr_1.4.5                  DBI_1.1.3                  
## [43] cachem_1.0.7                stringr_1.5.0              
## [45] zlibbioc_1.46.0             splines_4.3.0              
## [47] parallel_4.3.0              AnnotationDbi_1.62.0       
## [49] restfulr_0.0.15             BiocManager_1.30.20        
## [51] matrixStats_0.63.0          vctrs_0.6.2                
## [53] glmnet_4.1-7                Matrix_1.5-4               
## [55] jsonlite_1.8.4              bookdown_0.33              
## [57] hms_1.1.3                   bit64_4.0.5                
## [59] magick_2.7.4                foreach_1.5.2              
## [61] jquerylib_0.1.4             glue_1.6.2                 
## [63] codetools_0.2-19            stringi_1.7.12             
## [65] shape_1.4.6                 BiocIO_1.10.0              
## [67] tibble_3.2.1                pillar_1.9.0               
## [69] rappdirs_0.3.3              htmltools_0.5.5            
## [71] GenomeInfoDbData_1.2.10     R6_2.5.1                   
## [73] dbplyr_2.3.2                evaluate_0.20              
## [75] lattice_0.21-8              Biobase_2.60.0             
## [77] highr_0.10                  png_0.1-8                  
## [79] Rsamtools_2.16.0            memoise_2.0.1              
## [81] bslib_0.4.2                 Rcpp_1.0.10                
## [83] xfun_0.39                   MatrixGenerics_1.12.0      
## [85] pkgconfig_2.0.3