run_kallisto | R Documentation |
Runs the kallisto quant tool
run_kallisto(
input1 = NULL,
input2 = NULL,
index = NULL,
sample.name = NULL,
fusion = NULL,
strandedness = NULL,
out.dir = NULL,
threads = 10,
bootstrap = 100,
fragment.length = 300,
std.dev = 25,
parallel = FALSE,
cores = 4,
execute = TRUE,
kallisto = "",
version = FALSE
)
input1 |
List of the paths to files containing to the forward reads, required |
input2 |
List of the paths to files containing to the reverse reads, required for paired end sequence data |
index |
Path to the reference transcriptome kallisto index, required |
sample.name |
List of the sample names, required |
fusion |
Search for fusions used for Pizzly |
strandedness |
Strand spcific reads, values are "first" or "second" |
out.dir |
Name of the directory from the Kallisto output. If NULL, which is the default, a directory named "kallisto" is created in the current working directory. |
threads |
Number of threads for kallisto to use, default set to 10 |
bootstrap |
Number of bootstrap samples, default set to 100 |
fragment.length |
Estimated average insert fragment size, must be set for single end sequence data |
std.dev |
Estimated standard deviation of insert fragment size, must be set for single end sequence data |
parallel |
Run in parallel, default set to FALSE |
cores |
Number of cores/threads to use for parallel processing, default set to 4 |
execute |
Whether to execute the commands or not, default set to TRUE |
kallisto |
Path to the kallisto program, required |
version |
Returns the version number |
A list with the kallisto commands
## Not run:
trimmed_reads_dir <- "trimmed_reads"
mate1 <- list.files(path = trimmed_reads_dir, pattern = "*_R1_001.fastq$", full.names = TRUE)
mate2 <- list.files(path = trimmed_reads_dir, pattern = "*_R2_001.fastq$", full.names = TRUE)
sample_names <- unlist(lapply(strsplit(list.files(path = trimmed_reads_dir,
pattern = "*_R1_001.fastq$",
full.names = FALSE),"_"), `[[`, 1))
index <- "/export/buzz1/Genome/Homo_sapiens/Ensembl/GRCH38_p7/Sequence/Transcriptome/
KallistoIndex/GRCh38_p7.kall"
strandedness <- "second"
# Paired end
kallisto.cmds <- run_kallisto(input1 = mate1,
input2 = mate2,
index = transcriptome,
sample.name = sample.names,
strandedness = strandedness,
out.dir = kalisto.results.dir,
kallisto = "/software/kallisto_v0.45.1/kallisto")
# Single end
kallisto.cmds <- run_kallisto(input1 = mate1,
index = transcriptome,
sample.name = sample.names,
strandedness = strandedness,
fragment.length = 300,
std.dev = 25,
out.dir = kalisto.results.dir,
kallisto = "/software/kallisto_v0.45.1/kallisto")
## End(Not run)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.