library(Aerith)
library(ggplot2)
Stable isotope probing (SIP) experiments benefit from reliable theoretical spectra to interpret precursor and fragment ions. This vignette walks through the Aerith workflow for generating SIP-aware peptide spectra, emphasizing how each function supports isotope-resolved analysis. We start with unlabeled peptides and progress to labeled scenarios, outlining key parameters, interpretation tips, and best practices. Aerith integrates binomial-based isotope modeling, fragment intensity estimation, and visualization in a single toolkit, aligning with current SIP proteomics methodologies that rely on high accuracy precursor modeling and fragment annotation.
For unlabeled peptides, Aerith estimates precursor isotope distributions, calculates neutral masses, and visualizes theoretical spectra that reflect expected instrument observations.
The precursor_peak_calculator function returns the
isotopic distribution for a peptide sequence, while
calPepAtomCount reports elemental composition derived from
the amino acid content. The calPepPrecursorMass function
uses a binomial model to estimate the precursor monoisotopic mass given
a target isotope and its natural abundance. The sequence parameter
accepts standard one-letter amino acid codes, and the isotope argument
takes values such as "C13" or "N15". Adjust
the abundance parameter to match experimental labeling levels; use
natural abundance (for example 0.0107 for carbon) when no enrichment is
applied.
a <- precursor_peak_calculator("PEPTIDECCCC")
head(a, 5)
#> Mass Prob
#> 1 1439.483 0.40250131
#> 2 1440.485 0.27773163
#> 3 1441.484 0.18623515
#> 4 1442.485 0.08384669
#> 5 1443.484 0.03376080
calPepAtomCount("PEPTIDECCCC")
#> C H O N P S
#> 1 54 85 23 15 0 4
calPepPrecursorMass("PEPTIDECCCC", "C13", 0.0107)
#> [1] 1440.484
getPrecursorSpectra constructs theoretical mass spectra
for a peptide across specified charge states. The charges
argument accepts individual values or integer ranges, enabling quick
comparison of charge state envelopes. Interpreting the plot helps assess
isotope spacing and relative intensities, which should coincide with
observed MS1 signals when the peptide is unlabeled.
a <- getPrecursorSpectra("PEPTIDE", 2)
plot(a) + scale_x_continuous(breaks = seq(400, 405, by = 1)) + geom_linerange(linewidth = 0.2)
#> Scale for x is already present.
#> Adding another scale for x, which will replace the existing scale.
a <- getPrecursorSpectra("PEPTIDE", 1:2)
plot(a)
The resulting plots show the theoretical mass-to-charge distribution for the provided sequence. A narrow spacing with consistent intensity decay indicates a well-resolved isotope cluster typical for low charge state precursors.
BYion_peak_calculator_DIY enumerates b/y fragment ions
for the peptide, incorporating isotope substitution probabilities.
Provide the desired isotope label and its enrichment probability to
tailor the output to experimental conditions. For unlabeled peptides,
set the abundance to the natural isotope frequency to obtain expected
fragment masses.
df <- BYion_peak_calculator_DIY("PEPTIDE", "C13", 0.0107)
head(df, 5)
#> Mass Prob Kind
#> 1 147.0526 9.341439e-01 Y1
#> 2 148.0556 5.625045e-02 Y1
#> 3 149.0571 9.092087e-03 Y1
#> 4 150.0599 4.783419e-04 Y1
#> 5 151.0617 3.517561e-05 Y1
getSipBYionSpectra generates labeled or unlabeled
fragment spectra, while plotSipBYionLabel annotates
fragment peaks. The charge state range controls which fragments appear,
and the isotope count parameter governs the number of labeled atoms
considered. Use higher charge states when analyzing fragmentation data
collected in higher energy dissociation experiments.
a <- getSipBYionSpectra("KHRIPCDRK", "C13", 0.0107, 1:2, 2)
plot(a) + plotSipBYionLabel(a)
a <- getSipBYionSpectra("KHRIPCDRK", "C13", 0.0107, 1:2, 0)
plot(a) + plotSipBYionLabel(a)
a <- getSipBYionSpectra("KHRIPCDRK", "C13", 0.0107, 1:2, 2:3)
plot(a) + plotSipBYionLabel(a)
These plots highlight how isotope incorporation alters fragment mass distributions. Observe the position and intensity changes when varying isotope counts to anticipate how labeled fragments manifest in MS/MS spectra.
SIP experiments enrich specific isotopes within peptides, shifting precursor and fragment masses. Aerith accommodates enrichment levels from low to complete labeling, helping quantify incorporation levels and differentiate labeled populations.
precursor_peak_calculator_DIY accepts user-defined
isotope abundances, making it suitable for labeled designs. When
combined with calPepPrecursorMass, users can contrast
natural abundance and enriched scenarios. Selecting the isotope
(Atom) and enrichment (Prob) parameters should
reflect experimental labeling; for example, set Prob close
to 1 for near-complete enrichment or match pulse-labeling levels (e.g.,
0.55) for partial incorporation.
a <- precursor_peak_calculator_DIY("PEPTIDECCCC", "N15", 0.55)
head(a, 5)
#> Mass Prob
#> 1 1439.483 2.672840e-06
#> 2 1440.480 5.069828e-05
#> 3 1441.477 4.514761e-04
#> 4 1442.474 2.507844e-03
#> 5 1443.472 9.738476e-03
calPepPrecursorMass("PEPTIDECCCC", "C13", 0.55)
#> [1] 1469.581
The following calls illustrate how precursor mass estimates respond to different isotope species and enrichment levels. Interpret shifts in the resulting masses to gauge the labeling effect and to set accurate extraction windows when processing LC-MS/MS data.
calPepAtomCount("PEPTIDECCCC")
#> C H O N P S
#> 1 54 85 23 15 0 4
calPepPrecursorMass("PEPTIDECCCC", "C13", 0.0107)
#> [1] 1440.484
calPepPrecursorMass("PEPTIDECCCC", "C13", 0.5)
#> [1] 1466.571
calPepPrecursorMass("PEPTIDECCCC", "H2", 0.000115)
#> [1] 1440.484
calPepPrecursorMass("PEPTIDECCCC", "H2", 0.5)
#> [1] 1482.748
calPepPrecursorMass("PEPTIDECCCC", "N15", 0.00368)
#> [1] 1440.484
calPepPrecursorMass("PEPTIDECCCC", "N15", 0.5)
#> [1] 1447.463
calPepPrecursorMass("PEPTIDECCCC", "O18", 0.00205)
#> [1] 1440.484
calPepPrecursorMass("PEPTIDECCCC", "O18", 0.5)
#> [1] 1463.486
calPepPrecursorMass("PEPTIDECCCC", "S34", 0.0429)
#> [1] 1440.484
calPepPrecursorMass("PEPTIDECCCC", "S34", 0.5)
#> [1] 1444.473
calPepAtomCount("PEPTIDECCCC")
#> C H O N P S
#> 1 54 85 23 15 0 4
df <- precursor_peak_calculator_DIY("PEPTIDECCCC", "C13", 0.0107)
df$Mass[which.max(df$Prob)]
#> [1] 1439.483
df <- precursor_peak_calculator_DIY("PEPTIDECCCC", "C13", 0.5)
df$Mass[which.max(df$Prob)]
#> [1] 1466.571
df <- precursor_peak_calculator_DIY("PEPTIDECCCC", "O18", 0.00205)
df$Mass[which.max(df$Prob)]
#> [1] 1439.483
df <- precursor_peak_calculator_DIY("PEPTIDECCCC", "O18", 0.5)
df$Mass[which.max(df$Prob)]
#> [1] 1463.532
df <- precursor_peak_calculator_DIY("PEPTIDECCCC", "S34", 0.0429)
df$Mass[which.max(df$Prob)]
#> [1] 1439.483
df <- precursor_peak_calculator_DIY("PEPTIDECCCC", "S34", 0.5)
df$Mass[which.max(df$Prob)]
#> [1] 1443.476
getSipPrecursorSpectra visualizes isotope envelopes
under enrichment. Analyze how peak intensities redistribute as the
enrichment fraction increases. Wider distributions or shifted apex
masses signify successful labeling, enabling downstream quantitative
comparisons.
a <- getSipPrecursorSpectra("PEPTIDE", "C13", 0.55, 2)
plot(a) + scale_x_continuous(breaks = seq(400, 420, by = 1)) + geom_linerange(linewidth = 0.2)
#> Scale for x is already present.
#> Adding another scale for x, which will replace the existing scale.
a <- getSipPrecursorSpectra("PEPTIDE", "C13", 0.55, 1:2)
plot(a)
Expect the labeled spectra to shift toward higher m/z values relative to the unlabeled counterpart. Monitoring this displacement assists in verifying labeling efficiency across charge states.
When labeling impacts fragmentation, use
BYion_peak_calculator_DIY with enriched probabilities to
determine the dominant fragment masses. Filtering on Kind
enables investigation of particular ion series, which is useful when
tracking diagnostic fragment transitions.
df <- BYion_peak_calculator_DIY("PEPTIDE", "C13", 0.55)
head(df, 5)
#> Mass Prob Kind
#> 1 147.0526 0.01819015 Y1
#> 2 148.0560 0.11127364 Y1
#> 3 149.0593 0.27256130 Y1
#> 4 150.0627 0.33469720 Y1
#> 5 151.0660 0.20723904 Y1
df <- BYion_peak_calculator_DIY("PEPTIDECCCC", "C13", 0.0107)
df <- df[which(df$Kind == "Y6"), ]
df$Mass[which.max(df$Prob)]
#> [1] 902.2022
df <- BYion_peak_calculator_DIY("PEPTIDECCCC", "C13", 0.5)
df <- df[which(df$Kind == "Y6"), ]
df$Mass[which.max(df$Prob)]
#> [1] 917.2506
df <- BYion_peak_calculator_DIY("PEPTIDECCCC", "H2", 0.000115)
df <- df[which(df$Kind == "Y6"), ]
df$Mass[which.max(df$Prob)]
#> [1] 902.2022
df <- BYion_peak_calculator_DIY("PEPTIDECCCC", "H2", 0.5)
df <- df[which(df$Kind == "Y6"), ]
df$Mass[which.max(df$Prob)]
#> [1] 925.3425
df <- BYion_peak_calculator_DIY("PEPTIDECCCC", "N15", 0.00368)
df <- df[which(df$Kind == "Y6"), ]
df$Mass[which.max(df$Prob)]
#> [1] 902.2022
df <- BYion_peak_calculator_DIY("PEPTIDECCCC", "N15", 0.5)
df <- df[which(df$Kind == "Y6"), ]
df$Mass[which.max(df$Prob)]
#> [1] 908.1869
The fragment spectra for labeled peptides reveal how isotope
incorporation propagates across fragment ions. Examine the annotations
produced by plotSipBYionLabel to confirm which ion series
carry the label and to identify shifts that support peptide
identification and quantification.
a <- getSipBYionSpectra("KHRIPCDRK", "C13", 0.55, 1:2, 2)
plot(a) + plotSipBYionLabel(a)
a <- getSipBYionSpectra("KHRIPCDRK", "C13", 0.55, 1:2, 0)
plot(a) + plotSipBYionLabel(a)
a <- getSipBYionSpectra("KHRIPCDRK", "C13", 0.55, 1:2, 2:3)
plot(a) + plotSipBYionLabel(a)
Through these visualizations, Aerith showcases its advantage in modeling SIP-specific fragmentation patterns in line with current state-of-the-art approaches that couple theoretical spectra with high-resolution MS/MS data analysis.
sessionInfo()
#> R Under development (unstable) (2025-11-04 r88984)
#> Platform: aarch64-apple-darwin20
#> Running under: macOS Ventura 13.7.8
#>
#> Matrix products: default
#> BLAS: /System/Library/Frameworks/Accelerate.framework/Versions/A/Frameworks/vecLib.framework/Versions/A/libBLAS.dylib
#> LAPACK: /Library/Frameworks/R.framework/Versions/4.6-arm64/Resources/lib/libRlapack.dylib; LAPACK version 3.12.1
#>
#> locale:
#> [1] C
#>
#> time zone: America/New_York
#> tzcode source: internal
#>
#> attached base packages:
#> [1] stats graphics grDevices utils datasets methods base
#>
#> other attached packages:
#> [1] tidyr_1.3.1 ggplot2_4.0.1 stringr_1.6.0 dplyr_1.1.4 Aerith_0.99.11
#>
#> loaded via a namespace (and not attached):
#> [1] rlang_1.1.6 magrittr_2.0.4
#> [3] clue_0.3-66 matrixStats_1.5.0
#> [5] compiler_4.6.0 systemfonts_1.3.1
#> [7] vctrs_0.6.5 reshape2_1.4.5
#> [9] ProtGenerics_1.43.0 pkgconfig_2.0.3
#> [11] MetaboCoreUtils_1.19.1 crayon_1.5.3
#> [13] fastmap_1.2.0 XVector_0.51.0
#> [15] labeling_0.4.3 rmarkdown_2.30
#> [17] preprocessCore_1.73.0 ragg_1.5.0
#> [19] purrr_1.2.0 xfun_0.54
#> [21] MultiAssayExperiment_1.37.2 cachem_1.1.0
#> [23] jsonlite_2.0.0 DelayedArray_0.37.0
#> [25] BiocParallel_1.45.0 parallel_4.6.0
#> [27] cluster_2.1.8.1 R6_2.6.1
#> [29] bslib_0.9.0 stringi_1.8.7
#> [31] RColorBrewer_1.1-3 limma_3.67.0
#> [33] GenomicRanges_1.63.0 jquerylib_0.1.4
#> [35] Rcpp_1.1.0 Seqinfo_1.1.0
#> [37] SummarizedExperiment_1.41.0 iterators_1.0.14
#> [39] knitr_1.50 IRanges_2.45.0
#> [41] Matrix_1.7-4 igraph_2.2.1
#> [43] tidyselect_1.2.1 dichromat_2.0-0.1
#> [45] abind_1.4-8 yaml_2.3.11
#> [47] doParallel_1.0.17 codetools_0.2-20
#> [49] affy_1.89.0 lattice_0.22-7
#> [51] tibble_3.3.0 plyr_1.8.9
#> [53] Biobase_2.71.0 withr_3.0.2
#> [55] S7_0.2.1 evaluate_1.0.5
#> [57] Spectra_1.21.0 pillar_1.11.1
#> [59] affyio_1.81.0 BiocManager_1.30.27
#> [61] MatrixGenerics_1.23.0 foreach_1.5.2
#> [63] stats4_4.6.0 MSnbase_2.37.0
#> [65] MALDIquant_1.22.3 ncdf4_1.24
#> [67] generics_0.1.4 S4Vectors_0.49.0
#> [69] scales_1.4.0 glue_1.8.0
#> [71] lazyeval_0.2.2 tools_4.6.0
#> [73] mzID_1.49.0 data.table_1.17.8
#> [75] QFeatures_1.21.0 vsn_3.79.0
#> [77] mzR_2.45.0 fs_1.6.6
#> [79] XML_3.99-0.20 grid_4.6.0
#> [81] impute_1.85.0 MsCoreUtils_1.23.1
#> [83] PSMatch_1.15.0 cli_3.6.5
#> [85] textshaping_1.0.4 S4Arrays_1.11.1
#> [87] AnnotationFilter_1.35.0 pcaMethods_2.3.0
#> [89] gtable_0.3.6 sass_0.4.10
#> [91] digest_0.6.39 BiocGenerics_0.57.0
#> [93] SparseArray_1.11.7 ggrepel_0.9.6
#> [95] farver_2.1.2 htmltools_0.5.8.1
#> [97] lifecycle_1.0.4 statmod_1.5.1
#> [99] MASS_7.3-65