Mass spectrometry is a powerful tool in biomedical research. Advancements in label-free methods and MS instruments have enabled high-throughput proteomics profiling of biological samples with the increasingly deeper coverage of the proteome, but missing values are still ubiquitous in MS-based proteomics data. Single-cell MS-based proteomics has also started to become available, with increased frequency of missing values. The limpa package implements a pipeline for quantification and differential expression analysis of mass-spectrometry (MS) proteomics data, with probabilistic information recovery from missing values. The package accepts peptide-level data, usually including missing values, and produces complete protein quantifications without the need for imputation. A key feature is the ability to quantify an expression value for every protein and every sample, regardless of the number of quantified peptides for each protein or the number of missing values at the peptide level. Another key feature is the ability to propagate quantification uncertainty through to the differential expression analysis, using variance modeling and precision weights, increasing statistical power while avoiding false discoveries. limpa produces a linear model object suitable for downstream analysis with the limma package, allowing complex experimental designs and other downstream tasks such as the gene ontology or pathway analysis. The package name (“limpa”) is an acronym for “Linear Models for Proteomics Data”. limpa package version: 1.0.3.
limpa assumes that label-free proteomics profiling has been performed on a series of biological samples, and that peptide quantifications have been obtained for each sample (Luo et al 2023). limma reads the peptide or peptide-precursor data from popular software tools such as DIA-NN (https://github.com/vdemichev/DiaNN) (Demichev et al 2020), Spectronaut (https://staging.biognosys.com/software/spectronaut/) or MaxQuant (Tyanova et al 2016). In each case, the peptide quantifications are read into a wide-format (peptide by sample) matrix of log-intensities. The only other information required to start a limpa analysis is a vector specifying the protein or protein-group membership of each peptide.
Peptide quantifications typically include missing values for some samples, which cause considerable problems for downstream analyses. The missing values can’t just be ignored, because they occur more often at low intensities, and missing value imputation also introduces biases and problems. limpa uses statistical models developed by Li & Smyth (2023) to evaluate how much information can be recovered from the missing values, quantifying the overall missing value “mechanism” into a detection probability curve (DPC). limpa uses the DPC, together with a Bayesian model, to summarize the peptide quantifications into protein-level quantifications and to undertake limma-based differential expression analyses (Li 2024, Ritchie et al 2015).
If y.peptide is a matrix of peptide-level log2-intensities (including NAs), protein.id is a vector of protein IDs, and design is a design matrix, then the following code will quantify complete log2-expression for the proteins without missing values and will conduct a differential expression analysis defined by the design matrix.
library(limpa)
dpcfit <- dpc(y.peptide)
y.protein <- dpcQuant(y.peptide, protein.id, dpc=dpcfit)
fit <- dpcDE(y.protein, design)
fit <- eBayes(fit)
topTable(fit)Real proteomics datasets are too large to analyze in this vignette, but we can demonstrate a complete reproducible analysis using a small simulated data. First, generate the dataset:
> library(limpa)Loading required package: limma> set.seed(20241230)
> y.peptide <- simProteinDataSet()The dataset is stored as a limma EList object, with components E (log2-expression), genes (feature annotation) and targets (sample annotation).
The simulation function can generate any number of peptides or proteins but, by default, the dataset has 100 peptides belonging to 25 proteins and the samples are in two groups with \(n=5\) replicates in each group.
About 40% of the expression values are missing.
> dim(y.peptide)[1] 100  10> head(y.peptide$genes)             Protein DEStatus
Peptide001 Protein01    NotDE
Peptide002 Protein01    NotDE
Peptide003 Protein01    NotDE
Peptide004 Protein01    NotDE
Peptide005 Protein02    NotDE
Peptide006 Protein02    NotDE> table(y.peptide$targets$Group)
1 2 
5 5 > mean(is.na(y.peptide$E))[1] 0.389Next we estimate the intercept and slope of the detection probability curve, which relates the probability of detection to the underlying peptide expression level on the logit scale.
> dpcfit <- dpc(y.peptide)1 peptides are completely missing in all samples.> dpcfit$dpc  beta0   beta1 
-3.8185  0.7455 > plotDPC(dpcfit)Then we use the DPC to quantify the protein log2-expression values, using the DPC to represent the missing values.
> y.protein <- dpcQuant(y.peptide, "Protein", dpc=dpcfit)Estimating hyperparameters ...Quantifying proteins ...Proteins: 25 Peptides: 100There are no longer any missing values, and the samples now cluster into groups.
The plotMDSUsingSEs() function is similar to plotMDS() in the limma package, but takes account of the standard errors generated by dpcQuant().
> plotMDSUsingSEs(y.protein)> Group <- factor(y.peptide$targets$Group)
> Group.color <- Group
> levels(Group.color) <- c("blue","red")
> plotMDSUsingSEs(y.protein, pch=16, col=as.character(Group.color))Finally, we conduct a differential expression analysis using the limma package.
The dpcDE function calls limma’s voomaLmFit function, which was specially developed for use with limpa.
voomaLmFit computes precision weights, in a similar way to voom for RNA-seq, but instead of using count sizes it use the quantification precisions from dpcQuant.
The plot shows how dpcDE predicts the protein-wise variances from the quantification uncertainties and expression levels.
> design <- model.matrix(~Group)
> fit <- dpcDE(y.protein, design, plot=TRUE)> fit <- eBayes(fit)
> topTable(fit, coef=2)            Protein DEStatus NPeptides PropObs   logFC AveExpr       t
Protein23 Protein23       Up         4   0.950  0.9339   9.244  4.8083
Protein24 Protein24     Down         4   0.975 -0.8071   9.330 -4.0091
Protein22 Protein22       Up         4   1.000  0.6810   8.912  3.9458
Protein08 Protein08       Up         4   0.450  0.6313   4.728  2.3289
Protein13 Protein13    NotDE         4   0.650 -0.5068   5.976 -2.2626
Protein07 Protein07    NotDE         4   0.425  0.4689   4.285  1.7542
Protein11 Protein11    NotDE         4   0.650  0.3455   5.189  1.3688
Protein21 Protein21    NotDE         4   0.850  0.2861   8.486  1.4679
Protein10 Protein10    NotDE         4   0.375  0.2369   4.850  0.9797
Protein03 Protein03     Down         4   0.250 -0.2145   2.724 -0.8672
            P.Value adj.P.Val       B
Protein23 3.026e-06 7.566e-05  4.1967
Protein24 8.659e-05 9.230e-04  1.0706
Protein22 1.108e-04 9.230e-04  0.7351
Protein08 2.088e-02 1.238e-01 -3.6447
Protein13 2.475e-02 1.238e-01 -3.9470
Protein07 8.097e-02 3.374e-01 -4.7688
Protein11 1.726e-01 5.395e-01 -5.3940
Protein21 1.437e-01 5.134e-01 -5.5073
Protein10 3.284e-01 8.872e-01 -5.8837
Protein03 3.869e-01 8.872e-01 -5.9558This small dataset has five truly DE proteins. Four of the give are top-ranked in the DE results. The other DE protein is ranked 10th in the DE results and does not achieve statistical significant because it had only 17% detected observations and, hence, a high quantification uncertainty.
To view the log-expression values for the top DE protein, together with standard errors:
> plotProtein(y.protein, "Protein23", col=as.character(Group.color))For some applications, it is desirable to complete and impute the peptide-level intensities without summarizing them at the protein level. This is can done in limpa by
y.peptide.complete <- dpcImpute(y.peptide, dpc=dpcfit)> sessionInfo()R version 4.5.0 RC (2025-04-04 r88126)
Platform: x86_64-pc-linux-gnu
Running under: Ubuntu 24.04.2 LTS
Matrix products: default
BLAS:   /home/biocbuild/bbs-3.21-bioc/R/lib/libRblas.so 
LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0  LAPACK version 3.12.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] stats     graphics  grDevices utils     datasets  methods   base     
other attached packages:
[1] limpa_1.0.3      limma_3.64.1     BiocStyle_2.36.0
loaded via a namespace (and not attached):
 [1] cli_3.6.5           knitr_1.50          rlang_1.1.6        
 [4] magick_2.8.6        xfun_0.52           jsonlite_2.0.0     
 [7] data.table_1.17.4   statmod_1.5.0       htmltools_0.5.8.1  
