NCBSet-class: The NCBSet class

NCBSet-classR Documentation

The NCBSet class

Description

The NCBSet class is a container for a top-down proteomics experiment similar to the TopDownSet but instead of intensity values it just stores the information if a bond is covered by a N-terminal (encoded as 1), a C-terminal (encoded as 2) and/or bidirectional fragments (encoded as 3).

Usage

## S4 method for signature 'NCBSet'
bestConditions(
  object,
  n = ncol(object),
  minN = 0L,
  maximise = c("fragments", "bonds"),
  ...
)

## S4 method for signature 'NCBSet'
fragmentationMap(
  object,
  nCombinations = 10,
  cumCoverage = TRUE,
  maximise = c("fragments", "bonds"),
  labels = colnames(object),
  alphaIntensity = TRUE,
  ...
)

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

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

Arguments

object

NCBSet

n

integer, max number of combinations/iterations.

minN

integer, stop if there are less than minN additional fragments

maximise

character, optimisation targeting for the highest number of "fragments" (default) or "bonds".

nCombinations

integer, number of combinations to show (0 to avoid plotting them at all).

cumCoverage

logical, if TRUE (default) cumulative coverage of combinations is shown.

labels

character, overwrite x-axis labels.

alphaIntensity

logical, if TRUE (default) the alpha level of the color is used to show the colData(object)$AssignedIntensity. If FALSE the alpha is set to 1.

what

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

...

arguments passed to internal/other methods. added.

Value

An NCBSet object.

Methods (by generic)

  • bestConditions(NCBSet): Best combination of conditions.

    Finds the best combination of conditions for highest coverage of fragments or bonds. If there are two (or more conditions) that would add the same number of fragments/bonds the one with the highest assigned intensity is used. Use n to limit the number of iterations and combinations that should be returned. If minN is set at least minN fragments have to be added to the combinations. The function returns a 7-column matrix. The first column contains the index (Index) of the condition (column number). The next columns contain the newly added number of fragments or bonds (FragmentsAddedToCombination, BondsAddedToCombination), the fragments or bonds in the condition (FragmentsInCondition, BondsInCondition), and the cumulative coverage fragments or bonds (FragmentCoverage, BondCoverage).

  • fragmentationMap(NCBSet): Plot fragmentation map.

    Plots a fragmentation map of the Protein. Use nCombinations to add another plot with nCombinations combined conditions. If cumCoverage is TRUE (default) these combinations increase the coverage cumulatively.

  • summary(NCBSet): Summary statistics.

    Returns a matrix with some statistics: number of fragments, total/min/first quartile/median/mean/third quartile/maximum of intensity values.

Slots

rowViews

Biostrings::XStringViews, information about bonds (name, start, end, width, sequence), see Biostrings::XStringViews for details.

colData

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

assay

Matrix::dgCMatrix, coverage values of the bonds. The rows correspond to the bonds and the columns to the condition/run. It just stores values that are different from zero. If a bond is covered by an N-terminal fragment its encoded with 1, by an C-terminal fragmentl with 2 and by both fragment types/bidirectional by 3 respectively.

files

character, files that were imported.

processing

character, log messages.

Author(s)

Sebastian Gibb mail@sebastiangibb.de

See Also

  • An NCBSet is generated from an TopDownSet object.

  • Biostrings::XStringViews for the row view interface.

  • Matrix::dgCMatrix for technical details about the coverage storage.

Examples

## Example data
data(tds, package="topdownr")

## Aggregate technical replicates
tds <- aggregate(tds)

## Coercion into an NCBSet object
ncb <- as(tds, "NCBSet")

ncb

head(summary(ncb))

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

# Accessing colData
ncb$Mz

# Subsetting

# First 100 bonds
ncb[1:100]

# Just bond 152
ncb["bond152"]

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

# Plot fragmentation map
fragmentationMap(ncb)

sgibb/topdownr documentation built on Jan. 16, 2024, 12:14 a.m.