| Title: | Bioinformatics Library for Integrated Tools |
|---|---|
| Description: | An all-encompassing R toolkit designed to streamline the process of calling various bioinformatics software and then performing data analysis and visualization in R. With 'blit', users can easily integrate a wide array of bioinformatics command line tools into their workflows, leveraging the power of R for sophisticated data manipulation and graphical representation. |
| Authors: | Yun Peng [aut, cre] (ORCID: <https://orcid.org/0000-0003-2801-3332>), Shixiang Wang [aut] (ORCID: <https://orcid.org/0000-0001-9855-7357>), Jia Ding [ctb] (ORCID: <https://orcid.org/0009-0002-2530-0897>), Jennifer Lu [cph] (Author of the included scripts from Kraken2 and KrakenTools libraries), Li Song [cph] (Author of included scripts from TRUST4 library), X. Shirley Liu [cph] (Author of included scripts from TRUST4 library) |
| Maintainer: | Yun Peng <[email protected]> |
| License: | GPL (>= 3) |
| Version: | 0.2.0.9000 |
| Built: | 2026-05-16 08:30:36 UTC |
| Source: | https://github.com/WangLabCSU/blit |
The alleleCount program primarily exists to prevent code duplication
between some other projects, specifically AscatNGS and Battenberg.
allele_counter( hts_file, loci_file, ofile, ..., odir = getwd(), alleleCounter = NULL )allele_counter( hts_file, loci_file, ofile, ..., odir = getwd(), alleleCounter = NULL )
hts_file |
A string of path to sample HTS file. |
loci_file |
A string of path to loci file. |
ofile |
A string of path to the output file. |
... |
<dynamic dots> Additional arguments passed to |
odir |
A string of path to the output directory. |
alleleCounter |
A string of path to |
A command object.
Other commands:
bcftools(),
bedtools(),
bowtie2(),
bwa(),
cellranger(),
conda(),
fastp(),
fastq_pair(),
freebayes(),
gistic2(),
kraken2(),
kraken_tools(),
minimap2(),
perl(),
pyscenic(),
python(),
samtools(),
seqkit(),
snpEff(),
strelka(),
trust4(),
varscan()
micromamba
Manage Environment with micromamba
appmamba(...) install_appmamba(force = FALSE) uninstall_appmamba() appmamba_rc(edit = FALSE)appmamba(...) install_appmamba(force = FALSE) uninstall_appmamba() appmamba_rc(edit = FALSE)
... |
<dynamic dots> Additional arguments passed to
|
force |
A logical value indicating whether to reinstall |
edit |
A logical value indicating whether to open the config file for editing. |
appmamba(): blit utilizes micromamba to manage environments.
This function simply executes the specified micromamba command.
install_appmamba(): Install appmamba (micromamba).
uninstall_appmamba(): Remove appmamba (micromamba).
appmamba_rc(): Get the run commands config file of the micromamba.
install_appmamba() appmamba() appmamba("env", "list") # uninstall_appmamba() # Uninstall the `micromamba`install_appmamba() appmamba() appmamba("env", "list") # uninstall_appmamba() # Uninstall the `micromamba`
arg() is intended for user use, while arg0() is for developers and does
not perform argument validation.
arg(tag, value, indicator = FALSE, lgl2int = FALSE, format = "%s", sep = " ") arg0( tag, value, indicator = FALSE, lgl2int = FALSE, format = "%s", sep = " ", allow_null = FALSE, arg = caller_arg(value), call = caller_call() )arg(tag, value, indicator = FALSE, lgl2int = FALSE, format = "%s", sep = " ") arg0( tag, value, indicator = FALSE, lgl2int = FALSE, format = "%s", sep = " ", allow_null = FALSE, arg = caller_arg(value), call = caller_call() )
tag |
A string specifying argument tag, like "-i", "-o". |
value |
Value passed to the argument. |
indicator |
A logical value specifying whether value should be an
indicator of tag. If |
lgl2int |
A logical value indicates whether transfrom value |
format |
The format of the value, details see |
sep |
A character string used to separate |
allow_null |
A single logical value indicates whether |
arg |
An argument name as a string. This argument will be mentioned in error messages as the input that is at the origin of a problem. |
call |
The execution environment of a currently running function. |
A string.
BCFtools is a program for variant calling and manipulating files in the Variant Call Format (VCF) and its binary counterpart BCF. All commands work transparently with both VCFs and BCFs, both uncompressed and BGZF-compressed.
bcftools(subcmd = NULL, ..., bcftools = NULL)bcftools(subcmd = NULL, ..., bcftools = NULL)
subcmd |
Sub-Command of bcftools. Details see: |
... |
<dynamic dots> Additional arguments passed to |
bcftools |
A string of path to |
A command object.
Other commands:
allele_counter(),
bedtools(),
bowtie2(),
bwa(),
cellranger(),
conda(),
fastp(),
fastq_pair(),
freebayes(),
gistic2(),
kraken2(),
kraken_tools(),
minimap2(),
perl(),
pyscenic(),
python(),
samtools(),
seqkit(),
snpEff(),
strelka(),
trust4(),
varscan()
The bedtools is a powerful toolset for genome arithmetic.
bedtools(subcmd = NULL, ..., bedtools = NULL)bedtools(subcmd = NULL, ..., bedtools = NULL)
subcmd |
Sub-Command of bedtools. Details see: |
... |
<dynamic dots> Additional arguments passed to |
bedtools |
A string of path to |
A command object.
Other commands:
allele_counter(),
bcftools(),
bowtie2(),
bwa(),
cellranger(),
conda(),
fastp(),
fastq_pair(),
freebayes(),
gistic2(),
kraken2(),
kraken_tools(),
minimap2(),
perl(),
pyscenic(),
python(),
samtools(),
seqkit(),
snpEff(),
strelka(),
trust4(),
varscan()
Bowtie 2 is an ultrafast and memory-efficient tool for aligning sequencing reads to long reference sequences. It is particularly good at aligning reads of about 50 up to 100s or 1,000s of characters, and particularly good at aligning to relatively long (e.g. mammalian) genomes. Bowtie 2 indexes the genome with an FM Index to keep its memory footprint small: for the human genome, its memory footprint is typically around 3.2 GB. Bowtie 2 supports gapped, local, and paired-end alignment modes.
bowtie2(index, reads, ofile, ..., bowtie2 = NULL)bowtie2(index, reads, ofile, ..., bowtie2 = NULL)
index |
Path to bowtie2 index prefix (without file extensions). |
reads |
A character vector of FASTQ files used as input to bowtie2. |
ofile |
A string of path to the output sam file. |
... |
<dynamic dots> Additional arguments passed to |
bowtie2 |
A string of path to |
A command object.
Other commands:
allele_counter(),
bcftools(),
bedtools(),
bwa(),
cellranger(),
conda(),
fastp(),
fastq_pair(),
freebayes(),
gistic2(),
kraken2(),
kraken_tools(),
minimap2(),
perl(),
pyscenic(),
python(),
samtools(),
seqkit(),
snpEff(),
strelka(),
trust4(),
varscan()
BWA is a software package that aligns low-divergence sequences to a large reference genome, such as the human genome
bwa(subcmd = NULL, ..., bwa = NULL)bwa(subcmd = NULL, ..., bwa = NULL)
subcmd |
Sub-Command of BWA (e.g., "index", "mem"). |
... |
<dynamic dots> Additional arguments passed to |
bwa |
A string of path to |
A command object.
Other commands:
allele_counter(),
bcftools(),
bedtools(),
bowtie2(),
cellranger(),
conda(),
fastp(),
fastq_pair(),
freebayes(),
gistic2(),
kraken2(),
kraken_tools(),
minimap2(),
perl(),
pyscenic(),
python(),
samtools(),
seqkit(),
snpEff(),
strelka(),
trust4(),
varscan()
## Not run: # Index reference genome bwa("index", "-a", "bwtsw", "reference.fa") |> cmd_run() # Paired-end sequence alignment bwa("mem", "-t", "4", "reference.fa", "read1.fq", "read2.fq") |> cmd_run(stdout = "output.sam") # Single-end alignment (generate sai file) bwa("aln", "-t", "4", "reference.fa", "read.fq") |> cmd_run(stdout = "read.sai") ## End(Not run)## Not run: # Index reference genome bwa("index", "-a", "bwtsw", "reference.fa") |> cmd_run() # Paired-end sequence alignment bwa("mem", "-t", "4", "reference.fa", "read1.fq", "read2.fq") |> cmd_run(stdout = "output.sam") # Single-end alignment (generate sai file) bwa("aln", "-t", "4", "reference.fa", "read.fq") |> cmd_run(stdout = "read.sai") ## End(Not run)
Run cellranger
cellranger(subcmd = NULL, ..., cellranger = NULL)cellranger(subcmd = NULL, ..., cellranger = NULL)
subcmd |
Sub-Command of cellranger. |
... |
<dynamic dots> Additional arguments passed to |
cellranger |
A string of path to |
A command object.
https://www.10xgenomics.com/support/software/cell-ranger/latest
https://www.10xgenomics.com/support/software/cell-ranger/downloads#reference-downloads
Other commands:
allele_counter(),
bcftools(),
bedtools(),
bowtie2(),
bwa(),
conda(),
fastp(),
fastq_pair(),
freebayes(),
gistic2(),
kraken2(),
kraken_tools(),
minimap2(),
perl(),
pyscenic(),
python(),
samtools(),
seqkit(),
snpEff(),
strelka(),
trust4(),
varscan()
## Not run: fastq_dir # 10x raw fastq files directory genome_ref # Please download the transcriptome reference data cellranger( "count", sprintf("--fastqs=%s", fastq_dir), sprintf("--id=%s", basename(fastq_dir)), sprintf("--sample=%s", basename(fastq_dir)), sprintf("--localcores=%s", parallel::detectCores()), sprintf("--transcriptome=%s", genome_ref), sprintf("--chemistry=%s", shQuote("auto")), "--nosecondary" ) ## End(Not run)## Not run: fastq_dir # 10x raw fastq files directory genome_ref # Please download the transcriptome reference data cellranger( "count", sprintf("--fastqs=%s", fastq_dir), sprintf("--id=%s", basename(fastq_dir)), sprintf("--sample=%s", basename(fastq_dir)), sprintf("--localcores=%s", parallel::detectCores()), sprintf("--transcriptome=%s", genome_ref), sprintf("--chemistry=%s", shQuote("auto")), "--nosecondary" ) ## End(Not run)
conda-like environment prefix to the PATH environment variables Set
conda-like environment prefix to the
PATH environment variables.
cmd_conda(...)cmd_conda(...)
... |
Additional arguments passed to |
Schedule expressions to run
cmd_on_start(command, ...) cmd_on_exit(command, ...) cmd_on_fail(command, ...) cmd_on_succeed(command, ...)cmd_on_start(command, ...) cmd_on_exit(command, ...) cmd_on_fail(command, ...) cmd_on_succeed(command, ...)
command |
A |
... |
The expressions input will be captured with
|
cmd_on_start: The command object itself, with the start code updated.
cmd_on_exit: The command object itself, with the exit code updated.
cmd_on_fail: The command object itself, with the failure code updated.
cmd_on_succeed: The command object itself, with the successful code
updated.
cmd_on_start(): define the startup code of the command
cmd_on_exit(): define the exit code of the command
cmd_on_fail(): define the failure code of the command
cmd_on_succeed(): define the successful code of the command
Execute a list of commands
cmd_parallel( ..., stdouts = FALSE, stderrs = FALSE, stdins = NULL, stdout_callbacks = NULL, stderr_callbacks = NULL, timeouts = NULL, threads = NULL, verbose = TRUE )cmd_parallel( ..., stdouts = FALSE, stderrs = FALSE, stdins = NULL, stdout_callbacks = NULL, stderr_callbacks = NULL, timeouts = NULL, threads = NULL, verbose = TRUE )
... |
A list of |
stdouts, stderrs
|
Specifies how the output/error streams of the child process are handled. One of or a list of following values:
For When a single file path is specified, the stdout/stderr of all commands will be merged into this single file. |
stdins |
should the input be diverted? One of or a list of following values:
|
stdout_callbacks, stderr_callbacks
|
One of or a list of following values:
|
timeouts |
Timeout in seconds. Can be a single value or a list, specifying the maximum elapsed time for running the command in the separate process. |
threads |
Number of threads to use. |
verbose |
A single boolean value indicating whether the command execution should be verbose. |
A list of exit status invisiblely.
cmd_run: Run the command.
cmd_help: Print the help document for this command.
cmd_background: Run the command in the background. This function is
provided for completeness. Instead of using this function, we recommend
using cmd_parallel(), which can run multiple commands in the background
while ensuring that all processes are properly cleaned up when the process
exits.
cmd_run( command, stdout = TRUE, stderr = TRUE, stdin = TRUE, stdout_callback = NULL, stderr_callback = NULL, timeout = NULL, spinner = FALSE, verbose = TRUE ) cmd_help( command, stdout = TRUE, stderr = TRUE, stdout_callback = NULL, stderr_callback = NULL, verbose = TRUE ) cmd_background( command, stdout = FALSE, stderr = FALSE, stdin = NULL, verbose = TRUE )cmd_run( command, stdout = TRUE, stderr = TRUE, stdin = TRUE, stdout_callback = NULL, stderr_callback = NULL, timeout = NULL, spinner = FALSE, verbose = TRUE ) cmd_help( command, stdout = TRUE, stderr = TRUE, stdout_callback = NULL, stderr_callback = NULL, verbose = TRUE ) cmd_background( command, stdout = FALSE, stderr = FALSE, stdin = NULL, verbose = TRUE )
command |
A |
stdout, stderr
|
Specifies how the output/error streams of the child process are handled. Possible values include:
For For For When using a |
stdin |
should the input be diverted? Possible values include:
|
stdout_callback, stderr_callback
|
Possible values include:
|
timeout |
Timeout in seconds. This is a limit for the elapsed time running command in the separate process. |
spinner |
Whether to show a reassuring spinner while the process is running. |
verbose |
A single boolean value indicating whether the command execution should be verbose. |
cmd_run: Exit status invisiblely.
cmd_help: The input command invisiblely.
cmd_background: A process object.
Setup the context for the command
cmd_wd(command, wd = NULL) cmd_envvar(command, ..., action = "replace", sep = NULL) cmd_envpath(command, ..., action = "prefix", name = "PATH") cmd_condaenv(command, ..., root = NULL, action = "prefix")cmd_wd(command, wd = NULL) cmd_envvar(command, ..., action = "replace", sep = NULL) cmd_envpath(command, ..., action = "prefix", name = "PATH") cmd_condaenv(command, ..., root = NULL, action = "prefix")
command |
A |
wd |
A string or |
... |
<dynamic dots>:
|
action |
Should the new values |
sep |
A string to separate new and old value when |
name |
A string define the PATH environment variable name. You
can use this to define other |
root |
A string specifying the path to the conda root prefix. If not provided, the function searches for the root in the following order:
|
cmd_wd: The command object itself, with working directory updated.
cmd_envvar: The command object itself, with running environment
variable updated.
cmd_envpath: The command object itself, with running environment
variable specified in name updated.
cmd_condaenv: The command object itself, with running environment
variable PATH updated.
cmd_wd(): define the working directory.
cmd_envvar(): define the environment variables.
cmd_envpath(): define the PATH-like environment variables.
cmd_condaenv(): Set conda-like environment prefix to the PATH
environment variables.
Command is an R6 class used by developers to create new command. It should
not be used by end users.
new()
Create a new Command object.
Command$new(...)
...Additional argument passed into command.
build_command()
Build the command line
Command$build_command(help = FALSE, verbose = TRUE)
helpA boolean value indicating whether to build parameters for help document or not.
verboseA boolean value indicating whether the command execution should be verbose.
envirAn environment used to Execute command.
An atomic character combine the command and parameters.
get_on_start()
Get the command startup code
Command$get_on_start()
A list of quosures.
get_on_exit()
Get the command exit code
Command$get_on_exit()
A list of quosures.
get_on_fail()
Get the command failure code
Command$get_on_fail()
A list of quosures.
get_on_succeed()
Get the command succeessful code
Command$get_on_succeed()
A list of quosures.
print()
Build parameters to run command.
Command$print(indent = NULL)
indentA single integer number giving the space of indent.
The object itself.
clone()
The objects of this class are cloneable with this method.
Command$clone(deep = FALSE)
deepWhether to make a deep clone.
make_command
Run conda
conda(subcmd = NULL, ..., conda = NULL)conda(subcmd = NULL, ..., conda = NULL)
subcmd |
Sub-Command of conda. |
... |
<dynamic dots> Additional arguments passed to |
conda |
A string of path to |
A command object.
Other commands:
allele_counter(),
bcftools(),
bedtools(),
bowtie2(),
bwa(),
cellranger(),
fastp(),
fastq_pair(),
freebayes(),
gistic2(),
kraken2(),
kraken_tools(),
minimap2(),
perl(),
pyscenic(),
python(),
samtools(),
seqkit(),
snpEff(),
strelka(),
trust4(),
varscan()
Invoke a System Command
exec(cmd, ...)exec(cmd, ...)
cmd |
Command to be invoked, as a character string. |
... |
<dynamic dots> Additional arguments passed to |
A command object.
command collectionscmd_run(exec("echo", "$PATH"))cmd_run(exec("echo", "$PATH"))
The fastp is a tool designed to provide ultrafast all-in-one preprocessing
and quality control for FastQ data.
fastp(fq1, ofile1, ..., fq2 = NULL, ofile2 = NULL, fastp = NULL)fastp(fq1, ofile1, ..., fq2 = NULL, ofile2 = NULL, fastp = NULL)
fq1, fq2
|
A string of fastq file path. |
ofile1, ofile2
|
A string of path to the output fastq file. |
... |
<dynamic dots> Additional arguments passed to |
fastp |
A string of path to |
A command object.
Other commands:
allele_counter(),
bcftools(),
bedtools(),
bowtie2(),
bwa(),
cellranger(),
conda(),
fastq_pair(),
freebayes(),
gistic2(),
kraken2(),
kraken_tools(),
minimap2(),
perl(),
pyscenic(),
python(),
samtools(),
seqkit(),
snpEff(),
strelka(),
trust4(),
varscan()
Rewrite paired end fastq files to make sure that all reads have a mate and to separate out singletons.
Usually when you get paired end read files you have two files with a /1 sequence in one and a /2 sequence in the other (or a /f and /r or just two reads with the same ID). However, often when working with files from a third party source (e.g. the SRA) there are different numbers of reads in each file (because some reads fail QC). Spades, bowtie2 and other tools break because they demand paired end files have the same number of reads.
fastq_pair( fq1, fq2, ..., hash_table_size = NULL, max_hash_table_size = NULL, fastq_pair = NULL ) fastq_read_pair(fastq_files)fastq_pair( fq1, fq2, ..., hash_table_size = NULL, max_hash_table_size = NULL, fastq_pair = NULL ) fastq_read_pair(fastq_files)
fq1, fq2
|
A string of fastq file path. |
... |
<dynamic dots> Additional arguments passed to |
hash_table_size |
Size of hash table to use. |
max_hash_table_size |
Maximal hash table size to use. |
fastq_pair |
A string of path to |
fastq_files |
A character of the fastq file paths. |
A command object.
Other commands:
allele_counter(),
bcftools(),
bedtools(),
bowtie2(),
bwa(),
cellranger(),
conda(),
fastp(),
freebayes(),
gistic2(),
kraken2(),
kraken_tools(),
minimap2(),
perl(),
pyscenic(),
python(),
samtools(),
seqkit(),
snpEff(),
strelka(),
trust4(),
varscan()
freebayes is a Bayesian genetic variant detector designed to find small polymorphisms, specifically SNPs (single-nucleotide polymorphisms), indels (insertions and deletions), MNPs (multi-nucleotide polymorphisms), and complex events (composite insertion and substitution events) smaller than the length of a short-read sequencing alignment.
freebayes(ref, input, ofile, ..., freebayes = NULL)freebayes(ref, input, ofile, ..., freebayes = NULL)
ref |
A string of reference file path. |
input |
A string of path to the input bam file. |
ofile |
A string of path to the output vcf file. |
... |
<dynamic dots> Additional arguments passed to |
freebayes |
A string of path to |
A command object.
Other commands:
allele_counter(),
bcftools(),
bedtools(),
bowtie2(),
bwa(),
cellranger(),
conda(),
fastp(),
fastq_pair(),
gistic2(),
kraken2(),
kraken_tools(),
minimap2(),
perl(),
pyscenic(),
python(),
samtools(),
seqkit(),
snpEff(),
strelka(),
trust4(),
varscan()
The GISTIC module identifies regions of the genome that are significantly amplified or deleted across a set of samples. Each aberration is assigned a G-score that considers the amplitude of the aberration as well as the frequency of its occurrence across samples. False Discovery Rate q-values are then calculated for the aberrant regions, and regions with q-values below a user-defined threshold are considered significant. For each significant region, a "peak region" is identified, which is the part of the aberrant region with greatest amplitude and frequency of alteration. In addition, a "wide peak" is determined using a leave-one-out algorithm to allow for errors in the boundaries in a single sample. The "wide peak" boundaries are more robust for identifying the most likely gene targets in the region. Each significantly aberrant region is also tested to determine whether it results primarily from broad events (longer than half a chromosome arm), focal events, or significant levels of both. The GISTIC module reports the genomic locations and calculated q-values for the aberrant regions. It identifies the samples that exhibit each significant amplification or deletion, and it lists genes found in each "wide peak" region.
gistic2(seg, refgene, ..., odir = getwd(), gistic2 = NULL)gistic2(seg, refgene, ..., odir = getwd(), gistic2 = NULL)
seg |
A data.frame of segmented data. |
refgene |
Path to reference genome data input file (REQUIRED, see below for file description). |
... |
<dynamic dots> Additional arguments passed to |
odir |
A string of path to the output directory. |
gistic2 |
A string of path to |
A command object.
Other commands:
allele_counter(),
bcftools(),
bedtools(),
bowtie2(),
bwa(),
cellranger(),
conda(),
fastp(),
fastq_pair(),
freebayes(),
kraken2(),
kraken_tools(),
minimap2(),
perl(),
pyscenic(),
python(),
samtools(),
seqkit(),
snpEff(),
strelka(),
trust4(),
varscan()
These scripts are designed to help Kraken users with downstream analysis of Kraken results.
kraken_tools(script, ..., python = NULL)kraken_tools(script, ..., python = NULL)
script |
Name of the kraken2 script. One of
|
... |
<dynamic dots> Additional arguments passed to |
python |
A string of path to |
A command object.
Other commands:
allele_counter(),
bcftools(),
bedtools(),
bowtie2(),
bwa(),
cellranger(),
conda(),
fastp(),
fastq_pair(),
freebayes(),
gistic2(),
kraken2(),
minimap2(),
perl(),
pyscenic(),
python(),
samtools(),
seqkit(),
snpEff(),
strelka(),
trust4(),
varscan()
Kraken is a taxonomic sequence classifier that assigns taxonomic labels to DNA sequences. Kraken examines the k-mers within a query sequence and uses the information within those k-mers to query a database. That database maps k-mers to the lowest common ancestor (LCA) of all genomes known to contain a given k-mer.
kraken2( reads, ..., ofile = "kraken_output.txt", report = "kraken_report.txt", classified_out = NULL, unclassified_out = NULL, odir = getwd(), kraken2 = NULL )kraken2( reads, ..., ofile = "kraken_output.txt", report = "kraken_report.txt", classified_out = NULL, unclassified_out = NULL, odir = getwd(), kraken2 = NULL )
reads |
A character vector of FASTQ files used as input to Kraken2. Can be one file (single-end) or two files (paired-end). |
... |
<dynamic dots> Additional arguments passed to |
ofile |
A string of path to save kraken2 output. |
report |
A string of path to save kraken2 report. |
classified_out |
A string of path to save classified sequences, which should be a fastq file. |
unclassified_out |
A string of path to save unclassified sequences, which should be a fastq file. |
odir |
A string of path to the output directory. |
kraken2 |
A string of path to |
A command object.
Other commands:
allele_counter(),
bcftools(),
bedtools(),
bowtie2(),
bwa(),
cellranger(),
conda(),
fastp(),
fastq_pair(),
freebayes(),
gistic2(),
kraken_tools(),
minimap2(),
perl(),
pyscenic(),
python(),
samtools(),
seqkit(),
snpEff(),
strelka(),
trust4(),
varscan()
make_command is a helper function used by developers to create function for
a new Command object. It should not be used by end users.
make_command(name, fun, envir = parent.frame())make_command(name, fun, envir = parent.frame())
name |
A string of the function name. |
fun |
A function used to initialize the |
envir |
A environment used to bind the created function. |
A function.
Minimap2 is a versatile sequence alignment program that aligns DNA or mRNA sequences against a large reference database. Typical use cases include: (1) mapping PacBio or Oxford Nanopore genomic reads to the human genome; (2) finding overlaps between long reads with error rate up to ~15%; (3) splice-aware alignment of PacBio Iso-Seq or Nanopore cDNA or Direct RNA reads against a reference genome; (4) aligning Illumina single- or paired-end reads; (5) assembly-to-assembly alignment; (6) full-genome alignment between two closely related species with divergence below ~15%.
minimap2(method, ref, reads, ofile, ..., minimap2 = NULL)minimap2(method, ref, reads, ofile, ..., minimap2 = NULL)
method |
Mapping preset: "map-ont" (single-end long reads) or "sr" (paired-end short reads). |
ref |
Path to the reference genome FASTA file. |
reads |
A character vector of FASTQ files used as input to minimap2. |
ofile |
A string of path to the output SAM file. |
... |
<dynamic dots> Additional arguments passed to |
minimap2 |
A string of path to |
A command object.
Other commands:
allele_counter(),
bcftools(),
bedtools(),
bowtie2(),
bwa(),
cellranger(),
conda(),
fastp(),
fastq_pair(),
freebayes(),
gistic2(),
kraken2(),
kraken_tools(),
perl(),
pyscenic(),
python(),
samtools(),
seqkit(),
snpEff(),
strelka(),
trust4(),
varscan()
Perl is a highly capable, feature-rich programming language with over 36 years of development.
perl(..., perl = NULL)perl(..., perl = NULL)
... |
<dynamic dots> Additional arguments passed to |
perl |
A string of path to |
A command object.
Other commands:
allele_counter(),
bcftools(),
bedtools(),
bowtie2(),
bwa(),
cellranger(),
conda(),
fastp(),
fastq_pair(),
freebayes(),
gistic2(),
kraken2(),
kraken_tools(),
minimap2(),
pyscenic(),
python(),
samtools(),
seqkit(),
snpEff(),
strelka(),
trust4(),
varscan()
Run pyscenic
pyscenic(subcmd = NULL, ..., pyscenic = NULL)pyscenic(subcmd = NULL, ..., pyscenic = NULL)
subcmd |
Sub-Command of pyscenic. |
... |
<dynamic dots> Additional arguments passed to |
pyscenic |
A string of path to |
A command object.
Other commands:
allele_counter(),
bcftools(),
bedtools(),
bowtie2(),
bwa(),
cellranger(),
conda(),
fastp(),
fastq_pair(),
freebayes(),
gistic2(),
kraken2(),
kraken_tools(),
minimap2(),
perl(),
python(),
samtools(),
seqkit(),
snpEff(),
strelka(),
trust4(),
varscan()
Python is a programming language that lets you work quickly and integrate systems more effectively.
python(..., python = NULL)python(..., python = NULL)
... |
<dynamic dots> Additional arguments passed to |
python |
A string of path to |
A command object.
Other commands:
allele_counter(),
bcftools(),
bedtools(),
bowtie2(),
bwa(),
cellranger(),
conda(),
fastp(),
fastq_pair(),
freebayes(),
gistic2(),
kraken2(),
kraken_tools(),
minimap2(),
perl(),
pyscenic(),
samtools(),
seqkit(),
snpEff(),
strelka(),
trust4(),
varscan()
Python is a programming language that lets you work quickly and integrate systems more effectively.
samtools(subcmd = NULL, ..., samtools = NULL)samtools(subcmd = NULL, ..., samtools = NULL)
subcmd |
Sub-Command of samtools. Details see: |
... |
<dynamic dots> Additional arguments passed to |
samtools |
A string of path to |
A command object.
Other commands:
allele_counter(),
bcftools(),
bedtools(),
bowtie2(),
bwa(),
cellranger(),
conda(),
fastp(),
fastq_pair(),
freebayes(),
gistic2(),
kraken2(),
kraken_tools(),
minimap2(),
perl(),
pyscenic(),
python(),
seqkit(),
snpEff(),
strelka(),
trust4(),
varscan()
Run seqkit
seqkit(subcmd = NULL, ..., seqkit = NULL)seqkit(subcmd = NULL, ..., seqkit = NULL)
subcmd |
Sub-Command of seqkit. |
... |
<dynamic dots> Additional arguments passed to |
seqkit |
A string of path to |
A command object.
Other commands:
allele_counter(),
bcftools(),
bedtools(),
bowtie2(),
bwa(),
cellranger(),
conda(),
fastp(),
fastq_pair(),
freebayes(),
gistic2(),
kraken2(),
kraken_tools(),
minimap2(),
perl(),
pyscenic(),
python(),
samtools(),
snpEff(),
strelka(),
trust4(),
varscan()
Genetic variant annotation, and functional effect prediction toolbox. It annotates and predicts the effects of genetic variants on genes and proteins (such as amino acid changes).
snpEff(subcmd = NULL, ..., snpEff = NULL)snpEff(subcmd = NULL, ..., snpEff = NULL)
subcmd |
Sub-Command of snpEff. Details see: |
... |
<dynamic dots> Additional arguments passed to |
snpEff |
A string of path to |
A command object.
Other commands:
allele_counter(),
bcftools(),
bedtools(),
bowtie2(),
bwa(),
cellranger(),
conda(),
fastp(),
fastq_pair(),
freebayes(),
gistic2(),
kraken2(),
kraken_tools(),
minimap2(),
perl(),
pyscenic(),
python(),
samtools(),
seqkit(),
strelka(),
trust4(),
varscan()
Strelka is a fast and accurate small variant caller optimized for analysis of germline variation in small cohorts and somatic variation in tumor/normal sample pairs.
strelka(script, ...)strelka(script, ...)
script |
Name of the Strelka script. One of
|
... |
<dynamic dots> Additional arguments passed to |
strelka |
A string of path to |
A command object.
Other commands:
allele_counter(),
bcftools(),
bedtools(),
bowtie2(),
bwa(),
cellranger(),
conda(),
fastp(),
fastq_pair(),
freebayes(),
gistic2(),
kraken2(),
kraken_tools(),
minimap2(),
perl(),
pyscenic(),
python(),
samtools(),
seqkit(),
snpEff(),
trust4(),
varscan()
TRUST4: immune repertoire reconstruction from bulk and single-cell RNA-seq data
trust4( file1, ref_coordinate, ..., file2 = NULL, mode = NULL, ref_annot = NULL, ofile = NULL, odir = getwd(), trust4 = NULL ) trust4_imgt_annot( species = "Homo_sapien", ..., ofile = "IMGT+C.fa", odir = getwd(), perl = NULL ) trust4_gene_names(imgt_annot, ofile = "bcr_tcr_gene_name.txt", odir = getwd())trust4( file1, ref_coordinate, ..., file2 = NULL, mode = NULL, ref_annot = NULL, ofile = NULL, odir = getwd(), trust4 = NULL ) trust4_imgt_annot( species = "Homo_sapien", ..., ofile = "IMGT+C.fa", odir = getwd(), perl = NULL ) trust4_gene_names(imgt_annot, ofile = "bcr_tcr_gene_name.txt", odir = getwd())
file1 |
Path to bam file or fastq file. |
ref_coordinate |
Path to the fasta file coordinate and sequence of V/D/J/C genes. |
... |
|
file2 |
Path to the second paired-end read fastq file, only used for
|
mode |
One of "bam" or "fastq". If |
ref_annot |
Path to detailed V/D/J/C gene reference file, such as from IMGT database. (default: not used). (recommended). |
ofile |
|
odir |
A string of path to the output directory. |
trust4 |
A string of path to |
species |
Species to extract IMGT annotation, details see https://www.imgt.org//download/V-QUEST/IMGT_V-QUEST_reference_directory/. |
perl |
A string of path to |
imgt_annot |
Path of IMGT annotation file, created via
|
A command object.
Other commands:
allele_counter(),
bcftools(),
bedtools(),
bowtie2(),
bwa(),
cellranger(),
conda(),
fastp(),
fastq_pair(),
freebayes(),
gistic2(),
kraken2(),
kraken_tools(),
minimap2(),
perl(),
pyscenic(),
python(),
samtools(),
seqkit(),
snpEff(),
strelka(),
varscan()
VarScan is a platform-independent software tool developed at the Genome Institute at Washington University to detect variants in NGS data.
varscan(subcmd = NULL, ..., varscan = NULL)varscan(subcmd = NULL, ..., varscan = NULL)
subcmd |
Sub-Command of varscan. Details see: |
... |
<dynamic dots> Additional arguments passed to |
varscan |
A string of path to |
A command object.
Other commands:
allele_counter(),
bcftools(),
bedtools(),
bowtie2(),
bwa(),
cellranger(),
conda(),
fastp(),
fastq_pair(),
freebayes(),
gistic2(),
kraken2(),
kraken_tools(),
minimap2(),
perl(),
pyscenic(),
python(),
samtools(),
seqkit(),
snpEff(),
strelka(),
trust4()