omicRexposome
is an R package designed to work join with
rexposome
. The aim of omicRexposome
is to
perform analysis joining exposome and omic data with the goal to find
the relationship between a single or set of exposures (external
exposome) and the behavior of a gene, a group of CpGs, the level of a
protein, etc. Also to provide a series of tools to analyse exposome and
omic data using standard methods from Biocondcutor.
omicRexposome
is currently in development and not
available from CRAN nor Bioconductor. Anyway, the package can be
installed by using devtools R package and taking the source from
Bioinformatic Research Group in Epidemiology’s GitHub
repository.
This can be done by opening an R session and typing the following code:
User must take into account that this sentence do not install the packages’ dependencies.
Two different types of analyses can be done with
omicRexposome
:
Analysis | omicRexposome function |
---|---|
Association Study | association |
Integration Study | crossomics |
Both association and integration studies are based in objects of
class MultiDataSet
. A MultiDataSet
object is a
contained for multiple layers of sample information. Once the exposome
data and the omics data are encapsulated in a MultiDataSet
the object can be used for both association and integration studies.
The method association
requires a
MultiDataSet
object having to types of information: the
exposome data from an ExposomeSet
object and omic
information from objects of class ExpressionSet
,
MethylationSet
, SummarizedExperiment
or
others. ExposomeSet
objects are created with functions
read_exposome
and load_exposome
from
rexposome
R package (see next section Loading Exposome
Data) and encapsulates exposome data. The method
crossomics
expects a MultiDataSet
with any
number of different data-sets (at last two). Compared with
association
method, crossomics
do not requires
an ExposomeSet
.
In order to illustrate the capabilities of omicRexposome
and the exposome-omic analysis pipeline, we will use the data from
BRGdata
package. This package includes different omic-sets
including methylation, transcriptome and proteome data-sets and an
exposome-data set.
omicRexposome
and MultiDataSet
R packages
are loaded using the standard library command:
The association studies are performed using the method
association
. This method requires, at last four,
augments:
object
should be filled with a
MultiDataSet
object.formula
should be filled with an expression
containing the covariates used to adjust the model.expset
should be filled with the name that the
exposome-set receives in the MultiDataSet
object.omicset
should be filled with the name that
the omic-set receives in the MultiDataSet
object.The argument formula
should follow the pattern:
~sex+age
. The method association
will fill the
formula placing the exposures in the ExposomeSet
m between
~
and the covariates sex+age
.
association
implements the limma
pipeline
using lmFit
and eBayes
in the extraction
methods from MultiDataSet
. The method takes care of the
missing data in exposures, outcomes and omics data and locating and is
subsets both data-sets, exposome data and omic data, by common samples.
The argument method
allows to select the fitting method
used in lmFit
. By default it takes the value
"ls"
for least squares but it can also takes
"robust"
for robust regression.
The following subsections illustrates the usage of
association
with different types of omics data:
methylome, transcriptome and proteome.
First we get the exposome data from BRGdata
package that
we will use in the whole section.
## [1] "ExposomeSet"
## attr(,"package")
## [1] "rexposome"
The aim of this analysis is to perform an association test between
the gene expression levels and the exposures. So the first point is to
obtain the transcriptome data from the brgedata
package.
The association studies between exposures and transcriptome are done
in the same way that the ones with methylome. The method used is
association
, that takes as input an object of
MultiDataSet
class with both exposome and expression
data.
mds <- createMultiDataSet()
mds <- add_genexp(mds, brge_gexp)
mds <- add_exp(mds, brge_expo)
gexp <- association(mds, formula=~Sex+Age,
expset = "exposures", omicset = "expression")
We can have a look to the number of hits and the lambda score of each
analysis with the methods tableHits
and
tableLambda
, seen in the previous section.
## exposure hits lambda
## 1 BPA_p 19 0.9072377
## 2 BPA_t1 27 0.8807316
## 3 BPA_t3 56 0.9391129
## 4 Ben_p 19 0.8013466
## 5 Ben_t1 12 0.8234104
## 6 Ben_t2 9 0.8393350
## 7 Ben_t3 21 0.8301203
## 8 NO2_p 32 1.0281960
## 9 NO2_t1 16 0.7942881
## 10 NO2_t2 35 1.1482314
## 11 NO2_t3 31 0.8770931
## 12 PCB118 59 0.9308472
## 13 PCB138 38 1.0726221
## 14 PCB153 51 1.1743989
## 15 PCB180 17 0.9790750
Since most of all models have a lambda under one, we should consider
use Surrogate Variable Analysis. This can be done using the
same association
method but by setting the argument
sva
to "fast"
so the pipeline of
isva
and SmartSVA
R packages is applied. If
sva
is set to "slow"
the applied. pipeline is
the one from sva
R package.
gexp <- association(mds, formula=~Sex+Age,
expset = "exposures", omicset = "expression", sva = "fast")
We can re-check the results creating the same table than before:
## exposure hits lambda
## 1 BPA_p 50 0.9874143
## 2 BPA_t1 51 0.9433845
## 3 BPA_t3 61 0.9842234
## 4 Ben_p 76 1.0117743
## 5 Ben_t1 64 1.0115556
## 6 Ben_t2 71 1.0089886
## 7 Ben_t3 59 0.9968889
## 8 NO2_p 78 1.0116970
## 9 NO2_t1 67 1.0056894
## 10 NO2_t2 69 1.0209991
## 11 NO2_t3 49 0.9802463
## 12 PCB118 129 1.0518532
## 13 PCB138 67 1.0094398
## 14 PCB153 58 0.9925023
## 15 PCB180 67 0.9974135
The objects of class ResultSet
have a method called
plotAssociation
that allows to create QQ Plots (that are
another useful way to see if there are some inflation/deflation in the
P-Values).
gridExtra::grid.arrange(
plotAssociation(gexp, rid="Ben_p", type="qq") +
ggplot2::ggtitle("Transcriptome - Pb Association"),
plotAssociation(gexp, rid="BPA_p", type="qq") +
ggplot2::ggtitle("Transcriptome - THM Association"),
ncol=2
)
## Warning in ggplot2::geom_text(ggplot2::aes(x = -Inf, y = Inf, hjust = 0, : All aesthetics have length 1, but the data has 67528 rows.
## ℹ Please consider using `annotate()` or provide this layer with data containing
## a single row.
## All aesthetics have length 1, but the data has 67528 rows.
## ℹ Please consider using `annotate()` or provide this layer with data containing
## a single row.
Following this line, the same method plotAssociation
can
be used to create volcano plots.
The proteome data-set included in brgedata
has 47
proteins for 90 samples.
## ExpressionSet (storageMode: lockedEnvironment)
## assayData: 47 features, 90 samples
## element names: exprs
## protocolData: none
## phenoData
## sampleNames: x0001 x0002 ... x0090 (90 total)
## varLabels: age sex
## varMetadata: labelDescription
## featureData
## featureNames: Adiponectin_ok Alpha1AntitrypsinAAT_ok ...
## VitaminDBindingProte_ok (47 total)
## fvarLabels: chr start end
## fvarMetadata: labelDescription
## experimentData: use 'experimentData(object)'
## Annotation:
The association analysis between exposures and proteome is also done
using association
.
mds <- createMultiDataSet()
mds <- add_eset(mds, brge_prot, dataset.type ="proteome")
mds <- add_exp(mds, brge_expo)
prot <- association(mds, formula=~Sex+Age,
expset = "exposures", omicset = "proteome")
The tableHits
indicates that no association was found
between the 47 proteins and the exposures.
## exposure hits
## Ben_p Ben_p 0
## Ben_t1 Ben_t1 0
## Ben_t2 Ben_t2 0
## Ben_t3 Ben_t3 0
## BPA_p BPA_p 0
## BPA_t1 BPA_t1 0
## BPA_t3 BPA_t3 0
## NO2_p NO2_p 0
## NO2_t1 NO2_t1 1
## NO2_t2 NO2_t2 0
## NO2_t3 NO2_t3 0
## PCB118 PCB118 0
## PCB138 PCB138 0
## PCB153 PCB153 0
## PCB180 PCB180 0
This is also seen in the Manhattan plot for proteins that can be
obtained from plotAssociation
.
gridExtra::grid.arrange(
plotAssociation(prot, rid="Ben_p", type="protein") +
ggplot2::ggtitle("Proteome - Cd Association") +
ggplot2::geom_hline(yintercept = 1, color = "LightPink"),
plotAssociation(prot, rid="NO2_p", type="protein") +
ggplot2::ggtitle("Proteome - Cotinine Association") +
ggplot2::geom_hline(yintercept = 1, color = "LightPink"),
ncol=2
)
NOTE: A real Manhattan plot can be draw with
plot
method for ResultSet
objects by setting
the argument type
to "manhattan"
.
## used (Mb) gc trigger (Mb) max used (Mb)
## Ncells 8706510 465.0 15132181 808.2 15132181 808.2
## Vcells 23141672 176.6 347195268 2648.9 433993058 3311.2
omicRexposome
allows to study the relation between
exposures and omic-features from another perspective, different from the
association analyses. The integration analysis can be done, in
omicRexposome
using multi canonical correlation
analysis or using multiple co-inertia analysis. The first
methods is implemented in R package PMA
(CRAN) and the
second in omicade4
R package (Bioconductor). The two
methods are encapsulated in the crossomics
method.
The differences between association
and
crossomics
are that the first method test association
between two complete data-sets, by removing the samples having missing
values in any of the involved data-sets, and the second try to find
latent relationships between two or more sets.
Hence, we need to explore the missing data in the exposome data-set.
This can be done using the methods plotMissings
and
tableMissings
from rexposome
R package.
From the plot we can see that more of the exposures have up to 25% of missing values. Hence the first step in the integration analysis is to avoid missing values. so, we perform a fast imputation on the exposures side:
crossomics
function expects to obtain the different
data-sets in a single labelled-list, in the argument called
list
. The argument method
from
crossomics
function can be set to mcia
(for
multiple co-inertia analysis) or to mcca
(for
multi canonical correlation analysis).
The following code shows how to perform the integration of the
exposome and the proteome. The method crossomics
request a
MultiDataSet
object as input, containing the data-set to be
integrated.
mds <- createMultiDataSet()
mds <- add_genexp(mds, brge_gexp)
mds <- add_eset(mds, brge_prot, dataset.type = "proteome")
mds <- add_exp(mds, brge_expo)
cr_mcia <- crossomics(mds, method = "mcia", verbose = TRUE)
## 'svd' fail to convergence, 'eigen' used to perform singular value decomposition
## Object of class 'ResultSet'
## . created with: crossomics
## . sva:
## . method: mcia ( omicade4 )
## . #results: 1 ( error: 0 )
## . featureData: 3
## . expression: 67528x11
## . proteome: 47x3
## . exposures: 15x12
As can be seen, crossomics
returns an object of class
ResultSet
. In the integration process, the different
data-sets are subset by common samples. This is done taking advantage of
MultiDataSet
capabilities.
The same is done when method is set to mcca
.
We used an extra argument (permute
) into the previous
call to crossomics
using multi canonical correlation
analysis. This argument allows to set the internal argument
corresponding to permutations
and iterations
,
that are used to tune-up internal parameters.
When a ResultSet
is generated using
crossomics
the methods plotHits
,
plotLambda
and plotAssociation
can
NOT be used. But the plotIntegration
will
help us to understand what was done. This method allows to provide the
colors to be used on the plots:
The graphical representation of the results from a multiple co-inertia analysis is a composition of four different plots.
The first plot (first row, first column) is the samples space. It illustrates how the different data-sets are related in terms of intra-sample variability (each data-set has a different color). The second plot (first row, second column) shows the feature space. The features of each set are drawn on the same components so the relation between each data-set can be seen (the features are colored depending of the set were they belong).
The third plot (second row, first column) shows the inertia of each component. The two first plots are drawn on the first and second component. Finally, the fourth plot shows the behavior of the data-sets.
A radar plots is obtained when plotIntegration
is used
on a ResultSet
created though multi canonical
correlation analysis.
## Warning: `position_stack()` requires non-overlapping x intervals.
## `position_stack()` requires non-overlapping x intervals.
## `position_stack()` requires non-overlapping x intervals.
## `position_stack()` requires non-overlapping x intervals.
## Warning: ggrepel: 5 unlabeled data points (too many overlaps). Consider
## increasing max.overlaps
This plot shows the features of the three data-sets in the same 2D space.The relation between the features can be understood by proximity. This means that the features that clusters, or that are in the same quadrant are related and goes in a different direction than the features in the other quadrants.
## 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] stats graphics grDevices utils datasets methods base
##
## other attached packages:
## [1] rexposome_1.29.0 MultiDataSet_1.35.0 omicRexposome_1.29.0
## [4] Biobase_2.67.0 BiocGenerics_0.53.3 generics_0.1.3
## [7] BiocStyle_2.35.0
##
## loaded via a namespace (and not attached):
## [1] splines_4.4.2 norm_1.0-11.1
## [3] bitops_1.0-9 tibble_3.2.1
## [5] omicade4_1.47.0 XML_3.99-0.17
## [7] rpart_4.1.23 lifecycle_1.0.4
## [9] edgeR_4.5.0 lattice_0.22-6
## [11] MASS_7.3-61 flashClust_1.01-2
## [13] backports_1.5.0 magrittr_2.0.3
## [15] limma_3.63.2 Hmisc_5.2-0
## [17] sass_0.4.9 rmarkdown_2.29
## [19] jquerylib_0.1.4 yaml_2.3.10
## [21] minqa_1.2.8 DBI_1.2.3
## [23] buildtools_1.0.0 RColorBrewer_1.1-3
## [25] ade4_1.7-22 abind_1.4-8
## [27] zlibbioc_1.52.0 GenomicRanges_1.59.1
## [29] JADE_2.0-4 nnet_7.3-19
## [31] sandwich_3.1-1 sva_3.55.0
## [33] circlize_0.4.16 GenomeInfoDbData_1.2.13
## [35] IRanges_2.41.1 S4Vectors_0.45.2
## [37] ggrepel_0.9.6 maketools_1.3.1
## [39] genefilter_1.89.0 made4_1.81.0
## [41] SmartSVA_0.1.3 RSpectra_0.16-2
## [43] annotate_1.85.0 PMA_1.2-4
## [45] codetools_0.2-20 DelayedArray_0.33.2
## [47] DT_0.33 gmm_1.8
## [49] tidyselect_1.2.1 shape_1.4.6.1
## [51] farver_2.1.2 UCSC.utils_1.3.0
## [53] lme4_1.1-35.5 matrixStats_1.4.1
## [55] stats4_4.4.2 base64enc_0.1-3
## [57] jsonlite_1.8.9 Formula_1.2-5
## [59] survival_3.7-0 iterators_1.0.14
## [61] emmeans_1.10.5 foreach_1.5.2
## [63] tools_4.4.2 pryr_0.1.6
## [65] Rcpp_1.0.13-1 glue_1.8.0
## [67] gridExtra_2.3 SparseArray_1.7.2
## [69] xfun_0.49 mgcv_1.9-1
## [71] qvalue_2.39.0 MatrixGenerics_1.19.0
## [73] GenomeInfoDb_1.43.1 dplyr_1.1.4
## [75] withr_3.0.2 BiocManager_1.30.25
## [77] fastmap_1.2.0 boot_1.3-31
## [79] fansi_1.0.6 caTools_1.18.3
## [81] digest_0.6.37 R6_2.5.1
## [83] estimability_1.5.1 imputeLCMD_2.1
## [85] colorspace_2.1-1 isva_1.9
## [87] gtools_3.9.5 RSQLite_2.3.8
## [89] utf8_1.2.4 calibrate_1.7.7
## [91] data.table_1.16.2 httr_1.4.7
## [93] htmlwidgets_1.6.4 S4Arrays_1.7.1
## [95] scatterplot3d_0.3-44 pkgconfig_2.0.3
## [97] gtable_0.3.6 blob_1.2.4
## [99] impute_1.81.0 XVector_0.47.0
## [101] sys_3.4.3 htmltools_0.5.8.1
## [103] multcompView_0.1-10 clue_0.3-66
## [105] scales_1.3.0 tmvtnorm_1.6
## [107] leaps_3.2 png_0.1-8
## [109] corrplot_0.95 knitr_1.49
## [111] rstudioapi_0.17.1 reshape2_1.4.4
## [113] nloptr_2.1.1 checkmate_2.3.2
## [115] nlme_3.1-166 zoo_1.8-12
## [117] cachem_1.1.0 GlobalOptions_0.1.2
## [119] stringr_1.5.1 KernSmooth_2.23-24
## [121] parallel_4.4.2 foreign_0.8-87
## [123] AnnotationDbi_1.69.0 pillar_1.9.0
## [125] grid_4.4.2 fastICA_1.2-5.1
## [127] vctrs_0.6.5 gplots_3.2.0
## [129] pcaMethods_1.99.0 xtable_1.8-4
## [131] cluster_2.1.6 htmlTable_2.4.3
## [133] evaluate_1.0.1 qqman_0.1.9
## [135] mvtnorm_1.3-2 cli_3.6.3
## [137] locfit_1.5-9.10 compiler_4.4.2
## [139] rlang_1.1.4 crayon_1.5.3
## [141] labeling_0.4.3 plyr_1.8.9
## [143] stringi_1.8.4 lsr_0.5.2
## [145] BiocParallel_1.41.0 munsell_0.5.1
## [147] Biostrings_2.75.1 glmnet_4.1-8
## [149] Matrix_1.7-1 bit64_4.5.2
## [151] ggplot2_3.5.1 KEGGREST_1.47.0
## [153] statmod_1.5.0 FactoMineR_2.11
## [155] SummarizedExperiment_1.37.0 memoise_2.0.1
## [157] bslib_0.8.0 bit_4.5.0