autonomics is created to make cross-platform omics analysis intuitive and enjoyable. It contains import + preprocess + analyze + visualize one-liners for RNAseq, MS proteomics, SOMAscan and Metabolon platforms and a generic importer for data from any rectangular omics file. With a focus on not only automation but also customization, these importers have a flexible formula/ block/contrastdefs interface which allows to define any statistical formula, fit any general linear model, quantify any contrast, and use random effects or precision weights if required, building on top of the powerful limma platform for statistical analysis. It also offers exquisite support for analyzing complex designs such as the innovative contrastogram visualization to unravel and summarize the differential expression structure in complex designs. By autonomating repetitive tasks, generifying common steps, and intuifying complex designs, it makes cross-platform omics data analysis a fun experience. Try it and enjoy :-).
Autonomics offers import + preprocess + analyze + visualize one-liners for RNAseq, MS proteomics, SOMAscan and Metabolon platforms, as well a generic importer for data from any rectangular omics file. We will discuss each of these in more detail in the following sections, but all of them follow the following general steps:
~0+subgroupAutonomics provides ready-to-use importers for both count as well as BAM files, which read / preprocess / analyze RNAseq data. Specific to RNAseq is counting reads using Rsubread::featureCounts (for read_rna_bams)), normalizing for library composition (cpm) or gene length (tpm), estimating voom precision weights, and using these in linear modeling. Let’s illustrate both of these on example datasets:
Billing et al. (2019) studied the differentiation of embryonic (E00) into mesenchymal stem cells (E00.8, E01, E02, E05, E15, E30), with similar properties as bone-marrow mesenchymal stem cells (M00). Importing the RNAseq counts:
require(autonomics, quietly = TRUE)
## 
## Attaching package: 'autonomics'
## The following object is masked from 'package:stats':
## 
##     biplot
object <-  read_rnaseq_counts(file = download_data('billing19.rnacounts.txt'), 
                              pca = TRUE, fit = 'limma', plot = TRUE)
