readRibodata: Reads ribosomal and (optionally) rna data from alignment...

Description Usage Arguments Details Value Author(s) Examples

View source: R/readRibodata.R

Description

Reads BAM files, or flat text files (which may be compressed) containing strand, transcript name, start and sequence information for each alignment.

Usage

1
2
readRibodata(riboFiles, rnaFiles, columns = c(strand = 1, seqname = 2,
start = 3, sequence = 4), zeroIndexed = FALSE, header = FALSE, replicates, seqnames)

Arguments

riboFiles

Filenames of ribosomal alignments.

rnaFiles

Filenames of RNA alignments.

columns

Columns of alignment files containing strand, transcript (seqname) name, start of alignment, and sequence. Ignored for BAM files.

zeroIndexed

Are the alignments zero-indexed (i.e., the first base in a sequence is 0). Defaults to FALSE.

header

Does the (non-bam) alignment file have a header line? Defaults to FALSE.

replicates

Replicate information for the files.

seqnames

Transcript (seqname) names to be read into the object.

Details

If a filename ends in '.bam' or '.BAM' it will be assumed that the file is a BAM file and processed accordingly. Otherwise, it will be treated as a flat text file.

If using text files these should contain (as a minimum), the sequence of the read, the name of the sequence to which the read aligns, the strand to which it aligns, and the starting position of alignment. A Bowtie alignment (note that Bowtie, rather than Bowtie2, is recommended for short reads, which ribosome footprints are) using the option “–suppress 1,6,7,8” will generate this minimal data. It is by default assumed that the data are generated in this way, and the default columns specification for the default readRibodata function (see below) reflects this.

Value

riboData object.

Author(s)

Thomas J. Hardcastle

Examples

1
2
3
4
5
6
7
datadir <- system.file("extdata", package = "riboSeqR")
ribofiles <- paste(datadir, 
                   "/chlamy236_plus_deNovo_plusOnly_Index", c(17,3,5,7), sep = "")
rnafiles <- paste(datadir, 
                  "/chlamy236_plus_deNovo_plusOnly_Index", c(10,12,14,16), sep = "")

riboDat <- readRibodata(ribofiles, rnafiles, replicates = c("WT", "WT", "M", "M"))

Example output

Loading required package: GenomicRanges
Loading required package: stats4
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, sd, var, xtabs

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

    Filter, Find, Map, Position, Reduce, anyDuplicated, append,
    as.data.frame, basename, cbind, colMeans, colSums, colnames,
    dirname, 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,
    rowMeans, rowSums, rownames, sapply, setdiff, sort, table, tapply,
    union, unique, unsplit, which, which.max, which.min

Loading required package: S4Vectors

Attaching package: 'S4Vectors'

The following object is masked from 'package:base':

    expand.grid

Loading required package: IRanges
Loading required package: GenomeInfoDb
Loading required package: abind
Reading ribosomal files.......done!
Reading rna files.......done!

riboSeqR documentation built on Nov. 8, 2020, 8:23 p.m.