CalibraCurve 0.99.2
Targeted mass spectrometry based experiments (e.g. proteomics, metabolomics or lipidomics) are used for validation of biomarkers or for quantitative assays. During the development of these assays, calibration curves are a valuable tool to assess the upper and lower limit of quantification and the linearity of the obtained measurements. The resulting model may also allow to predict concentrations from measured intensities (absolute quantification).
CalibraCurve is a tool to generate and visualize calibration curves. The tool is already established as a KNIME workflow based on an R script (Kohl, Stepath, Bracht, Megger, Sitek, Marcus, and Eisenacher, 2020) and a service in de.NBI and ELIXIR Germany. The code base was re-worked into an R package without changing the underlying algorithm. Several improvements were made:
ggplot2
(Wickham, 2016) for nicer,
publication-ready visualizationCalibraCurve
To install the latest version of CalibraCurve
from Bioconductor use the
following code snippet:
if (!requireNamespace("BiocManager", quietly = TRUE)) {
install.packages("BiocManager")
}
BiocManager::install("CalibraCurve")
If you need help with CalibraCurve
, encountered a bug or have suggestions for
improvements leave a post in the
Bioconductor support site with the
CalibraCurve
tag or start an issue on our
github page.
CalibraCurve
We hope that CalibraCurve will be useful for your research. Please use the following information to cite the package and the overall approach. Thank you!
## Citation info
citation("CalibraCurve")
#> To cite package 'CalibraCurve' in publications use:
#>
#> Kohl M, Stepath M, Bracht T, Megger D, Sitek B, Marcus K, Eisenacher
#> M (2020). "CalibraCurve: A Tool for Calibration of Targeted MS-Based
#> Measurements." _Proteomics_, 1900143. doi:10.1002/pmic.201900143
#> <https://doi.org/10.1002/pmic.201900143>,
#> <https://analyticalsciencejournals.onlinelibrary.wiley.com/doi/10.1002/pmic.201900143>.
#>
#> A BibTeX entry for LaTeX users is
#>
#> @Article{,
#> title = {CalibraCurve: A Tool for Calibration of Targeted MS-Based Measurements},
#> author = {Michael Kohl and Markus Stepath and Thilo Bracht and Dominik A Megger and Barbara Sitek and Katrin Marcus and Martin Eisenacher},
#> year = {2020},
#> journal = {Proteomics},
#> volumne = {20},
#> pages = {1900143},
#> doi = {https://doi.org/10.1002/pmic.201900143},
#> url = {https://analyticalsciencejournals.onlinelibrary.wiley.com/doi/10.1002/pmic.201900143},
#> }
CalibraCurve
CalibraCurve’s main task is to calculate calibration curves for targeted mass spectrometry based data like proteomics, metabolomics or lipidomics. The experimental setup is the following: for the analyte of interest, samples with different known concentrations levels of this analyte are prepared. Ideally, these samples are measured in replicates, leading to intensities or areas as a measurement.
To calculate the calibration curve, two pieces of information are needed: The known concentration levels and the measured intensities for each sample.
There are two options to provide this information for CalibraCurve. First, a set of Excel files may be provided that contain at least a Concentration and a Measurement column each. Second, a SummarizedExperiment object may be provided. This object may contain measurements for different analytes, one analyte per row. Information about the concentration levels have to be provided in the colData part.
CalibraCurve contains test data from the msqc1
R package
(Kockmann, Trachsel, Panse, Wahlander, Selevsek, Grossmann, Wolski, and Schlapbach, 2016). More precisely, parts of the
dilution series from the msqc1_dil
object.
This data was filtered in the following way: - using only QTRAP, TSQVantage and QExactive instruments (PRM or SRM method) - use only data on the level of y-ions, no precursors - use only the heavy isotope version of the peptide
SummarizedExperiment
objects with the corresponding data for each peptide
sequence is given as .rds files.
Additionally, for the peptide sequence “GGPFSDSYR” the data is stored as .xlsx tables, separated by instrument and ion type.
A single .xlsx
, .csv
or .txt
file, containing at least a “Concentration”
and a “Measurement” column, can be imported using the readDataTable
function.
The number of the column belonging to the concentration levels has
to be defined as concCol
.
The number of the column belonging to the substance or analyte has
to be defined as measCol
.
The resulting object can then be directly used in the main function
CalibraCurve
to generate and visualize the Calibration curve.
Here, an example is given using a .xlsx
file:
library("CalibraCurve")
file <- system.file("extdata", "MSQC1_xlsx", "GGPFSDSYR_QExactive_y5.xlsx",
package = "CalibraCurve")
D <- readDataTable(dataPath = file,
fileType = "xlsx",
concCol = 16, # column "amount" containing concentrations
measCol = 12, # column "Area" containing measurements
naStrings = c("NA", "NaN", "Filtered", "#NV"),
sheet = 1)
print(head(D))
#> Concentration Measurement
#> 16 25 63532736
#> 17 25 72636312
#> 18 25 63185244
#> 13 50 142634912
#> 14 50 94896768
#> 15 50 67706488
The resulting object is a dataframe with two columns,
“Concentration” and “Measurement”.
The rows are sorted from lowest to highest concentration levels. This dataframe
can then be directly used in the main function
CalibraCurve
to generate and visualize the Calibration curve.
Alternatively, a whole folder of .xlsx
, .csv
or .txt
files can be
imported using the readMultipleTables
function. All files in the folder have
to follow the same structure, i.e. having the “Concentration” and “Measurement”
column in the same place.
Here is an example using a folder containing 9 .xlsx
files:
library("CalibraCurve")
folder <- system.file("extdata", "MSQC1_xlsx", package = "CalibraCurve")
D <- readMultipleTables(dataFolder = folder,
fileType = "xlsx",
concCol = 16,
measCol = 12)
print(D)
#> $GGPFSDSYR_QExactive_y5
#> Concentration Measurement
#> 16 25 63532736
#> 17 25 72636312
#> 18 25 63185244
#> 13 50 142634912
#> 14 50 94896768
#> 15 50 67706488
#> 10 200 584484352
#> 11 200 408839680
#> 12 200 337659296
#> 7 1000 2688073216
#> 8 1000 2024798720
#> 9 1000 1936675456
#> 4 2000 2988372480
#> 5 2000 2982280448
#> 6 2000 4745236480
#> 1 5000 10947928064
#> 2 5000 16187528192
#> 3 5000 13521389568
#>
#> $GGPFSDSYR_QExactive_y6
#> Concentration Measurement
#> 16 25 42505600
#> 17 25 49378136
#> 18 25 42739496
#> 13 50 94398144
#> 14 50 63172188
#> 15 50 44371872
#> 10 200 390347488
#> 11 200 270616320
#> 12 200 224925344
#> 7 1000 1801370496
#> 8 1000 1353777920
#> 9 1000 1304889088
#> 4 2000 1991971712
#> 5 2000 1988677376
#> 6 2000 3203979008
#> 1 5000 7434415616
#> 2 5000 10999006208
#> 3 5000 9158412288
#>
#> $GGPFSDSYR_QExactive_y7
#> Concentration Measurement
#> 16 25 4357065
#> 17 25 5220379
#> 18 25 4293157
#> 13 50 9211622
#> 14 50 6375533
#> 15 50 4412843
#> 10 200 40054380
#> 11 200 27407376
#> 12 200 12618111
#> 7 1000 178959088
#> 8 1000 137778240
#> 9 1000 129679664
#> 4 2000 190292000
#> 5 2000 196466864
#> 6 2000 331479904
#> 1 5000 723865344
#> 2 5000 1073278592
#> 3 5000 907921792
#>
#> $GGPFSDSYR_QTRAP_y5
#> Concentration Measurement
#> 16 25 636182
#> 17 25 624453
#> 18 25 656312
#> 13 50 1344918
#> 14 50 1391938
#> 15 50 1333177
#> 10 200 4956311
#> 11 200 5285020
#> 12 200 5365300
#> 7 1000 22445056
#> 8 1000 22948416
#> 9 1000 22717552
#> 4 2000 28129572
#> 5 2000 28504010
#> 6 2000 27675910
#> 1 5000 65790580
#> 2 5000 68766688
#> 3 5000 68546360
#>
#> $GGPFSDSYR_QTRAP_y6
#> Concentration Measurement
#> 16 25 158025
#> 17 25 156321
#> 18 25 167139
#> 13 50 304044
#> 14 50 318554
#> 15 50 317089
#> 10 200 1247797
#> 11 200 1283142
#> 12 200 1303153
#> 7 1000 5370834
#> 8 1000 5439073
#> 9 1000 5355637
#> 4 2000 6830844
#> 5 2000 6802424
#> 6 2000 6527015
#> 1 5000 16103898
#> 2 5000 17477076
#> 3 5000 17297898
#>
#> $GGPFSDSYR_QTRAP_y7
#> Concentration Measurement
#> 16 25 19337
#> 17 25 21675
#> 18 25 20447
#> 13 50 45252
#> 14 50 40547
#> 15 50 43635
#> 10 200 160839
#> 11 200 166228
#> 12 200 166188
#> 7 1000 708102
#> 8 1000 715868
#> 9 1000 701413
#> 4 2000 928816
#> 5 2000 914337
#> 6 2000 924336
#> 1 5000 2362928
#> 2 5000 2484053
#> 3 5000 2497078
#>
#> $GGPFSDSYR_TSQVantage_y5
#> Concentration Measurement
#> 16 25 726363
#> 17 25 776707
#> 18 25 638435
#> 13 50 938616
#> 14 50 1499013
#> 15 50 1646262
#> 10 200 5734246
#> 11 200 5328412
#> 12 200 5906098
#> 7 1000 31725796
#> 8 1000 34535764
#> 9 1000 32574730
#> 4 2000 59979412
#> 5 2000 59192548
#> 6 2000 55664044
#> 1 5000 156006784
#> 2 5000 158682160
#> 3 5000 157924304
#>
#> $GGPFSDSYR_TSQVantage_y6
#> Concentration Measurement
#> 16 25 374307
#> 17 25 398118
#> 18 25 349628
#> 13 50 475989
#> 14 50 830615
#> 15 50 894860
#> 10 200 3049708
#> 11 200 2837996
#> 12 200 3143589
#> 7 1000 16948022
#> 8 1000 18544714
#> 9 1000 17462780
#> 4 2000 32045960
#> 5 2000 31821972
#> 6 2000 29836298
#> 1 5000 82836872
#> 2 5000 83833328
#> 3 5000 82925984
#>
#> $GGPFSDSYR_TSQVantage_y7
#> Concentration Measurement
#> 16 25 37309
#> 17 25 51757
#> 18 25 38181
#> 13 50 56026
#> 14 50 89100
#> 15 50 98100
#> 10 200 334119
#> 11 200 318106
#> 12 200 343906
#> 7 1000 1833054
#> 8 1000 1976417
#> 9 1000 1892401
#> 4 2000 3511667
#> 5 2000 3492670
#> 6 2000 3200100
#> 1 5000 8909107
#> 2 5000 8990250
#> 3 5000 8934693
str(D)
#> List of 9
#> $ GGPFSDSYR_QExactive_y5 :'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 6.35e+07 7.26e+07 6.32e+07 1.43e+08 9.49e+07 ...
#> $ GGPFSDSYR_QExactive_y6 :'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 42505600 49378136 42739496 94398144 63172188 ...
#> $ GGPFSDSYR_QExactive_y7 :'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 4357065 5220379 4293157 9211622 6375533 ...
#> $ GGPFSDSYR_QTRAP_y5 :'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 636182 624453 656312 1344918 1391938 ...
#> $ GGPFSDSYR_QTRAP_y6 :'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 158025 156321 167139 304044 318554 ...
#> $ GGPFSDSYR_QTRAP_y7 :'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 19337 21675 20447 45252 40547 ...
#> $ GGPFSDSYR_TSQVantage_y5:'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 726363 776707 638435 938616 1499013 ...
#> $ GGPFSDSYR_TSQVantage_y6:'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 374307 398118 349628 475989 830615 ...
#> $ GGPFSDSYR_TSQVantage_y7:'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 37309 51757 38181 56026 89100 ...
Each file in the folder is imported and the resulting dataframes are combined as a list. Each list entry contains the data belonging to one of the input files, which here present different measurements of the same peptide that differ in their used fragment ion and machine.
A SummarizedExperiment object stored as an .rds file can be imported using the
function readDataSE
function.
The column name of the colData part belonging to the concentration levels have
to be defined as concColName
.
The column name of the rowData part belonging to the substance or analyte have
to be defined as substColName
.
Here, we will use the default options.
library("CalibraCurve")
file <- system.file("extdata", "MSQC1", "msqc1_dil_GGPFSDSYR.rds",
package = "CalibraCurve")
D <- readDataSE(file, concColName = "amount_fmol", substColName = "Substance")
print(D)
#> $QExactive_y5
#> Concentration Measurement
#> 106 25 63532736
#> 107 25 72636312
#> 108 25 63185244
#> 103 50 142634912
#> 104 50 94896768
#> 105 50 67706488
#> 100 200 584484352
#> 101 200 408839680
#> 102 200 337659296
#> 97 1000 2688073216
#> 98 1000 2024798720
#> 99 1000 1936675456
#> 94 2000 2988372480
#> 95 2000 2982280448
#> 96 2000 4745236480
#> 91 5000 10947928064
#> 92 5000 16187528192
#> 93 5000 13521389568
#>
#> $QExactive_y6
#> Concentration Measurement
#> 88 25 42505600
#> 89 25 49378136
#> 90 25 42739496
#> 85 50 94398144
#> 86 50 63172188
#> 87 50 44371872
#> 82 200 390347488
#> 83 200 270616320
#> 84 200 224925344
#> 79 1000 1801370496
#> 80 1000 1353777920
#> 81 1000 1304889088
#> 76 2000 1991971712
#> 77 2000 1988677376
#> 78 2000 3203979008
#> 73 5000 7434415616
#> 74 5000 10999006208
#> 75 5000 9158412288
#>
#> $QExactive_y7
#> Concentration Measurement
#> 70 25 4357065
#> 71 25 5220379
#> 72 25 4293157
#> 67 50 9211622
#> 68 50 6375533
#> 69 50 4412843
#> 64 200 40054380
#> 65 200 27407376
#> 66 200 12618111
#> 61 1000 178959088
#> 62 1000 137778240
#> 63 1000 129679664
#> 58 2000 190292000
#> 59 2000 196466864
#> 60 2000 331479904
#> 55 5000 723865344
#> 56 5000 1073278592
#> 57 5000 907921792
#>
#> $QTRAP_y5
#> Concentration Measurement
#> 52 25 636182
#> 53 25 624453
#> 54 25 656312
#> 49 50 1344918
#> 50 50 1391938
#> 51 50 1333177
#> 46 200 4956311
#> 47 200 5285020
#> 48 200 5365300
#> 43 1000 22445056
#> 44 1000 22948416
#> 45 1000 22717552
#> 40 2000 28129572
#> 41 2000 28504010
#> 42 2000 27675910
#> 37 5000 65790580
#> 38 5000 68766688
#> 39 5000 68546360
#>
#> $QTRAP_y6
#> Concentration Measurement
#> 34 25 158025
#> 35 25 156321
#> 36 25 167139
#> 31 50 304044
#> 32 50 318554
#> 33 50 317089
#> 28 200 1247797
#> 29 200 1283142
#> 30 200 1303153
#> 25 1000 5370834
#> 26 1000 5439073
#> 27 1000 5355637
#> 22 2000 6830844
#> 23 2000 6802424
#> 24 2000 6527015
#> 19 5000 16103898
#> 20 5000 17477076
#> 21 5000 17297898
#>
#> $QTRAP_y7
#> Concentration Measurement
#> 16 25 19337
#> 17 25 21675
#> 18 25 20447
#> 13 50 45252
#> 14 50 40547
#> 15 50 43635
#> 10 200 160839
#> 11 200 166228
#> 12 200 166188
#> 7 1000 708102
#> 8 1000 715868
#> 9 1000 701413
#> 4 2000 928816
#> 5 2000 914337
#> 6 2000 924336
#> 1 5000 2362928
#> 2 5000 2484053
#> 3 5000 2497078
#>
#> $TSQVantage_y5
#> Concentration Measurement
#> 160 25 726363
#> 161 25 776707
#> 162 25 638435
#> 157 50 938616
#> 158 50 1499013
#> 159 50 1646262
#> 154 200 5734246
#> 155 200 5328412
#> 156 200 5906098
#> 151 1000 31725796
#> 152 1000 34535764
#> 153 1000 32574730
#> 148 2000 59979412
#> 149 2000 59192548
#> 150 2000 55664044
#> 145 5000 156006784
#> 146 5000 158682160
#> 147 5000 157924304
#>
#> $TSQVantage_y6
#> Concentration Measurement
#> 142 25 374307
#> 143 25 398118
#> 144 25 349628
#> 139 50 475989
#> 140 50 830615
#> 141 50 894860
#> 136 200 3049708
#> 137 200 2837996
#> 138 200 3143589
#> 133 1000 16948022
#> 134 1000 18544714
#> 135 1000 17462780
#> 130 2000 32045960
#> 131 2000 31821972
#> 132 2000 29836298
#> 127 5000 82836872
#> 128 5000 83833328
#> 129 5000 82925984
#>
#> $TSQVantage_y7
#> Concentration Measurement
#> 124 25 37309
#> 125 25 51757
#> 126 25 38181
#> 121 50 56026
#> 122 50 89100
#> 123 50 98100
#> 118 200 334119
#> 119 200 318106
#> 120 200 343906
#> 115 1000 1833054
#> 116 1000 1976417
#> 117 1000 1892401
#> 112 2000 3511667
#> 113 2000 3492670
#> 114 2000 3200100
#> 109 5000 8909107
#> 110 5000 8990250
#> 111 5000 8934693
str(D)
#> List of 9
#> $ QExactive_y5 :'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 6.35e+07 7.26e+07 6.32e+07 1.43e+08 9.49e+07 ...
#> $ QExactive_y6 :'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 42505600 49378136 42739496 94398144 63172188 ...
#> $ QExactive_y7 :'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 4357065 5220379 4293157 9211622 6375533 ...
#> $ QTRAP_y5 :'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 636182 624453 656312 1344918 1391938 ...
#> $ QTRAP_y6 :'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 158025 156321 167139 304044 318554 ...
#> $ QTRAP_y7 :'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 19337 21675 20447 45252 40547 ...
#> $ TSQVantage_y5:'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 726363 776707 638435 938616 1499013 ...
#> $ TSQVantage_y6:'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 374307 398118 349628 475989 830615 ...
#> $ TSQVantage_y7:'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 37309 51757 38181 56026 89100 ...
The resulting object is a list of dataframes, each with two columns, “Concentration” and “Measurement”. The rows are sorted from lowest to highest concentration levels.
This list of dataframes can then be directly used in the main function
CalibraCurve
to generate and visualize the Calibration curve.
Alternatively, also a SummarizedExperiment object, that was already imported
in R, can be read in directly by readDataSE
:
library("CalibraCurve")
file <- system.file("extdata", "MSQC1", "msqc1_dil_GGPFSDSYR.rds",
package = "CalibraCurve")
SE_data <- readRDS(file)
D <- readDataSE(rawDataSE = SE_data,
concColName = "amount_fmol",
substColName = "Substance")
print(D)
#> $QExactive_y5
#> Concentration Measurement
#> 106 25 63532736
#> 107 25 72636312
#> 108 25 63185244
#> 103 50 142634912
#> 104 50 94896768
#> 105 50 67706488
#> 100 200 584484352
#> 101 200 408839680
#> 102 200 337659296
#> 97 1000 2688073216
#> 98 1000 2024798720
#> 99 1000 1936675456
#> 94 2000 2988372480
#> 95 2000 2982280448
#> 96 2000 4745236480
#> 91 5000 10947928064
#> 92 5000 16187528192
#> 93 5000 13521389568
#>
#> $QExactive_y6
#> Concentration Measurement
#> 88 25 42505600
#> 89 25 49378136
#> 90 25 42739496
#> 85 50 94398144
#> 86 50 63172188
#> 87 50 44371872
#> 82 200 390347488
#> 83 200 270616320
#> 84 200 224925344
#> 79 1000 1801370496
#> 80 1000 1353777920
#> 81 1000 1304889088
#> 76 2000 1991971712
#> 77 2000 1988677376
#> 78 2000 3203979008
#> 73 5000 7434415616
#> 74 5000 10999006208
#> 75 5000 9158412288
#>
#> $QExactive_y7
#> Concentration Measurement
#> 70 25 4357065
#> 71 25 5220379
#> 72 25 4293157
#> 67 50 9211622
#> 68 50 6375533
#> 69 50 4412843
#> 64 200 40054380
#> 65 200 27407376
#> 66 200 12618111
#> 61 1000 178959088
#> 62 1000 137778240
#> 63 1000 129679664
#> 58 2000 190292000
#> 59 2000 196466864
#> 60 2000 331479904
#> 55 5000 723865344
#> 56 5000 1073278592
#> 57 5000 907921792
#>
#> $QTRAP_y5
#> Concentration Measurement
#> 52 25 636182
#> 53 25 624453
#> 54 25 656312
#> 49 50 1344918
#> 50 50 1391938
#> 51 50 1333177
#> 46 200 4956311
#> 47 200 5285020
#> 48 200 5365300
#> 43 1000 22445056
#> 44 1000 22948416
#> 45 1000 22717552
#> 40 2000 28129572
#> 41 2000 28504010
#> 42 2000 27675910
#> 37 5000 65790580
#> 38 5000 68766688
#> 39 5000 68546360
#>
#> $QTRAP_y6
#> Concentration Measurement
#> 34 25 158025
#> 35 25 156321
#> 36 25 167139
#> 31 50 304044
#> 32 50 318554
#> 33 50 317089
#> 28 200 1247797
#> 29 200 1283142
#> 30 200 1303153
#> 25 1000 5370834
#> 26 1000 5439073
#> 27 1000 5355637
#> 22 2000 6830844
#> 23 2000 6802424
#> 24 2000 6527015
#> 19 5000 16103898
#> 20 5000 17477076
#> 21 5000 17297898
#>
#> $QTRAP_y7
#> Concentration Measurement
#> 16 25 19337
#> 17 25 21675
#> 18 25 20447
#> 13 50 45252
#> 14 50 40547
#> 15 50 43635
#> 10 200 160839
#> 11 200 166228
#> 12 200 166188
#> 7 1000 708102
#> 8 1000 715868
#> 9 1000 701413
#> 4 2000 928816
#> 5 2000 914337
#> 6 2000 924336
#> 1 5000 2362928
#> 2 5000 2484053
#> 3 5000 2497078
#>
#> $TSQVantage_y5
#> Concentration Measurement
#> 160 25 726363
#> 161 25 776707
#> 162 25 638435
#> 157 50 938616
#> 158 50 1499013
#> 159 50 1646262
#> 154 200 5734246
#> 155 200 5328412
#> 156 200 5906098
#> 151 1000 31725796
#> 152 1000 34535764
#> 153 1000 32574730
#> 148 2000 59979412
#> 149 2000 59192548
#> 150 2000 55664044
#> 145 5000 156006784
#> 146 5000 158682160
#> 147 5000 157924304
#>
#> $TSQVantage_y6
#> Concentration Measurement
#> 142 25 374307
#> 143 25 398118
#> 144 25 349628
#> 139 50 475989
#> 140 50 830615
#> 141 50 894860
#> 136 200 3049708
#> 137 200 2837996
#> 138 200 3143589
#> 133 1000 16948022
#> 134 1000 18544714
#> 135 1000 17462780
#> 130 2000 32045960
#> 131 2000 31821972
#> 132 2000 29836298
#> 127 5000 82836872
#> 128 5000 83833328
#> 129 5000 82925984
#>
#> $TSQVantage_y7
#> Concentration Measurement
#> 124 25 37309
#> 125 25 51757
#> 126 25 38181
#> 121 50 56026
#> 122 50 89100
#> 123 50 98100
#> 118 200 334119
#> 119 200 318106
#> 120 200 343906
#> 115 1000 1833054
#> 116 1000 1976417
#> 117 1000 1892401
#> 112 2000 3511667
#> 113 2000 3492670
#> 114 2000 3200100
#> 109 5000 8909107
#> 110 5000 8990250
#> 111 5000 8934693
str(D)
#> List of 9
#> $ QExactive_y5 :'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 6.35e+07 7.26e+07 6.32e+07 1.43e+08 9.49e+07 ...
#> $ QExactive_y6 :'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 42505600 49378136 42739496 94398144 63172188 ...
#> $ QExactive_y7 :'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 4357065 5220379 4293157 9211622 6375533 ...
#> $ QTRAP_y5 :'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 636182 624453 656312 1344918 1391938 ...
#> $ QTRAP_y6 :'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 158025 156321 167139 304044 318554 ...
#> $ QTRAP_y7 :'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 19337 21675 20447 45252 40547 ...
#> $ TSQVantage_y5:'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 726363 776707 638435 938616 1499013 ...
#> $ TSQVantage_y6:'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 374307 398118 349628 475989 830615 ...
#> $ TSQVantage_y7:'data.frame': 18 obs. of 2 variables:
#> ..$ Concentration: num [1:18] 25 25 25 50 50 50 200 200 200 1000 ...
#> ..$ Measurement : num [1:18] 37309 51757 38181 56026 89100 ...
As soon as the data have been imported, they can directly be
used inside the CalibraCurve
function. CalibraCurve has a lot of options,
e.g. thresholds or customizing options for the graphics. All of these options
have sensible default values, so you can just try out CalibraCurve without
thinking about these. However for optimal results or publication-ready
visualizations some adjustments may be necessary.
Details and examples about customizing the visualization are given in the
vignette “Customizing the visualizations of CalibraCurve”.
RES <- CalibraCurve(D)
CalibraCurve outputs several tables and two kinds of plots: the calibration
curve and the response factor plots. The calibration curves are stort in the
list entry plot_CC_list
of the result object.
By default, each calibration curve is plotted separately.
Please be aware that by default the plots are not exported (you have to specify
an output_path
for that).
As an example, we plot here the calibration curve for the fourth analyte:
print(RES$plot_CC_list[[4]])
#> Warning in ggplot2::scale_y_continuous(trans = "log10"): log-10 transformation
#> introduced infinite values.
The calibration curve is plotted as a red line, both axes are shown on logarithmic scale. The grey background shows the linear range (here between XY and 1000 fmol), the curve was only fitted to data points inside this range. The grey data point outside of the linear range visibly do not fit to the curve.
Additionally, we can now look at the corresponding response factor plot:
print(RES$plot_RF_list[[4]])
For each data point the so called response factor is calculated. For the data point \(i\) it is defined as (intensity - y intercept)/concentration.
The nominator describes the signal produced by the analyte (intensity minus the intersect) and the denominator the concentration or amount of the analyte. Ideally the mean response factors for each concentration level (shown as the larger dots) would lie in a straight line and not vary much for the different concentration levels. In orange, thresholds based on 80% and 120% of the overall mean response factor (within the linear range) are plotted. Values within these borders are coloured in green, the ones outside in pink. For all concentration levels within the final linear range the response factors lie within these borders, serving as a further quality check of the results and linearity of the response.
For each analyte, two result tables are produced to give inside into the steps of the algorithm. The first one gives an overview of the results on the level of concentration levels:
print(RES$RES[[4]]$result_table_conc_levels)
#> substance concentration mean_measurement estimated_measurement
#> 1 QTRAP_y5 25 638982.3 661747
#> 2 QTRAP_y5 50 1356677.7 1275920
#> 3 QTRAP_y5 200 5202210.3 4960959
#> 4 QTRAP_y5 1000 22703674.7 24614497
#> 5 QTRAP_y5 2000 28103164.0 49181421
#> 6 QTRAP_y5 5000 67701209.3 122882190
#> removed_while_cleaning CV CV_within_thres preliminary_linear_range
#> 1 FALSE 2.521674 TRUE TRUE
#> 2 FALSE 2.292034 TRUE TRUE
#> 3 FALSE 4.165634 TRUE TRUE
#> 4 FALSE 1.109806 TRUE TRUE
#> 5 FALSE 1.475567 TRUE TRUE
#> 6 FALSE 2.449464 TRUE TRUE
#> mean_percentage_bias SD_percentage_bias CV_percentage_bias
#> 1 3.706559 2.623536 70.78090
#> 2 6.574496 2.531495 38.50477
#> 3 4.973153 4.303218 86.52898
#> 4 7.778030 1.025634 13.18630
#> 5 NA NA NA
#> 6 NA NA NA
#> mean_response_factor RF_within_thres final_linear_range
#> 1 23656.34 TRUE TRUE
#> 2 26182.07 TRUE TRUE
#> 3 25773.18 TRUE TRUE
#> 4 22656.10 TRUE TRUE
#> 5 14027.80 FALSE FALSE
#> 6 13530.73 FALSE FALSE
There is one row for each concentration level. It contains the mean measurement, the estimated measurement by the linear model and further information on the different steps of the algorithm (e.g. if the concentration level was removed during the cleaning step or if it was resent in the preliminary linear range). The last column shows which concentration levels ended up in the final linear range.
For this specific analyte, no concentration level was removed during the cleaning step and all levels fulfilled the CV criterion and where included in the preliminary linear range. The two highest concentration levels were subsequently removed during the calculation of the final linear range.
The second table offers a more detailed view for each individual replicate:
print(RES$RES[[4]]$result_table_obs)
#> substance concentration measurement removed_while_cleaning percentage_bias
#> 1 QTRAP_y5 25 636182 FALSE 4.1625105
#> 2 QTRAP_y5 25 624453 FALSE 4.1625105
#> 3 QTRAP_y5 25 656312 FALSE 4.1625105
#> 4 QTRAP_y5 50 1344918 FALSE 6.0722327
#> 5 QTRAP_y5 50 1391938 FALSE 6.0722327
#> 6 QTRAP_y5 50 1333177 FALSE 6.0722327
#> 7 QTRAP_y5 200 4956311 FALSE 0.8849328
#> 8 QTRAP_y5 200 5285020 FALSE 0.8849328
#> 9 QTRAP_y5 200 5365300 FALSE 0.8849328
#> 10 QTRAP_y5 1000 22445056 FALSE 5.6171380
#> 11 QTRAP_y5 1000 22948416 FALSE 5.6171380
#> 12 QTRAP_y5 1000 22717552 FALSE 5.6171380
#> 13 QTRAP_y5 2000 28129572 FALSE NA
#> 14 QTRAP_y5 2000 28504010 FALSE NA
#> 15 QTRAP_y5 2000 27675910 FALSE NA
#> 16 QTRAP_y5 5000 65790580 FALSE NA
#> 17 QTRAP_y5 5000 68766688 FALSE NA
#> 18 QTRAP_y5 5000 68546360 FALSE NA
#> response_factor RF_within_thres final_linear_range
#> 1 23544.32 TRUE TRUE
#> 2 23075.16 TRUE TRUE
#> 3 24349.52 TRUE TRUE
#> 4 25946.88 TRUE TRUE
#> 5 26887.28 TRUE TRUE
#> 6 25712.06 TRUE TRUE
#> 7 24543.69 TRUE TRUE
#> 8 26187.23 TRUE TRUE
#> 9 26588.63 TRUE TRUE
#> 10 22397.48 TRUE TRUE
#> 11 22900.84 TRUE TRUE
#> 12 22669.98 TRUE TRUE
#> 13 14041.00 FALSE FALSE
#> 14 14228.22 FALSE FALSE
#> 15 13814.17 FALSE FALSE
#> 16 13148.60 FALSE FALSE
#> 17 13743.82 FALSE FALSE
#> 18 13699.76 FALSE FALSE
This table offers a more detailed view, e.g. it can be seen that every single replicate within the final linear range also fulfills the response factore criterion.
If multiple analytes are analyzed, a further summary table is generated:
print(RES$summary_tab)
#> substance intercept coeff r2 LLOQ ULOQ
#> QExactive_y5 QExactive_y5 10373596.693 2110785.5772 0.9009888 25 1000
#> QExactive_y6 QExactive_y6 7332829.885 1405113.0256 0.8976782 25 1000
#> QExactive_y7 QExactive_y7 1036751.715 133030.2906 0.8404977 25 1000
#> QTRAP_y5 QTRAP_y5 47573.936 24566.9233 0.9912075 25 1000
#> QTRAP_y6 QTRAP_y6 16843.745 5837.8402 0.9918128 25 1000
#> QTRAP_y7 QTRAP_y7 1965.934 770.7869 0.9885624 25 1000
#> TSQVantage_y5 TSQVantage_y5 -72408.091 30413.7556 0.9824233 25 5000
#> TSQVantage_y6 TSQVantage_y6 -43283.577 16247.5431 0.9808731 25 5000
#> TSQVantage_y7 TSQVantage_y7 -2770.665 1758.9092 0.9789053 25 5000
#> eq
#> QExactive_y5 y = 1e+07 + 2.1e+06 * x (R2 = 0.901)
#> QExactive_y6 y = 7.3e+06 + 1.4e+06 * x (R2 = 0.898)
#> QExactive_y7 y = 1e+06 + 1.3e+05 * x (R2 = 0.84)
#> QTRAP_y5 y = 4.8e+04 + 2.5e+04 * x (R2 = 0.991)
#> QTRAP_y6 y = 1.7e+04 + 5.8e+03 * x (R2 = 0.992)
#> QTRAP_y7 y = 2e+03 + 7.7e+02 * x (R2 = 0.989)
#> TSQVantage_y5 y = -7.2e+04 + 3e+04 * x (R2 = 0.982)
#> TSQVantage_y6 y = -4.3e+04 + 1.6e+04 * x (R2 = 0.981)
#> TSQVantage_y7 y = -2.8e+03 + 1.8e+03 * x (R2 = 0.979)
This table contains one row per analyte and contains more detailed information on the linear model (coefficients and R squared) and the lower and upper limit of quantification. In this case it can be seen that the QExactive leads to the worst R squared value overall, while the TSQCVantage lead to a wider linear range.
The resulting calibration curves can be used to predict the analytes’
concentration based on measured intensities. You need to define a vector with
measured intensities and
feed it into the predictConcentration
function together with the Calibracurve
result object for this particular analyte. Imagine we have measured three new
samples, with an intensity of 1e6, 1e7 and 1e8, respectively.
newdata <- c(1000000, #1e6
10000000, #1e7
100000000) # 1e8
CC_res = RES$RES[[4]]
predictConcentration(CC_res = RES$RES[[4]], newdata = newdata)
#> Warning in predictConcentration(CC_res = RES$RES[[4]], newdata = newdata): At least one predicted concentration is outside the final
#> linear range. Results may be unreliable.
#> intensity predConc linear_range
#> 1 1e+06 38.76863 TRUE
#> 2 1e+07 405.11488 TRUE
#> 3 1e+08 4068.57728 FALSE
In the table, the predicted concentrations based on the linear model are given. Please note the warning message: The predicted concentration for the intensity 1e8 falls outside of the linear range for this analyte. Therefore this value may be unreliable, like is shown by the warning message.
#> ─ Session info ───────────────────────────────────────────────────────────────────────────────────────────────────────
#> setting value
#> version R version 4.5.1 Patched (2025-08-23 r88802)
#> os Ubuntu 24.04.3 LTS
#> system x86_64, linux-gnu
#> ui X11
#> language (EN)
#> collate C
#> ctype en_US.UTF-8
#> tz America/New_York
#> date 2025-10-16
#> pandoc 2.7.3 @ /usr/bin/ (via rmarkdown)
#> quarto 1.7.32 @ /usr/local/bin/quarto
#>
#> ─ Packages ───────────────────────────────────────────────────────────────────────────────────────────────────────────
#> package * version date (UTC) lib source
#> abind 1.4-8 2024-09-12 [2] CRAN (R 4.5.1)
#> backports 1.5.0 2024-05-23 [2] CRAN (R 4.5.1)
#> bibtex 0.5.1 2023-01-26 [2] CRAN (R 4.5.1)
#> Biobase 2.69.1 2025-10-16 [2] Bioconductor 3.22 (R 4.5.1)
#> BiocGenerics 0.55.3 2025-10-16 [2] Bioconductor 3.22 (R 4.5.1)
#> BiocManager 1.30.26 2025-06-05 [2] CRAN (R 4.5.1)
#> BiocStyle * 2.37.1 2025-10-16 [2] Bioconductor 3.22 (R 4.5.1)
#> bookdown 0.45 2025-10-03 [2] CRAN (R 4.5.1)
#> bslib 0.9.0 2025-01-30 [2] CRAN (R 4.5.1)
#> cachem 1.1.0 2024-05-16 [2] CRAN (R 4.5.1)
#> CalibraCurve * 0.99.2 2025-10-16 [1] Bioconductor 3.22 (R 4.5.1)
#> checkmate 2.3.3 2025-08-18 [2] CRAN (R 4.5.1)
#> cli 3.6.5 2025-04-23 [2] CRAN (R 4.5.1)
#> crayon 1.5.3 2024-06-20 [2] CRAN (R 4.5.1)
#> DelayedArray 0.35.3 2025-10-16 [2] Bioconductor 3.22 (R 4.5.1)
#> dichromat 2.0-0.1 2022-05-02 [2] CRAN (R 4.5.1)
#> digest 0.6.37 2024-08-19 [2] CRAN (R 4.5.1)
#> dplyr 1.1.4 2023-11-17 [2] CRAN (R 4.5.1)
#> evaluate 1.0.5 2025-08-27 [2] CRAN (R 4.5.1)
#> farver 2.1.2 2024-05-13 [2] CRAN (R 4.5.1)
#> fastmap 1.2.0 2024-05-15 [2] CRAN (R 4.5.1)
#> generics 0.1.4 2025-05-09 [2] CRAN (R 4.5.1)
#> GenomicRanges 1.61.5 2025-10-16 [2] Bioconductor 3.22 (R 4.5.1)
#> ggplot2 4.0.0 2025-09-11 [2] CRAN (R 4.5.1)
#> glue 1.8.0 2024-09-30 [2] CRAN (R 4.5.1)
#> gtable 0.3.6 2024-10-25 [2] CRAN (R 4.5.1)
#> htmltools 0.5.8.1 2024-04-04 [2] CRAN (R 4.5.1)
#> httr 1.4.7 2023-08-15 [2] CRAN (R 4.5.1)
#> IRanges 2.43.5 2025-10-16 [2] Bioconductor 3.22 (R 4.5.1)
#> jquerylib 0.1.4 2021-04-26 [2] CRAN (R 4.5.1)
#> jsonlite 2.0.0 2025-03-27 [2] CRAN (R 4.5.1)
#> knitr 1.50 2025-03-16 [2] CRAN (R 4.5.1)
#> labeling 0.4.3 2023-08-29 [2] CRAN (R 4.5.1)
#> lattice 0.22-7 2025-04-02 [3] CRAN (R 4.5.1)
#> lifecycle 1.0.4 2023-11-07 [2] CRAN (R 4.5.1)
#> lubridate 1.9.4 2024-12-08 [2] CRAN (R 4.5.1)
#> magrittr 2.0.4 2025-09-12 [2] CRAN (R 4.5.1)
#> Matrix 1.7-4 2025-08-28 [3] CRAN (R 4.5.1)
#> MatrixGenerics 1.21.0 2025-10-16 [2] Bioconductor 3.22 (R 4.5.1)
#> matrixStats 1.5.0 2025-01-07 [2] CRAN (R 4.5.1)
#> openxlsx 4.2.8 2025-01-25 [2] CRAN (R 4.5.1)
#> pillar 1.11.1 2025-09-17 [2] CRAN (R 4.5.1)
#> pkgconfig 2.0.3 2019-09-22 [2] CRAN (R 4.5.1)
#> plyr 1.8.9 2023-10-02 [2] CRAN (R 4.5.1)
#> purrr 1.1.0 2025-07-10 [2] CRAN (R 4.5.1)
#> R6 2.6.1 2025-02-15 [2] CRAN (R 4.5.1)
#> RColorBrewer 1.1-3 2022-04-03 [2] CRAN (R 4.5.1)
#> Rcpp 1.1.0 2025-07-02 [2] CRAN (R 4.5.1)
#> RefManageR * 1.4.0 2022-09-30 [2] CRAN (R 4.5.1)
#> rlang 1.1.6 2025-04-11 [2] CRAN (R 4.5.1)
#> rmarkdown 2.30 2025-09-28 [2] CRAN (R 4.5.1)
#> S4Arrays 1.9.1 2025-10-16 [2] Bioconductor 3.22 (R 4.5.1)
#> S4Vectors 0.47.4 2025-10-16 [2] Bioconductor 3.22 (R 4.5.1)
#> S7 0.2.0 2024-11-07 [2] CRAN (R 4.5.1)
#> sass 0.4.10 2025-04-11 [2] CRAN (R 4.5.1)
#> scales 1.4.0 2025-04-24 [2] CRAN (R 4.5.1)
#> Seqinfo 0.99.2 2025-10-16 [2] Bioconductor 3.22 (R 4.5.1)
#> sessioninfo * 1.2.3 2025-02-05 [2] CRAN (R 4.5.1)
#> SparseArray 1.9.1 2025-10-16 [2] Bioconductor 3.22 (R 4.5.1)
#> stringi 1.8.7 2025-03-27 [2] CRAN (R 4.5.1)
#> stringr 1.5.2 2025-09-08 [2] CRAN (R 4.5.1)
#> SummarizedExperiment 1.39.2 2025-10-16 [2] Bioconductor 3.22 (R 4.5.1)
#> tibble 3.3.0 2025-06-08 [2] CRAN (R 4.5.1)
#> tidyr 1.3.1 2024-01-24 [2] CRAN (R 4.5.1)
#> tidyselect 1.2.1 2024-03-11 [2] CRAN (R 4.5.1)
#> timechange 0.3.0 2024-01-18 [2] CRAN (R 4.5.1)
#> vctrs 0.6.5 2023-12-01 [2] CRAN (R 4.5.1)
#> withr 3.0.2 2024-10-28 [2] CRAN (R 4.5.1)
#> xfun 0.53 2025-08-19 [2] CRAN (R 4.5.1)
#> xml2 1.4.0 2025-08-20 [2] CRAN (R 4.5.1)
#> XVector 0.49.1 2025-10-16 [2] Bioconductor 3.22 (R 4.5.1)
#> yaml 2.3.10 2024-07-26 [2] CRAN (R 4.5.1)
#> zip 2.3.3 2025-05-13 [2] CRAN (R 4.5.1)
#>
#> [1] /tmp/Rtmptowtjq/Rinst3ca11811b3798
#> [2] /home/biocbuild/bbs-3.22-bioc/R/site-library
#> [3] /home/biocbuild/bbs-3.22-bioc/R/library
#> * ── Packages attached to the search path.
#>
#> ──────────────────────────────────────────────────────────────────────────────────────────────────────────────────────
This vignette was generated using BiocStyle (Oleś, 2025) with knitr (Xie, 2025) and rmarkdown (Allaire, Xie, Dervieux, McPherson, Luraschi, Ushey, Atkins, Wickham, Cheng, Chang, and Iannone, 2025) running behind the scenes.
Citations made with RefManageR (McLean, 2017).
Allaire, J., Y. Xie, et al. (2025). rmarkdown: Dynamic Documents for R. R package version 2.30. URL: https://github.com/rstudio/rmarkdown.
Kockmann, T., C. Trachsel, et al. (2016). “Targeted proteomics coming of age - SRM, PRM, and DIA performance evaluated from a core facility perspective”. In: PROTEOMICS. DOI: 10.1002/pmic.201500502. URL: http://onlinelibrary.wiley.com/doi/10.1002/pmic.201500502/full.
Kohl, M., M. Stepath, et al. (2020). “CalibraCurve: A Tool for Calibration of Targeted MS-Based Measurements”. In: Proteomics, p. 1900143. DOI: https://doi.org/10.1002/pmic.201900143. URL: https://analyticalsciencejournals.onlinelibrary.wiley.com/doi/10.1002/pmic.201900143.
McLean, M. W. (2017). “RefManageR: Import and Manage BibTeX and BibLaTeX References in R”. In: The Journal of Open Source Software. DOI: 10.21105/joss.00338.
Oleś, A. (2025). BiocStyle: Standard styles for vignettes and other Bioconductor documents. R package version 2.37.1. DOI: 10.18129/B9.bioc.BiocStyle. URL: https://bioconductor.org/packages/BiocStyle.
Wickham, H. (2016). ggplot2: Elegant Graphics for Data Analysis. Springer-Verlag New York. ISBN: 978-3-319-24277-4. URL: https://ggplot2.tidyverse.org.
Xie, Y. (2025). knitr: A General-Purpose Package for Dynamic Report Generation in R. R package version 1.50. URL: https://yihui.org/knitr/.