##  Read ~/.cache/autonomics/aecd361c012f2_billing19.rnacounts.txt
##      Infer subgroup from sample_ids
##      Infer subgroup from sample_ids
##      Preprocess
##          Keep 22855/56269 features: count >= 10 (~0 + subgroup)
##          pseudocount 0.5
##          cpm:    tmm scale libsizes
##              cpm
##          voom: ~0 + subgroup
##          log2
## Warning: `guides(<scale> = FALSE)` is deprecated. Please use `guides(<scale> =
## "none")` instead.
##      Add PCA
##      lmFit(~0 + subgroup, weights = assays(object)$weights)
##              contrast ndown  nup
##          1: E00.8-E00     8   36
##          2: E01-E00.8     0    0
##          3:   E02-E01     0   31
##          4:   E05-E02   124   18
##          5:   E15-E05  3179 4165
##          6:   E30-E15  2104 1198
##          7:   M00-E30  1579 1523The sfile/sfileby arguments in both importers allow to merge in sample data from a separate file by a specified column. This is especially useful in clinical datasets, such as the GSE161731 Covid-19 dataset, for which both a counts file as well as a sample file can be downloaded from GEO:
basedir <- file.path(tempdir(), 'datasets')
dir.create(basedir, showWarnings = FALSE, recursive = TRUE)
if (!dir.exists(file.path(basedir, 'GSE161731'))){
    GEOquery::getGEOSuppFiles("GSE161731", baseDir=basedir)
}
## Setting options('download.file.method.GEOquery'='auto')
## Setting options('GEOquery.inmemory.gpl'=FALSE)
##                                                                     size isdir
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_counts.csv.gz       8347405 FALSE
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_counts_key.csv.gz      2443 FALSE
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_key.csv.gz             2398 FALSE
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_xpr_nlcpm.csv.gz   16839511 FALSE
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_xpr_tpm_geo.txt.gz 23859414 FALSE
##                                                                 mode
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_counts.csv.gz       644
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_counts_key.csv.gz   644
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_key.csv.gz          644
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_xpr_nlcpm.csv.gz    644
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_xpr_tpm_geo.txt.gz  644
##                                                                               mtime
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_counts.csv.gz      2022-11-01 16:51:36
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_counts_key.csv.gz  2022-11-01 16:51:36
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_key.csv.gz         2022-11-01 16:51:37
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_xpr_nlcpm.csv.gz   2022-11-01 16:51:41
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_xpr_tpm_geo.txt.gz 2022-11-01 16:51:48
##                                                                               ctime
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_counts.csv.gz      2022-11-01 16:51:36
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_counts_key.csv.gz  2022-11-01 16:51:36
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_key.csv.gz         2022-11-01 16:51:37
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_xpr_nlcpm.csv.gz   2022-11-01 16:51:41
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_xpr_tpm_geo.txt.gz 2022-11-01 16:51:48
##                                                                               atime
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_counts.csv.gz      2022-11-01 16:51:33
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_counts_key.csv.gz  2022-11-01 16:51:36
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_key.csv.gz         2022-11-01 16:51:36
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_xpr_nlcpm.csv.gz   2022-11-01 16:51:37
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_xpr_tpm_geo.txt.gz 2022-11-01 16:51:41
##                                                                  uid  gid
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_counts.csv.gz      1002 1002
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_counts_key.csv.gz  1002 1002
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_key.csv.gz         1002 1002
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_xpr_nlcpm.csv.gz   1002 1002
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_xpr_tpm_geo.txt.gz 1002 1002
##                                                                     uname
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_counts.csv.gz      biocbuild
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_counts_key.csv.gz  biocbuild
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_key.csv.gz         biocbuild
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_xpr_nlcpm.csv.gz   biocbuild
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_xpr_tpm_geo.txt.gz biocbuild
##                                                                    grname
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_counts.csv.gz      biocbuild
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_counts_key.csv.gz  biocbuild
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_key.csv.gz         biocbuild
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_xpr_nlcpm.csv.gz   biocbuild
## /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_xpr_tpm_geo.txt.gz biocbuild
object <- read_rnaseq_counts(
  file  = file.path(basedir, 'GSE161731/GSE161731_counts.csv.gz'), 
  sfile = file.path(basedir, 'GSE161731/GSE161731_counts_key.csv.gz'), 
  sfileby = 'rna_id', 
  subgroupvar = 'gender', 
  pca = TRUE, 
  fit = 'limma', 
  plot = TRUE)
##  Read /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_counts.csv.gz
##      Merge sdata: /tmp/RtmplnPYan/datasets/GSE161731/GSE161731_counts_key.csv.gz
##      Preprocess
##          Keep 18974/60675 features: count >= 10 (~0 + gender)
##          pseudocount 0.5
##          cpm:    tmm scale libsizes
##              cpm
##          voom: ~0 + gender
##          log2
## Warning: `guides(<scale> = FALSE)` is deprecated. Please use `guides(<scale> =
## "none")` instead.
##      Add PCA
##      lmFit(~0 + gender, weights = assays(object)$weights)
##                contrast ndown nup
##          1: Male-Female    62  61The block argument in both importers allows to model e.g. subject_id as a random effect, which is very common in clinical data. Repeating the GEO Covid-19 analysis with subject_id as a block effect. Outcomment the line below, which takes a few minutes to complete.
#object <- read_rnaseq_counts(
#  file  = '~/autonomicscache/datasets/GSE161731/GSE161731_counts.csv.gz', 
#  sfile = '~/autonomicscache/datasets/GSE161731/GSE161731_counts_key.csv.gz', 
#  sfileby = 'rna_id', 
#  subgroupvar = 'gender', 
#  block = 'subject_id')Billing et al. (2016) compared three types of stem cells: embryonic (E), embryonic differentiated into mesenchymal (EM), and bone-marrow mesenchymal (BM). Importing a downsized version of the RNAseq BAM files (with only a subset of all reads):
# not run to avoid issues with R CMD CHECK 
if (requireNamespace('Rsubread')){
  object <- read_rnaseq_bams(
                dir    = download_data('billing16.bam.zip'), 
                paired = TRUE, 
                genome = 'hg38', 
                pca    = TRUE, 
                fit   = 'limma', 
                plot   = TRUE)
}Proper preprocessing leads to increased power:
  file <- download_data('billing19.rnacounts.txt')
