detectPeak: Detection and quantification of peaks for a ptrSet object.

detectPeakR Documentation

Detection and quantification of peaks for a ptrSet object.

Description

The detectPeak function detects peaks on the average total spectrum around nominal masses, for all files present in ptrSet which have not already been processed. The temporal evolution of each peak is then evaluated by using a two-dimensional penalized spline regression. Finally, the expiration points (if defined in the ptrSet) are averaged, and a t-test is performed between expiration and ambient air. The peakList can be accessed with the getPeakList function, which returns the information about the detected peaks in each file as a list of ExpressionSet objects.The peak detection steps within each file are as follows:
for each nominal mass:

  • correction of the calibration drift

  • peak detection on the average spectrum

  • estimation of temporal evolution

  • t-test between expiration and ambient air

Usage

detectPeak(
  x,
  ppm = 130,
  minIntensity = 10,
  minIntensityRate = 0.01,
  mzNominal = NULL,
  resolutionRange = NULL,
  fctFit = NULL,
  smoothPenalty = NULL,
  parallelize = FALSE,
  nbCores = 2,
  saving = TRUE,
  saveDir = getParameters(x)$saveDir,
  ...
)

## S4 method for signature 'ptrRaw'
detectPeak(
  x,
  ppm = 130,
  minIntensity = 10,
  minIntensityRate = 0.01,
  mzNominal = NULL,
  resolutionRange = NULL,
  fctFit = NULL,
  smoothPenalty = NULL,
  timeLimit,
  knots = NULL,
  mzPrimaryIon = 21.022,
  ...
)

## S4 method for signature 'ptrSet'
detectPeak(
  x,
  ppm = 130,
  minIntensity = 10,
  minIntensityRate = 0.01,
  mzNominal = NULL,
  resolutionRange = NULL,
  fctFit = NULL,
  smoothPenalty = 0,
  parallelize = FALSE,
  nbCores = 2,
  saving = TRUE,
  saveDir = getParameters(x)$saveDir,
  ...
)

Arguments

x

a ptrSet object

ppm

minimum distance in ppm between two peaks

minIntensity

minimum intensity for peak detection. The final threshold for peak detection will be: max(minIntensity, threshold noise ). The threshold noise corresponds to max(max(noise around the nominal mass), minIntensityRate * max(intensity within the nominal mass)

minIntensityRate

Fraction of the maximum intensity to be used for noise thresholding

mzNominal

nominal masses at which peaks will be detected; if NULL, all nominal masses of the mass axis

resolutionRange

vector with the minimum, average, and maximum resolution of the PTR instrument. If NULL, the values are estimated by using the calibration peaks.

fctFit

function for the peak quantification: should be sech2 or averagePeak. If NULL, the best function is selected by using the calibration peaks

smoothPenalty

second order penalty coefficient of the p-spline used for two-dimensional regression. If NULL, the coefficient is estimated by generalized cross validation (GCV) criteria

parallelize

Boolean. If TRUE, loops over files are parallelized

nbCores

number of cluster to use for parallel computation.

saving

boolean. If TRUE, the object will be saved in saveDir with the setName parameter of the createPtrSet function

saveDir

The directory where the ptrSet object will be saved as .RData. If NULL, nothing will be saved

...

may be used to pass parameters to the processFileTemporal function

timeLimit

index time of the expiration limits and background. Should be provided by timeLimits function

knots

numeric vector corresponding to the knot values, which used for the two dimensional regression for each file. Should be provided by defineKnots function

mzPrimaryIon

the exact mass of the primary ion isotope

Value

an S4 ptrSet object, that contains the input ptrSet completed with the peakLists.

References

Muller et al 2014, Holzinger et al 2015, Marx and Eilers 1992

Examples


## For ptrRaw object
library(ptairData)
filePath <- system.file('extdata/exhaledAir/ind1', 'ind1-1.h5', 
package = 'ptairData')
raw <- readRaw(filePath,mzCalibRef=c(21.022,59.049),calib=TRUE)
timeLimit<-timeLimits(raw,fracMaxTIC=0.7)
knots<-defineKnots(object = raw,timeLimit=timeLimit)
raw <- detectPeak(raw, timeLimit=timeLimit, mzNominal = c(21,59),
smoothPenalty=0,knots=knots,resolutionRange=c(2000,5000,8000))

## For a ptrSet object
library(ptairData)
directory <- system.file("extdata/exhaledAir",  package = "ptairData")
exhaledPtrset<-createPtrSet(dir=directory,setName="exhaledPtrset",
mzCalibRef=c(21.022,59.049),
fracMaxTIC=0.9,saveDir= NULL)
exhaledPtrset  <- detectPeak(exhaledPtrset)
peakListEset<-getPeakList(exhaledPtrset)
Biobase::fData(peakListEset[[1]])
Biobase::exprs(peakListEset[[1]])

camilleroquencourt/ptairMS documentation built on April 24, 2024, 9:03 p.m.