Rbowtie.Rd
The following functions can be used to call the bowtie and bowtie-build binaries.
We recommend to use the QuasR package instead of using bowtie
and bowtie_build
directly. QuasR provides a simpler interface than Rbowtie and covers the whole analysis workflow of typical ultra-high throughput sequencing experiments, starting from the raw sequence reads, over pre-processing and alignment, up to quantification.
bowtie_build(references, outdir, ..., prefix = "index", force = FALSE,
strict = TRUE, execute = TRUE)
bowtie(sequences, index, ..., type = c("single", "paired", "crossbow"),
outfile, force = FALSE, strict = TRUE, execute = TRUE)
bowtie_build_usage()
bowtie_usage()
bowtie_version()
Character vector. The path to the files containing the references for which to build a bowtie index.
Character scalar. The path to the output directory in
which to store the bowtie index. If the directory already exists,
the function function will cast an error, unless
force==TRUE
.
Character scalar. The prefix to use for the bowtie index files.
If type
is either single
or
crossbow
, a character vector of filenames if additional
argument c==FALSE
, otherwise a vector of read sequences. If
type
is paired
, a list of filnames ore sequences of
length 2, where the first list item corresponds to the fist mate
pair sequences, and the second list item to the second mate pair
sequences.
Character scalar. The path to the bowtie index and prefix
to align against, in the form </path/to/index>/<prefix>
.
Character scalar, one in c("single", "paired",
"crossbow")
. If single
, the input sequences are
interpreted as single reads. If paired
, they are supposed
to be mate pair reads and if crossbow
, they are considered
to be Crossbow-style reads.
Character scalar. A path to a files used for the alignment output. If missing, the alignments will be returned as a regular R character vector.
Logical. Force overwriting of outdir
or
outfile
.
Logical. Turn off strict checking of input arguments.
Logical scalar. Whether to execute the assembled shell
command. If FALSE
, return a string with the command.
Additional arguments to be passed on to the binaries. See below for details.
All additional arguments in ...
are interpreted as additional
parameters to be passed on to the binaries. For flags, those are
supposed to be logicals (e.g., quiet=TRUE
will be translated
into --q
, q=TRUE
in -q
, and so on). Parameters
with additional input are supposed to be character or numeric vectors,
where the individual vector elements are collapsed into a single
comma-separated string (e.g., k=2
is translated into k
2
, bmax=100
into --bmax 100
, 3=letters[c(1,2,3)]
into
-3 a,b,c
, and so on). Note that some arguments to the bowtie
binary will be ignored if they are already handled as explicit
function arguments. See the output of bowtie_usage()
and
bowtie_build_usage()
for details about available parameters.
The output generated by calling the binaries. For bowtie_build
this is typically a report of the index generation, for
bowtie
this can be a vector of aligments (if outfile
is
missing), otherwise an empty character scalar.
bowtie_usage()
and bowtie_build_usage()
return the usage
information for the respective binaries.
bowtie_version()
return the bowtie versions information.
Langmead B, Trapnell C, Pop M, Salzberg SL. Ultrafast and memory-efficient alignment of short DNA sequences to the human genome. Genome Biology 10:R25.
Langmead B, Schatz M, Lin J, Pop M, Salzberg SL. Searching for SNPs with cloud computing. Genome Biology 10:R134.
td <- tempdir()
## Building a bowtie index
refs <- dir(system.file(package="Rbowtie", "samples", "refs"),
full=TRUE)
tmp <- bowtie_build(references=refs, outdir=file.path(td, "index"),
force=TRUE)
head(tmp)
#> [1] "Settings:"
#> [2] " Output files: \"/var/folders/24/8k48jl6d249_n_qfxwsl6xvm0000gn/T//Rtmp4wtgdn/index/index.*.ebwt\""
#> [3] " Line rate: 6 (line is 64 bytes)"
#> [4] " Lines per side: 1 (side is 64 bytes)"
#> [5] " Offset rate: 5 (one in 32)"
#> [6] " FTable chars: 10"
dir(file.path(td, "index"))
#> [1] "index.1.ebwt" "index.2.ebwt" "index.3.ebwt" "index.4.ebwt"
#> [5] "index.rev.1.ebwt" "index.rev.2.ebwt"
tmp2 <- bowtie_build(references=refs, outdir=file.path(td,"indexColor"),
force=TRUE, C=TRUE)
#> Warning: running command ''/Users/runner/work/_temp/Library/Rbowtie/bowtie-build' -C '/Users/runner/work/_temp/Library/Rbowtie/samples/refs/chr1.fa','/Users/runner/work/_temp/Library/Rbowtie/samples/refs/chr2.fa','/Users/runner/work/_temp/Library/Rbowtie/samples/refs/chr3.fa' '/var/folders/24/8k48jl6d249_n_qfxwsl6xvm0000gn/T//Rtmp4wtgdn/indexColor/index'' had status 1
dir(file.path(td, "indexColor"))
#> character(0)
head(tmp2)
#> character(0)
## Alignments
reads <- system.file(package="Rbowtie", "samples", "reads", "reads.fastq")
tmp <- bowtie(sequences=reads, index=file.path(td, "index", "index"))
tmp
#> [1] "HWUSI-EAS1513_0012:6:48:5769:946#0/1\t+\tchr1\t818\tTGGAGTTCATGTGAACAGCATATAATTAGGTCTTGCTTTTATTCTAATATGATATATTCTATTAGTTGATATATTTAAACTGCCACAACTTAATGTAATTA\tBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBB\t0\t"
#> [2] "HWUSI-EAS1513_0012:6:48:6908:952#0/1\t+\tchr2\t1132\tAACATAGTGAAGAAACCTCATAGAACCAGAGAAGTTGTGTGTTTCCATGCAGCTAAAATCTTGCTTAATCTGAGAAGCTCAGTTCTATGAAGCCTTACAAC\tBIGGFMONROYYYYY__________ILJJL_____WMQQQPTPPV__b___QQ____W__BBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBB\t0\t"
#> [3] "HWUSI-EAS1513_0012:6:48:8070:953#0/1\t+\tchr1\t7542\tGTCTGTCTAGTGGAGGGGAAGCCGGTGGAGCGATCGTGTGGAAGAGGAAGCGGGATCTGCAGGAGTTACTGCAGAGTATAAAGAAAATCAATGTTGGGGTA\tBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBB\t0\t"
bowtie(sequences=reads, index=file.path(td, "index", "index"),
outfile=file.path(td, "alignments.txt"), best=TRUE, force=TRUE)
readLines(file.path(td, "alignments.txt"))
#> [1] "HWUSI-EAS1513_0012:6:48:5769:946#0/1\t+\tchr1\t818\tTGGAGTTCATGTGAACAGCATATAATTAGGTCTTGCTTTTATTCTAATATGATATATTCTATTAGTTGATATATTTAAACTGCCACAACTTAATGTAATTA\tBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBB\t0\t"
#> [2] "HWUSI-EAS1513_0012:6:48:6908:952#0/1\t+\tchr2\t1132\tAACATAGTGAAGAAACCTCATAGAACCAGAGAAGTTGTGTGTTTCCATGCAGCTAAAATCTTGCTTAATCTGAGAAGCTCAGTTCTATGAAGCCTTACAAC\tBIGGFMONROYYYYY__________ILJJL_____WMQQQPTPPV__b___QQ____W__BBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBB\t0\t"
#> [3] "HWUSI-EAS1513_0012:6:48:8070:953#0/1\t+\tchr1\t7542\tGTCTGTCTAGTGGAGGGGAAGCCGGTGGAGCGATCGTGTGGAAGAGGAAGCGGGATCTGCAGGAGTTACTGCAGAGTATAAAGAAAATCAATGTTGGGGTA\tBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBB\t0\t"
bowtie(sequences=list("TGGGTGGGGTATTCTAGAAATTTCTATTAATCCT",
"TCTGTTCAAGTCAGATGGTCACCAATCTGAAGAC"),
index=file.path(td, "index", "index"), type="paired", c=TRUE)