XStringQuality-class: PhredQuality, SolexaQuality and IlluminaQuality objects

Description Usage Arguments Details Alphabet and encoding Author(s) See Also Examples

Description

Objects for storing string quality measures.

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
## Constructors:
PhredQuality(x)
SolexaQuality(x)
IlluminaQuality(x)

## alphabet and encoding
## S4 method for signature 'XStringQuality'
alphabet(x)
## S4 method for signature 'XStringQuality'
encoding(x)

Arguments

x

Either a character vector, BString, BStringSet, integer vector, or number vector of error probabilities.

Details

PhredQuality objects store characters that are interpreted as [0 - 99] quality measures by subtracting 33 from their ASCII decimal representation (e.g. ! = 0, " = 1, \# = 2, ...). Quality measures q encode probabilities as -10 * log10(p).

SolexaQuality objects store characters that are interpreted as [-5 - 99] quality measures by subtracting 64 from their ASCII decimal representation (e.g. ; = -5, < = -4, = = -3, ...). Quality measures q encode probabilities as -10 * (log10(p) - log10(1 - p)).

IlluminaQuality objects store characters that are interpreted as [0 - 99] quality measures by subtracting 64 from their ASCII decimal representation (e.g. @ = 0, A = 1, B = 2, ...). Quality measures q encode probabilities as -10 * log10(p)

Alphabet and encoding

In the code snippets below, x is an XStringQuality object.

alphabet(x): Valid letters in this quality score; not all letters are encountered in actual sequencing runs.

encoding(x): Map between letters and their corresponding integer encoding. Use as.integer and as.numeric to coerce objects to their integer and probability representations.

Author(s)

P. Aboyoun

See Also

pairwiseAlignment, PairwiseAlignments-class, DNAString-class, BStringSet-class

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
PhredQuality(0:40)
SolexaQuality(0:40)
IlluminaQuality(0:40)

pq <- PhredQuality(c("*+,-./", "0123456789:;"))
qs <- as(pq, "IntegerList")  # quality scores
qs
as(qs, "PhredQuality")
p <- as(pq, "NumericList")  # probabilities
as(p, "PhredQuality")

PhredQuality(seq(1e-4,0.5,length=10))
SolexaQuality(seq(1e-4,0.5,length=10))
IlluminaQuality(seq(1e-4,0.5,length=10))

x <- SolexaQuality(BStringSet(c(a="@ABC", b="abcd")))
as(x, "IntegerList")  # quality scores
as(x, "NumericList")  # probabilities
as.matrix(x)          # quality scores

Biostrings documentation built on Nov. 8, 2020, 11:12 p.m.