ChromBackendSpectra: Chromatographic Data Backend for Spectra Objects

View source: R/ChromBackendSpectra.R

ChromBackendSpectraR Documentation

Chromatographic Data Backend for Spectra Objects

Description

The ChromBackendSpectra class extends ChromBackendMemory, inheriting all its slots and methods while providing additional functionality for summarizing chromatographic data from Spectra::Spectra() objects.

It can be initialized with a Spectra object, which is stored in the spectra slot of the backend. Users can also provide a data.frame containing chromatographic metadata, stored in chromData. This metadata filters the Spectra object and generates peaksData. If chromData is not provided, a default data.frame is created from the Spectra data. An "rtmin", "rtmax", "mzmin", and "mzmax" column will be created by condensing the Spectra data corresponding to each unique combination of the factorize.by variables.

The dataOrigin core chromatogram variable should reflect the dataOrigin of the Spectra object. The factorize.by parameter defines the variables for grouping Spectra data into chromatographic data. The default is c("msLevel", "dataOrigin"), which will define separate chromatograms for each combination of msLevel and dataOrigin. These variables must be in both Spectra and chromData (if provided).

The summarize.method parameter defines how spectral data intensity is summarized:

  • "sum": Sums intensity to create a Total Ion Chromatogram (TIC).

  • "max": Takes max intensity for a Base Peak Chromatogram (BPC).

If chromData or its factorization columns are modified, the factorize() method must be called to update chromSpectraIndex.

Usage

ChromBackendSpectra()

## S4 method for signature 'ChromBackendSpectra'
backendInitialize(
  object,
  spectra = Spectra::Spectra(),
  factorize.by = c("msLevel", "dataOrigin"),
  summarize.method = c("sum", "max"),
  chromData = fillCoreChromVariables(),
  ...
)

chromSpectraIndex(object)

Arguments

object

A ChromBackendSpectra object.

spectra

A Spectra object.

factorize.by

A character vector of variables for grouping Spectra data into chromatographic data. Default: c("msLevel", "dataOrigin"). If chromData is provided, it must contain these columns.

summarize.method

A character string specifying intensity summary: "sum" (default) or "max".

chromData

A data.frame with chromatographic data for use in backendInitialize(). If missing, a default is generated. Columns like rtmin, rtmax, mzmin, and mzmax must be provided and not contain NA values. Use -Inf/Inf for unspecified values. The "dataOrigin" column must match the Spectra object's "dataOrigin".

...

Additional parameters.

Details

No peaksData is stored until the user calls a function that generates it (e.g., rtime(), peaksData(), intensity()). The peaksData slot replacement is unsupported — modifications are temporary to optimize memory. The inMemory slot indicates this with TRUE.

ChromBackendSpectra should reuse ChromBackendMemory methods whenever possible to keep implementations simple.

Note

ensure that it returns a factor

Author(s)

Philippine Louail, Johannes Rainer.

Examples

library(Spectra)
library(MsBackendMetaboLights)

## Get Spectra data from MetaboLights
be <- backendInitialize(MsBackendMetaboLights(),
    mtblsId = "MTBLS39",
    filePattern = c("63B.cdf")
)
s <- Spectra(be)

## Initialize ChromBackendSpectra
be_empty <- new("ChromBackendSpectra")
be <- backendInitialize(be_empty, s)

## replace the msLevel data
msLevel(be) <- c(1L, 2L, 3L)

## re-factorize the data
be <- factorize(be)

## Create BPC : we summarize the intensity present in the Spectra object
## by the maximum value, thus creating a Base Peak Chromatogram.
be <- backendInitialize(be_empty, s, summarize.method = "max")

## Can now see the details of this bpc by looking at the chromData of our
## object
chromData(be)

## Another possibilities is to create eics from the Spectra object.
## Here we create an EIC with a specific m/z and retention time window.
df <- data.frame(mzmin = 100.01, mzmax = 100.02 , rtmin = 50, rtmax = 100)
be <- backendInitialize(be_empty, s, summarize.method = "sum")
chromData(be) <- cbind(chromData(be), df)

## now when we call the peaksData function, we will get the intensity
## of the spectra object that are in the m/z and retention time window
## defined in the chromData.
peaksData(be)


rformassspectrometry/Chromatograms documentation built on April 12, 2025, 6:53 p.m.