srduplicated | R Documentation |
These generics order, rank, sort, and find duplicates in short read
objects, including fastq-encoded qualities. srorder
,
srrank
and srsort
differ from the default functions
rank
, order
and sort
in that sorting is based on
an internally-defined order rather than, e.g., the order implied by
LC_COLLATE
.
srorder(x, ...)
srrank(x, ...)
srsort(x, ...)
srduplicated(x, ...)
x |
The object to be sorted, ranked, ordered, or to have duplicates identified; see the examples below for objects for which methods are defined. |
... |
Additional arguments available for use by methods; usually ignored. |
Unlike sort
and friends, the implementation does not preserve
order of duplicated elements. Like duplicated
, one element in
each set of duplicates is marked as FALSE
.
srrank
settles ties using the “min” criterion described
in rank
, i.e., identical elements are ranked equal to the
rank of the first occurrence of the sorted element.
The following methods are defined, in addition to methods described in class-specific documentation:
signature(x = "XStringSet")
:
signature(x = "XStringSet")
:
signature(x = "XStringSet")
:
Apply srorder
, srrank
, srsort
,
srduplicated
to
XStringSet
objects such
as those returned by sread
.
signature(x = "ShortRead")
:
signature(x = "ShortRead")
:
signature(x = "ShortRead")
:
Apply srorder
, srrank
, srsort
,
srduplicated
to
XStringSet
objects to
the sread
component of ShortRead
and
derived objects.
The functions return the following values:
srorder |
An integer vector the same length as |
srrank |
An integer vector the same length as |
srsort |
An instance of |
srduplicated |
A logical vector the same length as |
Martin Morgan <mtmorgan@fhcrc.org>
showMethods("srsort")
showMethods("srorder")
showMethods("srduplicated")
sp <- SolexaPath(system.file('extdata', package='ShortRead'))
rfq <- readFastq(analysisPath(sp), pattern="s_1_sequence.txt")
sum(srduplicated(sread(rfq)))
srsort(sread(rfq))
srsort(quality(rfq))
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.