This vignette describes how to visualize quantitative mass spectrometry data contained in a QFeatures object. This vignette is distributed under a CC BY-SA license.
QFeatures 1.14.2
To demonstrate the data visualization of QFeatures
, we first perform
a quick processing of the hlpsms
example data. We load the data and
read it as a QFeautres
object. See the processing
vignette
for more details about data processing with QFeatures
.
library("QFeatures")
data(hlpsms)
hl <- readQFeatures(hlpsms, quantCols = 1:10, name = "psms")
We then aggregate the psms to peptides, and the peptodes to proteins.
hl <- aggregateFeatures(hl, "psms", "Sequence", name = "peptides", fun = colMeans)
## Your row data contain missing values. Please read the relevant
## section(s) in the aggregateFeatures manual page regarding the effects
## of missing values on data aggregation.
hl <- aggregateFeatures(hl, "peptides", "ProteinGroupAccessions", name = "proteins", fun = colMeans)
We also add the TMT tags that were used to multiplex the samples. The
data is added to the colData
of the QFeatures
object and will
allow us to demonstrate how to plot data from the colData
.
hl$tag <- c("126", "127N", "127C", "128N", "128C", "129N", "129C",
"130N", "130C", "131")
The dataset is now ready for data exploration.
QFeatures
hierarchyQFeatures
objects can contain several assays as the data goes through
the processing workflow. The plot
function provides an overview of
all the assays present in the dataset, showing also the hierarchical
relationships between the assays as determined by the AssayLinks
.
plot(hl)
This plot is rather simple with only three assays, but some processing
workflows may involve more steps. The feat3
example data illustrates
the different possible relationships: one parent to one child, multiple
parents to one child and one parent to multiple children.
data("feat3")
plot(feat3)
Note that some datasets may contain many assays, for instance because
the MS experiment consists of hundreds of batches. This can lead to an
overcrowded plot. Therefore, you can also explore this hierarchy of
assays through an interactive plot, supported by the plotly
package
(Sievert (2020)). You can use the viewer panel to zoom in and out and
navigate across the tree(s).
plot(hl, interactive = TRUE)
The quantitative data is retrieved using assay()
, the feature
metadata is retrieved using rowData()
on the assay of interest, and
the sample metadata is retrieved using colData()
. Once retrieved,
the data can be supplied to the base R data exploration tools. Here
are some examples:
proteins
assay.plot(assay(hl, "proteins")[1, ])
.n
from the
protein rowData
.hist(rowData(hl)[["proteins"]]$.n)
tag
from the
colData
.table(hl$tag)
##
## 126 127C 127N 128C 128N 129C 129N 130C 130N 131
## 1 1 1 1 1 1 1 1 1 1
ggplot2
ggplot2
is a powerful tool for data visualization in R
and is part
of the tidyverse
package ecosystem (Wickham et al. (2019)). It produces
elegant and publication-ready plots in a few lines of code. ggplot2
can be used to explore QFeatures
object, similarly to the base
functions shown above. Note that ggplot2
expects data.frame
or
tibble
objects whereas the quantitative data in QFeatures
are
encoded as matrix
(or matrix-like objects, see
?SummarizedExperiment
) and the rowData
and colData
are encoded
as DataFrame
. This is easily circumvented by converting those
objects to data.frame
s or tibble
s. See here how we reproduce the
plot above using ggplot2
.
library("ggplot2")
df <- data.frame(rowData(hl)[["proteins"]])
ggplot(df) +
aes(x = .n) +
geom_histogram()
We refer the reader to the ggplot2
package website for more information
about the wide variety of functions that the package offers and for
tutorials and cheatsheets.
Another useful package for quantitative data exploration is
ComplexHeatmap
(Gu, Eils, and Schlesner (2016)). It is part of the Bioconductor project
(Gentleman et al. (2004)) and facilitates visualization of matrix objects as
heatmap. See here an example where we plot the protein data.
library(ComplexHeatmap)
Heatmap(matrix = assay(hl, "proteins"),
show_row_names = FALSE)
ComplexHeatmap
also allows to add row and/or column annotations.
Let’s add the predicted protein location as row annotation.
ha <- rowAnnotation(markers = rowData(hl)[["proteins"]]$markers)
Heatmap(matrix = assay(hl, "proteins"),
show_row_names = FALSE,
left_annotation = ha)
More advanced usage of ComplexHeatmap
is described in the package
reference book.
In this section, we show how to combine in a single table
different pieces of information available in a QFeatures
object,
that are quantitation data, feature metadata and sample metadata. The
QFeatures
package provides the longFormat
function that converts a
QFeatures
object into a long table. Long tables are very useful when
using ggplot2
for data visualization. For instance, suppose we want
to visualize the distribution of protein quantitation (present in the
proteins
assay) with respect to the different acquisition tags
(present in the colData
) for each predicted cell location separately
(present in the rowData
of the assays). Furthermore, we link the
quantitation values coming from the same protein using lines. This can
all be plotted at once in a few lines of code.
lf <- longFormat(hl[, , "proteins"],
rowvars = "markers",
colvars = "tag")
## Warning: 'experiments' dropped; see 'drops()'
## harmonizing input:
## removing 20 sampleMap rows not in names(experiments)
ggplot(data.frame(lf)) +
aes(x = tag,
y = value,
group = rowname) +
geom_line() +
facet_wrap(~ markers, scales = "free_y", ncol = 3)
longFormat
allows to retrieve and combine all available data from a
Qfeatures
object. We here demonstrate the ease to combine different
pieces that could highlight sample specific and/or feature specific
effects on data quantitation.
Finally, a simply shiny
app allows to explore and visualise the
respective assays of a QFeatures
object.
display(hl)
A dropdown menu in the side bar allows the user to select an assay of interest, which can then be visualised as a heatmap (figure 1), as a quantitative table (figure 2) or a row data table (figure 3).
## R version 4.4.1 (2024-06-14)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 22.04.4 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.19-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] grid stats4 stats graphics grDevices utils datasets
## [8] methods base
##
## other attached packages:
## [1] ComplexHeatmap_2.20.0 gplots_3.1.3.1
## [3] dplyr_1.1.4 ggplot2_3.5.1
## [5] QFeatures_1.14.2 MultiAssayExperiment_1.30.2
## [7] SummarizedExperiment_1.34.0 Biobase_2.64.0
## [9] GenomicRanges_1.56.1 GenomeInfoDb_1.40.1
## [11] IRanges_2.38.1 S4Vectors_0.42.1
## [13] BiocGenerics_0.50.0 MatrixGenerics_1.16.0
## [15] matrixStats_1.3.0 BiocStyle_2.32.1
##
## loaded via a namespace (and not attached):
## [1] bitops_1.0-7 rlang_1.1.4 magrittr_2.0.3
## [4] clue_0.3-65 GetoptLong_1.0.5 compiler_4.4.1
## [7] png_0.1-8 vctrs_0.6.5 reshape2_1.4.4
## [10] stringr_1.5.1 ProtGenerics_1.36.0 shape_1.4.6.1
## [13] pkgconfig_2.0.3 crayon_1.5.3 fastmap_1.2.0
## [16] magick_2.8.3 XVector_0.44.0 labeling_0.4.3
## [19] caTools_1.18.2 utf8_1.2.4 rmarkdown_2.27
## [22] UCSC.utils_1.0.0 tinytex_0.51 purrr_1.0.2
## [25] xfun_0.45 zlibbioc_1.50.0 cachem_1.1.0
## [28] jsonlite_1.8.8 highr_0.11 DelayedArray_0.30.1
## [31] parallel_4.4.1 cluster_2.1.6 R6_2.5.1
## [34] bslib_0.7.0 stringi_1.8.4 RColorBrewer_1.1-3
## [37] limma_3.60.3 jquerylib_0.1.4 Rcpp_1.0.12
## [40] bookdown_0.40 iterators_1.0.14 knitr_1.48
## [43] BiocBaseUtils_1.6.0 Matrix_1.7-0 igraph_2.0.3
## [46] tidyselect_1.2.1 abind_1.4-5 yaml_2.3.9
## [49] doParallel_1.0.17 codetools_0.2-20 lattice_0.22-6
## [52] tibble_3.2.1 plyr_1.8.9 withr_3.0.0
## [55] evaluate_0.24.0 circlize_0.4.16 pillar_1.9.0
## [58] BiocManager_1.30.23 KernSmooth_2.23-24 foreach_1.5.2
## [61] generics_0.1.3 munsell_0.5.1 scales_1.3.0
## [64] gtools_3.9.5 glue_1.7.0 lazyeval_0.2.2
## [67] tools_4.4.1 Cairo_1.6-2 tidyr_1.3.1
## [70] MsCoreUtils_1.16.0 msdata_0.44.0 colorspace_2.1-0
## [73] GenomeInfoDbData_1.2.12 cli_3.6.3 fansi_1.0.6
## [76] S4Arrays_1.4.1 AnnotationFilter_1.28.0 gtable_0.3.5
## [79] sass_0.4.9 digest_0.6.36 SparseArray_1.4.8
## [82] rjson_0.2.21 farver_2.1.2 htmltools_0.5.8.1
## [85] lifecycle_1.0.4 httr_1.4.7 GlobalOptions_0.1.2
## [88] statmod_1.5.0 MASS_7.3-61
Gentleman, Robert C., Vincent J. Carey, Douglas M. Bates, Ben Bolstad, Marcel Dettling, Sandrine Dudoit, Byron Ellis, et al. 2004. “Bioconductor: Open Software Development for Computational Biology and Bioinformatics.” Genome Biol 5 (10): –80. https://doi.org/10.1186/gb-2004-5-10-r80.
Gu, Zuguang, Roland Eils, and Matthias Schlesner. 2016. “Complex Heatmaps Reveal Patterns and Correlations in Multidimensional Genomic Data.” Bioinformatics 32 (18): 2847–9.
Sievert, Carson. 2020. Interactive Web-Based Data Visualization with R, Plotly, and Shiny. Chapman; Hall/CRC. https://plotly-r.com.
Wickham, Hadley, Mara Averick, Jennifer Bryan, Winston Chang, Lucy McGowan, Romain François, Garrett Grolemund, et al. 2019. “Welcome to the Tidyverse.” J. Open Source Softw. 4 (43): 1686.