combi package: vignette

Introduction

This package implements a novel data integration model for sample-wise integration of different views. It accounts for compositionality and employs a non-parametric mean-variance trend for sequence count data. The resulting model can be conveniently plotted to allow for explorative visualization of variability shared over different views.

Installation

The package can be installed and loaded using the following commands:

library(BiocManager)
BiocManager::install("combi", update = FALSE)
library(devtools)
install_github("CenterForStatistics-UGent/combi")
suppressPackageStartupMessages(library(combi))
cat("combi package version", as.character(packageVersion("combi")), "\n")
## combi package version 1.12.1
data(Zhang)

Unconstrained integration

For an unconstrained ordination, a named list of datasets with overlapping samples must be supplied. The datasets can currently be supplied as a raw data matrix (with features in the columns), or as a phyloseq, SummarizedExperiment or ExpressionSet object. In addition, information on the required distribution (“quasi” for quasi-likelihood fitting, “gaussian” for normal data) and compositional nature (TRUE/FALSE) should be supplied

microMetaboInt = combi(list(microbiome = zhangMicrobio, metabolomics = zhangMetabo),
    distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE)

One can print basic infor about the ordination

microMetaboInt
## Unconstrained combi ordination of 2 dimensions on 2 views with 42 samples.
## Views and number of features were:
##  microbiome: 130
##  metabolomics: 174
##  Importance parameters of dimensions 1 to 2 are 117 and 44.9

A simple plot function is available for the result, for samples and shapes, a data frame should also be supplied

plot(microMetaboInt)

plot(microMetaboInt, samDf = zhangMetavars, samCol = "ABX")

By default, only the most important features (furthest away from the origin) are shown. To show all features, one can resort to point cloud plots or density plots as follows:

plot(microMetaboInt, samDf = zhangMetavars, samCol = "ABX", featurePlot = "points")

plot(microMetaboInt, samDf = zhangMetavars, samCol = "ABX", featurePlot = "density")

The drawback is that now no feature labels are shown.

Adding projections

As an aid to interpretation of compositional views, links between features can be plotted and projected onto samples by providing their names or approximate coordinates

# First define the plot, and return the coordinates
mmPlot = plot(microMetaboInt, samDf = zhangMetavars, samCol = "ABX", returnCoords = TRUE,
    featNum = 10)
# Providing feature names, and sample coordinates, but any combination is
# allowed
addLink(mmPlot, links = cbind("Staphylococcus_819c11", "OTU929ffc"), Views = 1, samples = c(0,
    1))

Coordinates

Finally, one can extract the coordinates for use in third-party software

coords = extractCoords(microMetaboInt, Dim = c(1, 2))

Constrained integration

For a constrained ordination also a data frame of sample variables should be supplied

microMetaboIntConstr = combi(list(microbiome = zhangMicrobio, metabolomics = zhangMetabo),
    distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE,
    covariates = zhangMetavars)

Also here we can get a quick overview

microMetaboIntConstr
## Constrained combi ordination of 2 dimensions on 2 views with 42 samples.
## Views and number of features were:
##  microbiome: 130
##  metabolomics: 174
##  Number of sample variables included was 4,
## for which 6 parameters were estimated per dimension.
##  Importance parameters of dimensions 1 to 2 are 34.2 and 21.4

and plot the ordination

plot(microMetaboIntConstr, samDf = zhangMetavars, samCol = "ABX")

Diagnostics

Convergence of the iterative algorithm can be assessed as follows:

convPlot(microMetaboInt)

Influence of the different views can be investigated through

inflPlot(microMetaboInt, samples = 1:20, plotType = "boxplot")

FAQ

Why are not all my samples shown in the constrained ordination?

Confusion often arises as to why less distinct sample dots are shown than there are samples in the constrained ordination. This occurs when few, categorical constraining variables are supplied as below.

# Linear with only 2 variables
microMetaboIntConstr2Vars = combi(list(microbiome = zhangMicrobio, metabolomics = zhangMetabo),
    distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE,
    covariates = zhangMetavars[, c("Sex", "ABX")])

Every constrained sample score is a linear combination of constraining variables, and in this case there are only 2x2=4 distinct sample scores possible, leading to the sample dots from the same gender and treatment to be plotted on top of each other.

plot(microMetaboIntConstr2Vars, samDf = zhangMetavars, samCol = "ABX")

In general, it is best to include all measured sample variables in a constrained analysis, and let the combi-algorithm find out which ones are the most important drivers of variability.

The combi function crashes, what should I do

