Contents

The following content is descibed in more detail in Egloff et al. (2018) (under review NMETH-A35040).

library(NestLink)
stopifnot(require(specL))

1 Input Pool Frequency

aa_pool_x8 <- c(rep('A', 12), rep('S', 0), rep('T', 12), rep('N', 12),
    rep('Q', 12), rep('D', 8),  rep('E', 0), rep('V', 12), rep('L', 0),
    rep('F', 0), rep('Y', 8), rep('W', 0), rep('G', 12), rep('P', 12))

aa_pool_1_2_9_10 <- c(rep('A', 8), rep('S', 7), rep('T', 7), rep('N', 6),
    rep('Q', 6), rep('D', 8), rep('E', 8), rep('V', 9), rep('L', 6),
    rep('F', 5), rep('Y', 9), rep('W', 6),  rep('G', 15), rep('P', 0))

aa_pool_3_8 <- c(rep('A', 5), rep('S', 4), rep('T', 5), rep('N', 2),
    rep('Q', 2), rep('D', 8), rep('E', 8), rep('V', 7), rep('L', 5),
    rep('F', 4), rep('Y', 6), rep('W', 4),  rep('G', 12), rep('P', 28))

2 Sanity Check

table(aa_pool_x8)
## aa_pool_x8
##  A  D  G  N  P  Q  T  V  Y 
## 12  8 12 12 12 12 12 12  8
length(aa_pool_x8)
## [1] 100
table(aa_pool_1_2_9_10)
## aa_pool_1_2_9_10
##  A  D  E  F  G  L  N  Q  S  T  V  W  Y 
##  8  8  8  5 15  6  6  6  7  7  9  6  9
length(aa_pool_1_2_9_10)
## [1] 100
table(aa_pool_3_8)
## aa_pool_3_8
##  A  D  E  F  G  L  N  P  Q  S  T  V  W  Y 
##  5  8  8  4 12  5  2 28  2  4  5  7  4  6
length(aa_pool_3_8)
## [1] 100

3 Compose Peptides

3.1 GPGXXXXXXXX(VR|VSR|VFGIR|VSGER) peptide

replicate(10, compose_GPGx8cTerm(pool=aa_pool_x8))
##  [1] "GPGVGVATTVTVSR"   "GPGAQGNYTDTVFGIR" "GPGGGPTQNNDVSR"   "GPGTYPNPVQDVFR"  
##  [5] "GPGAAYAGYGVVFGIR" "GPGTQVATYGYVFR"   "GPGTAYNDVPTVFGIR" "GPGAYTQNNQAVSR"  
##  [9] "GPGGPNQPYNYVSR"   "GPGVNPPTPYQVFR"

3.2 GPYYXXXXXXYYR peptide

compose_GPx10R(aa_pool_1_2_9_10, aa_pool_3_8)
## [1] "GPFEGPVEYPGQR"

4 Generate peptides

set.seed(2)
(sample.size <- 3E+04)
## [1] 30000
peptides.GPGx8cTerm <- replicate(sample.size, compose_GPGx8cTerm(pool=aa_pool_x8))
peptides.GPx10R <- replicate(sample.size, compose_GPx10R(aa_pool_1_2_9_10, aa_pool_3_8))
# write.table(peptides.GPGx8cTerm, file='/tmp/pp.txt')

5 Peptide mass

5.1 Compute peptide mass

library(protViz)
(smp.peptide <- compose_GPGx8cTerm(aa_pool_x8))
## [1] "GPGPDDTDTYGVFR"
parentIonMass(smp.peptide)
## [1] 1496.665
pim.GPGx8cTerm <- unlist(lapply(peptides.GPGx8cTerm, function(x){parentIonMass(x)}))
pim.GPx10R <- unlist(lapply(peptides.GPx10R, function(x){parentIonMass(x)}))
pim.iRT <-  unlist(lapply(as.character(iRTpeptides$peptide), function(x){parentIonMass(x)}))

5.2 Draw parent ion mass histogram