[10] tinytex_0.57        sass_0.4.10         rmarkdown_2.29     
[13] evaluate_1.0.3      jquerylib_0.1.4     fastmap_1.2.0      
[16] yaml_2.3.10         lifecycle_1.0.4     bookdown_0.43      
[19] BiocManager_1.30.25 compiler_4.5.0      Rcpp_1.0.14        
[22] digest_0.6.37       R6_2.6.1            magrittr_2.0.3     
[25] bslib_0.9.0         tools_4.5.0         cachem_1.1.0       The limpa pipeline starts with a matrix of peptide precursor intensities (rows for peptides and columns for samples) and a character vector of protein IDs.
The input data can be conveniently supplied as a limma EList object, but a plain numeric matrix containing the log-intensities is also acceptable.
Non-detected peptides should be entered as NA.
limpa includes the functions readDIANN() and readSpectronaut(), which read peptide or peptide-precursor level data output by the popular software tools DIA-NN or Spectronaut respectively.
For example,
> y.peptide <- readDIANN("Report.tsv")reads the Report.tsv file written by DIA-NN from the current working directory.
We are aware that DIA-NN no longer exports its main report in the .tsv format, but in the .parquet format from version 2.0.
To process DIA-NN data in the .parquet format, use
> y.peptide <- readDIANN("Report.parquet", format="parquet")For data searched with match-between-run (MBR), as suggested by Thierry Nordmann (Max Planck Institute of Biochemistry), we recommend to filter the observations by the following
> y.peptide <- readDIANN("report.tsv",
+              q.columns = c("Q.Value","Lib.Q.Value","Lib.PG.Q.Value"),
+              q.cutoffs = 0.01)If searched without MBR, we recommend
> y.peptide <- readDIANN("report.tsv",
+              q.columns = c("Q.Value","Global.Q.Value","Global.PG.Q.Value"),
+              q.cutoffs = 0.01)After reading in the data, it is common to filter out non-proteotypic peptides by
> y.peptide <- filterNonProteotypicPeptides(y.peptide)and to filter out compound protein groups (protein groups mapped to two or more protein IDs) by
> y.peptide <- filterCompoundProteins(y.peptide)These filtering steps are not required by limpa but help with downstream interpretation of the results.
Protein groups can also be filtered by number of peptides, typically to remove proteins with only one detected precursor:
> y.peptide <- filterSingletonPeptides(y.peptide, min.n.peptides = 2)This step is entirely optional.
We generally recommend that users keep all proteins in order to retain maximum information.
The dpcQuant() function will still quantify complete data even for proteins with just one peptide.
Note that peptides do not need to be filtered based on the proportion of detected or missing values. limpa can process peptides correctly even if the number of detections is small.
Missing values for some peptides in some samples has complicated the analysis of MS proteomics data.
Peptides with very low expression values are frequently not detected, but peptides at high expression levels may also be undetected for a variety of reasons that are not completely understood or easily predictable, for example ambiguity of their elution profile with that of other peptides.
If y is the true expression level of a particular peptide in a particular sample (on the log2 scale), then limpa assumes that the probability of detection is given by
\[P(D | y) = F(\beta_0 + \beta_1 y)\]
where \(D\) indicates detection, \(\beta_0\) and \(\beta_1\) are the intercept and slope of the DPC and \(F\) is the logistic function, given by plogis in R.
This probability relationship is called the detection probability curve (DPC) in limpa.
The slope \(\beta_1\) measures how dependent the missing value process is on the underlying expression level.
A slope of zero would means completely random missing values, while very large slopes correspond to left censoring.
The DPC allows limpa to recover information in a probabilistic manner from the missing values.
The larger the slope, the more information there is to recover.
We typically find \(\beta_1\) values between about 0.7 and 1 to be representative of real MS data.
The DPC is difficult to estimate because y is only observed for detected peptides, and the detected values are a biased representation of the complete values that in principle might have been observed had the missing value mechanism not operated.
limpa uses a mathematical exponential tilting argument to represent the DPC in terms of observed values only, which provides a means to estimate the DPC from real data.
The DPC slope \(\beta_1\) is nevertheless often under-estimated if the variability of each peptide is large.
limpa uses the DPC, together with a Bayesian model, to estimate the expression level of each protein in each sample. It fits an additive model \[\mu_{ij} = \gamma_i + \delta_j\] to the peptide log-intensities for each protein, where \(\mu_{ij}\) is the expected log-expression of peptide \(j\) in sample \(i\), \(\gamma_i\) is the log-expression level of the protein in sample \(i\) and \(\delta_j\) is the baseline effect for peptide \(j\). A sum-to-zero constraint is applied to the peptide effects \(\delta_j\) so that the protein expression \(\gamma_i\) represents the average log-expression of the peptides in sample \(i\). The log-likelihood consists of squared residuals for each observed peptide value and the log probability of being missing for each non-detected peptide value. A multivariate normal prior is also applied to the protein log-expression values, where the prior is estimated from the global data for all proteins. DPC-Quant maximizes the log-posterior with respect to the \(\gamma_i\)s and \(\delta_j\)s, and the final \(\gamma_i\)s become the protein quantifications. DPC-Quant also returns the posterior standard error with which each log-expression value is estimated.
Finally, the protein log2-expression values and associated uncertainties are passed to the voomaLmFitWithImputation function, which computes precision weights for each observation and fits proteinwise linear models.
voomaLmFitWithImputation combines features from the voomaLmFit and voomLmFit functions in the limma and edgeR packages.
It uses both protein expression and the quantification standard errors to predict the protein-wise variances and, hence, to construct precision weights for downstream linear modelling.
This allows the uncertainty associated with missing values imputation to be propagated through to the differential expression analysis.
voomaLmFitWithImputation also gives special consideration to instances where all the expression values for a particular protein are imputed for one or more treatment conditions, ensuring robust differential expression analyses even for proteins with only a small proportion of detected peptides.
limpa’s dpcDE function is a wrapper function, passing the appropriate standard errors from dpcQuant to voomaLmFitWithImputation.
The limpa package is fully compatible with limma analysis pipelines, allowing any complex experimental design and other downstream tasks such as the gene ontology or pathway analysis.
The dpcDE() function accepts any argument that voomaLmFitWithImputation (or voomaLmFit) does.
For example, dpcDE(y.protein, design, sample.weights=TRUE) can be used to downweight outlier samples.
Or dpc(y.protein, design, block=subject) could be used to model the correlation between repeated observations on the same subject.
Any questions about limpa can be sent to the Bioconductor support forum.
Li M, Cobbold SA, Smyth GK (2025). Quantification and differential analysis of mass spectrometry proteomics data with probabilistic recovery of information from missing values. bioRxiv 2025/651125. doi:10.1101/2025.04.28.651125
The limpa project was supported by Chan Zuckerberg Initiative EOSS grant 2021-237445, by Melbourne Research and CSL Translational Data Science Scholarships to ML, and by NHMRC Investigator Grant 2025645 to GS.
Demichev V, Messner CB, Vernardis SI, Lilley KS, Ralser M (2020). DIA-NN: neural networks and interference correction enable deep proteome coverage in high throughput. Nature Methods 17(1), 41-44.
Li M, Smyth GK (2023). Neither random nor censored: estimating intensity-dependent probabilities for missing values in label-free proteomics. Bioinformatics 39(5), btad200. doi:10.1093/bioinformatics/btad200
Li M, Cobbold SA, Smyth GK (2025). Quantification and differential analysis of mass spectrometry proteomics data with probabilistic recovery of information from missing values. bioRxiv 2025/651125. doi:10.1101/2025.04.28.651125
Li M (2024). Linear Models and Empirical Bayes Methods for Mass Spectrometry-based Proteomics Data. PhD Thesis, University of Melbourne. http://hdl.handle.net/11343/351600
Lou R, Cao Y, Li S, Lang X, Li Y, Zhang Y, Shui W (2023). Benchmarking commonly used software suites and analysis workflows for DIA proteomics and phosphoproteomics. Nature Communications 14(1), 94.
Ritchie ME, Phipson B, Wu D, Hu Y, Law CW, Shi W, Smyth GK (2015). limma powers differential expression analyses for RNA-sequencing and microarray studies. Nucleic Acids Research 43, e47. doi:10.1093/nar/gkv007
Tyanova S, Temu T, Cox J (2016). The MaxQuant computational platform for mass spectrometry-based shotgun proteomics. Nature Protocols 11, 2301-2319.