Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Read PASEF DDA MS2 precursor information #20

Open
wants to merge 6 commits into
base: main
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from 4 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
8 changes: 6 additions & 2 deletions DESCRIPTION
Original file line number Diff line number Diff line change
@@ -1,6 +1,6 @@
Package: MsBackendTimsTof
Title: Mass spectrometry Data Backend for Bruker TimsTOF Files
Version: 0.1.2
Version: 0.2.0
Authors@R:
c(person(given = "Johannes", family = "Rainer",
email = "[email protected]",
Expand All @@ -10,6 +10,10 @@ Authors@R:
email = "[email protected]",
role = "aut",
comment = c(ORCID = "0000-0001-9438-6909")),
person(given = "Roger", family = "Gine Bertomeu",
email = "[email protected]",
comment = c(ORCID = "0000-0003-0288-9619"),
role = "aut"),
jorainer marked this conversation as resolved.
Show resolved Hide resolved
person(given = "Steffen", family = "Neumann",
email = "[email protected]",
role = "ctb",
Expand Down Expand Up @@ -48,4 +52,4 @@ BugReports: https://github.com/RforMassSpectrometry/MsBackendTimsTof/issues
URL: https://github.com/RforMassSpectrometry/MsBackendTimsTof
biocViews: Infrastructure, MassSpectrometry, Metabolomics, DataImport
Roxygen: list(markdown=TRUE)
RoxygenNote: 7.1.2
RoxygenNote: 7.2.2
10 changes: 10 additions & 0 deletions NAMESPACE
Original file line number Diff line number Diff line change
Expand Up @@ -4,8 +4,10 @@ export(MsBackendTimsTof)
exportClasses(MsBackendTimsTof)
exportMethods(length)
importFrom(BiocParallel,bplapply)
importFrom(BiocParallel,bpmapply)
importFrom(BiocParallel,bpparam)
importFrom(IRanges,NumericList)
importFrom(MsCoreUtils,between)
importFrom(MsCoreUtils,i2index)
importFrom(MsCoreUtils,rbindFill)
importFrom(S4Vectors,DataFrame)
Expand All @@ -19,11 +21,19 @@ importFrom(opentimsr,OpenTIMS)
importFrom(opentimsr,opentims_set_threads)
importFrom(opentimsr,query)
importFrom(opentimsr,setup_bruker_so)
importFrom(opentimsr,table2df)
importFrom(stats,setNames)
importFrom(utils,capture.output)
importMethodsFrom(Spectra,backendInitialize)
importMethodsFrom(Spectra,collisionEnergy)
importMethodsFrom(Spectra,dataStorage)
importMethodsFrom(Spectra,isolationWindowLowerMz)
importMethodsFrom(Spectra,isolationWindowTargetMz)
importMethodsFrom(Spectra,isolationWindowUpperMz)
importMethodsFrom(Spectra,msLevel)
importMethodsFrom(Spectra,peaksVariables)
importMethodsFrom(Spectra,precursorCharge)
importMethodsFrom(Spectra,precursorIntensity)
importMethodsFrom(Spectra,precursorMz)
importMethodsFrom(Spectra,spectraData)
importMethodsFrom(Spectra,spectraVariables)
68 changes: 68 additions & 0 deletions R/MsBackendTimsTof-functions.R
Original file line number Diff line number Diff line change
Expand Up @@ -288,6 +288,11 @@ MsBackendTimsTof <- function() {
# can we assume that tms@all_coulmns is the same for all the TimsTOF?
.TIMSTOF_COLUMNS <- c("mz", "intensity", "tof", "inv_ion_mobility")

.TIMSTOF_MS2_COLUMNS <- c("precursorMz", "precursorCharge",
"precursorIntensity", "collisionEnergy",
"isolationWindowLowerMz", "isolationWindowTargetMz",
"isolationWindowUpperMz")

#' @importFrom methods as
#'
#' @importFrom S4Vectors DataFrame
Expand Down Expand Up @@ -349,6 +354,19 @@ MsBackendTimsTof <- function() {
res[["dataStorage"]] <- dataStorage(x)
core_cols <- core_cols[core_cols != "dataStorage"]
}
if (any(core_cols %in% .TIMSTOF_MS2_COLUMNS)) {
res_ms2_data <- .calculate_core_ms2_information(x)
for (col in core_cols[core_cols %in% .TIMSTOF_MS2_COLUMNS]) {
if (col == "precursorCharge"){
res[[col]] <- as.integer(
res_ms2_data[, which(.TIMSTOF_MS2_COLUMNS == col)]
)
} else {
res[[col]] <- res_ms2_data[, which(.TIMSTOF_MS2_COLUMNS == col)]
}
core_cols <- core_cols[core_cols != col]
}
}
if (length(core_cols)) {
res[core_cols] <- lapply(coreSpectraVariables()[core_cols],
function(z, n) rep(as(NA, z), n), length(x))
Expand All @@ -357,3 +375,53 @@ MsBackendTimsTof <- function() {
as(res, "DataFrame")
else DataFrame()
}

#' @importFrom BiocParallel bpmapply bpparam
.calculate_core_ms2_information <- function(x){
if (is.null(names(x@fileNames))){
return(matrix(numeric(), nrow = 0, ncol = length(.TIMSTOF_MS2_COLUMNS)))
}
tbls <- bpmapply(FUN = .do_calculate_core_ms2_information,
x = names(x@fileNames),
indices = lapply(split(x@indices[,1:2], f = x@indices[,3]),
matrix, ncol = 2),
SIMPLIFY = FALSE,
BPPARAM = bpparam())
output <- do.call(rbind, tbls)
output[order(order(x@indices[,3])), ]
}

#' @importFrom opentimsr OpenTIMS CloseTIMS table2df
#' @importFrom MsCoreUtils between
.do_calculate_core_ms2_information <- function(x, indices){
tms <- opentimsr::OpenTIMS(x)
on.exit(opentimsr::CloseTIMS(tms))
if (missing(indices)) {
indices <- unique(.query_tims(tms, frames = tms@frames$Id,
columns = c("frame", "scan")))
}
output <- matrix(NA_real_,
nrow = nrow(indices),
ncol = length(.TIMSTOF_MS2_COLUMNS))
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

maybe also add the column names with

ncol = length(.TIMSTOF_MS2_COLUMNS),
dimnames = list(NULL, .TIMSTOF_MS2_COLUMNS))

tbl <- opentimsr::table2df(tms, c("PASEFFrameMsMsInfo", "Precursors"))
jorainer marked this conversation as resolved.
Show resolved Hide resolved
for (i in seq(nrow(tbl[[1]]))) {
row <- tbl[[1]][i, ]
prec <- tbl[[2]][row$Precursor,]
target_rows <- which(indices[,1] == row$Frame &
MsCoreUtils::between(indices[,2],
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

The MsCoreUtils:: should not be needed here, since you are importing between from MsCoreUtils.

c(row$ScanNumBegin,
row$ScanNumEnd)))
if (length(target_rows)) {
output[target_rows, ] <-
rep(c(prec$MonoisotopicMz, ## precursorMz
prec$Charge, ## precursorCharge
prec$Intensity, ## precursorIntensity
row$CollisionEnergy, ## collisionEnergy
row$IsolationMz - row$IsolationWidth, ## isolationWindowLowerMz
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Hm, I'm wondering how the isolation width is defined - just to be sure that it should not be isolationMz - isolationWidth/2 ... can you please check?

row$IsolationMz, ## isolationWindowTargetMz
row$IsolationMz + row$IsolationWidth), ## isolationWindowUpperMz
each = length(target_rows))
}
}
output
}
56 changes: 56 additions & 0 deletions R/MsBackendTimsTof.R
Original file line number Diff line number Diff line change
Expand Up @@ -288,6 +288,62 @@ setMethod("msLevel", "MsBackendTimsTof", function(object, ...) {
.get_msLevel(object)
})

#' @importMethodsFrom Spectra precursorMz
#'
#' @rdname MsBackendTimsTof
setMethod("precursorMz", "MsBackendTimsTof", function(object, ...) {
jorainer marked this conversation as resolved.
Show resolved Hide resolved
tbl <- .calculate_core_ms2_information(object)
tbl[, .TIMSTOF_MS2_COLUMNS == "precursorMz"]
})

#' @importMethodsFrom Spectra precursorCharge
#'
#' @rdname MsBackendTimsTof
setMethod("precursorCharge", "MsBackendTimsTof", function(object, ...) {
tbl <- .calculate_core_ms2_information(object)
tbl[, .TIMSTOF_MS2_COLUMNS == "precursorCharge"]
jorainer marked this conversation as resolved.
Show resolved Hide resolved
})

#' @importMethodsFrom Spectra precursorIntensity
#'
#' @rdname MsBackendTimsTof
setMethod("precursorIntensity", "MsBackendTimsTof", function(object, ...) {
tbl <- .calculate_core_ms2_information(object)
tbl[, .TIMSTOF_MS2_COLUMNS == "precursorIntensity"]
})

#' @importMethodsFrom Spectra collisionEnergy
#'
#' @rdname MsBackendTimsTof
setMethod("collisionEnergy", "MsBackendTimsTof", function(object, ...) {
tbl <- .calculate_core_ms2_information(object)
tbl[, .TIMSTOF_MS2_COLUMNS == "collisionEnergy"]
})

#' @importMethodsFrom Spectra isolationWindowLowerMz
#'
#' @rdname MsBackendTimsTof
setMethod("isolationWindowLowerMz", "MsBackendTimsTof", function(object, ...) {
tbl <- .calculate_core_ms2_information(object)
tbl[, .TIMSTOF_MS2_COLUMNS == "isolationWindowLowerMz"]
})

#' @importMethodsFrom Spectra isolationWindowTargetMz
#'
#' @rdname MsBackendTimsTof
setMethod("isolationWindowTargetMz", "MsBackendTimsTof", function(object, ...) {
tbl <- .calculate_core_ms2_information(object)
tbl[, .TIMSTOF_MS2_COLUMNS == "isolationWindowTargetMz"]
})

#' @importMethodsFrom Spectra isolationWindowUpperMz
#'
#' @rdname MsBackendTimsTof
setMethod("isolationWindowUpperMz", "MsBackendTimsTof", function(object, ...) {
tbl <- .calculate_core_ms2_information(object)
tbl[, .TIMSTOF_MS2_COLUMNS == "isolationWindowUpperMz"]
})

#' @rdname MsBackendTimsTof
setMethod("$", "MsBackendTimsTof", function(x, name) {
if (!any(spectraVariables(x) == name))
Expand Down
21 changes: 21 additions & 0 deletions man/MsBackendTimsTof.Rd

Some generated files are not rendered by default. Learn more about how customized files appear on GitHub.

15 changes: 15 additions & 0 deletions tests/testthat/test_MsBackendTimsTof-functions.R
Original file line number Diff line number Diff line change
Expand Up @@ -145,3 +145,18 @@ test_that(".inv_ion_mobility works", {
res_2 <- .inv_ion_mobility(be_2)
expect_equal(res_2, res[idx])
})

test_that("calculate_core_ms2_information works", {
tbl <- .calculate_core_ms2_information(be)
expect_true(is.numeric(tbl))
expect_equal(ncol(tbl), length(.TIMSTOF_MS2_COLUMNS))
expect_equal(nrow(tbl), nrow(be@indices))
expect_true(!all(is.na(tbl)))

#Check that it also works with subsets of scans
be_subset <- be[c(120:130, 17000:17010)] #From both file 1 and 2
tbl <- .calculate_core_ms2_information(be_subset)
expect_equal(ncol(tbl), length(.TIMSTOF_MS2_COLUMNS))
expect_equal(nrow(tbl), nrow(be_subset@indices))
})

42 changes: 42 additions & 0 deletions tests/testthat/test_MsBackendTimsTof.R
Original file line number Diff line number Diff line change
Expand Up @@ -184,6 +184,48 @@ test_that("msLevel,MsBackendTimsTof works", {
expect_identical(res, match(.get_frame_columns(be, "MsMsType"), c(0L, 8L)))
})

test_that("precursorMz,MsBackendTimsTof works", {
expect_identical(precursorMz(MsBackendTimsTof()), numeric(0))
res <- precursorMz(be)
expect_identical(length(res), length(be))
})

test_that("precursorCharge,MsBackendTimsTof works", {
expect_identical(precursorCharge(MsBackendTimsTof()), numeric(0))
res <- precursorCharge(be)
expect_identical(length(res), length(be))
})

test_that("precursorIntensity,MsBackendTimsTof works", {
expect_identical(precursorIntensity(MsBackendTimsTof()), numeric(0))
res <- precursorIntensity(be)
expect_identical(length(res), length(be))
})

test_that("collisionEnergy,MsBackendTimsTof works", {
expect_identical(collisionEnergy(MsBackendTimsTof()), numeric(0))
res <- collisionEnergy(be)
expect_identical(length(res), length(be))
})

test_that("isolationWindowLowerMz,MsBackendTimsTof works", {
expect_identical(isolationWindowLowerMz(MsBackendTimsTof()), numeric(0))
res <- isolationWindowLowerMz(be)
expect_identical(length(res), length(be))
})

test_that("isolationWindowTargetMz,MsBackendTimsTof works", {
expect_identical(isolationWindowTargetMz(MsBackendTimsTof()), numeric(0))
res <- isolationWindowTargetMz(be)
expect_identical(length(res), length(be))
})

test_that("isolationWindowUpperMz,MsBackendTimsTof works", {
expect_identical(isolationWindowUpperMz(MsBackendTimsTof()), numeric(0))
res <- isolationWindowUpperMz(be)
expect_identical(length(res), length(be))
})

test_that("$,MsBackendTimsTof works", {
res_all <- spectraData(be)

Expand Down