Skip to content

Commit fbfd6c9

Browse files
authored
Merge pull request #303 from jorainer/devel
feat: extract electronBeamEnergy MS term
2 parents 0368dbf + b90fb90 commit fbfd6c9

File tree

4 files changed

+15
-3
lines changed

4 files changed

+15
-3
lines changed

DESCRIPTION

+1-1
Original file line numberDiff line numberDiff line change
@@ -2,7 +2,7 @@ Package: mzR
22
Type: Package
33
Title: parser for netCDF, mzXML and mzML and mzIdentML files
44
(mass spectrometry data)
5-
Version: 2.41.1
5+
Version: 2.41.2
66
Author: Bernd Fischer, Steffen Neumann, Laurent Gatto, Qiang Kou, Johannes Rainer
77
Authors@R: c(
88
person("Steffen", "Neumann", email="[email protected]", role=c("aut","cre")),

NEWS

+4
Original file line numberDiff line numberDiff line change
@@ -1,3 +1,7 @@
1+
CHANGES IN VERSION 2.41.2
2+
-------------------------
3+
o Report also electron beam energy (MS:1003410) for EAD data.
4+
15
CHANGES IN VERSION 2.41.1
26
-------------------------
37
o Fix compilation error with stricter compiler checks

inst/unitTests/test_pwiz.R

+2
Original file line numberDiff line numberDiff line change
@@ -20,6 +20,8 @@ test_mzXML <- function() {
2020
checkTrue(any(colnames(hdr) == "scanWindowUpperLimit"))
2121
checkTrue(all(!is.na(hdr$scanWindowLowerLimit)))
2222
checkTrue(all(!is.na(hdr$scanWindowUpperLimit)))
23+
checkTrue(any(colnames(hdr) == "electronBeamEnergy"))
24+
checkTrue(all(is.na(hdr$electronBeamEnergy)))
2325
hdr <- header(mzxml,1)
2426
checkTrue(is.list(hdr))
2527
hdr <- header(mzxml, 2:3)

src/RcppPwiz.cpp

+8-2
Original file line numberDiff line numberDiff line change
@@ -200,7 +200,8 @@ Rcpp::DataFrame RcppPwiz::getScanHeaderInfo (Rcpp::IntegerVector whichScan) {
200200
Rcpp::NumericVector isolationWindowUpperOffset(N_scans);
201201
Rcpp::NumericVector scanWindowLowerLimit(N_scans);
202202
Rcpp::NumericVector scanWindowUpperLimit(N_scans);
203-
203+
Rcpp::NumericVector electronBeamEnergy(N_scans);
204+
204205
for (size_t i = 0; i < N_scans; i++) {
205206
int current_scan = whichScan[i];
206207
size_t current_index = static_cast<size_t>(current_scan - 1);
@@ -262,6 +263,8 @@ Rcpp::DataFrame RcppPwiz::getScanHeaderInfo (Rcpp::IntegerVector whichScan) {
262263
const Precursor& precursor = sp->precursors[0];
263264
// collisionEnergy
264265
collisionEnergy[i] = precursor.activation.cvParam(MS_collision_energy).valueAs<double>();
266+
// electronBeamEnergy - for EAD
267+
electronBeamEnergy[i] = precursor.activation.cvParam(MS_electron_beam_energy).value.empty() ? NA_REAL : precursor.activation.cvParam(MS_electron_beam_energy).valueAs<double>();
265268
// precursorScanNum
266269
size_t precursorIndex = slp->find(precursor.spectrumID);
267270
if (precursorIndex < N) {
@@ -288,6 +291,7 @@ Rcpp::DataFrame RcppPwiz::getScanHeaderInfo (Rcpp::IntegerVector whichScan) {
288291
}
289292
} else {
290293
collisionEnergy[i] = NA_REAL;
294+
electronBeamEnergy[i] = NA_REAL;
291295
precursorScanNum[i] = NA_INTEGER;
292296
precursorMZ[i] = NA_REAL;
293297
precursorCharge[i] = NA_INTEGER;
@@ -302,7 +306,7 @@ Rcpp::DataFrame RcppPwiz::getScanHeaderInfo (Rcpp::IntegerVector whichScan) {
302306
}
303307
}
304308

305-
Rcpp::List header(31);
309+
Rcpp::List header(32);
306310
std::vector<std::string> names;
307311
size_t i = 0;
308312
names.push_back("seqNum");
@@ -367,6 +371,8 @@ Rcpp::DataFrame RcppPwiz::getScanHeaderInfo (Rcpp::IntegerVector whichScan) {
367371
header[i++] = Rcpp::wrap(scanWindowLowerLimit);
368372
names.push_back("scanWindowUpperLimit");
369373
header[i++] = Rcpp::wrap(scanWindowUpperLimit);
374+
names.push_back("electronBeamEnergy");
375+
header[i++] = Rcpp::wrap(electronBeamEnergy);
370376
header.attr("names") = names;
371377

372378
return header;

0 commit comments

Comments
 (0)