Title: | Direct Access to Orbitrap Data and Beyond |
---|---|
Description: | This package wraps the functionality of the Thermo Fisher Scientic RawFileReader .NET 8.0 assembly. Within the R environment, spectra and chromatograms are represented by S3 objects. The package provides basic functions to download and install the required third-party libraries. The package is developed, tested, and used at the Functional Genomics Center Zurich, Switzerland. |
Authors: | Christian Panse [aut, cre] , Leonardo Schwarz [ctb] , Tobias Kockmann [aut] |
Maintainer: | Christian Panse <[email protected]> |
License: | GPL-3 |
Version: | 1.15.8 |
Built: | 2024-11-13 03:23:51 UTC |
Source: | https://github.com/bioc/rawrr |
Benchmark spectra per second
.benchmark(rawfile)
.benchmark(rawfile)
rawfile |
the name of the raw file containing the mass spectrometry data from the Thermo Fisher Scientific instrument. |
Christian Panse 2024-11-05
sampleFilePath() |> rawrr:::.benchmark() -> S plot (count / runTimeInSec ~ count, log='xy', data=S, sub = paste0("Overall runtime took ", round(sum(S$runTimeInSec), 3), " seconds."), xlab = 'number of random generated scan ids', main = "benchmark spectra per second")
sampleFilePath() |> rawrr:::.benchmark() -> S plot (count / runTimeInSec ~ count, log='xy', data=S, sub = paste0("Overall runtime took ", round(sum(S$runTimeInSec), 3), " seconds."), xlab = 'number of random generated scan ids', main = "benchmark spectra per second")
Download and install the Thermo Fisher Scientific .NET 8.0 nupkgs
.downloadNupkgs(sourceUrl = .thermofisherlsmsUrl(), force = TRUE)
.downloadNupkgs(sourceUrl = .thermofisherlsmsUrl(), force = TRUE)
sourceUrl |
url of nupkgs. |
force |
if |
An (invisible) vector of integer code, 0 for success and non-zero for failure. For the "wget" and "curl" methods this is the status code returned by the external program.
Christian Panse <[email protected]>, 2021, 2024
URL for Thermo Fisher .NET assemblies
.thermofisherlsmsUrl()
.thermofisherlsmsUrl()
an URL
deriving area under the curve (AUC)
auc.rawrrChromatogram(x)
auc.rawrrChromatogram(x)
x |
an rawrrChromatogram object contains |
A numeric value.
Base peak of a spectrum
basePeak(x)
basePeak(x)
x |
A rawrrSpectrum object |
A double vector of length two. The first component is the base peak position (m/z). The second component is the base peak intensity.
S <- readSpectrum(rawfile = sampleFilePath(), 1) basePeak(S[[1]])
S <- readSpectrum(rawfile = sampleFilePath(), 1) basePeak(S[[1]])
rawrr.exe
console application.builds rawrr.exe
file from C# source code requiring
.NET SDK. The console application rawrr.exe
is used by the package's reader functions through a system2 call
or a textConnection.
To use this function, ensure that the local RawFileReader NuGet packages
are added to the NuGet source list. You can accomplish this by
downloading the necessary packages with
rawrr:::.downloadNupkgs()
and subsequently running
rawrr:::.addNupkgSource()
.
buildRawrrExe()
buildRawrrExe()
the return value of the system2 command.
Christian Panse <[email protected]>, 2021, 2024
https://www.mono-project.com/docs/advanced/assemblies-and-the-gac/, 2020
https://docs.microsoft.com/en-us/dotnet/csharp/language-reference/compiler-options/advanced
Retrieve dependent scan(s) of a scan listed in scan index
dependentScan(x, scanNumber)
dependentScan(x, scanNumber)
x |
A scan index returned by |
scanNumber |
The scan number that should be inspected for dependent scans. |
The scan number of the dependent scan(s).
Idx <- readIndex(rawfile = sampleFilePath()) dependentScan(Idx, scanNumber = 1)
Idx <- readIndex(rawfile = sampleFilePath()) dependentScan(Idx, scanNumber = 1)
Is FAIMS Voltage on?
faimsVoltageOn(x)
faimsVoltageOn(x)
x |
A rawrrSpectrum object |
A boolean
S <- readSpectrum(rawfile = sampleFilePath(), 1:10) try(faimsVoltageOn(S[[1]]))
S <- readSpectrum(rawfile = sampleFilePath(), 1:10) try(faimsVoltageOn(S[[1]]))
determine scan numbers which match a specified filter
filter(rawfile, filter = "ms", precision = 10, tmpdir = tempdir())
filter(rawfile, filter = "ms", precision = 10, tmpdir = tempdir())
rawfile |
the name of the raw file containing the mass spectrometry data from the Thermo Fisher Scientific instrument. |
filter |
scan filter string, e.g., |
precision |
mass precision, default is 10. |
tmpdir |
defines the directory used to store temporary data generated by
the .NET assembly |
a vecntor of integer values.
deprecated installRawFileReaderDLLs
installRawFileReaderDLLs()
installRawFileReaderDLLs()
rawrr
assemblydownloads and installs the rawrr.exe
.NET assembly in
the directory provided by rawrrAssemblyPath()
.
installRawrrExe( sourceUrl = "https://fgcz-ms.uzh.ch/~cpanse/rawrr/dotnet/", force = FALSE, ... )
installRawrrExe( sourceUrl = "https://fgcz-ms.uzh.ch/~cpanse/rawrr/dotnet/", force = FALSE, ... )
sourceUrl |
url of |
force |
if |
... |
other parameter for |
The console application rawrr
is used by the package's
reader functions through a system2 call.
An integer code, 0 for success and non-zero for failure. For the "wget" and "curl" methods this is the status code returned by the external program.
rawrrChromatogram
Function to check if an object is an instance of class rawrrChromatogram
is.rawrrChromatogram(x)
is.rawrrChromatogram(x)
x |
x any R object to be tested. |
TRUE
or FALSE
Tobias Kockmann, 2020.
rawfile <- sampleFilePath() C <- readChromatogram(rawfile, mass = 445.1181, tol = 10) is.rawrrChromatogram(C[[1]])
rawfile <- sampleFilePath() C <- readChromatogram(rawfile, mass = 445.1181, tol = 10) is.rawrrChromatogram(C[[1]])
rawrrSpectrum
Function to check if an object is an instance of class rawrrSpectrum
is.rawrrSpectrum(x)
is.rawrrSpectrum(x)
x |
any R object to be tested. |
TRUE
or FALSE
S <- rawrr::sampleFilePath() |> rawrr::readSpectrum(scan = 1:10) rawrr::is.rawrrSpectrum(S[[1]])
S <- rawrr::sampleFilePath() |> rawrr::readSpectrum(scan = 1:10) rawrr::is.rawrrSpectrum(S[[1]])
rawrrSpectrumSet
Function to check if an object is an instance of class rawrrSpectrumSet
is.rawrrSpectrumSet(x)
is.rawrrSpectrumSet(x)
x |
any R object to be tested. |
TRUE
or FALSE
rawrr::sampleFilePath() |> rawrr::readSpectrum(scan = 1:10) |> rawrr::is.rawrrSpectrumSet()
rawrr::sampleFilePath() |> rawrr::readSpectrum(scan = 1:10) |> rawrr::is.rawrrSpectrumSet()
Make accessor function for key value pair returned by RawFileReader
makeAccessor(key, returnType = "integer")
makeAccessor(key, returnType = "integer")
key |
An object name found in instance of class |
returnType |
The type used for casting of values |
This function factory creates accessor functions for class
rawrrSpectrum
.
An accessor function
Tobias Kockmann, 2020
S <- readSpectrum(rawfile = sampleFilePath(), 1:10) maxIonTime <- makeAccessor(key = "Max. Ion Time (ms):", returnType = "double") maxIonTime(S[[1]])
S <- readSpectrum(rawfile = sampleFilePath(), 1:10) maxIonTime <- makeAccessor(key = "Max. Ion Time (ms):", returnType = "double") maxIonTime(S[[1]])
Acquisition/scan range of spectrum
massRange(x)
massRange(x)
x |
A rawrrSpectrum object |
A double vector of length two. The first component is the start m/z, the second is the stop m/z value used by the detector during data acquisition. Also referred to as scan range.
S <- readSpectrum(rawfile = sampleFilePath(), 1) massRange(S[[1]])
S <- readSpectrum(rawfile = sampleFilePath(), 1) massRange(S[[1]])
Retrieve master scan of scan listed in scan index
masterScan(x, scanNumber)
masterScan(x, scanNumber)
x |
A scan index returned by |
scanNumber |
The scan number that should be inspected for the presence of a master scan. |
Returns the scan number of the master scan or NA if no master scan exists.
rawrr::sampleFilePath() |> rawrr::readIndex() |> rawrr::masterScan(scanNumber = 1)
rawrr::sampleFilePath() |> rawrr::readIndex() |> rawrr::masterScan(scanNumber = 1)
rawrrSpectrum
Developer function.
new_rawrrSpectrum( scan = numeric(), massRange = numeric(), scanType = character(), StartTime = numeric(), centroidStream = logical(), mZ = numeric(), intensity = numeric() )
new_rawrrSpectrum( scan = numeric(), massRange = numeric(), scanType = character(), StartTime = numeric(), centroidStream = logical(), mZ = numeric(), intensity = numeric() )
scan |
scan number |
massRange |
Mass range covered by spectrum |
scanType |
Character string describing the scan type. |
StartTime |
Retention time in minutes |
centroidStream |
Logical indicating if centroided data is available |
mZ |
m/z values |
intensity |
Intensity values |
Object of class rawrrSpectrum
Tobias Kockmann, 2020.
rawrrChromatogram
objectsPlot rawrrChromatogram
objects
## S3 method for class 'rawrrChromatogram' plot(x, legend = TRUE, ...)
## S3 method for class 'rawrrChromatogram' plot(x, legend = TRUE, ...)
x |
A |
legend |
Should legend be printed? |
... |
Passes additional arguments. |
This function creates a plot.
Tobias Kockmann, 2020.
rawfile <- sampleFilePath() C <- readChromatogram(rawfile, mass = 445.1181, tol = 10) plot(C[[1]])
rawfile <- sampleFilePath() C <- readChromatogram(rawfile, mass = 445.1181, tol = 10) plot(C[[1]])
rawrrChromatogramSet
objectsPlot rawrrChromatogramSet
objects
## S3 method for class 'rawrrChromatogramSet' plot(x, diagnostic = FALSE, ...)
## S3 method for class 'rawrrChromatogramSet' plot(x, diagnostic = FALSE, ...)
x |
A |
diagnostic |
Show diagnostic legend? |
... |
Passes additional arguments. |
This function creates a plot.
Tobias Kockmann, 2020.
rawrrSpectrum
Plot method for objects of class rawrrSpectrum
.
## S3 method for class 'rawrrSpectrum' plot( x, relative = TRUE, centroid = FALSE, SN = FALSE, legend = TRUE, diagnostic = FALSE, ... )
## S3 method for class 'rawrrSpectrum' plot( x, relative = TRUE, centroid = FALSE, SN = FALSE, legend = TRUE, diagnostic = FALSE, ... )
x |
an object of class |
relative |
If set to |
centroid |
Should centroided data be used for plotting? |
SN |
Should Signal/Noise be used for plotting? |
legend |
Should legend be printed? |
diagnostic |
Should this option be applied? The default is |
... |
function passes arbitrary additional arguments. |
plot.rawrrSpectrum
is a low level function that calls
base::plot
for plotting rawrrSpectrum
objects. It passes all
additional arguments to plot()
Is usually called by method dispatch.
This function creates a plot.
Tobias Kockmann, 2020.
Print method imitate the look and feel of Thermo Fisher Scientific FreeStyle's output
## S3 method for class 'rawrrSpectrum' print(x, ...)
## S3 method for class 'rawrrSpectrum' print(x, ...)
x |
an |
... |
Arguments to be passed to methods. |
This function creates a print message.
Christian Panse and Tobias Kockmann, 2020.
Derives the path where all .NET assemblies are stored.
rawrrAssemblyPath()
rawrrAssemblyPath()
path
installRawrrExe
and buildRawrrExe
rawrrAssemblyPath()
rawrrAssemblyPath()
rawrrSpectrum
objectsHigh-level constructor for instances of class
rawrrSpectrum
, also named helper function. Currently, mainly to support
testing and for demonstration.
rawrrSpectrum(sim = "TESTPEPTIDE")
rawrrSpectrum(sim = "TESTPEPTIDE")
sim |
Either |
Function returns a validated rawrrSpectrum
object
Tobias Kockmann, 2020.
plot(rawrrSpectrum(sim = "TESTPEPTIDE")) rawrrSpectrum(sim = "example_1")
plot(rawrrSpectrum(sim = "TESTPEPTIDE")) rawrrSpectrum(sim = "example_1")
Extracts chromatographic data from a raw file.
readChromatogram( rawfile, mass = NULL, tol = 10, filter = "ms", type = "xic", tmpdir = tempdir() )
readChromatogram( rawfile, mass = NULL, tol = 10, filter = "ms", type = "xic", tmpdir = tempdir() )
rawfile |
the name of the raw file containing the mass spectrometry data from the Thermo Fisher Scientific instrument. |
mass |
a vector of mass values iff |
tol |
mass tolerance in ppm iff |
filter |
defines the scan filter, default is |
type |
|
tmpdir |
defines the directory used to store temporary data generated by
the .NET assembly |
Chromatograms come in different flavors but are always signal intensity values as a function of time. Signal intensities can be point estimates from scanning detectors or plain intensities from non-scanning detectors, e.g., UV trace. Scanning detector (mass analyzers) point estimates can be defined in different ways by, for instance, summing all signals of a given spectrum (total ion chromatogram or TIC), or by extracting signal around an expected value (extracted ion chromatogram = XIC), or by using the maximum signal contained in a spectrum (base peak chromatogram = BPC). On top, chromatograms can be computed from pre-filtered lists of scans. A total ion chromatogram (TIC), for instance, is typically generated by iterating over all MS1-level scans.
chromatogram object(s) containing of a vector of times
and a
corresponding vector of intensities
.
Christian Trachsel, Tobias Kockmann and Christian Panse <[email protected]> 2018, 2019, 2020.
Automated quality control sample 1 (autoQC01) analyzed across different Thermo Scientific mass spectrometers, MSV000086542.
The Thermo Fisher Scientific RawFileReader C# code snippets https://planetorbitrap.com/rawfilereader.
https://massive.ucsd.edu/ProteoSAFe/dataset.jsp?accession=MSV000086542
# Example 1: not meaningful but proof-of-concept (rawfile <- rawrr::sampleFilePath()) rawrr::readChromatogram(rawfile, mass=c(669.8381, 726.8357), tol=1000) |> plot() rawrr::readChromatogram(rawfile, type='bpc') |> plot() rawrr::readChromatogram(rawfile, type='tic') |> plot() # Example 2: extract iRT peptides if (require(ExperimentHub) & require(protViz)){ iRTpeptide <- c("LGGNEQVTR", "YILAGVENSK", "GTFIIDPGGVIR", "GTFIIDPAAVIR", "GAGSSEPVTGLDAK", "TPVISGGPYEYR", "VEATFGVDESNAK", "TPVITGAPYEYR", "DGLDAASYYAPVR", "ADVTPADFSEWSK", "LFLQFGAQGSPFLK") # fetch via ExperimentHub library(ExperimentHub) eh <- ExperimentHub::ExperimentHub() EH4547 <- normalizePath(eh[["EH4547"]]) (rawfile <- paste0(EH4547, ".raw")) if (!file.exists(rawfile)){ file.link(EH4547, rawfile) } op <- par(mfrow=c(2,1)) readChromatogram(rawfile, type='bpc') |> plot() readChromatogram(rawfile, type='tic') |> plot() par(op) # derive [2H+] ions ((protViz::parentIonMass(iRTpeptide) + 1.008) / 2) |> readChromatogram(rawfile=rawfile) |> plot() }
# Example 1: not meaningful but proof-of-concept (rawfile <- rawrr::sampleFilePath()) rawrr::readChromatogram(rawfile, mass=c(669.8381, 726.8357), tol=1000) |> plot() rawrr::readChromatogram(rawfile, type='bpc') |> plot() rawrr::readChromatogram(rawfile, type='tic') |> plot() # Example 2: extract iRT peptides if (require(ExperimentHub) & require(protViz)){ iRTpeptide <- c("LGGNEQVTR", "YILAGVENSK", "GTFIIDPGGVIR", "GTFIIDPAAVIR", "GAGSSEPVTGLDAK", "TPVISGGPYEYR", "VEATFGVDESNAK", "TPVITGAPYEYR", "DGLDAASYYAPVR", "ADVTPADFSEWSK", "LFLQFGAQGSPFLK") # fetch via ExperimentHub library(ExperimentHub) eh <- ExperimentHub::ExperimentHub() EH4547 <- normalizePath(eh[["EH4547"]]) (rawfile <- paste0(EH4547, ".raw")) if (!file.exists(rawfile)){ file.link(EH4547, rawfile) } op <- par(mfrow=c(2,1)) readChromatogram(rawfile, type='bpc') |> plot() readChromatogram(rawfile, type='tic') |> plot() par(op) # derive [2H+] ions ((protViz::parentIonMass(iRTpeptide) + 1.008) / 2) |> readChromatogram(rawfile=rawfile) |> plot() }
This function extracts the meta information from a given raw file.
readFileHeader(rawfile)
readFileHeader(rawfile)
rawfile |
the name of the raw file containing the mass spectrometry data from the Thermo Fisher Scientific instrument. |
A list object containing the following entries: RAW file version, Creation date, Operator, Number of instruments, Description, Instrument model, Instrument name, Serial number, Software version, Firmware version, Units, Mass resolution, Number of scans, Number of ms2 scans, Scan range, Time range, Mass range, Scan filter (first scan), Scan filter (last scan), Total number of filters, Sample name, Sample id, Sample type, Sample comment, Sample vial, Sample volume, Sample injection volume, Sample row number, Sample dilution factor, or Sample barcode.
Tobias Kockmann and Christian Panse 2018, 2019, 2020.
Thermo Fisher Scientific's NewRawfileReader C# code snippets https://planetorbitrap.com/rawfilereader.
rawrr::sampleFilePath() |> readFileHeader()
rawrr::sampleFilePath() |> readFileHeader()
Read scan index
readIndex(rawfile)
readIndex(rawfile)
rawfile |
the name of the raw file containing the mass spectrometry data from the Thermo Fisher Scientific instrument. |
returns a data.frame
with the column names
scan, scanType, StartTime, precursorMass, MSOrder, charge, masterScan, and
dependencyType of all spectra.
Tobias Kockmann and Christian Panse <[email protected]>, 2020, 2021
Idx <- rawrr::sampleFilePath() |> rawrr::readIndex() table(Idx$scanType) plot(Idx$StartTime, Idx$precursorMass, col=as.factor(Idx$charge), pch=16) table(Idx$MSOrder)
Idx <- rawrr::sampleFilePath() |> rawrr::readIndex() table(Idx$scanType) plot(Idx$StartTime, Idx$precursorMass, col=as.factor(Idx$charge), pch=16) table(Idx$MSOrder)
The function derives spectra of a given raw file and a given vector of scan numbers.
readSpectrum( rawfile, scan = NULL, tmpdir = tempdir(), validate = FALSE, mode = "" )
readSpectrum( rawfile, scan = NULL, tmpdir = tempdir(), validate = FALSE, mode = "" )
rawfile |
the name of the raw file containing the mass spectrometry data from the Thermo Fisher Scientific instrument. |
scan |
a vector of requested scan numbers. |
tmpdir |
defines the directory used to store temporary data generated by
the .NET assembly |
validate |
boolean default is |
mode |
if |
All mass spectra are recorded by scanning detectors (mass analyzers)
that log signal intensities for ranges of mass to charge ratios (m/z), also
referred to as position. These recordings can be of continuous nature,
so-called profile data (p), or appear centroided (c) in case discrete
information (tuples of position and intensity values) are sufficient.
This heavily compacted data structure is often called a peak list.
In addition to signal intensities, a peak list can also cover additional
peak attributes like peak resolution (R), charge (z), or local noise
estimates. In short, the additional attributes further described the nature
of the original profile signal or help to group peak lists with respect to
their molecular nature or processing history. A well-known example is the
assignment of peaks to peak groups that constitute isotope patterns
(M, M+1, M+2, ...).
The names of objects encapsulated within rawrrSpectrum
instances are
keys returned by the Thermo Fisher Scientific New RawFileReader API and the
corresponding values become data parts of the objects, typically vectors.
a nested list of rawrrSpectrum
objects containing more than 50
values of scan information, e.g., the charge state, two vectors containing
the mZ and its corresponding intensity values or the AGC information,
mass calibration, ion optics ...
Tobias Kockmann and Christian Panse <[email protected]> 2018, 2019, 2020, 2021
C# code snippets of the NewRawfileReader library https://planetorbitrap.com/rawfilereader.
Universal Spectrum Explorer: https://www.proteomicsdb.org/use/ doi:10.1021/acs.jproteome.1c00096
https://massive.ucsd.edu/ProteoSAFe/dataset.jsp?accession=MSV000086542
# Example 1 S <- rawrr::sampleFilePath() |> rawrr::readSpectrum(scan = 1:9) S[[1]] names(S[[1]]) plot(S[[1]]) # Example 2 - find best peptide spectrum match using the |> pipe operator # fetch via ExperimentHub if (require(ExperimentHub) & require(protViz)){ eh <- ExperimentHub::ExperimentHub() EH4547 <- normalizePath(eh[["EH4547"]]) (rawfile <- paste0(EH4547, ".raw")) if (!file.exists(rawfile)){ file.link(EH4547, rawfile) } GAG <- "GAGSSEPVTGLDAK" .bestPeptideSpectrumMatch <- function(rawfile, sequence="GAGSSEPVTGLDAK"){ readIndex(rawfile) |> subset(abs((1.008 + (protViz::parentIonMass(sequence) - 1.008) / 2) - precursorMass) < 0.001, select = scan) |> unlist() |> readSpectrum(rawfile = rawfile) |> lapply(function(x) { y <- protViz::psm(sequence = GAG, spec=x, plot=FALSE); y$scan <- x$scan; y }) |> lapply(FUN= function(x){ score <- sum(abs(x$mZ.Da.error) < 0.01); cbind(scan=x$scan, score=score) }) |> (function(x) as.data.frame(Reduce(rbind, x)))() |> subset(score > 0) |> (function(x) x[order(x$score, decreasing = TRUE), 'scan'])() |> head(1) } start_time <- Sys.time() bestMatch <- .bestPeptideSpectrumMatch(rawfile, GAG) |> rawrr::readSpectrum(rawfile=rawfile) |> lapply(function(x) protViz::peakplot(peptideSequence = GAG, x)) end_time <- Sys.time() end_time - start_time # Example 3 # using proteomicsdb \doi{10.1101/2020.09.08.287557} # through https://www.proteomicsdb.org/use/ .UniversalSpectrumExplorer <- function(x, sequence){ m <- protViz::psm( sequence, x) cat(paste(x$mZ[m$idx], "\t", x$intensity[m$idx]), sep = "\n") } rawrr::readSpectrum(rawfile=rawfile, 11091) |> lapply(function(x).UniversalSpectrumExplorer(x, sequence = GAG)) }
# Example 1 S <- rawrr::sampleFilePath() |> rawrr::readSpectrum(scan = 1:9) S[[1]] names(S[[1]]) plot(S[[1]]) # Example 2 - find best peptide spectrum match using the |> pipe operator # fetch via ExperimentHub if (require(ExperimentHub) & require(protViz)){ eh <- ExperimentHub::ExperimentHub() EH4547 <- normalizePath(eh[["EH4547"]]) (rawfile <- paste0(EH4547, ".raw")) if (!file.exists(rawfile)){ file.link(EH4547, rawfile) } GAG <- "GAGSSEPVTGLDAK" .bestPeptideSpectrumMatch <- function(rawfile, sequence="GAGSSEPVTGLDAK"){ readIndex(rawfile) |> subset(abs((1.008 + (protViz::parentIonMass(sequence) - 1.008) / 2) - precursorMass) < 0.001, select = scan) |> unlist() |> readSpectrum(rawfile = rawfile) |> lapply(function(x) { y <- protViz::psm(sequence = GAG, spec=x, plot=FALSE); y$scan <- x$scan; y }) |> lapply(FUN= function(x){ score <- sum(abs(x$mZ.Da.error) < 0.01); cbind(scan=x$scan, score=score) }) |> (function(x) as.data.frame(Reduce(rbind, x)))() |> subset(score > 0) |> (function(x) x[order(x$score, decreasing = TRUE), 'scan'])() |> head(1) } start_time <- Sys.time() bestMatch <- .bestPeptideSpectrumMatch(rawfile, GAG) |> rawrr::readSpectrum(rawfile=rawfile) |> lapply(function(x) protViz::peakplot(peptideSequence = GAG, x)) end_time <- Sys.time() end_time - start_time # Example 3 # using proteomicsdb \doi{10.1101/2020.09.08.287557} # through https://www.proteomicsdb.org/use/ .UniversalSpectrumExplorer <- function(x, sequence){ m <- protViz::psm( sequence, x) cat(paste(x$mZ[m$idx], "\t", x$intensity[m$idx]), sep = "\n") } rawrr::readSpectrum(rawfile=rawfile, 11091) |> lapply(function(x).UniversalSpectrumExplorer(x, sequence = GAG)) }
Read and extract scan trailer from TFS raw files.
readTrailer(rawfile, label = NULL)
readTrailer(rawfile, label = NULL)
rawfile |
the name of the raw file containing the mass spectrometry data from the Thermo Fisher Scientific instrument. |
label |
if NULL; the function scans for all available labels. |
an vector of trailers or values of a given trailer. Of note, the values are usually returned as a character.
rawrr::sampleFilePath() |> rawrr:::readTrailer() rawrr::sampleFilePath() |> rawrr:::readTrailer("AGC:") |> head()
rawrr::sampleFilePath() |> rawrr:::readTrailer() rawrr::sampleFilePath() |> rawrr:::readTrailer("AGC:") |> head()
sample.raw
BLOBThe binary example file sample.raw, shipped with the package, contains 574 Fourier-transformed Orbitrap spectra (FTMS) recorded on a Thermo Fisher Scientific Q Exactive HF-X. The mass spectrometer was operated in line with a nano electrospray source (NSI) in positive mode (+). All spectra were written to disk after applying centroiding (c) and lock mass correction.
sampleFilePath()
sampleFilePath()
Thermo Fisher Scientific Q Exactive HF-X raw file
of size 1.5M bytes and checksum
MD5 (sample.raw) = fe67058456c79af7442316c474d20e96
.
Additional raw data for
demonstration and extended testing is available through
MSV000086542
and the
tartare package.
Lions love raw meat!
file path of the sample.raw location.
Tobias Kockmann, 2018, 2019.
Bioconductor tartare package.
Automated quality control sample 1 (autoQC01) analyzed across different Thermo Fisher Scientific mass spectrometers, MSV000086542.
sampleFilePath()
sampleFilePath()
rawrrSpectrum
objectsAccessor function for scan number of rawrrSpectrum
objects
scanNumber(x)
scanNumber(x)
x |
A rawrrSpectrum object |
This accessor function returns the scan number of a mass spectrum
stored as rawrrSpectrum
object. Scan numbers are equal to the scan index
running from 1 to
with
being the last scan of a raw file.
The scan number of type integer
S <- readSpectrum(rawfile = sampleFilePath(), 1:10) scanNumber(S[[1]])
S <- readSpectrum(rawfile = sampleFilePath(), 1:10) scanNumber(S[[1]])
Text summary of chromatogram
## S3 method for class 'rawrrChromatogram' summary(object, ...)
## S3 method for class 'rawrrChromatogram' summary(object, ...)
object |
A |
... |
Function passes additional arguments. |
A rawrrChromatogram
object
C <- readChromatogram(rawfile = sampleFilePath(), mass = c(445.1181, 519.1367)) summary(C[[1]]) summary(C[[2]])
C <- readChromatogram(rawfile = sampleFilePath(), mass = c(445.1181, 519.1367)) summary(C[[1]]) summary(C[[2]])
Basic summary function
## S3 method for class 'rawrrSpectrum' summary(object, ...)
## S3 method for class 'rawrrSpectrum' summary(object, ...)
object |
an |
... |
Arguments to be passed to methods. |
This function creates a print message.
Christian Panse and Tobias Kockmann, 2020.
Total ion current of a spectrum
tic(x)
tic(x)
x |
A rawrrSpectrum object |
A double vector of length one.
S <- readSpectrum(rawfile = sampleFilePath(), 1) tic(S[[1]])
S <- readSpectrum(rawfile = sampleFilePath(), 1) tic(S[[1]])
Checks the validity of an readIndex
returned object.
validate_rawrrIndex(x)
validate_rawrrIndex(x)
x |
object to be validated. |
Validated data.frame
of readIndex
object
Tobias Kockmann and Christian Panse, 2020-12-09.
rawrr::sampleFilePath() |> rawrr::readIndex() |> rawrr::validate_rawrrIndex()
rawrr::sampleFilePath() |> rawrr::readIndex() |> rawrr::validate_rawrrIndex()
Checks the validity of rawrrSpectrum
object attributes.
validate_rawrrSpectrum(x)
validate_rawrrSpectrum(x)
x |
object to be validated. |
Validated rawrrSpectrum
object
Tobias Kockmann and Christian Panse, 2020.
S <- rawrr::sampleFilePath() |> rawrr::readSpectrum(scan = 1:9, validate=TRUE)
S <- rawrr::sampleFilePath() |> rawrr::readSpectrum(scan = 1:9, validate=TRUE)