QFeatures 1.0.0
QFeaturesWe are going to use a subset of the CPTAC study 6 containing
conditions A and B (Paulovich et al. 2010). The peptide-level data, as
processed by MaxQuant (Cox and Mann 2008) is available in the msdata
package:
basename(f <- msdata::quant(pattern = "cptac", full.names = TRUE))## [1] "cptac_a_b_peptides.txt"From the names of the columns, we see that the quantitative columns,
starting with "Intensity." (note the space!) are at positions 56 to
61.
names(read.delim(f))##  [1] "Sequence"                 "N.term.cleavage.window"  
##  [3] "C.term.cleavage.window"   "Amino.acid.before"       
##  [5] "First.amino.acid"         "Second.amino.acid"       
##  [7] "Second.last.amino.acid"   "Last.amino.acid"         
##  [9] "Amino.acid.after"         "A.Count"                 
## [11] "R.Count"                  "N.Count"                 
## [13] "D.Count"                  "C.Count"                 
## [15] "Q.Count"                  "E.Count"                 
## [17] "G.Count"                  "H.Count"                 
## [19] "I.Count"                  "L.Count"                 
## [21] "K.Count"                  "M.Count"                 
## [23] "F.Count"                  "P.Count"                 
## [25] "S.Count"                  "T.Count"                 
## [27] "W.Count"                  "Y.Count"                 
## [29] "V.Count"                  "U.Count"                 
## [31] "Length"                   "Missed.cleavages"        
## [33] "Mass"                     "Proteins"                
## [35] "Leading.razor.protein"    "Start.position"          
## [37] "End.position"             "Unique..Groups."         
## [39] "Unique..Proteins."        "Charges"                 
## [41] "PEP"                      "Score"                   
## [43] "Identification.type.6A_7" "Identification.type.6A_8"
## [45] "Identification.type.6A_9" "Identification.type.6B_7"
## [47] "Identification.type.6B_8" "Identification.type.6B_9"
## [49] "Experiment.6A_7"          "Experiment.6A_8"         
## [51] "Experiment.6A_9"          "Experiment.6B_7"         
## [53] "Experiment.6B_8"          "Experiment.6B_9"         
## [55] "Intensity"                "Intensity.6A_7"          
## [57] "Intensity.6A_8"           "Intensity.6A_9"          
## [59] "Intensity.6B_7"           "Intensity.6B_8"          
## [61] "Intensity.6B_9"           "Reverse"                 
## [63] "Potential.contaminant"    "id"                      
## [65] "Protein.group.IDs"        "Mod..peptide.IDs"        
## [67] "Evidence.IDs"             "MS.MS.IDs"               
## [69] "Best.MS.MS"               "Oxidation..M..site.IDs"  
## [71] "MS.MS.Count"(i <- grep("Intensity\\.", names(read.delim(f))))## [1] 56 57 58 59 60 61We now read these data using the readQFeatures function. The peptide
level expression data will be imported into R as an instance of class
QFeatures named cptac with an assay named peptides. We also use
the fnames argument to set the row-names of the peptides assay to
the peptide sequences.
library("QFeatures")
cptac <- readQFeatures(f, ecol = i, sep = "\t", name = "peptides", fnames = "Sequence")Below we update the sample (column) annotations to encode the two groups, 6A and 6B, and the original sample numbers.
cptac$group <- rep(c("6A", "6B"), each = 3)
cptac$sample <- rep(7:9, 2)
colData(cptac)## DataFrame with 6 rows and 2 columns
##                      group    sample
##                <character> <integer>
## Intensity.6A_7          6A         7
## Intensity.6A_8          6A         8
## Intensity.6A_9          6A         9
## Intensity.6B_7          6B         7
## Intensity.6B_8          6B         8
## Intensity.6B_9          6B         9filterFeatures(cptac, ~ Reverse == "")## An instance of class QFeatures containing 1 assays:
##  [1] peptides: SummarizedExperiment with 11436 rows and 6 columnsfilterFeatures(cptac, ~ Potential.contaminant == "")## An instance of class QFeatures containing 1 assays:
##  [1] peptides: SummarizedExperiment with 11385 rows and 6 columnslibrary("magrittr")
cptac <- cptac %>%
    filterFeatures(~ Reverse == "") %>%
    filterFeatures(~ Potential.contaminant == "")The spreadsheet that was read above contained numerous variables that are returned by MaxQuant, but not necessarily necessary in the frame of a downstream statistical analysis.
