#alignment task bwa { String name String read1 String read2 String output_dir String ref command <<< if [ ! -d ${output_dir}/alignment ];then mkdir ${output_dir}/alignment fi bwa mem -R '@RG\tID:group_n\tLB:library_n\tPL:BGI\tPU:unit1\tSM:${name}' -M -t 10 ${ref} ${read1} ${read2} | \ samtools view -@ 10 -bh -o - | samtools sort -@ 10 -o ${output_dir}/alignment/${name}.sorted.bam samtools index ${output_dir}/alignment/${name}.sorted.bam >>> output { String sorted_bam = "${output_dir}/alignment/${name}.sorted.bam" } } #remove PCR duplicates task markduplicates_genecore { String name String ref String sorted_bam String output_dir command <<< if [ ! -d ${output_dir}/alignment ];then mkdir ${output_dir}/alignment fi gencore -i ${sorted_bam} \ -o ${output_dir}/alignment/${name}.rmdup.bam \ -r ${ref} \ -u UMI \ -j ${output_dir}/alignment/${name}_gencore.json \ -h ${output_dir}/alignment/${name}_gencore.html samtools index ${output_dir}/alignment/${name}.rmdup.bam >>> output { String rmdup_bam = "${output_dir}/alignment/${name}.rmdup.bam" String genecore_json = "${output_dir}/alignment/${name}_gencore.json" } } task markduplicates_picard { String name String ref String sorted_bam String output_dir command <<< if [ ! -d ${output_dir}/alignment ];then mkdir ${output_dir}/alignment fi java -XX:+UseParallelGC -XX:ParallelGCThreads=8 -Xmx12G -jar $PICARD MarkDuplicates \ I=${sorted_bam} \ O=${output_dir}/alignment/${name}.rmdup.bam \ CREATE_INDEX=true \ M=${output_dir}/alignment/${name}.rmdup.metrics.txt \ R=${ref} >>> output { String rmdup_bam = "${output_dir}/alignment/${name}.rmdup.bam" } } # generater mpileup file task generater_mpileup { String name String rmdup_bam String ref String bed String output_dir command <<< samtools mpileup -Bq 20 -Q 20 -f ${ref} -l ${bed} ${rmdup_bam} -o ${output_dir}/alignment/${name}.pileup >>> output { String pileup = "${output_dir}/alignment/${name}.pileup" } } workflow alignment { Boolean run=true String tumor String tumor_r1 String tumor_r2 String? normal String? normal_r1 String? normal_r2 Boolean umi String ref String bed String output_dir if (run) { # 单样本 if (!defined(normal)) { call bwa as bwa_tumor { input: name=tumor, ref=ref, output_dir=output_dir, read1=tumor_r1, read2=tumor_r2 } if (umi) { call markduplicates_genecore as tumor_markduplicates_genecore { input: name=tumor, ref=ref, output_dir=output_dir, sorted_bam=bwa_tumor.sorted_bam, } } if (!umi) { call markduplicates_picard as tumor_markduplicates_picard { input: name=tumor, ref=ref, output_dir=output_dir, sorted_bam=bwa_tumor.sorted_bam, } } } # 双样本 if (defined(normal)) { call bwa as bwa_tumor_control { input: name=tumor, ref=ref, output_dir=output_dir, read1=tumor_r1, read2=tumor_r2 } call bwa as bwa_normal_control { input: name=normal, ref=ref, output_dir=output_dir, read1=normal_r1, read2=normal_r2 } call markduplicates_picard as normal_markduplicates_picard_control { input: name=normal, ref=ref, output_dir=output_dir, sorted_bam=bwa_tumor_control.sorted_bam, } if (umi) { call markduplicates_genecore as tumor_markduplicates_genecore_control { input: name=tumor, ref=ref, output_dir=output_dir, sorted_bam=bwa_normal_control.sorted_bam, } } if (!umi) { call markduplicates_picard as tumor_markduplicates_picard_control { input: name=tumor, ref=ref, output_dir=output_dir, sorted_bam=bwa_normal_control.sorted_bam, } } } } output { String tumor_sorted_bam = "${output_dir}/alignment/${tumor}.sorted.bam" String tumor_rmdup_bam = "${output_dir}/alignment/${tumor}.rmdup.bam" String tumor_pileup = "${output_dir}/alignment/${tumor}.pileup" String normal_sorted_bam = "${output_dir}/alignment/${normal}.sorted.bam" String normal_rmdup_bam = "${output_dir}/alignment/${normal}.rmdup.bam" String normal_pileup = "${output_dir}/alignment/${normal}.pileup" String genecore_json = "${output_dir}/alignment/${tumor}_gencore.json" } }