# log2counts
  object <- read_rnaseq_counts(file, 
       cpm = FALSE, voom = FALSE, fit = 'limma', verbose = FALSE, plot = FALSE)
##      Infer subgroup from sample_ids
  colSums(summarize_fit(object, 'limma')[, -1])
## ndown   nup 
##   456    61
# log2cpm
  object <- read_rnaseq_counts(file,
       cpm = TRUE,  voom = FALSE, fit = 'limma', verbose = FALSE, plot=FALSE)
##      Infer subgroup from sample_ids
  colSums(summarize_fit(object, 'limma')[, -1])
## ndown   nup 
##  5458  5107
# log2cpm + voom
  object <- read_rnaseq_counts(file,  # log2 cpm + voom
       cpm = TRUE,  voom = TRUE,  fit = 'limma', verbose = FALSE, plot=FALSE)
##      Infer subgroup from sample_ids
  colSums(summarize_fit(object, 'limma')[, -1])
## ndown   nup 
##  6994  6971A popular approach in (DDA) MS proteomics data analysis is to use MaxQuant to produce proteinGroups.txt and phospho (STY)Sites.txt files. These can be imported using read_proteingroups / read_phosphosites, which :
An LFQ intensities example is the Fukuda et al. (2020) dataset, which compares the proteome of 30d zebrafish juveniles versus adults using label-free quantification. In the volcano plot measured values are shown with circles, imputed values as triangles.
object <- read_proteingroups(file = download_data('fukuda20.proteingroups.txt'), 
                             pca = TRUE, fit = 'limma', plot = TRUE)
##  Read ~/.cache/autonomics/1c5c8970b3f7cd_fukuda20.proteingroups.txt
##      Infer subgroup from sample_ids
##      Replace 0->NA for 15654/41112 values (in 3998/6852 features of 6/6 samples)
##      Log2 transform
##  Filter features
##      Retain 6735/6852 features: ~Reverse != "+"
##      Retain 6692/6735 features: contaminant != '+'
##      Retain 6679/6692 features: non-zero, non-NA, and non-NaN for some sample
##      Filter 4534/6679 features: expr > 0 for at least two samples in some subgroup
##  Transform exprs
##      Impute systematic nondetects for 762/4534 features in 6/6 samples
## Warning: `guides(<scale> = FALSE)` is deprecated. Please use `guides(<scale> =
## "none")` instead.
##      Add PCA
##      lmFit(~0 + subgroup)
##                 contrast ndown nup
##          1: Adult-X30dpt  1620 876Normalized ratios were used in the proteomics profiling of the earlier described Billing et al. (2019) study, which examined the differentiation of embryonic stem cells (E00) into mesenchymal stem cells (E01,E02, E05, E15, E30), and compared these to bone-marrow derived mesenchymal stem cells (M00). The proteomics profiling experiment was performed using MS1 labeling: a light (L) labeled internal standard was created by mixing all samples in equal amounts, the subsequent samples were either medium (M) or heavy (H) labeled. Not all label combinations were of interest, and the select_subgroups argument allows to specifies which subgroup combinations to retain:
object <- read_proteingroups(
              file = download_data('billing19.proteingroups.txt'),
              select_subgroups = c('E00_STD', 'E01_STD', 'E02_STD',
                                   'E05_STD', 'E15_STD', 'E30_STD', 'M00_STD'), 
              pca = TRUE, fit = 'limma', plot = TRUE)
