Version: | 1.9.3 |
Date: | 2024-10-02 |
Title: | Reads Mass Spectrometry Data in Bruker *flex Format |
Depends: | R (≥ 3.3.0) |
Suggests: | testthat |
Description: | Reads data files acquired by Bruker Daltonics' matrix-assisted laser desorption/ionization-time-of-flight mass spectrometer of the *flex series. |
License: | GPL (≥ 3) |
URL: | https://strimmerlab.github.io/software/maldiquant/, https://github.com/sgibb/readBrukerFlexData/ |
BugReports: | https://github.com/sgibb/readBrukerFlexData/issues/ |
LazyLoad: | yes |
RoxygenNote: | 7.3.1 |
Encoding: | UTF-8 |
NeedsCompilation: | no |
Packaged: | 2024-10-02 21:01:19 UTC; sebastian |
Author: | Sebastian Gibb |
Maintainer: | Sebastian Gibb <mail@sebastiangibb.de> |
Repository: | CRAN |
Date/Publication: | 2024-10-02 21:20:02 UTC |
The readBrukerFlexData Package
Description
The readBrukerFlexData package reads data files acquired by MALDI-TOF MS on
Bruker Daltonics machines of the *flex series. (autoflex,
microflex, ultraflex).
The package was developed without any knowledge nor even support by
Bruker Daltonics.
All trademarks are owned by or licensed to Bruker Daltonics.
Author(s)
Sebastian Gibb mail@sebastiangibb.de
References
https://github.com/sgibb/readBrukerFlexData
See Also
Useful links:
Report bugs at https://github.com/sgibb/readBrukerFlexData/issues/
Change precision.
Description
This function converts double values to double values in a given precision (only correctly working for cut a higher precision to a lower one; e.g. IEEE 754 double precision to IEEE 754 single precision)
Usage
.changePrecision(x, size)
Arguments
x |
|
size |
|
Value
double
.
See Also
Cubic calibration
Description
Only basic support (not 100% identical results) for Bruker Daltonics' cubic calibration (encoded in '"##$NTBCal"' as '"V1.0CTOF2CalibrationConstants"'). This is an internal function and should normally not used by the user.
Usage
.ctof2calibration(tof, d)
Arguments
tof |
|
d |
|
Value
double
, mz values.
Note
Bruker Daltonics doesn't explain how the cubic calibration works.
All formula are results of “trial and error”.
That is why mass calculated by .ctof2calibration
differs little from original mass.
See also https://github.com/sgibb/readBrukerFlexData/issues/3.
Author(s)
Alan Race, Samuel Granjeaud, Sebastian Gibb
See Also
https://github.com/sgibb/readBrukerFlexData/issues/3
Converts double to single precision.
Description
This function simulates the conversion of floating point numbers from double
precision (64bit, R: double()
, C: double
) to single precision
(32bit, R: none
, C: float
). It follows IEEE 754 standard.
Usage
.double2singlePrecision(x)
Arguments
x |
|
Details
The same could be done in C by using casts:
double precision32(double value) { float x=value; return (double)x; }
Value
double
(in single precision).
References
IEEE 754 standard: https://ieeexplore.ieee.org/xpl/articleDetails.jsp?arnumber=30711&filter=AND(p_Publication_Number:2355)
See Also
Examples
## load library
library("readBrukerFlexData")
## show more details
oldDigits <- options()$digits
options(digits=22)
## a test number
num <- 1/3
num
readBrukerFlexData:::.double2singlePrecision(num)
## reset digits option
options(digits=oldDigits)
Extracts High Precision Calibration coefficients.
Description
This is a helper function to extract coefficients and constants from metaData$hpcStr.
Usage
.extractHPCConstants(hpcStr)
Arguments
hpcStr |
metaData$hpcStr, which store coefficents hpcStr <- "V1.0CHPCData Order 10 vCoeff V1.0VectorDouble 11 -0.48579224953906053 0.0009361303203700988 -6.92711401708155e-008 -1.0992953299897006e-009 1.1718229914003113e-012 -5.392578762547374e-016 9.0176664604755316e-020 1.9704001597871883e-023 -1.1794161284667635e-026 2.0351573912658823e-030 -1.2617853301428769e-034 c2 -0.046701600316874939 c0 237.64781433281422 minMass 736.50266799999997 maxMass 3698.6377320000001 bUse 1 endCHPCData" |
Value
a list
:
hpcConstants$coefficients
: double, vector of coefficentshpcConstants$calibrationConstant0
: c0hpcConstants$calibrationConstant2
: c2
See Also
Extracts values from acqu file.
Description
This function is a helper function to extract double values from acqu file.
Usage
.grepAcquDoubleValue(patternStr, srcStr)
Arguments
patternStr |
|
srcStr |
|
Value
double
vector of the value given in patternStr
See Also
Extracts values from acqu file.
Description
This function is a helper function to extract values from acqu file.
Usage
.grepAcquValue(patternStr, srcStr)
Arguments
patternStr |
|
srcStr |
|
Value
character
vector of the value given in patternStr
See Also
.grepAcquDoubleValue
,
.readAcquFile
High Precision Calibration
Description
Only basic support (not 100% identical results) for Bruker Daltonics' HPC.
HPC stands for High Precision Calibration.
This is an internal function and should normally not used by the user.
Usage
.hpc(mass, minMass, maxMass, hpcCoefficients)
Arguments
mass |
|
minMass |
|
maxMass |
|
hpcCoefficients |
|
Details
Bruker Daltonics doesn't explain how HPC works. All formula are results of
“trial and error”. That is why mass calculated by .hpc
differs little from original HPC mass.
(In example file 214 of 24860 mass are incorrect;
deviations: min: 6.103515625e-05, max: 0.02935791015625.)
In the manual of mass spectrometry instruments of Bruker Daltonics machines
the *flex series you can find an article about HPC principles:
Gobom, J. and Mueller, M. and Egelhofer V. and Theiss, D. and
Lehrach, H. and Nordhoff, E. (2002)
“A Calibration Method That Simplifies and Improves Accurate
Determination of Peptide Molecular mass by MALDI-TOF MS.”,
Anal Chem 74: 3915-3923
https://pubmed.ncbi.nlm.nih.gov/12175185/
Value
A vector of HPC corrected mass (double
).
Note
Please note that .hpc is not correct! You have been warned.
References
Gobom, J. and Mueller, M. and Egelhofer V. and Theiss, D. and
Lehrach, H. and Nordhoff, E. (2002)
“A Calibration Method That Simplifies and Improves Accurate
Determination of Peptide Molecular mass by MALDI-TOF MS.”,
Anal Chem 74: 3915-3923
https://pubmed.ncbi.nlm.nih.gov/12175185/
See Also
readBrukerFlexDir
,
readBrukerFlexFile
,
.double2singlePrecision
Examples
## load library
library("readBrukerFlexData")
## get examples directory
exampleDirectory <- system.file("Examples", package="readBrukerFlexData")
## read example spectra
## please note: filterZeroIntensities=TRUE is used for compatibility reasons.
## You should NOT use it!
noHpcSpec <- readBrukerFlexFile(file.path(exampleDirectory,
"hpc/fid/0_A20/1/1SRef/fid"), filterZeroIntensities=TRUE, useHpc=FALSE)
hpcSpec <- readBrukerFlexFile(file.path(exampleDirectory,
"hpc/fid/0_A20/1/1SRef/fid"), filterZeroIntensities=TRUE)
## plot spectrum
plot(noHpcSpec$spectrum$mass, noHpcSpec$spectrum$intensity, type="l",
col="red", xlim=c(1296, 1300))
lines(hpcSpec$spectrum$mass, hpcSpec$spectrum$intensity, type="l",
col="green", xlim=c(1296, 1300))
legend(x="topright", legend=c("no hpc", "hpc"), col=c("red", "green"), lwd=1)
## show difference between .hpc and original HPC
## load mzXML generated by Bruker Daltonics CompassXport 1.3.5
## you could do it like this:
#library("readMzXmlData")
#cpSpecHpcMzXml <- readMzXmlFile(file.path(exampleDirectory,
# "hpc/mzXML/hpc.mzXML"))
## or easily use:
data(cpSpecHpcMzXml)
## reduce R double precision to single precision because our CompassXport 1.3.5
## supports only mzXML with precision=32 (only for compatibility reasons)
noHpcSpec$spectrum$mass32 <-
readBrukerFlexData:::.double2singlePrecision(noHpcSpec$spectrum$mass)
hpcSpec$spectrum$mass32 <-
readBrukerFlexData:::.double2singlePrecision(hpcSpec$spectrum$mass)
## calculate deviance
d <- noHpcSpec$spectrum$mass32-cpSpecHpcMzXml$spectrum$mass
dHPC <- hpcSpec$spectrum$mass32-cpSpecHpcMzXml$spectrum$mass
## a little summary
cat("without .hpc:\n",
"not matching: ", length(cpSpecHpcMzXml$spectrum$mass[d!=0]), " of ",
length(cpSpecHpcMzXml$spectrum$mass), "; range: ",
range(abs(d[d!=0])), "\nwith .hpc:\n",
"not matching: ", length(cpSpecHpcMzXml$spectrum$mass[dHPC!=0]), " of ",
length(cpSpecHpcMzXml$spectrum$mass), "; range: ",
range(abs(d[dHPC!=0])), "\n")
##
## doing things manually
##
hpcMass <- readBrukerFlexData:::.hpc(mass=noHpcSpec$spectrum$mass,
minMass=noHpcSpec$metaData$hpc$limits["minMass"],
maxMass=noHpcSpec$metaData$hpc$limits["maxMass"],
hpcCoefficients=noHpcSpec$metaData$hpc$coefficients)
Reads an acqu file.
Description
This function reads constants, calibrations values and a lot of more from acqu files.
Usage
.readAcquFile(fidFile, verbose = FALSE)
Arguments
fidFile |
|
verbose |
|
Value
a list
with metadata
See Also
Reads binary fid file.
Description
This function reads a binary fid file. A fid file contains intensities for all measured time points.
Usage
.readFidFile(fidFile, nIntensities, endian = "little")
Arguments
fidFile |
|
nIntensities |
number of data entries (total count; get from acqu file) |
endian |
endianness of the fid file (‘little’ or ‘big’; default: ‘little’) |
Value
A vector of intensity values.
See Also
.readAcquFile
,
readBrukerFlexFile
Creates sample name.
Description
This function guesses the name of current spot from filename.
Usage
.sampleName(fidFile)
Arguments
fidFile |
|
Value
character
, sample name: e.g ‘Pankreas_HB_L_061019_A10’
Note
WARNING: if the 4th upper directory hasn't an unique name you will get equal names for your spot list
e.g. the following will create a list with 4 elements but
only 2 unique spot names (2-100kDa:0_A1 and 2-100kDa:0_B1)
./Run1/2-100kDa/0_A1/1/1SLin/fid
./Run1/2-100kDa/0_B1/1/1SLin/fid
./Run2/2-100kDa/0_A1/1/1SLin/fid
./Run2/2-100kDa/0_B1/1/1SLin/fid
See Also
.readAcquFile
,
readBrukerFlexFile
Calculates mass from time-of-flight values.
Description
This function calculates mass from time of flight values
based on the following article:
Titulaer, M. K. and Siccama, I. and Dekker, J. L. and
van Rijswijk, A. L. and Heeren R. M. and Sillevis Smitt, P. A. and
Luider, T. M. (2006)
“A database application for pre-processing, storage and comparison
of mass spectra derived from patients and controls”,
BMC Bioinformatics, 7: 403
https://pubmed.ncbi.nlm.nih.gov/16953879/
Usage
.tof2mass(tof, c1, c2, c3)
Arguments
tof |
|
c1 |
metaData$calibrationConstants[1] |
c2 |
metaData$calibrationConstants[2] |
c3 |
metaData$calibrationConstants[3] |
Details
Arguments are imported from metadata (acqu-file).
Value
double
, calculated mass
References
Titulaer, M. K. and Siccama, I. and Dekker, J. L. and
van Rijswijk, A. L. and Heeren R. M. and Sillevis Smitt, P. A. and
Luider, T. M. (2006)
“A database application for pre-processing, storage and comparison
of mass spectra derived from patients and controls”,
BMC Bioinformatics, 7: 403
https://pubmed.ncbi.nlm.nih.gov/16953879/
Mass spectrum generated by Bruker Daltonics CompassXport
Description
This dataset was generated by Bruker Daltonics CompassXport and imported by
readMzXmlFile
to R. It is only needed for
comparison between Bruker Daltonics' HPC and
.hpc
.
Format
A list containing a mass and an intensity vector.
Source
Examples/hpc/mzXML/hpc.mzXML
See Also
Removed functions in package readBrukerFlexData
Description
These functions are defunct and no longer available.
Details
- mqReadBrukerFlex:
use
importBrukerFlex
instead.
Reads recursively mass spectrometry data in Bruker Daltonics XMASS format.
Description
This function leads recursively all mass spectrometry data in Bruker Daltonics XMASS format in a specified directory.
Usage
readBrukerFlexDir(
brukerFlexDir,
removeCalibrationScans = TRUE,
removeMetaData = FALSE,
useHpc = TRUE,
useSpectraNames = TRUE,
filterZeroIntensities = FALSE,
verbose = FALSE
)
Arguments
brukerFlexDir |
|
removeCalibrationScans |
|
removeMetaData |
|
useHpc |
|
useSpectraNames |
|
filterZeroIntensities |
|
verbose |
|
Details
See readBrukerFlexFile
.
Value
A list
of spectra.
[[1]]$spectrum$mass
: A vector of calculated mass.[[1]]$spectrum$intensity
: A vector of intensity values.[[1]]$metaData
: A list of metaData depending on read spectrum.
See Also
importBrukerFlex
,
readBrukerFlexFile
,
.hpc
Examples
## load library
library("readBrukerFlexData")
## get examples directory
exampleDirectory <- system.file("Examples", package="readBrukerFlexData")
## read example spectra
spec <- readBrukerFlexDir(file.path(exampleDirectory,
"2010_05_19_Gibb_C8_A1"))
## plot spectra
plot(spec[[1]]$spectrum$mass, spec[[1]]$spectrum$intensity, type="n")
l <- length(spec)
legendStr <- character(l)
for (i in seq(along=spec)) {
lines(spec[[i]]$spectrum$mass, spec[[i]]$spectrum$intensity, type="l",
col=rainbow(l)[i])
legendStr[i] <- spec[[i]]$metaData$fullName
}
## draw legend
legend(x="topright", legend=legendStr, col=rainbow(l), lwd=1)
Reads mass spectrometry data in Bruker Daltonics XMASS format.
Description
This function reads mass spectrometry data in Bruker Daltonics XMASS format used by Bruker Daltonics mass spectrometer of *flex series (autoflex, microflex, ultraflex).
Usage
readBrukerFlexFile(
fidFile,
removeMetaData = FALSE,
useHpc = TRUE,
filterZeroIntensities = FALSE,
keepNegativeIntensities = FALSE,
verbose = FALSE
)
Arguments
fidFile |
|
removeMetaData |
|
useHpc |
|
filterZeroIntensities |
|
keepNegativeIntensities |
|
verbose |
|
Details
readBrukerFlexFile
has to import the following data to calculating
mass from acqu file:
acqu-value | becomes metaData | description |
$BYTORDA | metaData$byteOrder | endianness of fid file |
$TD | metaData$number | total number of measured time periods |
$DELAY | metaData$timeDelay | first measured intensity after metaData$timeDelay ns |
$DW | metaData$timeDelta | ns between measured time periods |
$ML1 | metaData$calibrationConstants[1] | mass calibration constant |
$ML2 | metaData$calibrationConstants[2] | mass calibration constant |
$ML3 | metaData$calibrationConstants[3] | mass calibration constant |
If High Precision Calibration (HPC) is used, readBrukerFlexFile
needs:
acqu-value | becomes metaData | description |
$HPClBHi | metaData$hpc$limits[“maxMass”] | upper mass threshold |
$HPClBLo | metaData$hpc$limits[“minMass”] | lower mass threshold |
$HPClOrd | metaData$hpc$order | polynomial order |
$HPClUse | metaData$hpc$use | maybe using of HPC? (seems to be always “yes” in our test data) |
$HPCStr | metaData$hpc$coefficients | polynomial coefficients in a string |
readBrukerFlexFile
tries also to import [optional]:
acqu-value | becomes metaData | description |
DATATYPE | metaData$dataType | e.g CONTINUOUS MASS SPECTRUM |
SPECTROMETER/DATASYSTEM | metaData$dataSystem | e.g. Bruker Flex Series |
.SPECTROMETER TYPE | metaData$spectrometerType | e.g. TOF |
.INLET | metaData$inlet | DIRECT |
.IONIZATION MODE | metaData$ionizationMode | e.g. LD+ |
$DATE | metaData$date | same as $AQ_DATE but often only "0" |
$ACQMETH | metaData$acquisitionMethod | path to method file |
$AQ_DATE | metaData$acquisitionDate | acquisition date |
$AQ_mod | metaData$acquisitionMode | acquisition mode |
$AQOP_m | metaData$acquisitionOperatorMode, metaData$tofMode | LINEAR / REFLECTOR |
$ATTEN | metaData$laserAttenuation | laser beam attenuation |
$CMT[1:4] | metaData$comments | comments |
$DEFLON | metaData$deflection | deflection ON/OFF |
$DIGTYP | metaData$digitizerType | type of digitizer |
$DPCAL1 | metaData$deflectionPulserCal1 | deflection pulser cal 1 |
$DPMASS | metaData$deflectionPulserMass | deflection pulser mass |
$FCVer | metaData$flexControlVersion | Version of Bruker Daltonics FlexControl software |
$ID_raw | metaData$id | spectrum id |
$INSTRUM | metaData$instrument | e.g. AUTOFLEX |
$InstrID | metaData$instrumentId | ID of mass spectrometer |
$InstTyp | metaData$instrumentType | instrument type |
$Lift1 | metaData$lift[1] | LIFT constant? |
$Lift2 | metaData$lift[2] | LIFT constant? |
$Masserr | metaData$massError | initial mass error in ppm |
$NoSHOTS | metaData$laserShots | number of applied laser shots |
$PATCHNO | metaData$patch | sample postion on target |
$PATH | metaData$path | original file path (on Bruker *flex series controller PC) |
$REPHZ | metaData$laserRepetition | laser repetition rate in Hz |
$SPOTNO | metaData$spot | same as $PATCHNO
(in older files often empty, that's why readBrukerFlexFile
uses $PATHNO instead) |
$SPType | metaData$spectrumType | e.g. TOF |
$TgIDS | metaData$target$id | target ids |
$TgCount | metaData$target$count | number of measurements with this target |
$TgSer | metaData$target$serialNumber | target serial number |
$TgTyp | metaData$target$typeNumber | target type number |
$TLift | metaData$tlift | LIFT constant? |
import from file path:
value | becomes metaData | description |
full current path to fid file | metaData$file | path on local machine |
sample name | metaData$sampleName | - |
filterZeroIntensities
: Change default value is not recommended!
If TRUE
all intensities equal zero are removed.
This parameter exists only to be compatible to
Bruker Daltonics CompassXport's mzXML export function. For details see:
‘Release Notes for CompassXport 3.0.3’,
cap. 6 ‘Filtering of Zero Intensities’:
“Bruker Daltonics' Acquisition Software will compress Analysis raw
data. To save on operation time and to keep export file sizes small,
CompassXport 3.0.3 will filter out zero (0.0) intensities
when exporting to mzXML or mzData ...”
keepNegativeIntensities
: Change default value is
not recommended! If TRUE
negative intensity values are not
replaced by zero. This parameter exists only to be compatible to
Bruker Daltonics CompassXport.
Value
A list
of spectra and metadata.
spectrum$mass
: A vector of calculated mass.spectrum$tof
: A vector of time-of-flight data.spectrum$intensity
: A vector of intensity values.metaData
: A list of metaData depending on read spectrum.
See Also
https://github.com/sgibb/readBrukerFlexData/wiki,
importBrukerFlex
,
readBrukerFlexDir
,
.hpc
Examples
## load library
library("readBrukerFlexData")
## get examples directory
exampleDirectory <- system.file("Examples", package="readBrukerFlexData")
## read example spectrum
spec <- readBrukerFlexFile(file.path(exampleDirectory,
"2010_05_19_Gibb_C8_A1/0_A1/1/1SLin/fid"))
## print metaData
print(spec$metaData)
## plot spectrum
plot(spec$spectrum$mass, spec$spectrum$intensity, type="l", col="red")