View source: R/utils_imports.R
fimport | R Documentation |
Wraps around ORFik file format loaders and rtracklayer::import and tries to speed up loading with the use of data.table. Supports gzip, gz, bgz compression formats. Also safer chromosome naming with the argument chrStyle
fimport(path, chrStyle = NULL, param = NULL, strandMode = 0)
path |
a character path to file (1 or 2 files), or data.table with 2 colums(forward&reverse) or a GRanges/Galignment/GAlignmentPairs object etc. If it is ranged object it will presume to be already loaded, so will return the object as it is, updating the seqlevelsStyle if given. |
chrStyle |
a GRanges object, TxDb, FaFile,
, a |
param |
By default (i.e. |
strandMode |
numeric, default 0. Only used for paired end bam files. One of (0: strand = *, 1: first read of pair is +, 2: first read of pair is -). See ?strandMode. Note: Sets default to 0 instead of 1, as readGAlignmentPairs uses 1. This is to guarantee hits, but will also make mismatches of overlapping transcripts in opposite directions. |
NOTE: For wig/bigWig files you can send in 2 files, so that it automatically merges forward and reverse stranded objects. You can also just send 1 wig/bigWig file, it will then have "*" as strand.
a GAlignments
/GRanges
object,
depending on input.
Other utils:
bedToGR()
,
convertToOneBasedRanges()
,
export.bed12()
,
export.bigWig()
,
export.fstwig()
,
export.wiggle()
,
findFa()
,
fread.bed()
,
optimizeReads()
,
readBam()
,
readBigWig()
,
readWig()
bam_file <- system.file("extdata/Danio_rerio_sample", "ribo-seq.bam", package = "ORFik")
fimport(bam_file)
# Certain chromosome naming
fimport(bam_file, "NCBI")
# Paired end bam strandMode 1:
fimport(bam_file, strandMode = 1)
# (will have no effect in this case, since it is not paired end)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.