if (!requireNamespace("BiocManager", quietly = TRUE))
install.packages("BiocManager")
BiocManager::install("RcwlPipelines")
The development version is also available to download from Github.
BiocManager::install("hubentu/RcwlPipelines")
library(RcwlPipelines)
The R scripts to build the CWL tools and pipelines based on the Rcwl
package are stored in the “tools” and “pipelines” folder respectively.
The function cwlTools
can be used to collect the available scripts.
The cachePath
can be your existing cache directory or a new folder.
tools <- cwlTools(cachePath = tempdir())
tools
#> class: BiocFileCache
#> bfccache: /tmp/RtmpalrYxr
#> bfccount: 36
#> For more information see: bfcinfo() or bfcquery()
The full paths can be pulled from the “fpath” column.
library(dplyr)
bfcinfo(tools) %>% select(rname, fpath)
#> # A tibble: 36 x 2
#> rname fpath
#> <chr> <chr>
#> 1 HTseq /tmp/RtmpbYdKVx/Rinst6deb5d0870f6/RcwlPipelines/tools/HTseq.R
#> 2 STAR /tmp/RtmpbYdKVx/Rinst6deb5d0870f6/RcwlPipelines/tools/STAR.R
#> 3 blastn /tmp/RtmpbYdKVx/Rinst6deb5d0870f6/RcwlPipelines/tools/blastn…
#> 4 bowtie2_ali… /tmp/RtmpbYdKVx/Rinst6deb5d0870f6/RcwlPipelines/tools/bowtie…
#> 5 bowtie2_bui… /tmp/RtmpbYdKVx/Rinst6deb5d0870f6/RcwlPipelines/tools/bowtie…
#> 6 bwa /tmp/RtmpbYdKVx/Rinst6deb5d0870f6/RcwlPipelines/tools/bwa.R
#> 7 bwa_index /tmp/RtmpbYdKVx/Rinst6deb5d0870f6/RcwlPipelines/tools/bwa_in…
#> 8 cutadapt /tmp/RtmpbYdKVx/Rinst6deb5d0870f6/RcwlPipelines/tools/cutada…
#> 9 fastQC /tmp/RtmpbYdKVx/Rinst6deb5d0870f6/RcwlPipelines/tools/fastQC…
#> 10 featureCoun… /tmp/RtmpbYdKVx/Rinst6deb5d0870f6/RcwlPipelines/tools/featur…
#> # … with 26 more rows
We can develop a pipline by utilizing the available tools. For
example, a simple alignment pipelines with mapping and marking
duplicates can be built from the tools
.
First, we load the required tools, bwa, samtools and picard markduplicates.
scripts <- bfcinfo(tools) %>%
filter(rname %in% c("bwa",
"samtools_samTobam",
"samtools_sortBam",
"samtools_index",
"markdup")) %>%
pull(rpath)
invisible(sapply(scripts, source))
Next, we define the input parameters.
p1 <- InputParam(id = "threads", type = "int")
p2 <- InputParam(id = "RG", type = "string")
p3 <- InputParam(id = "Ref", type = "string")
p4 <- InputParam(id = "FQ1", type = "File")
p5 <- InputParam(id = "FQ2", type = "File?")
Then we define the pipeline steps, from raw fastqs to duplicates marked alignments.
## bwa
s1 <- Step(id = "bwa", run = bwa,
In = list(threads = "threads",
RG = "RG",
Ref = "Ref",
FQ1 = "FQ1",
FQ2 = "FQ2"))
## sam to bam
s2 <- Step(id = "sam2bam", run = sam2bam,
In = list(sam = "bwa/sam"))
## sort bam
s3 <- Step(id = "sortBam", run = sortBam,
In = list(bam = "sam2bam/bam"))
## mark duplicates
s4 <- Step(id = "markdup", run = markdup,
In = list(ibam = "sortBam/sbam",
obam = list(
valueFrom="$(inputs.ibam.nameroot).mdup.bam"),
matrix = list(
valueFrom="$(inputs.ibam.nameroot).markdup.txt")))
## index bam
s5 <- Step(id = "idxBam", run = samtools_index,
In = list(bam = "markdup/mBam"))
Last, we define the outputs and connect the steps to a new pipeline.
req1 <- list(class = "StepInputExpressionRequirement")
req2 <- list(class = "InlineJavascriptRequirement")
## outputs
o1 <- OutputParam(id = "Bam", type = "File", outputSource = "markdup/mBam")
o2 <- OutputParam(id = "Idx", type = "File", outputSource = "idxBam/idx")
## stepParam
Align <- cwlStepParam(requirements = list(req1, req2),
inputs = InputParamList(p1, p2, p3, p4, p5),
outputs = OutputParamList(o1, o2))
## build pipeline
Align <- Align + s1 + s2 + s3 + s4 + s5
The pipeline is ready for use. We can plot the pipeline with
plotCWL
from the Rcwl
package.
plotCWL(Align)
There are mainly 4 pipelines are collected in this package. Here is a brief introduction to these pipelines. More pipelines and tools are expected to be included in the future.
The pipeline can be used to preprocess DNA sequences in fastq format. It can take paired fastqs, read groups from multiple batches as input.
data(alignMerge)
inputs(alignMerge)
#> inputs:
#> idBam (string):
#> RG (string[]):
#> threads (int):
#> Ref (File):
#> FQ1s (File[]):
#> FQ2s (File[]):
The pipeline includes two steps and several jobs will be run in each step.
bwaAlign
: bwa alignment by read groups.runs(runs(alignMerge)[[1]])
#> List of length 4
#> names(4): bwa sam2bam sortBam idxBam
bwa
: To align fastqs and read groups to reference genome with bwa
.sam2bam
: To convert the alignments in “sam” format to “bam”
format with samtools
.sortBam
: To sort the “bam” file by coordinates with samtools
.idxBam
: To index “bam” file with samtools
.mergeBamDup
: Merge by samples and markduplicates.runs(runs(alignMerge)[[2]])
#> List of length 4
#> names(4): mergeBam markdup samtools_index samtools_flagstat
mergeBam
: To merge bam files from multiple batches with picard
.markdup
: To mark duplicates with picard
.samtools_index
: To index bam file with samtools
.samtools_flagstat
: To summarize flags in bam with samtools
.The final bam files after duplicates marked, bam index, duplicates matrix, and flag statistics summary will be in the output folder.
outputs(alignMerge)
#> outputs:
#> oBam:
#> type: File
#> outputSource: mergeBamDup/oBam
#> matrix:
#> type: File
#> outputSource: mergeBamDup/matrix
#> Idx:
#> type: File
#> outputSource: mergeBamDup/Idx
#> stat:
#> type: File
#> outputSource: mergeBamDup/stat
Here you can find an example to run the pipeline.
https://hubentu.github.io/others/Rcwl/application.html#dnaseq-alignment-pipeline
The pipeline was built with reads quality summary, STAR
alignment,
quantification by featureCounts
and RSeQC
quality control. Here
are the inputs.
data(rnaseq_Sf)
inputs(rnaseq_Sf)
#> inputs:
#> in_seqfiles (File[]):
#> in_prefix (string):
#> in_genomeDir (Directory):
#> in_GTFfile (File):
#> in_runThreadN (int): 1
The pipeline includes 6 steps.
fastqc
: To run quality summary for raw fastqs with fastqc
.STAR
: To align fastqs with STAR
.samtools_index
: To index aligned bam file.samtools_flagstat
: To summary alignment flags.featureCounts
: To quantify gene abundances.RSeQC
: Several steps included.gtfToGenePred
: To convert GTF annotation to “genePred” format.genePredToBed
: To convert “genePred” annotation to “bed” format.r_distribution
: To run reads distribution over genome features.gCoverage
: To summarize read coverage over gene body.The outputs and logs from alignment, quantification and QC steps are
collected together to the output folder. A final QC report could be
generated by multiqc
, which is also available in the data package.
An example to run the pipeline.
https://hubentu.github.io/others/Rcwl/application.html#rnaseq-pipeline
The Multi-Center Mutation Calling in Multiple Cancers project (MC3)
pipeline was developed by the TCGA group, which was implemented with
CWL and available at https://github.com/OpenGenomics/mc3. Two major
steps, variant calling and merging VCFs to MAF, was imported to the
mc3
pipeline in this package.
The steps of the pipeline was built on the CWL files from its github repository, which were also contained in the package. Thereforce, we need the load the pipleine by sourcing it from the script.
The pipeline requires a paired of tumor/normal BAM files and genomic annotation files as inputs.
bfcinfo(tools) %>% filter(rname == "mc3") %>% pull(rpath) %>% source
inputs(mc3)
#> inputs:
#> tumorID (string):
#> normalID (string):
#> tumor (File):
#> normal (File):
#> bed_file (File?):
#> centromere (File):
#> cosmic (File):
#> dbsnp (File):
#> refFasta (File):
#> vepData (Directory):
Two steps are included.
call_variants
: To call variants by 7 pipelines.covert
: To merge VCFs and convert to MAFThe merged VCF and converted MAF files will be collected to the output folder.
outputs(mc3)
#> outputs:
#> outmaf:
#> type: File
#> outputSource: convert/outmaf
#> outvcf:
#> type: File
#> outputSource: convert/vepvcf
Here is an examples to run the pipeline.
https://hubentu.github.io/others/Rcwl/application.html#mc3-somatic-variant-calling-pipeline
The GATK4 best practice pipeline for germline variant calling was
implemented with Workflow Description Language (WDL). We wrapped the
WDL pipeline into 3 steps with Rcwl
. The details of the pipeline can
be find here:
https://software.broadinstitute.org/gatk/best-practices/workflow?id=11145
GAlign
GATK alignment.The fastqs, sample information and customized json files for WDL are
required as inputs. Multiple steps will run in this step, including
bwa
alignment, mark duplicates and base quality recalibration. GATK
ready BAM files will be collected to the output directory.
hapCall
HaplotypeCaller.The GATK ready BAM and customized json files are inputs in this step. The local paths of GATK bundle files are required to be modified in your json file. A “gVCF” files will be generated.
jdCall
Joint variant discoveryThis step will combine the “gVCF” files and then call germline variants in all samples. The paths of the local bundle files are also required to be changed in the json template file. The final VCF file of germline variants will be collected.
An example to run the pipeline.
https://hubentu.github.io/others/Rcwl/application.html#gatk4-germline-variant-calling-pipeline
sessionInfo()
#> R version 3.6.0 (2019-04-26)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 18.04.2 LTS
#>
#> Matrix products: default
#> BLAS: /home/biocbuild/bbs-3.9-bioc/R/lib/libRblas.so
#> LAPACK: /home/biocbuild/bbs-3.9-bioc/R/lib/libRlapack.so
#>
#> locale:
#> [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
#> [3] LC_TIME=en_US.UTF-8 LC_COLLATE=C
#> [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
#> [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
#> [9] LC_ADDRESS=C LC_TELEPHONE=C
#> [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
#>
#> attached base packages:
#> [1] parallel stats4 stats graphics grDevices utils datasets
#> [8] methods base
#>
#> other attached packages:
#> [1] dplyr_0.8.1 RcwlPipelines_1.0.4 BiocFileCache_1.8.0
#> [4] dbplyr_1.4.2 Rcwl_1.0.4 S4Vectors_0.22.0
#> [7] BiocGenerics_0.30.0 yaml_2.2.0 BiocStyle_2.12.0
#>
#> loaded via a namespace (and not attached):
#> [1] httr_1.4.0 viridis_0.5.1 tidyr_0.8.3
#> [4] bit64_0.9-7 jsonlite_1.6 viridisLite_0.3.0
#> [7] shiny_1.3.2 assertthat_0.2.1 BiocManager_1.30.4
#> [10] blob_1.1.1 base64url_1.4 progress_1.2.2
#> [13] pillar_1.4.1 RSQLite_2.1.1 backports_1.1.4
#> [16] glue_1.3.1 downloader_0.4 digest_0.6.19
#> [19] RColorBrewer_1.1-2 promises_1.0.1 checkmate_1.9.3
#> [22] colorspace_1.4-1 htmltools_0.3.6 httpuv_1.5.1
#> [25] XML_3.98-1.20 pkgconfig_2.0.2 bookdown_0.11
#> [28] DiagrammeR_1.0.1 purrr_0.3.2 xtable_1.8-4
#> [31] scales_1.0.0 brew_1.0-6 later_0.8.0
#> [34] BiocParallel_1.18.0 tibble_2.1.3 ggplot2_3.2.0
#> [37] influenceR_0.1.0 withr_2.1.2 lazyeval_0.2.2
#> [40] cli_1.1.0 rgexf_0.15.3 magrittr_1.5
#> [43] crayon_1.3.4 mime_0.7 memoise_1.1.0
#> [46] evaluate_0.14 fansi_0.4.0 Rook_1.1-1
#> [49] tools_3.6.0 data.table_1.12.2 prettyunits_1.0.2
#> [52] hms_0.4.2 stringr_1.4.0 munsell_0.5.0
#> [55] compiler_3.6.0 rlang_0.3.4 grid_3.6.0
#> [58] rstudioapi_0.10 rappdirs_0.3.1 htmlwidgets_1.3
#> [61] visNetwork_2.0.7 igraph_1.2.4.1 rmarkdown_1.13
#> [64] gtable_0.3.0 curl_3.3 DBI_1.0.0
#> [67] R6_2.4.0 gridExtra_2.3 knitr_1.23
#> [70] zeallot_0.1.0 utf8_1.1.4 bit_1.1-14
#> [73] readr_1.3.1 stringi_1.4.3 Rcpp_1.0.1
#> [76] vctrs_0.1.0 batchtools_0.9.11 tidyselect_0.2.5
#> [79] xfun_0.7