rowDataNames(cptac)## CharacterList of length 1
## [["peptides"]] Sequence N.term.cleavage.window ... MS.MS.CountThe only ones that we will be needing below are the peptides sequences
and the protein identifiers. Below, we store these variables of
interest and filter them using the selectRowData function.
rowvars <- c("Sequence", "Proteins", "Leading.razor.protein")
cptac <- selectRowData(cptac, rowvars)
rowDataNames(cptac)## CharacterList of length 1
## [["peptides"]] Sequence Proteins Leading.razor.proteinMissing values can be very numerous in certain proteomics experiments
and need to be dealt with carefully. The first step is to assess their
presence across samples and features. But before being able to do so,
we need to replace 0 by NA, given that MaxQuant encodes missing data
with a 0 using the zeroIsNA function.
cptac <- zeroIsNA(cptac, i = seq_along(cptac))
nNA(cptac, i = seq_along(cptac))## $nNA
## DataFrame with 1 row and 3 columns
##         assay       nNA       pNA
##   <character> <integer> <numeric>
## 1    peptides     30609   44.9194
## 
## $nNArows
## DataFrame with 11357 rows and 4 columns
##             assay          name       nNA       pNA
##       <character>   <character> <integer> <numeric>
## 1        peptides AAAAGAGGAG...         4   66.6667
## 2        peptides     AAAALAGGK         0    0.0000
## 3        peptides    AAAALAGGKK         0    0.0000
## 4        peptides AAADALSDLE...         0    0.0000
## 5        peptides AAADALSDLE...         0    0.0000
## ...           ...           ...       ...       ...
## 11353    peptides YYSIYDLGNN...         6  100.0000
## 11354    peptides YYTFNGPNYN...         3   50.0000
## 11355    peptides    YYTITEVATR         4   66.6667
## 11356    peptides YYTVFDRDNN...         6  100.0000
## 11357    peptides YYTVFDRDNN...         6  100.0000
## 
## $nNAcols
## DataFrame with 6 rows and 4 columns
##         assay          name       nNA       pNA
##   <character>   <character> <integer> <numeric>
## 1    peptides Intensity....      4669   41.1112
## 2    peptides Intensity....      5388   47.4421
## 3    peptides Intensity....      5224   45.9981
## 4    peptides Intensity....      4651   40.9527
## 5    peptides Intensity....      5470   48.1641
## 6    peptides Intensity....      5207   45.8484The output of the nNA function tells us that
In this dataset, we have such a high number of peptides without any
data because the 6 samples are a subset of a larger dataset, and these
peptides happened to be absent in groups A and B. Below, we use
filterNA to remove all the peptides that contain one or more missing
values by using pNA = 0 (which also is the default value).
cptac <- filterNA(cptac, i = seq_along(cptac), pNA = 0)
cptac## An instance of class QFeatures containing 1 assays:
##  [1] peptides: SummarizedExperiment with 4051 rows and 6 columnsI we wanted to keep peptides that have up to 90% of missing values,
corresponsing in this case to those that have only one value (i.e 5/6
percent of missing values), we could have set pNA to 0.9.
The impute method can be used to perform missing value imputation
using a variety of imputation methods. The method takes an instance of
class QFeatures (or a SummarizedExperiment) as input, an a
character naming the desired method (see ?impute for the complete
list with details) and returns a new instance of class QFeatures (or
SummarizedExperiment) with imputed data.
As described in more details in (Lazar et al. 2016), there are two types of mechanisms resulting in missing values in LC/MSMS experiments.
Missing values resulting from absence of detection of a feature, despite ions being present at detectable concentrations. For example in the case of ion suppression or as a result from the stochastic, data-dependent nature of the MS acquisition method. These missing value are expected to be randomly distributed in the data and are defined as missing at random (MAR) or missing completely at random (MCAR).
Biologically relevant missing values, resulting from the absence of the low abundance of ions (below the limit of detection of the instrument). These missing values are not expected to be randomly distributed in the data and are defined as missing not at random (MNAR).
MAR and MCAR values can be reasonably well tackled by many imputation methods. MNAR data, however, requires some knowledge about the underlying mechanism that generates the missing data, to be able to attempt data imputation. MNAR features should ideally be imputed with a left-censor (for example using a deterministic or probabilistic minimum value) method. Conversely, it is recommended to use hot deck methods (for example nearest neighbour, maximum likelihood, etc) when data are missing at random.
Figure 1: Mixed imputation method
Black cells represent presence of quantitation values and light grey corresponds to missing data. The two groups of interest are depicted in green and blue along the heatmap columns. Two classes of proteins are annotated on the left: yellow are proteins with randomly occurring missing values (if any) while proteins in brown are candidates for non-random missing value imputation.
It is anticipated that the identification of both classes of missing values will depend on various factors, such as feature intensities and experimental design. Below, we use perform mixed imputation, applying nearest neighbour imputation on the 654 features that are assumed to contain randomly distributed missing values (if any) (yellow on figure 1) and a deterministic minimum value imputation on the 35 proteins that display a non-random pattern of missing values (brown on figure 1).
When analysing continuous data using parametric methods (such as t-test or linear models), it is often necessary to log-transform the data. The figure below (left) show that how our data is mainly composed of small values with a long tail of larger ones, which is a typical pattern of quantitative omics data.
Below, we use the logTransform function to log2-transform our
data. This time, instead of overwriting the peptides assay, we are
going to create a new one to contain the log2-transformed data.
cptac <- addAssay(cptac,
                  logTransform(cptac[[1]]),
                  name = "peptides_log")