##  Read ~/.cache/autonomics/1c5c891d6b1d02_billing19.proteingroups.txt
##      Infer subgroup from sample_ids
##          Retain 21/33 samples: ~subgroup %in% select_subgroups
##      Replace NaN->NA for 77204/205338 values (in 5423/9778 features of 21/21 samples)
##      Log2 transform
##  Filter features
##      Retain 9566/9778 features: ~Reverse != "+"
##      Retain 9470/9566 features: contaminant != '+'
##      Retain 7864/9470 features: non-zero, non-NA, and non-NaN for some sample
##      Filter 7094/7864 features: expr != 0 for at least two samples in some subgroup
##  Transform exprs
## Warning: `guides(<scale> = FALSE)` is deprecated. Please use `guides(<scale> =
## "none")` instead.
##      Add PCA
##      lmFit(~0 + subgroup)
##                    contrast ndown  nup
##          1: E01_STD-E00_STD     1    2
##          2: E02_STD-E01_STD     1    2
##          3: E05_STD-E02_STD     4    6
##          4: E15_STD-E05_STD  1350 1454
##          5: E30_STD-E15_STD   821  405
##          6: M00_STD-E30_STD  1705  579The phosphosites can be read in a similar fashion:
object <- read_phosphosites(
              file = download_data('billing19.phosphosites.txt'), 
              proteinfile = download_data('billing19.proteingroups.txt'), 
              select_subgroups = c('E00_STD', 'E01_STD', 'E02_STD',
                                   'E05_STD', 'E15_STD', 'E30_STD', 'M00_STD'), 
              pca = TRUE, fit = 'limma', plot = TRUE)
##  Read ~/.cache/autonomics/1c5c891d6b1d02_billing19.proteingroups.txt
##      Infer subgroup from sample_ids
##          Retain 21/33 samples: ~subgroup %in% select_subgroups
##      Replace NaN->NA for 77204/205338 values (in 5423/9778 features of 21/21 samples)
##      Log2 transform
##  Read ~/.cache/autonomics/1c5c8911a67ae8_billing19.phosphosites.txt
##      Infer subgroup from sample_ids
##          Retain 21/33 samples: ~subgroup %in% select_subgroups
##      Replace NaN->NA for 11144/296226 values (in 1290/14106 features of 21/21 samples)
##      Log2 transform
##  Filter features
##      Retain 13968/14106 features: ~Reverse != "+"
##      Retain 13826/13968 features: contaminant != '+'
##      Retain 9691/13826 features: non-zero, non-NA, and non-NaN for some sample
##      Filter 7072/9691 features: expr != 0 for at least two samples in some subgroup
##      Retain 6034/7072 features: ~`Localization prob` >= min_localization_prob
##  Add occupancies(phospho) = values(phospho) - values(proteins)
##  Transform exprs
## Warning: `guides(<scale> = FALSE)` is deprecated. Please use `guides(<scale> =
## "none")` instead.
##      Add PCA
##      lmFit(~0 + subgroup)
##                    contrast ndown  nup
##          1: E01_STD-E00_STD     0    4
##          2: E02_STD-E01_STD     1    3
##          3: E05_STD-E02_STD    31   20
##          4: E15_STD-E05_STD  1000 1085
##          5: E30_STD-E15_STD   649  477
##          6: M00_STD-E30_STD  1338  224Metabolon delivers metabolomics results in the form of an xlsx file, with values of interest likely in the (second) OrigScale sheet. Features are laid out in rows, samples in columns, and additional feature/sample data are self-contained in the file. Metabolon data can be read/analyzed with the autonomics function read_metabolon, as illustrated below on the dataset of Halama and and Atkin (2018), who investigates the effects of hypoglycemia on a number of subjects (SUB) at four different time points (t0, t1, t2, t3):
object <- read_metabolon(
        file  = download_data('atkin18.metabolon.xlsx'),  subgroupvar = 'SET',
        pca = TRUE, fit = 'limma', block = 'SUB', plot = TRUE)