The combi method works by iteratively solving systems of non-linear equations through Newton-Raphson. This can lead to numerical instability, with errors like “infinite values returned by jacobian” in the nleqslv function. This class of problems has no general solution, but is mainly a matter of trial and error. Following two things can be tried:

  1. Tweaking the prevCutOff and minFraction parameters to remove more sparse features.
  2. Tweaking the initPower parameter, which will lead to different starting values, and hopefully a solution path that is numerically more stable.
# Linear with only 2 variables
microMetaboTweak = combi(list(microbiome = zhangMicrobio, metabolomics = zhangMetabo),
    distributions = c("quasi", "gaussian"), compositional = c(TRUE, FALSE), logTransformGaussian = FALSE,
    initPower = 1.5, minFraction = 0.25, prevCutOff = 0.8)

Session info

This vignette was generated with following version of R:

sessionInfo()
## R version 4.3.1 (2023-06-16)
## 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] stats     graphics  grDevices utils     datasets  methods   base     
## 
## other attached packages:
## [1] combi_1.12.1 DBI_1.1.3   
## 
## loaded via a namespace (and not attached):
##   [1] ade4_1.7-22                 tidyselect_1.2.0           
##   [3] farver_2.1.1                dplyr_1.1.2                
##   [5] Biostrings_2.68.1           bitops_1.0-7               
##   [7] fastmap_1.1.1               RCurl_1.98-1.12            
##   [9] phyloseq_1.44.0             digest_0.6.33              
##  [11] lifecycle_1.0.3             cluster_2.1.4              
##  [13] survival_3.5-5              magrittr_2.0.3             
##  [15] compiler_4.3.1              rlang_1.1.1                
##  [17] sass_0.4.7                  tools_4.3.1                
##  [19] igraph_1.5.0.1              utf8_1.2.3                 
##  [21] yaml_2.3.7                  data.table_1.14.8          
##  [23] knitr_1.43                  BB_2019.10-1               
##  [25] labeling_0.4.2              S4Arrays_1.0.5             
##  [27] DelayedArray_0.26.7         alabama_2022.4-1           
##  [29] plyr_1.8.8                  abind_1.4-5                
##  [31] withr_2.5.0                 nleqslv_3.3.4              
##  [33] numDeriv_2016.8-1.1         BiocGenerics_0.46.0        
##  [35] grid_4.3.1                  stats4_4.3.1               
##  [37] fansi_1.0.4                 multtest_2.56.0            
##  [39] biomformat_1.28.0           colorspace_2.1-0           
##  [41] Rhdf5lib_1.22.0             ggplot2_3.4.2              
##  [43] scales_1.2.1                iterators_1.0.14           
##  [45] MASS_7.3-60                 isoband_0.2.7              
##  [47] SummarizedExperiment_1.30.2 cli_3.6.1                  
##  [49] vegan_2.6-4                 rmarkdown_2.23             
##  [51] crayon_1.5.2                generics_0.1.3             
##  [53] reshape2_1.4.4              ape_5.7-1                  
##  [55] cachem_1.0.8                rhdf5_2.44.0               
##  [57] stringr_1.5.0               zlibbioc_1.46.0            
##  [59] splines_4.3.1               parallel_4.3.1             
##  [61] formatR_1.14                XVector_0.40.0             
##  [63] matrixStats_1.0.0           vctrs_0.6.3                
##  [65] Matrix_1.6-0                jsonlite_1.8.7             
##  [67] SparseM_1.81                IRanges_2.34.1             
##  [69] S4Vectors_0.38.1            tensor_1.5                 
##  [71] foreach_1.5.2               limma_3.56.2               
##  [73] jquerylib_0.1.4             glue_1.6.2                 
##  [75] codetools_0.2-19            stringi_1.7.12             
##  [77] gtable_0.3.3                GenomeInfoDb_1.36.1        
##  [79] GenomicRanges_1.52.0        quadprog_1.5-8             
##  [81] munsell_0.5.0               tibble_3.2.1               
##  [83] pillar_1.9.0                cobs_1.3-5                 
##  [85] htmltools_0.5.5             quantreg_5.96              
##  [87] rhdf5filters_1.12.1         GenomeInfoDbData_1.2.10    
##  [89] R6_2.5.1                    evaluate_0.21              
##  [91] lattice_0.21-8              Biobase_2.60.0             
##  [93] highr_0.10                  bslib_0.5.0                
##  [95] MatrixModels_0.5-2          Rcpp_1.0.11                
##  [97] permute_0.9-7               nlme_3.1-162               
##  [99] mgcv_1.9-0                  xfun_0.39                  
## [101] MatrixGenerics_1.12.3       pkgconfig_2.0.3