(pim.min <- min(pim.GPGx8cTerm, pim.GPx10R))
## [1] 1037.512
(pim.max <- max(pim.GPGx8cTerm, pim.GPx10R))
## [1] 1890.877
(pim.breaks <- seq(round(pim.min - 1) , round(pim.max + 1) , length=75))
##  [1] 1037.000 1048.554 1060.108 1071.662 1083.216 1094.770 1106.324 1117.878
##  [9] 1129.432 1140.986 1152.541 1164.095 1175.649 1187.203 1198.757 1210.311
## [17] 1221.865 1233.419 1244.973 1256.527 1268.081 1279.635 1291.189 1302.743
## [25] 1314.297 1325.851 1337.405 1348.959 1360.514 1372.068 1383.622 1395.176
## [33] 1406.730 1418.284 1429.838 1441.392 1452.946 1464.500 1476.054 1487.608
## [41] 1499.162 1510.716 1522.270 1533.824 1545.378 1556.932 1568.486 1580.041
## [49] 1591.595 1603.149 1614.703 1626.257 1637.811 1649.365 1660.919 1672.473
## [57] 1684.027 1695.581 1707.135 1718.689 1730.243 1741.797 1753.351 1764.905
## [65] 1776.459 1788.014 1799.568 1811.122 1822.676 1834.230 1845.784 1857.338
## [73] 1868.892 1880.446 1892.000
hist(pim.GPGx8cTerm, breaks=pim.breaks, probability = TRUE, 
     col='#1111AAAA', xlab='peptide mass [Dalton]', ylim=c(0, 0.006))
hist(pim.GPx10R, breaks=pim.breaks,
     probability = TRUE, add=TRUE, col='#11AA1188')
abline(v=pim.iRT, col='grey')
legend("topleft", c('GPGx8cTerm', 'GPx10R', 'iRT'), 
     fill=c('#1111AAAA', '#11AA1133', 'grey'))

6 Hydrophobicity

6.1 Compute Hydrophobicity value using SSRC

the SSRC model, see Krokhin et al. (2004), is implemented as ssrc function in protViz.

For a sanity check we apply the ssrc function to a real world LC-MS run peptideStd consits of a digest of the FETUIN_BOVINE protein (400 amol) shipped with specL Panse et al. (2015).

library(specL)
ssrc <- sapply(peptideStd, function(x){ssrc(x$peptideSequence)})
rt <- unlist(lapply(peptideStd, function(x){x$rt}))
plot(ssrc, rt); abline(ssrc.lm <- lm(rt ~ ssrc), col='red'); 
legend("topleft", paste("spearman", round(cor(ssrc, rt, method='spearman'),2)))

here we apply ssrc to the simulated flycodes and iRT peptides Escher et al. (2012).

hyd.GPGx8cTerm <- ssrc(peptides.GPGx8cTerm)
hyd.GPx10R <- ssrc(peptides.GPx10R)
hyd.iRT <- ssrc(as.character(iRTpeptides$peptide))

(hyd.min <- min(hyd.GPGx8cTerm, hyd.GPx10R))
## [1] -7.63055
(hyd.max <- max(hyd.GPGx8cTerm, hyd.GPx10R))
## [1] 65.12112
hyd.breaks <- seq(round(hyd.min - 1) , round(hyd.max + 1) , length=75)

6.2 Draw hydrophobicity histogram

hist(hyd.GPGx8cTerm, breaks = hyd.breaks, probability = TRUE, 
     col='#1111AAAA', xlab='hydrophobicity', 
     ylim=c(0, 0.06),
     main='Histogram')
hist(hyd.GPx10R, breaks = hyd.breaks, probability = TRUE, add=TRUE, col='#11AA1188')
  abline(v=hyd.iRT, col='grey')
legend("topleft", c('GPGx8cTerm', 'GPx10R', 'iRT'),  fill=c('#1111AAAA', '#11AA1133', 'grey'))

7 Quality Control (QC)

7.1 QC of composed peptides

7.1.1 Input