##  Read: ~/.cache/autonomics/1c5c89372776af_atkin18.metabolon.xlsx
##      Log2 transform
## Warning: `guides(<scale> = FALSE)` is deprecated. Please use `guides(<scale> =
## "none")` instead.
##      Add PCA
##      lmFit(~0 + SET | SUB)
##              dupcor `SUB`
##             contrast ndown nup
##          1:    t1-t0   326  20
##          2:    t2-t1    39  28
##          3:    t3-t2    58 359Somascan results are returned in a txt file with extension .adat. Features are laid out in columns and samples in rows. Rectangles with feature/sample data are also contained in the file. The autonomics function read_somascan reads/analyzes SOMAscan files, additionally filtering samples/features for quality and type. This is illustrated on the above described dataset from Halama and and Atkin (2018), who investigated the effects of hypoglycemia on a number of subjects (SUB) at four different time points (t0, t1, t2, t3):
object <- read_somascan(
            file = download_data('atkin18.somascan.adat'), 
            pca = TRUE, fit = 'limma', block = 'Subject_ID', plot = TRUE)
##  Read: ~/.cache/autonomics/1c5c89650c86ed_atkin18.somascan.adat
##      Retain 1125/1129 features: ~ColCheck %in% ~feature_quality
##      Log2 transform
## Warning: `guides(<scale> = FALSE)` is deprecated. Please use `guides(<scale> =
## "none")` instead.
##      Add PCA
##      lmFit(~0 + SampleGroup | Subject_ID)
##              dupcor `Subject_ID`
##             contrast ndown nup
##          1:    t1-t0     1   2
##          2:    t2-t1     2   2
##          3:    t3-t2    65  47sessionInfo()
## R version 4.2.1 (2022-06-23)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.5 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.16-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.16-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] stats     graphics  grDevices utils     datasets  methods   base     
## 
## other attached packages:
## [1] autonomics_1.6.0 BiocStyle_2.26.0
## 
## loaded via a namespace (and not attached):
##   [1] assertive.base_0.0-9        colorspace_2.0-3           
##   [3] ellipsis_0.3.2              XVector_0.38.0             
##   [5] GenomicRanges_1.50.0        farver_2.1.1               
##   [7] assertive.sets_0.0-3        MultiAssayExperiment_1.24.0
##   [9] ggrepel_0.9.1               bit64_4.0.5                
##  [11] fansi_1.0.3                 xml2_1.3.3                 
##  [13] assertive.data.uk_0.0-2     codetools_0.2-18           
##  [15] R.methodsS3_1.8.2           cachem_1.0.6               
##  [17] knitr_1.40                  jsonlite_1.8.3             
##  [19] assertive_0.3-6             assertive.data.us_0.0-2    
##  [21] rematch_1.0.1               dbplyr_2.2.1               
##  [23] R.oo_1.25.0                 BiocManager_1.30.19        
##  [25] readr_2.1.3                 compiler_4.2.1             
##  [27] httr_1.4.4                  assertthat_0.2.1           
##  [29] Matrix_1.5-1                fastmap_1.1.0              
##  [31] limma_3.54.0                cli_3.4.1                  
##  [33] htmltools_0.5.3             tools_4.2.1                
##  [35] gtable_0.3.1                glue_1.6.2                 
##  [37] GenomeInfoDbData_1.2.9      dplyr_1.0.10               
##  [39] rappdirs_0.3.3              Rcpp_1.0.9                 
##  [41] Biobase_2.58.0              cellranger_1.1.0           
##  [43] jquerylib_0.1.4             vctrs_0.5.0                
##  [45] assertive.files_0.0-2       assertive.datetimes_0.0-3  
##  [47] assertive.models_0.0-2      xfun_0.34                  
##  [49] stringr_1.4.1               lifecycle_1.0.3            
##  [51] statmod_1.4.37              edgeR_3.40.0               
##  [53] zlibbioc_1.44.0             scales_1.2.1               
##  [55] pcaMethods_1.90.0           hms_1.1.2                  
##  [57] MatrixGenerics_1.10.0       parallel_4.2.1             
##  [59] SummarizedExperiment_1.28.0 GEOquery_2.66.0            
##  [61] assertive.matrices_0.0-2    assertive.strings_0.0-3    
##  [63] yaml_2.3.6                  curl_4.3.3                 
##  [65] memoise_2.0.1               gridExtra_2.3              
##  [67] ggplot2_3.3.6               sass_0.4.2                 
##  [69] stringi_1.7.8               RSQLite_2.2.18             
##  [71] highr_0.9                   S4Vectors_0.36.0           
##  [73] BiocGenerics_0.44.0         filelock_1.0.2             
##  [75] GenomeInfoDb_1.34.0         rlang_1.0.6                
##  [77] pkgconfig_2.0.3             matrixStats_0.62.0         
##  [79] bitops_1.0-7                evaluate_0.17              
##  [81] lattice_0.20-45             assertive.data_0.0-3       
##  [83] purrr_0.3.5                 labeling_0.4.2             
##  [85] assertive.properties_0.0-5  bit_4.0.4                  
##  [87] tidyselect_1.2.0            assertive.code_0.0-3       
##  [89] magrittr_2.0.3              bookdown_0.29              
##  [91] R6_2.5.1                    IRanges_2.32.0             
##  [93] magick_2.7.3                generics_0.1.3             
##  [95] DelayedArray_0.24.0         DBI_1.1.3                  
##  [97] pillar_1.8.1                withr_2.5.0                
##  [99] assertive.numbers_0.0-2     abind_1.4-5                
## [101] RCurl_1.98-1.9              tibble_3.1.8               
## [103] assertive.types_0.0-3       utf8_1.2.2                 
## [105] BiocFileCache_2.6.0         tzdb_0.3.0                 
## [107] rmarkdown_2.17              locfit_1.5-9.6             
## [109] grid_4.2.1                  readxl_1.4.1               
## [111] data.table_1.14.4           blob_1.2.3                 
## [113] digest_0.6.30               tidyr_1.2.1                
## [115] R.utils_2.12.1              stats4_4.2.1               
## [117] munsell_0.5.0               assertive.reflection_0.0-5 
## [119] bslib_0.4.0A M Billing, S S Dib, A M Bhagwat, I T da Silva, R D Drummond, S Hayat, R Al-Mismar, H Ben-Hamidane, N Goswami, K Engholm-Keller, M R Larsen, K Suhre, A Rafii, J Graummann (2019). Mol Cell Proteomics. 18(10):1950-1966. doi: 10.1074/mcp.RA119.001356.
A M Billing, H Ben Hamidane, S S Dib, R J Cotton, A M Bhagwat, P Kumar, S Hayat, N A Yousri, N Goswami, K Suhre, A Rafii, J Graumann (2016). Comprehensive transcriptomics and proteomics characterization of human mesenchymal stem cells reveals source specific cellular markers. Sci Rep. 9;6:21507. doi: 10.1038/srep21507.
R Fukuda, R Marin-Juez, H El-Sammak, A Beisaw, R Ramadass, C Kuenne, S Guenther, A Konzer, A M Bhagwat, J Graumann, D YR Stainier (2020). EMBO Rep. 21(8): e49752. doi: 10.15252/embr.201949752
A Halama, M Kulinski, S Dib, S B Zaghlool, K S Siveen, A Iskandarani, J Zierer 3, K S Prabhu, N J Satheesh, A M Bhagwat, S Uddin, G Kastenmüller, O Elemento, S S Gross, K Suhre (2018). Accelerated lipid catabolism and autophagy are cancer survival mechanisms under inhibited glutaminolysis. Cancer Lett., 430:133-147. doi:10.1016/j.canlet.2018.05.017
A Halama, H Kahal, A M Bhagwat, J Zierer, T Sathyapalan, J Graumann, K Suhre, S L Atkin (2018). Metabolic and proteomics signatures of hypoglycaemia in type 2 diabetes. Diabetes, obesity and metabolism, https://doi.org/10.1111/dom.13602