cptac## An instance of class QFeatures containing 2 assays:
##  [1] peptides: SummarizedExperiment with 4051 rows and 6 columns 
##  [2] peptides_log: SummarizedExperiment with 4051 rows and 6 columnspar(mfrow = c(1, 2))
limma::plotDensities(assay(cptac[[1]]))
limma::plotDensities(assay(cptac[[2]]))
Figure 2: Quantitative data in its original scale (left) and log2-transformed (right)
Assays in QFeatures objects can be normalised with the normalize
function. The type of normalisation is defined by the method
argument; below, we use quantiles normalisation, store the normalised
data into a new experiment, and visualise the resulting data.
cptac <- addAssay(cptac,
                  normalize(cptac[["peptides_log"]], method = "quantiles"),
                  name = "peptides_norm")## Loading required namespace: preprocessCorecptac## An instance of class QFeatures containing 3 assays:
##  [1] peptides: SummarizedExperiment with 4051 rows and 6 columns 
##  [2] peptides_log: SummarizedExperiment with 4051 rows and 6 columns 
##  [3] peptides_norm: SummarizedExperiment with 4051 rows and 6 columnspar(mfrow = c(1, 2))
limma::plotDensities(assay(cptac[["peptides_log"]]))
limma::plotDensities(assay(cptac[["peptides_norm"]]))
Figure 3: Distribution of log2 peptide intensities before (left) and after (right) quantiles normalisation
At this stage, it is possible to directly use the peptide-level intensities to perform a statistical analysis (Goeminne, Gevaert, and Clement 2016), or aggregate the peptide-level data into protein intensities, and perform the differential expression analysis at the protein level.
To aggregate feature data, we can use the aggregateFeatures function
that takes the following inputs:
QFeatures instance that contains the peptide
quantitation data - "cptac" in our example;i: the name or index of the assay that contains the
(normalised) peptide quantitation data - "peptides_norm" in our
case;fcol: the feature variable (in the assay above) to be used to
define what peptides to aggregate - "Proteins" here, given that we
want to aggregate all peptides that belong to one protein (group);name: the name of the new aggregates assay - "proteins" in this case;fun, the function that will compute this
aggregation - we will be using the default value, namely
robustSummary (Sticker et al. 2019).cptac <- aggregateFeatures(cptac, i = "peptides_norm", fcol = "Proteins", name = "proteins")
cptac## An instance of class QFeatures containing 4 assays:
##  [1] peptides: SummarizedExperiment with 4051 rows and 6 columns 
##  [2] peptides_log: SummarizedExperiment with 4051 rows and 6 columns 
##  [3] peptides_norm: SummarizedExperiment with 4051 rows and 6 columns 
##  [4] proteins: SummarizedExperiment with 1125 rows and 6 columnsWe obtain a final 1125 quantified proteins in the new proteins
assay. Below, we display the quantitation data for the first 6
proteins and their respective variables. The latter shown that number
of peptides that were using during the aggregation step (.n column).
head(assay(cptac[["proteins"]]))##                                      Intensity.6A_7 Intensity.6A_8
## P00918ups|CAH2_HUMAN_UPS                   17.30380       16.89828
## P01008ups|ANT3_HUMAN_UPS;CON__P41361       16.89434       15.98417
## P01127ups|PDGFB_HUMAN_UPS                  16.52022       16.81015
## P02144ups|MYG_HUMAN_UPS                    16.88925       16.46800
## P02753ups|RETBP_HUMAN_UPS                  17.84540       16.65529
## P02787ups|TRFE_HUMAN_UPS                   16.82634       16.88629
##                                      Intensity.6A_9 Intensity.6B_7
## P00918ups|CAH2_HUMAN_UPS                   16.63440       18.28533
## P01008ups|ANT3_HUMAN_UPS;CON__P41361       16.30053       16.82910
## P01127ups|PDGFB_HUMAN_UPS                  16.83638       18.20670
## P02144ups|MYG_HUMAN_UPS                    17.27137       17.87518
## P02753ups|RETBP_HUMAN_UPS                  16.49888       18.39589
## P02787ups|TRFE_HUMAN_UPS                   16.34263       18.14808
##                                      Intensity.6B_8 Intensity.6B_9
## P00918ups|CAH2_HUMAN_UPS                   18.57310       18.47074
## P01008ups|ANT3_HUMAN_UPS;CON__P41361       16.72590       16.44177
## P01127ups|PDGFB_HUMAN_UPS                  18.79476       17.14791
## P02144ups|MYG_HUMAN_UPS                    18.58260       18.29267
## P02753ups|RETBP_HUMAN_UPS                  17.73285       18.13718
## P02787ups|TRFE_HUMAN_UPS                   18.52184       18.15342rowData(cptac[["proteins"]])## DataFrame with 1125 rows and 3 columns
##                                           Proteins Leading.razor.protein
##                                        <character>           <character>
## P00918ups|CAH2_HUMAN_UPS             P00918ups|...         P00918ups|...
## P01008ups|ANT3_HUMAN_UPS;CON__P41361 P01008ups|...         P01008ups|...
## P01127ups|PDGFB_HUMAN_UPS            P01127ups|...         P01127ups|...
## P02144ups|MYG_HUMAN_UPS              P02144ups|...         P02144ups|...
## P02753ups|RETBP_HUMAN_UPS            P02753ups|...         P02753ups|...
## ...                                            ...                   ...
## sp|Q99207|NOP14_YEAST                sp|Q99207|...         sp|Q99207|...
## sp|Q99216|PNO1_YEAST                 sp|Q99216|...         sp|Q99216|...
## sp|Q99257|MEX67_YEAST                sp|Q99257|...         sp|Q99257|...
## sp|Q99258|RIB3_YEAST                 sp|Q99258|...         sp|Q99258|...
## sp|Q99383|HRP1_YEAST                 sp|Q99383|...         sp|Q99383|...
##                                             .n
##                                      <integer>
## P00918ups|CAH2_HUMAN_UPS                     1
## P01008ups|ANT3_HUMAN_UPS;CON__P41361         1
## P01127ups|PDGFB_HUMAN_UPS                    1
## P02144ups|MYG_HUMAN_UPS                      1
## P02753ups|RETBP_HUMAN_UPS                    2
## ...                                        ...
## sp|Q99207|NOP14_YEAST                        1
## sp|Q99216|PNO1_YEAST                         1
## sp|Q99257|MEX67_YEAST                        2
## sp|Q99258|RIB3_YEAST                         2
## sp|Q99383|HRP1_YEAST                         2We can get a quick overview of this .n variable by computing the
table below, that shows us that we have 405 proteins that are based on
a single peptides, 230 that are based on two, 119 that are based on
three, … and a single protein that is the results of aggregating 44
peptides.
table(rowData(cptac[["proteins"]])$.n)## 
##   1   2   3   4   5   6   7   8   9  10  11  12  13  14  15  16  17  18  19  20 
## 405 230 119  84  64  53  37  29  24  24  13   9   4   3   3   7   3   1   1   1 
##  21  22  23  24  25  30  31  33  44 
##   1   2   2   1   1   1   1   1   1Let’s choose P02787ups|TRFE_HUMAN_UPS and visualise its expression
pattern in the 2 groups at the protein and peptide level.
library("ggplot2")
library("dplyr")
longFormat(cptac["P02787ups|TRFE_HUMAN_UPS", ]) %>%
    as.data.frame() %>%
    mutate(group = ifelse(grepl("A", colname), "A", "B")) %>%
    mutate(sample = sub("Intensity\\.", "", colname)) %>%
    ggplot(aes(x = sample, y = value, colour = rowname, shape = group)) +
    geom_point() +
    facet_grid(~ assay)## harmonizing input:
