Uses minimap2 to align sequences agains a reference databse.
Uses options '-ax splice -t 12 -k14 –secondary=no fa_file
fq_in
'
Usage
minimap2_align(
config,
fa_file,
fq_in,
annot,
outdir,
minimap2 = NA,
k8 = NA,
samtools = NA,
prefix = NULL,
threads = 1
)
Arguments
- config
Parsed list of FLAMES config file
- fa_file
Path to the fasta file used as a reference database for alignment
- fq_in
File path to the fastq file used as a query sequence file
- annot
Genome annotation file used to create junction bed files
- outdir
Output folder
- minimap2
Path to minimap2 binary
- k8
Path to the k8 Javascript shell binary
- samtools
path to the samtools binary, required for large datasets since
Rsamtools
does not supportCSI
indexing- prefix
String, the prefix (e.g. sample name) for the outputted BAM file
- threads
Integer, threads for minimap2 to use, see minimap2 documentation for details, FLAMES will try to detect cores if this parameter is not provided.
Examples
temp_path <- tempfile()
bfc <- BiocFileCache::BiocFileCache(temp_path, ask = FALSE)
file_url <- 'https://raw.githubusercontent.com/OliverVoogd/FLAMESData/master/data'
fastq1 <- bfc[[names(BiocFileCache::bfcadd(bfc, 'Fastq1', paste(file_url, 'fastq/sample1.fastq.gz', sep = '/')))]]
genome_fa <- bfc[[names(BiocFileCache::bfcadd(bfc, 'genome.fa', paste(file_url, 'SIRV_isoforms_multi-fasta_170612a.fasta', sep = '/')))]]
annotation <- bfc[[names(BiocFileCache::bfcadd(bfc, 'annot.gtf', paste(file_url, 'SIRV_isoforms_multi-fasta-annotation_C_170612a.gtf', sep = '/')))]]
outdir <- tempfile()
dir.create(outdir)
minimap2_align(
config = jsonlite::fromJSON(
system.file("extdata", "config_sclr_nanopore_3end.json", package = 'FLAMES')
),
fa_file = genome_fa,
fq_in = fastq1,
annot = annotation,
outdir = outdir
)
#> 02:13:33 AM Fri Oct 25 2024 minimap2_align
#> total mapped primary
#> /tmp/RtmpoS7Kzz/file1f6d66ff3a3f/align2genome.bam 2503 1772 2500
#> secondary
#> /tmp/RtmpoS7Kzz/file1f6d66ff3a3f/align2genome.bam 0