round(table(aa_pool_x8)/length(aa_pool_x8), 2)
## aa_pool_x8
##    A    D    G    N    P    Q    T    V    Y 
## 0.12 0.08 0.12 0.12 0.12 0.12 0.12 0.12 0.08

7.1.2 Output

peptide2aa <- function(seq, from=4, to=4+8){
  unlist(lapply(seq, function(x){strsplit(substr(x, from, to), '')}))
}
peptides.GPGx8cTerm.aa <- peptide2aa(peptides.GPGx8cTerm)
round(table(peptides.GPGx8cTerm.aa)/length(peptides.GPGx8cTerm.aa), 2)
## peptides.GPGx8cTerm.aa
##    A    D    G    N    P    Q    T    V    Y 
## 0.11 0.07 0.11 0.11 0.11 0.11 0.11 0.22 0.07
peptides.GPx10R.aa <- peptide2aa(peptides.GPx10R, from=3, to=12)
round(table(peptides.GPx10R.aa)/length(peptides.GPx10R.aa), 2)
## peptides.GPx10R.aa
##    A    D    E    F    G    L    N    P    Q    S    T    V    W    Y 
## 0.06 0.08 0.08 0.04 0.13 0.05 0.04 0.17 0.04 0.05 0.06 0.08 0.05 0.07

7.2 Count GP patterns

sample.size 
## [1] 30000
length(grep('^GP(.*)GP(.*)R$', peptides.GPGx8cTerm))
## [1] 6319
length(grep('^GP(.*)GP(.*)R$', peptides.GPx10R))
## [1] 5959

7.3 Compute AA frequency table

count the peptides having the same AA composition

sample.size 
## [1] 30000
table(table(tt<-unlist(lapply(peptides.GPGx8cTerm, 
  function(x){paste(sort(unlist(strsplit(x, ''))), collapse='')}))))
## 
##    1    2    3    4    5    6    7    8    9   10   11   12   13   14   16   17 
## 9541 3606 1607  792  427  204  104   50   34   20    6    5    6    2    1    1
# write.table(tt, file='GPGx8cTerm.txt')
table(table(unlist(lapply(peptides.GPx10R, 
  function(x){paste(sort(unlist(strsplit(x, ''))), collapse='')}))))
## 
##     1     2     3     4     5 
## 24844  2104   265    32     5

the NestLink function plot_in_silico_LCMS_map graphs the LC-MS maps.

par(mfrow=c(2, 2))
h <- NestLink:::.plot_in_silico_LCMS_map(peptides.GPGx8cTerm, main='GPGx8cTerm')
h <- NestLink:::.plot_in_silico_LCMS_map(peptides.GPx10R, main='GPx10R')

8 Session info

Here is the output of the sessionInfo() commmand.

