View source: R/functions-TopDownSet.R
readTopDownFiles | R Documentation |
It creates an TopDownSet object and is its only constructor.
readTopDownFiles(
path,
pattern = ".*",
type = c("a", "b", "c", "x", "y", "z"),
modifications = c("Carbamidomethyl", "Acetyl", "Met-loss"),
customModifications = data.frame(),
adducts = data.frame(),
neutralLoss = PSMatch::defaultNeutralLoss(),
sequenceOrder = c("original", "random", "inverse"),
tolerance = 5e-06,
redundantIonMatch = c("remove", "closest"),
redundantFragmentMatch = c("remove", "closest"),
dropNonInformativeColumns = TRUE,
sampleColumns = c("Mz", "AgcTarget", "EtdReagentTarget", "EtdActivation",
"CidActivation", "HcdActivation", "UvpdActivation"),
conditions = "ScanDescription",
verbose = interactive()
)
path |
|
pattern |
|
type |
|
modifications |
|
customModifications |
|
adducts |
|
neutralLoss |
|
sequenceOrder |
|
tolerance |
|
redundantIonMatch |
|
redundantFragmentMatch |
|
dropNonInformativeColumns |
logical, should columns with just one identical value across all runs be removed? |
sampleColumns |
|
conditions |
|
verbose |
|
readTopDownFiles
reads and processes all top-down files, namely:
.fasta
(protein sequence)
.mzML
(spectra)
.experiments.csv
(method/fragmentation conditions)
.txt
(scan header information)
customModifications
: additional to the provided unimod modifications
available through the modifications
argument customModifications
allow to
apply user-definied modifications to the output of
PSMatch::calculateFragments()
.
The customModifications
argument takes a
data.frame
with the mass
to add, the name
of the modification, the
location (could be the position of the amino acid or "N-term"/"C-term"),
whether the modification is always seen (variable=FALSE
) or both, the
modified and unmodified amino acid are present (variable=TRUE
), e.g.
for Activation (which is available via modification="Acetyl"
)
data.frame(mass=42.010565, name="Acetyl", location="N-term", variable=FALSE)
or variable one (that could be present or not):
data.frame(mass=365.132, name="Custom", location=10, variable=TRUE)
If the customModifications
data.frame
contains multiple columns the
modifications are applied from row one to the last row one each time.
adducts
: Thermo's Xtract
allows some mistakes in deisotoping, mostly it
allows +/- C13-C12
and +/- H+
.
The adducts
argument takes a
data.frame
with the mass
to add, the name
that should assign to these
new fragments and an information to
whom the modification should be
applied, e.g. for H+
on z
,
data.frame(mass=1.008, name="zpH", to="z")
.
Please note: The adducts
are added to the output of
PSMatch::calculateFragments()
.
That has some limitations, e.g.
neutral loss calculation could not be done in
topdownr-package.
If neutral loss should be applied on adducts you have to create
additional rows, e.g.:
data.frame(mass=c(1.008, 1.008), name=c("cpH", "cpH_"), to=c("c", "c_"))
.
A TopDownSet
object.
PSMatch::calculateFragments()
,
PSMatch::defaultNeutralLoss()
if (require("topdownrdata")) {
# add H+ to z and no neutral loss of water
tds <- readTopDownFiles(
topdownrdata::topDownDataPath("myoglobin"),
## Use an artifical pattern to load just the fasta
## file and files from m/z == 1211, ETD reagent
## target 1e6 and first replicate to keep runtime
## of the example short
pattern=".*fasta.gz$|1211_.*1e6_1",
adducts=data.frame(mass=1.008, name="zpH", to="z"),
neutralLoss=PSMatch::defaultNeutralLoss(
disableWaterLoss=c("Cterm", "D", "E", "S", "T")),
tolerance=25e-6
)
}
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.