##   removing 12 sampleMap rows not in names(experiments)
## harmonizing input:
##   removing 12 sampleMap rows not in names(experiments)
Figure 4: Expression intensities for the protein P02787ups|TRFE_HUMAN_UPS (right, green) and its peptides (left) in groups A (circles) and B (triangles)
## R version 4.0.3 (2020-10-10)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 18.04.5 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.12-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.12-bioc/R/lib/libRlapack.so
## 
## 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       
## 
## attached base packages:
## [1] parallel  stats4    stats     graphics  grDevices utils     datasets 
## [8] methods   base     
## 
## other attached packages:
##  [1] gplots_3.1.0                magrittr_1.5               
##  [3] dplyr_1.0.2                 ggplot2_3.3.2              
##  [5] QFeatures_1.0.0             MultiAssayExperiment_1.16.0
##  [7] SummarizedExperiment_1.20.0 Biobase_2.50.0             
##  [9] GenomicRanges_1.42.0        GenomeInfoDb_1.26.0        
## [11] IRanges_2.24.0              S4Vectors_0.28.0           
## [13] BiocGenerics_0.36.0         MatrixGenerics_1.2.0       
## [15] matrixStats_0.57.0          BiocStyle_2.18.0           
## 
## loaded via a namespace (and not attached):
##  [1] Rcpp_1.0.5              msdata_0.29.0           lattice_0.20-41        
##  [4] gtools_3.8.2            digest_0.6.27           R6_2.4.1               
##  [7] evaluate_0.14           highr_0.8               pillar_1.4.6           
## [10] zlibbioc_1.36.0         rlang_0.4.8             lazyeval_0.2.2         
## [13] magick_2.5.0            Matrix_1.2-18           preprocessCore_1.52.0  
## [16] rmarkdown_2.5           labeling_0.4.2          stringr_1.4.0          
## [19] ProtGenerics_1.22.0     RCurl_1.98-1.2          munsell_0.5.0          
## [22] DelayedArray_0.16.0     compiler_4.0.3          xfun_0.18              
## [25] pkgconfig_2.0.3         htmltools_0.5.0         tidyselect_1.1.0       
## [28] tibble_3.0.4            GenomeInfoDbData_1.2.4  bookdown_0.21          
## [31] crayon_1.3.4            withr_2.3.0             MASS_7.3-53            
## [34] bitops_1.0-6            grid_4.0.3              gtable_0.3.0           
## [37] lifecycle_0.2.0         AnnotationFilter_1.14.0 MsCoreUtils_1.2.0      
## [40] scales_1.1.1            KernSmooth_2.23-17      stringi_1.5.3          
## [43] farver_2.0.3            XVector_0.30.0          limma_3.46.0           
## [46] ellipsis_0.3.1          generics_0.0.2          vctrs_0.3.4            
## [49] tools_4.0.3             glue_1.4.2              purrr_0.3.4            
## [52] yaml_2.2.1              colorspace_1.4-1        BiocManager_1.30.10    
## [55] caTools_1.18.0          knitr_1.30Cox, J, and M Mann. 2008. “MaxQuant Enables High Peptide Identification Rates, Individualized P.p.b.-range Mass Accuracies and Proteome-Wide Protein Quantification.” Nat Biotechnol 26 (12):1367–72. https://doi.org/10.1038/nbt.1511.
Goeminne, L J, K Gevaert, and L Clement. 2016. “Peptide-Level Robust Ridge Regression Improves Estimation, Sensitivity, and Specificity in Data-Dependent Quantitative Label-Free Shotgun Proteomics.” Mol Cell Proteomics 15 (2):657–68. https://doi.org/10.1074/mcp.M115.055897.
Lazar, C, L Gatto, M Ferro, C Bruley, and T Burger. 2016. “Accounting for the Multiple Natures of Missing Values in Label-Free Quantitative Proteomics Data Sets to Compare Imputation Strategies.” J Proteome Res 15 (4):1116–25. https://doi.org/10.1021/acs.jproteome.5b00981.
Paulovich, Amanda G, Dean Billheimer, Amy-Joan L Ham, Lorenzo Vega-Montoto, Paul A Rudnick, David L Tabb, Pei Wang, et al. 2010. “Interlaboratory Study Characterizing a Yeast Performance Standard for Benchmarking LC-MS Platform Performance.” Mol. Cell. Proteomics 9 (2):242–54.
Sticker, Adriaan, Ludger Goeminne, Lennart Martens, and Lieven Clement. 2019. “Robust Summarization and Inference in Proteome-Wide Label-Free Quantification.” bioRxiv. Cold Spring Harbor Laboratory. https://doi.org/10.1101/668863.