Last data update: 2014.03.03

R: 'SffHeader'
SffHeader-classR Documentation

SffHeader

Description

Class SffHeader contains meta-data stored in the header of the SFF files read in.

Objects from this class are the result of readSffHeader, or from the result of readSff. The resulting object will contain a header slot which is a list. If multiple sff files were processed by either of the above functions, this list will contain meta-data about each of these files in corresponding positions in the list.

Meta-data included is defined in the SFF file specifications and include:

filename:

The name of the file that was read in.

magic_number:

779314790, which encodes the string ".sff"

version:

Version number

index_offset:

An optional field which indicates the position of a read index within the file.

index_length:

An optional field which indicates the length of a read index within the file.

number_of_reads:

Stores the number of reads in the file.

header_length:

The number of bytes required by header fields.

key_length:

The length of the key sequence used for these reads.

number_of_flows_per_read:

The number of flows carried out during sequencing of the reads.

flowgram_format_code:

Indicates the format of the flowgram encoding. Currently "1" is the only valid value.

flow_chars:

Indicates the nucleotide bases used for each sequencing flow.

key_sequence:

Nucleotide sequence use for these reads.

Slots

header:

Object of class "list", containing data frames or lists of data frames summarizing a description of the SFF files.

Methods

header

signature(object = "SffHeader"): access the header slot of object, returning a list object.

Author(s)

Matt Settles <msettles@uidaho.edu>

See Also

SffReads,SffReadsQ.

Examples

showClass("SffHeader")

## The readSffHeader can be used to extract header information from one or more sff files:
sffFiles = c(system.file("extdata", "SmallTorrentTest.sff", package = "rSFFreader"), 
             system.file("extdata","Small454Test.sff",package="rSFFreader"))
header <- readSffHeader(sffFiles)
header
header(header)[[1]]$number_of_reads
header(header)[[2]]$number_of_reads

## Header information is also retrieved when using readSff:
sff <- readSff(sffFiles)

## Number of flows for dataset 1 and 2:
header(sff)[[1]]$number_of_flows_per_read
header(sff)[[2]]$number_of_flows_per_read

Results


R version 3.3.1 (2016-06-21) -- "Bug in Your Hair"
Copyright (C) 2016 The R Foundation for Statistical Computing
Platform: x86_64-pc-linux-gnu (64-bit)

R is free software and comes with ABSOLUTELY NO WARRANTY.
You are welcome to redistribute it under certain conditions.
Type 'license()' or 'licence()' for distribution details.

R is a collaborative project with many contributors.
Type 'contributors()' for more information and
'citation()' on how to cite R or R packages in publications.

Type 'demo()' for some demos, 'help()' for on-line help, or
'help.start()' for an HTML browser interface to help.
Type 'q()' to quit R.

> library(rSFFreader)
Loading required package: ShortRead
Loading required package: BiocGenerics
Loading required package: parallel

Attaching package: 'BiocGenerics'

The following objects are masked from 'package:parallel':

    clusterApply, clusterApplyLB, clusterCall, clusterEvalQ,
    clusterExport, clusterMap, parApply, parCapply, parLapply,
    parLapplyLB, parRapply, parSapply, parSapplyLB

The following objects are masked from 'package:stats':

    IQR, mad, xtabs

The following objects are masked from 'package:base':

    Filter, Find, Map, Position, Reduce, anyDuplicated, append,
    as.data.frame, cbind, colnames, do.call, duplicated, eval, evalq,
    get, grep, grepl, intersect, is.unsorted, lapply, lengths, mapply,
    match, mget, order, paste, pmax, pmax.int, pmin, pmin.int, rank,
    rbind, rownames, sapply, setdiff, sort, table, tapply, union,
    unique, unsplit

Loading required package: BiocParallel
Loading required package: Biostrings
Loading required package: S4Vectors
Loading required package: stats4

Attaching package: 'S4Vectors'

The following objects are masked from 'package:base':

    colMeans, colSums, expand.grid, rowMeans, rowSums

Loading required package: IRanges
Loading required package: XVector
Loading required package: Rsamtools
Loading required package: GenomeInfoDb
Loading required package: GenomicRanges
Loading required package: GenomicAlignments
Loading required package: SummarizedExperiment
Loading required package: Biobase
Welcome to Bioconductor

    Vignettes contain introductory material; view with
    'browseVignettes()'. To cite Bioconductor, see
    'citation("Biobase")', and for packages 'citation("pkgname")'.

> png(filename="/home/ddbj/snapshot/RGM3/R_BC/result/rSFFreader/SffHeader-class.Rd_%03d_medium.png", width=480, height=480)
> ### Name: SffHeader-class
> ### Title: 'SffHeader'
> ### Aliases: class:SffHeader SffHeader-class SffHeader
> ###   show,SffHeader-method detail,SffHeader-method header,SffHeader-method
> ###   header
> ### Keywords: classes
> 
> ### ** Examples
> 
> showClass("SffHeader")
Class "SffHeader" [package "rSFFreader"]

Slots:
             
Name:  header
Class:   list

Known Subclasses: 
Class "SffReads", directly
Class "SffReadsQ", by class "SffReads", distance 2
> 
> ## The readSffHeader can be used to extract header information from one or more sff files:
> sffFiles = c(system.file("extdata", "SmallTorrentTest.sff", package = "rSFFreader"), 
+              system.file("extdata","Small454Test.sff",package="rSFFreader"))
> header <- readSffHeader(sffFiles)
reading header for sff file:/home/ddbj/local/lib64/R/library/rSFFreader/extdata/SmallTorrentTest.sff
reading header for sff file:/home/ddbj/local/lib64/R/library/rSFFreader/extdata/Small454Test.sff
> header
class: SffHeader
> header(header)[[1]]$number_of_reads
[1] 1000
> header(header)[[2]]$number_of_reads
[1] 1000
> 
> ## Header information is also retrieved when using readSff:
> sff <- readSff(sffFiles)
Total number of reads to be read: 2000
reading header for sff file:/home/ddbj/local/lib64/R/library/rSFFreader/extdata/SmallTorrentTest.sff
reading header for sff file:/home/ddbj/local/lib64/R/library/rSFFreader/extdata/Small454Test.sff
reading file:/home/ddbj/local/lib64/R/library/rSFFreader/extdata/SmallTorrentTest.sff
reading file:/home/ddbj/local/lib64/R/library/rSFFreader/extdata/Small454Test.sff
> 
> ## Number of flows for dataset 1 and 2:
> header(sff)[[1]]$number_of_flows_per_read
[1] 640
> header(sff)[[2]]$number_of_flows_per_read
[1] 1600
> 
> 
> 
> 
> 
> 
> dev.off()
null device 
          1 
>