Skip to content

Commit

Permalink
feat: add filterMsLevel method for MsExperiment and XcmsExperiment
Browse files Browse the repository at this point in the history
- Add the `filterMsLevel` methods for `MsExperiment` and `XcmsExperiment`.
  • Loading branch information
jorainer committed Nov 27, 2023
1 parent d6b0198 commit b397afe
Show file tree
Hide file tree
Showing 7 changed files with 70 additions and 2 deletions.
8 changes: 8 additions & 0 deletions R/MsExperiment.R
Original file line number Diff line number Diff line change
Expand Up @@ -20,6 +20,14 @@ setMethod("filterMz", "MsExperiment",
filterMzRange(object, mz, msLevel.)
})

#' @rdname XcmsExperiment
setMethod("filterMsLevel", "MsExperiment",
function(object, msLevel. = uniqueMsLevels(object)) {
message("Filter spectra")
.mse_filter_spectra(object, filterMsLevel,
msLevel. = msLevel.)
})

#' @rdname XcmsExperiment
setMethod("uniqueMsLevels", "MsExperiment", function(object) {
uniqueMsLevels(spectra(object))
Expand Down
19 changes: 19 additions & 0 deletions R/XcmsExperiment.R
Original file line number Diff line number Diff line change
Expand Up @@ -83,6 +83,9 @@
#' `"isolationWindowUpperMz"` (columns in `chromPeakData`) do not contain
#' the user provided `mz`.
#'
#' - `filterMsLevel`: filter the data of the `XcmsExperiment` or `MsExperiment`
#' to keep only data of the MS level(s) specified with parameter `msLevel.`.
#'
#' - `filterMz`, `filterMzRange`: filter the spectra within an
#' `XcmsExperiment` or `MsExperiment` to the specified m/z range (parameter
#' `mz`). For `XcmsExperiment` also identified chromatographic peaks and
Expand Down Expand Up @@ -482,6 +485,8 @@
#' length equal to the numer of rows of the parameters `mz` and `rt`
#' defining the m/z and rt regions from which the chromatograms should
#' be created. Defaults to `msLevel = 1L`.
#' for `filterMsLevel`: `integer` defining the MS level(s) to which the
#' data should be subset.
#'
#' @param mz For `chromPeaks` and `featureDefinitions`: `numeric(2)` optionally
#' defining the m/z range for which chromatographic peaks or feature
Expand Down Expand Up @@ -816,6 +821,20 @@ setMethod(
callNextMethod(object = object, mz = mz, msLevel. = msLevel.)
})

#' @rdname XcmsExperiment
setMethod(
"filterMsLevel", "XcmsExperiment",
function(object, msLevel. = uniqueMsLevels(object)) {
if (!length(msLevel.))
return(object)
if (hasChromPeaks(object)) {
keep <- chromPeakData(object)$ms_level %in% msLevel.
object <- .filter_chrom_peaks(object, idx = base::which(keep))
}
callNextMethod(object = object, msLevel. = msLevel.)
})


################################################################################
## chromatographic peaks
################################################################################
Expand Down
2 changes: 2 additions & 0 deletions R/methods-XCMSnExp.R
Original file line number Diff line number Diff line change
Expand Up @@ -1460,6 +1460,8 @@ setMethod("smooth", "XCMSnExp", function(x, method = c("SavitzkyGolay",
setAs(from = "XCMSnExp", to = "xcmsSet", def = .XCMSnExp2xcmsSet)

#' @rdname XcmsExperiment
#'
#' @name XcmsExperiment
setAs(from = "XcmsExperiment", to = "xcmsSet", def = .XCMSnExp2xcmsSet)

#' @rdname XCMSnExp-peak-grouping-results
Expand Down
7 changes: 7 additions & 0 deletions inst/NEWS
Original file line number Diff line number Diff line change
@@ -1,3 +1,10 @@
Changes in version 4.1.x
----------------------

- Add `filterMsLevel` method for `MsExperiment` and `XcmsExperiment`.
- Ensure chunk-wise processing of Spectra (introduced with version 1.13.2) is
disabled when xcms is using its own chunk-wise processing.

Changes in version 4.1.1
----------------------

Expand Down
14 changes: 12 additions & 2 deletions man/XcmsExperiment.Rd

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

12 changes: 12 additions & 0 deletions tests/testthat/test_MsExperiment.R
Original file line number Diff line number Diff line change
Expand Up @@ -16,6 +16,18 @@ test_that("filterRt,MsExperiment works", {
expect_equal(rtime(spectra(res)), rtime(spectra(mse)))
})

test_that("filterMsLevel,MsExperiment works", {
res <- filterMsLevel(mse, msLevel = 1)
expect_equal(rtime(res), rtime(mse))

res <- filterMsLevel(mse, msLevel = integer())
expect_equal(rtime(res), rtime(mse))

res <- filterMsLevel(mse, msLevel = 2L)
expect_equal(rtime(res), numeric())
expect_equal(length(res), 3L)
})

test_that("filterFile,MsExperiment works", {
res <- filterFile(mse)
expect_s4_class(res, "MsExperiment")
Expand Down
10 changes: 10 additions & 0 deletions tests/testthat/test_XcmsExperiment.R
Original file line number Diff line number Diff line change
Expand Up @@ -194,6 +194,16 @@ test_that("filterRt,XcmsExperiment works", {
seq_len(nrow(chromPeaks(res))))
})

test_that("filterMsLevel,XcmsExperiment works", {
res <- filterMsLevel(xmse, c(1L, 2L))
expect_equal(rtime(res), rtime(xmse))
expect_equal(chromPeaks(res), chromPeaks(xmse))

res <- filterMsLevel(xmse, msLevel = 2L)
expect_equal(rtime(res), numeric())
expect_true(nrow(chromPeaks(res)) == 0L)
})

test_that("filterFile,XcmsExperiment works", {
res <- filterFile(xmse)
expect_s4_class(res, "XcmsExperiment")
Expand Down

0 comments on commit b397afe

Please sign in to comment.