## R Under development (unstable) (2022-10-25 r83175)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 22.04.1 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       
## 
## attached base packages:
## [1] stats4    stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
##  [1] specL_1.33.0                seqinr_4.2-16              
##  [3] RSQLite_2.2.18              DBI_1.1.3                  
##  [5] knitr_1.40                  scales_1.2.1               
##  [7] ggplot2_3.3.6               NestLink_1.15.0            
##  [9] ShortRead_1.57.0            GenomicAlignments_1.35.0   
## [11] SummarizedExperiment_1.29.0 Biobase_2.59.0             
## [13] MatrixGenerics_1.11.0       matrixStats_0.62.0         
## [15] Rsamtools_2.15.0            GenomicRanges_1.51.0       
## [17] BiocParallel_1.33.0         protViz_0.7.3              
## [19] gplots_3.1.3                Biostrings_2.67.0          
## [21] GenomeInfoDb_1.35.0         XVector_0.39.0             
## [23] IRanges_2.33.0              S4Vectors_0.37.0           
## [25] ExperimentHub_2.7.0         AnnotationHub_3.7.0        
## [27] BiocFileCache_2.7.0         dbplyr_2.2.1               
## [29] BiocGenerics_0.45.0         BiocStyle_2.27.0           
## 
## loaded via a namespace (and not attached):
##  [1] bitops_1.0-7                  deldir_1.0-6                 
##  [3] rlang_1.0.6                   magrittr_2.0.3               
##  [5] ade4_1.7-20                   compiler_4.3.0               
##  [7] mgcv_1.8-41                   png_0.1-7                    
##  [9] vctrs_0.5.0                   stringr_1.4.1                
## [11] pkgconfig_2.0.3               crayon_1.5.2                 
## [13] fastmap_1.1.0                 magick_2.7.3                 
## [15] ellipsis_0.3.2                labeling_0.4.2               
## [17] caTools_1.18.2                utf8_1.2.2                   
## [19] promises_1.2.0.1              rmarkdown_2.17               
## [21] purrr_0.3.5                   bit_4.0.4                    
## [23] xfun_0.34                     zlibbioc_1.45.0              
## [25] cachem_1.0.6                  jsonlite_1.8.3               
## [27] blob_1.2.3                    highr_0.9                    
## [29] later_1.3.0                   DelayedArray_0.25.0          
## [31] interactiveDisplayBase_1.37.0 jpeg_0.1-9                   
## [33] parallel_4.3.0                R6_2.5.1                     
## [35] bslib_0.4.0                   stringi_1.7.8                
## [37] RColorBrewer_1.1-3            jquerylib_0.1.4              
## [39] Rcpp_1.0.9                    bookdown_0.29                
## [41] assertthat_0.2.1              splines_4.3.0                
## [43] httpuv_1.6.6                  Matrix_1.5-1                 
## [45] tidyselect_1.2.0              yaml_2.3.6                   
## [47] codetools_0.2-18              hwriter_1.3.2.1              
## [49] curl_4.3.3                    lattice_0.20-45              
## [51] tibble_3.1.8                  withr_2.5.0                  
## [53] shiny_1.7.3                   KEGGREST_1.39.0              
## [55] evaluate_0.17                 pillar_1.8.1                 
## [57] BiocManager_1.30.19           filelock_1.0.2               
## [59] KernSmooth_2.23-20            generics_0.1.3               
## [61] RCurl_1.98-1.9                BiocVersion_3.17.0           
## [63] munsell_0.5.0                 gtools_3.9.3                 
## [65] xtable_1.8-4                  glue_1.6.2                   
## [67] tools_4.3.0                   interp_1.1-3                 
## [69] grid_4.3.0                    latticeExtra_0.6-30          
## [71] colorspace_2.0-3              AnnotationDbi_1.61.0         
## [73] nlme_3.1-160                  GenomeInfoDbData_1.2.9       
## [75] cli_3.4.1                     rappdirs_0.3.3               
## [77] fansi_1.0.3                   dplyr_1.0.10                 
## [79] gtable_0.3.1                  sass_0.4.2                   
## [81] digest_0.6.30                 farver_2.1.1                 
## [83] memoise_2.0.1                 htmltools_0.5.3              
## [85] lifecycle_1.0.3               httr_1.4.4                   
## [87] mime_0.12                     MASS_7.3-58.1                
## [89] bit64_4.0.5

References

Escher, C., L. Reiter, B. MacLean, R. Ossola, F. Herzog, J. Chilton, M. J. MacCoss, and O. Rinner. 2012. “Using iRT, a normalized retention time for more targeted measurement of peptides.” Proteomics 12 (8): 1111–21.

Krokhin, O. V., R. Craig, V. Spicer, W. Ens, K. G. Standing, R. C. Beavis, and J. A. Wilkins. 2004. “An improved model for prediction of retention times of tryptic peptides in ion pair reversed-phase HPLC: its application to protein peptide mapping by off-line HPLC-MALDI MS.” Mol. Cell Proteomics 3 (9): 908–19.

Panse, C., C. Trachsel, J. Grossmann, and R. Schlapbach. 2015. “specL–an R/Bioconductor package to prepare peptide spectrum matches for use in targeted proteomics.” Bioinformatics 31 (13): 2228–31.