TopDownSet-class: The TopDownSet class

Description Usage Arguments Details Value Methods (by generic) Slots Coercion Author(s) See Also Examples

Description

The TopDownSet class is a container for a whole top-down proteomics experiment.

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
## S4 method for signature 'TopDownSet'
aggregate(x, by = x$Sample, ...)

## S4 method for signature 'TopDownSet,TopDownSet'
combine(x, y)

## S4 method for signature 'TopDownSet'
filterCv(object, threshold, by = object$Sample, ...)

## S4 method for signature 'TopDownSet'
filterInjectionTime(
  object,
  maxDeviation = log2(3),
  keepTopN = 2,
  by = object$Sample,
  ...
)

## S4 method for signature 'TopDownSet'
filterIntensity(object, threshold, relative = TRUE, ...)

## S4 method for signature 'TopDownSet'
filterNonReplicatedFragments(object, minN = 2, by = object$Sample, ...)

## S4 method for signature 'TopDownSet'
normalize(object, method = "TIC", ...)

## S4 method for signature 'TopDownSet,missing'
plot(x, y, ..., verbose = interactive())

## S4 method for signature 'TopDownSet'
show(object)

## S4 method for signature 'TopDownSet'
summary(object, what = c("conditions", "fragments"), ...)

Arguments

x, object

TopDownSet

by

list, grouping variable, in general it refers to technical

y

missing, not used.

threshold

double, threshold variable.

maxDeviation

double, maximal allowed deviation in the log2 injection time in ms in comparison to the median ion injection time.

keepTopN

integer, how many technical replicates should be kept?

relative

logical, if relative is TRUE all fragments with intensity below threshold * max(intensity) per fragment are removed, otherwise all fragments below threshold are removed.

minN

numeric, if less than minN of a fragment are found across technical replicates it is removed.

method

character, normalisation method, currently just "TIC" for Total Ion Current normalisation of the scans/conditions (column-wise normalisation) is supported.

verbose

logical, verbose output?

what

character, specifies whether "conditions" (columns; default) or "fragments" (rows) should be summarized.

...

arguments passed to internal/other methods. replicates (that's why the default is the "Sample" column in colData).

Details

See vignette("analysis", package="topdownr") for a detailed example how to work with TopDownSet objects.

Value

An TopDownSet object.

Methods (by generic)

Slots

rowViews

FragmentViews, information about fragments (name, type, sequence, mass, charge), see FragmentViews for details.

colData

S4Vectors::DataFrame, information about the MS2 experiments and the fragmentation conditions.

assay

Matrix::dgCMatrix, intensity values of the fragments. The rows correspond to the fragments and the columns to the condition/run. It just stores values that are different from zero.

files

character, files that were imported.

tolerance

double, tolerance in ppm that were used for matching the experimental mz values to the theoretical fragments.

redundantMatching

character, matching strategies for redundant ion/fragment matches. See redundantIonMatch and redundantFragmentMatch in readTopDownFiles() for details.

processing

character, log messages.

Coercion

'as(object, "MSnSet"): Coerce an TopDownSet object into an MSnbase::MSnSet object.

'as(object, "NCBSet"): Coerce an TopDownSet object into an NCBSet object.

Author(s)

Sebastian Gibb mail@sebastiangibb.de

See Also

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
## Example data
data(tds, package="topdownr")

tds

head(summary(tds))

# Accessing slots
rowViews(tds)
colData(tds)
head(assayData(tds))

# Accessing colData
tds$Mz
tds$FilterString

# Subsetting

# First 100 fragments
tds[1:100]

# All c fragments
tds["c"]

# Just c 152
tds["c152"]

# Condition 1 to 10
tds[, 1:10]

# Filtering
# Filter all intensities that don't have at least 10 % of the highest
# intensity per fragment.
tds <- filterIntensity(tds, threshold=0.1)

# Filter all conditions with a CV above 30 % (across technical replicates)
tds <- filterCv(tds, threshold=30)

# Filter all conditions with a large deviation in injection time
tds <- filterInjectionTime(tds, maxDeviation=log2(3), keepTopN=2)

# Filter all conditions where fragments don't replicate
tds <- filterNonReplicatedFragments(tds)

# Normalise by TIC
tds <- normalize(tds)

# Aggregate technical replicates
tds <- aggregate(tds)

head(summary(tds))

# Coercion
as(tds, "NCBSet")

if (require("MSnbase")) {
    as(tds, "MSnSet")
}
## Not run: 
# plot a single condition
# pseudo-code (replace topdownset with your object)
plot(topdownset[,1])

# plot the whole object
pdf("topdown-spectra.pdf", paper="a4r", width=12)
# pseudo-code (replace topdownset with your object)
plot(topdownset)
dev.off()

## End(Not run)

topdownr documentation built on Nov. 8, 2020, 8:10 p.m.