-
Notifications
You must be signed in to change notification settings - Fork 705
/
rnaseq.nf
executable file
·872 lines (784 loc) · 34.9 KB
/
rnaseq.nf
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
/*
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
VALIDATE INPUTS
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
*/
def valid_params = [
aligners : ['star_salmon', 'star_rsem', 'hisat2'],
pseudoaligners : ['salmon'],
rseqc_modules : ['bam_stat', 'inner_distance', 'infer_experiment', 'junction_annotation', 'junction_saturation', 'read_distribution', 'read_duplication', 'tin']
]
def summary_params = NfcoreSchema.paramsSummaryMap(workflow, params)
// Validate input parameters
WorkflowRnaseq.initialise(params, log, valid_params)
// Check input path parameters to see if they exist
checkPathParamList = [
params.input, params.multiqc_config,
params.fasta, params.transcript_fasta, params.additional_fasta,
params.gtf, params.gff, params.gene_bed,
params.ribo_database_manifest, params.splicesites,
params.star_index, params.hisat2_index, params.rsem_index, params.salmon_index
]
for (param in checkPathParamList) { if (param) { file(param, checkIfExists: true) } }
// Check mandatory parameters
if (params.input) { ch_input = file(params.input) } else { exit 1, 'Input samplesheet not specified!' }
// Check rRNA databases for sortmerna
if (params.remove_ribo_rna) {
ch_ribo_db = file(params.ribo_database_manifest, checkIfExists: true)
if (ch_ribo_db.isEmpty()) {exit 1, "File provided with --ribo_database_manifest is empty: ${ch_ribo_db.getName()}!"}
}
// Check if file with list of fastas is provided when running BBSplit
if (!params.skip_bbsplit && !params.bbsplit_index && params.bbsplit_fasta_list) {
ch_bbsplit_fasta_list = file(params.bbsplit_fasta_list, checkIfExists: true)
if (ch_bbsplit_fasta_list.isEmpty()) {exit 1, "File provided with --bbsplit_fasta_list is empty: ${ch_bbsplit_fasta_list.getName()}!"}
}
// Check alignment parameters
def prepareToolIndices = []
if (!params.skip_bbsplit) { prepareToolIndices << 'bbsplit' }
if (!params.skip_alignment) { prepareToolIndices << params.aligner }
if (params.pseudo_aligner) { prepareToolIndices << params.pseudo_aligner }
// Get RSeqC modules to run
def rseqc_modules = params.rseqc_modules ? params.rseqc_modules.split(',').collect{ it.trim().toLowerCase() } : []
if (params.bam_csi_index) {
for (rseqc_module in ['read_distribution', 'inner_distance', 'tin']) {
if (rseqc_modules.contains(rseqc_module)) {
rseqc_modules.remove(rseqc_module)
}
}
}
// Stage dummy file to be used as an optional input where required
ch_dummy_file = file("$projectDir/assets/dummy_file.txt", checkIfExists: true)
// Check if an AWS iGenome has been provided to use the appropriate version of STAR
def is_aws_igenome = false
if (params.fasta && params.gtf) {
if ((file(params.fasta).getName() - '.gz' == 'genome.fa') && (file(params.gtf).getName() - '.gz' == 'genes.gtf')) {
is_aws_igenome = true
}
}
/*
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
CONFIG FILES
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
*/
ch_multiqc_config = Channel.fromPath("$projectDir/assets/multiqc_config.yml", checkIfExists: true)
ch_multiqc_custom_config = params.multiqc_config ? Channel.fromPath( params.multiqc_config, checkIfExists: true ) : Channel.empty()
ch_multiqc_logo = params.multiqc_logo ? Channel.fromPath( params.multiqc_logo, checkIfExists: true ) : Channel.empty()
ch_multiqc_custom_methods_description = params.multiqc_methods_description ? file(params.multiqc_methods_description, checkIfExists: true) : file("$projectDir/assets/methods_description_template.yml", checkIfExists: true)
// Header files for MultiQC
ch_pca_header_multiqc = file("$projectDir/assets/multiqc/deseq2_pca_header.txt", checkIfExists: true)
ch_clustering_header_multiqc = file("$projectDir/assets/multiqc/deseq2_clustering_header.txt", checkIfExists: true)
ch_biotypes_header_multiqc = file("$projectDir/assets/multiqc/biotypes_header.txt", checkIfExists: true)
/*
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
IMPORT LOCAL MODULES/SUBWORKFLOWS
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
*/
//
// MODULE: Loaded from modules/local/
//
include { BEDTOOLS_GENOMECOV } from '../modules/local/bedtools_genomecov'
include { DESEQ2_QC as DESEQ2_QC_STAR_SALMON } from '../modules/local/deseq2_qc'
include { DESEQ2_QC as DESEQ2_QC_RSEM } from '../modules/local/deseq2_qc'
include { DESEQ2_QC as DESEQ2_QC_SALMON } from '../modules/local/deseq2_qc'
include { DUPRADAR } from '../modules/local/dupradar'
include { MULTIQC } from '../modules/local/multiqc'
include { MULTIQC_CUSTOM_BIOTYPE } from '../modules/local/multiqc_custom_biotype'
include { UMITOOLS_PREPAREFORRSEM as UMITOOLS_PREPAREFORSALMON } from '../modules/local/umitools_prepareforrsem.nf'
include { MULTIQC_TSV_FROM_LIST as MULTIQC_TSV_FAIL_MAPPED } from '../modules/local/multiqc_tsv_from_list'
include { MULTIQC_TSV_FROM_LIST as MULTIQC_TSV_FAIL_TRIMMED } from '../modules/local/multiqc_tsv_from_list'
include { MULTIQC_TSV_FROM_LIST as MULTIQC_TSV_STRAND_CHECK } from '../modules/local/multiqc_tsv_from_list'
//
// SUBWORKFLOW: Consisting of a mix of local and nf-core/modules
//
include { INPUT_CHECK } from '../subworkflows/local/input_check'
include { PREPARE_GENOME } from '../subworkflows/local/prepare_genome'
include { ALIGN_STAR } from '../subworkflows/local/align_star'
include { QUANTIFY_RSEM } from '../subworkflows/local/quantify_rsem'
include { QUANTIFY_SALMON as QUANTIFY_STAR_SALMON } from '../subworkflows/local/quantify_salmon'
include { QUANTIFY_SALMON as QUANTIFY_SALMON } from '../subworkflows/local/quantify_salmon'
/*
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
IMPORT NF-CORE MODULES/SUBWORKFLOWS
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
*/
//
// MODULE: Installed directly from nf-core/modules
//
include { CAT_FASTQ } from '../modules/nf-core/cat/fastq/main'
include { BBMAP_BBSPLIT } from '../modules/nf-core/bbmap/bbsplit/main'
include { SAMTOOLS_SORT } from '../modules/nf-core/samtools/sort/main'
include { PRESEQ_LCEXTRAP } from '../modules/nf-core/preseq/lcextrap/main'
include { QUALIMAP_RNASEQ } from '../modules/nf-core/qualimap/rnaseq/main'
include { SORTMERNA } from '../modules/nf-core/sortmerna/main'
include { STRINGTIE_STRINGTIE } from '../modules/nf-core/stringtie/stringtie/main'
include { SUBREAD_FEATURECOUNTS } from '../modules/nf-core/subread/featurecounts/main'
include { CUSTOM_DUMPSOFTWAREVERSIONS } from '../modules/nf-core/custom/dumpsoftwareversions/main'
//
// SUBWORKFLOW: Consisting entirely of nf-core/modules
//
include { FASTQ_SUBSAMPLE_FQ_SALMON } from '../subworkflows/nf-core/fastq_subsample_fq_salmon/main'
include { FASTQ_FASTQC_UMITOOLS_TRIMGALORE } from '../subworkflows/nf-core/fastq_fastqc_umitools_trimgalore/main'
include { FASTQ_ALIGN_HISAT2 } from '../subworkflows/nf-core/fastq_align_hisat2/main'
include { BAM_SORT_STATS_SAMTOOLS } from '../subworkflows/nf-core/bam_sort_stats_samtools/main'
include { BAM_MARKDUPLICATES_PICARD } from '../subworkflows/nf-core/bam_markduplicates_picard/main'
include { BAM_RSEQC } from '../subworkflows/nf-core/bam_rseqc/main'
include { BAM_DEDUP_STATS_SAMTOOLS_UMITOOLS as BAM_DEDUP_STATS_SAMTOOLS_UMITOOLS_GENOME } from '../subworkflows/nf-core/bam_dedup_stats_samtools_umitools/main'
include { BAM_DEDUP_STATS_SAMTOOLS_UMITOOLS as BAM_DEDUP_STATS_SAMTOOLS_UMITOOLS_TRANSCRIPTOME } from '../subworkflows/nf-core/bam_dedup_stats_samtools_umitools/main'
include { BEDGRAPH_BEDCLIP_BEDGRAPHTOBIGWIG as BEDGRAPH_BEDCLIP_BEDGRAPHTOBIGWIG_FORWARD } from '../subworkflows/nf-core/bedgraph_bedclip_bedgraphtobigwig/main'
include { BEDGRAPH_BEDCLIP_BEDGRAPHTOBIGWIG as BEDGRAPH_BEDCLIP_BEDGRAPHTOBIGWIG_REVERSE } from '../subworkflows/nf-core/bedgraph_bedclip_bedgraphtobigwig/main'
/*
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
RUN MAIN WORKFLOW
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
*/
// Info required for completion email and summary
def multiqc_report = []
def pass_percent_mapped = [:]
def fail_percent_mapped = [:]
workflow RNASEQ {
ch_versions = Channel.empty()
//
// SUBWORKFLOW: Uncompress and prepare reference genome files
//
def biotype = params.gencode ? "gene_type" : params.featurecounts_group_type
PREPARE_GENOME (
prepareToolIndices,
biotype,
is_aws_igenome
)
ch_versions = ch_versions.mix(PREPARE_GENOME.out.versions)
// Check if contigs in genome fasta file > 512 Mbp
if (!params.skip_alignment && !params.bam_csi_index) {
PREPARE_GENOME
.out
.fai
.map { WorkflowRnaseq.checkMaxContigSize(it, log) }
}
//
// SUBWORKFLOW: Read in samplesheet, validate and stage input files
//
INPUT_CHECK (
ch_input
)
.reads
.map {
meta, fastq ->
new_id = meta.id - ~/_T\d+/
[ meta + [id: new_id], fastq ]
}
.groupTuple()
.branch {
meta, fastq ->
single : fastq.size() == 1
return [ meta, fastq.flatten() ]
multiple: fastq.size() > 1
return [ meta, fastq.flatten() ]
}
.set { ch_fastq }
ch_versions = ch_versions.mix(INPUT_CHECK.out.versions)
//
// MODULE: Concatenate FastQ files from same sample if required
//
CAT_FASTQ (
ch_fastq.multiple
)
.reads
.mix(ch_fastq.single)
.set { ch_cat_fastq }
ch_versions = ch_versions.mix(CAT_FASTQ.out.versions.first().ifEmpty(null))
// Branch FastQ channels if 'auto' specified to infer strandedness
ch_cat_fastq
.branch {
meta, fastq ->
auto_strand : meta.strandedness == 'auto'
return [ meta, fastq ]
known_strand: meta.strandedness != 'auto'
return [ meta, fastq ]
}
.set { ch_strand_fastq }
//
// SUBWORKFLOW: Sub-sample FastQ files and pseudo-align with Salmon to auto-infer strandedness
//
FASTQ_SUBSAMPLE_FQ_SALMON (
ch_strand_fastq.auto_strand,
PREPARE_GENOME.out.fasta,
PREPARE_GENOME.out.transcript_fasta,
PREPARE_GENOME.out.gtf,
PREPARE_GENOME.out.salmon_index,
!params.salmon_index && !('salmon' in prepareToolIndices)
)
ch_versions = ch_versions.mix(FASTQ_SUBSAMPLE_FQ_SALMON.out.versions)
FASTQ_SUBSAMPLE_FQ_SALMON
.out
.json_info
.join(ch_strand_fastq.auto_strand)
.map { meta, json, reads ->
return [ meta + [ strandedness: WorkflowRnaseq.getSalmonInferredStrandedness(json) ], reads ]
}
.mix(ch_strand_fastq.known_strand)
.set { ch_strand_inferred_fastq }
//
// SUBWORKFLOW: Read QC, extract UMI and trim adapters
//
FASTQ_FASTQC_UMITOOLS_TRIMGALORE (
ch_strand_inferred_fastq,
params.skip_fastqc || params.skip_qc,
params.with_umi,
params.skip_umi_extract,
params.skip_trimming,
params.umi_discard_read
)
ch_versions = ch_versions.mix(FASTQ_FASTQC_UMITOOLS_TRIMGALORE.out.versions)
//
// Filter channels to get samples that passed minimum trimmed read count
//
ch_fail_trimming_multiqc = Channel.empty()
ch_filtered_reads = FASTQ_FASTQC_UMITOOLS_TRIMGALORE.out.reads
if (!params.skip_trimming) {
ch_filtered_reads
.join(FASTQ_FASTQC_UMITOOLS_TRIMGALORE.out.trim_log, remainder: true)
.map {
meta, reads, trim_log ->
if (trim_log) {
if (!meta.single_end) {
trim_log = trim_log[-1]
}
num_reads = WorkflowRnaseq.getTrimGaloreReadsAfterFiltering(trim_log)
[ meta, reads, num_reads ]
} else {
[ meta, reads, params.min_trimmed_reads + 1 ]
}
}
.set { ch_num_trimmed_reads }
ch_num_trimmed_reads
.map { meta, reads, num_reads -> if (num_reads > params.min_trimmed_reads) [ meta, reads ] }
.set { ch_filtered_reads }
ch_num_trimmed_reads
.map {
meta, reads, num_reads ->
if (num_reads <= params.min_trimmed_reads) {
return [ "$meta.id\t$num_reads" ]
}
}
.set { ch_num_trimmed_reads }
MULTIQC_TSV_FAIL_TRIMMED (
ch_num_trimmed_reads.collect(),
["Sample", "Reads after trimming"],
'fail_trimmed_samples'
)
.set { ch_fail_trimming_multiqc }
}
//
// MODULE: Remove genome contaminant reads
//
if (!params.skip_bbsplit) {
BBMAP_BBSPLIT (
ch_filtered_reads,
PREPARE_GENOME.out.bbsplit_index,
[],
[ [], [] ],
false
)
.primary_fastq
.set { ch_filtered_reads }
ch_versions = ch_versions.mix(BBMAP_BBSPLIT.out.versions.first())
}
//
// MODULE: Remove ribosomal RNA reads
//
ch_sortmerna_multiqc = Channel.empty()
if (params.remove_ribo_rna) {
ch_sortmerna_fastas = Channel.from(ch_ribo_db.readLines()).map { row -> file(row, checkIfExists: true) }.collect()
SORTMERNA (
ch_filtered_reads,
ch_sortmerna_fastas
)
.reads
.set { ch_filtered_reads }
ch_sortmerna_multiqc = SORTMERNA.out.log
ch_versions = ch_versions.mix(SORTMERNA.out.versions.first())
}
//
// SUBWORKFLOW: Alignment with STAR and gene/transcript quantification with Salmon
//
ch_genome_bam = Channel.empty()
ch_genome_bam_index = Channel.empty()
ch_samtools_stats = Channel.empty()
ch_samtools_flagstat = Channel.empty()
ch_samtools_idxstats = Channel.empty()
ch_star_multiqc = Channel.empty()
ch_aligner_pca_multiqc = Channel.empty()
ch_aligner_clustering_multiqc = Channel.empty()
if (!params.skip_alignment && params.aligner == 'star_salmon') {
ALIGN_STAR (
ch_filtered_reads,
PREPARE_GENOME.out.star_index,
PREPARE_GENOME.out.gtf,
params.star_ignore_sjdbgtf,
'',
params.seq_center ?: '',
is_aws_igenome,
PREPARE_GENOME.out.fasta
)
ch_genome_bam = ALIGN_STAR.out.bam
ch_genome_bam_index = ALIGN_STAR.out.bai
ch_transcriptome_bam = ALIGN_STAR.out.bam_transcript
ch_samtools_stats = ALIGN_STAR.out.stats
ch_samtools_flagstat = ALIGN_STAR.out.flagstat
ch_samtools_idxstats = ALIGN_STAR.out.idxstats
ch_star_multiqc = ALIGN_STAR.out.log_final
if (params.bam_csi_index) {
ch_genome_bam_index = ALIGN_STAR.out.csi
}
ch_versions = ch_versions.mix(ALIGN_STAR.out.versions)
//
// SUBWORKFLOW: Remove duplicate reads from BAM file based on UMIs
//
if (params.with_umi) {
// Deduplicate genome BAM file before downstream analysis
BAM_DEDUP_STATS_SAMTOOLS_UMITOOLS_GENOME (
ch_genome_bam.join(ch_genome_bam_index, by: [0]),
params.umitools_dedup_stats
)
ch_genome_bam = BAM_DEDUP_STATS_SAMTOOLS_UMITOOLS_GENOME.out.bam
ch_genome_bam_index = BAM_DEDUP_STATS_SAMTOOLS_UMITOOLS_GENOME.out.bai
ch_samtools_stats = BAM_DEDUP_STATS_SAMTOOLS_UMITOOLS_GENOME.out.stats
ch_samtools_flagstat = BAM_DEDUP_STATS_SAMTOOLS_UMITOOLS_GENOME.out.flagstat
ch_samtools_idxstats = BAM_DEDUP_STATS_SAMTOOLS_UMITOOLS_GENOME.out.idxstats
if (params.bam_csi_index) {
ch_genome_bam_index = BAM_DEDUP_STATS_SAMTOOLS_UMITOOLS_GENOME.out.csi
}
ch_versions = ch_versions.mix(BAM_DEDUP_STATS_SAMTOOLS_UMITOOLS_GENOME.out.versions)
// Co-ordinate sort, index and run stats on transcriptome BAM
BAM_SORT_STATS_SAMTOOLS (
ch_transcriptome_bam,
PREPARE_GENOME.out.fasta
)
ch_transcriptome_sorted_bam = BAM_SORT_STATS_SAMTOOLS.out.bam
ch_transcriptome_sorted_bai = BAM_SORT_STATS_SAMTOOLS.out.bai
// Deduplicate transcriptome BAM file before read counting with Salmon
BAM_DEDUP_STATS_SAMTOOLS_UMITOOLS_TRANSCRIPTOME (
ch_transcriptome_sorted_bam.join(ch_transcriptome_sorted_bai, by: [0]),
params.umitools_dedup_stats
)
// Name sort BAM before passing to Salmon
SAMTOOLS_SORT (
BAM_DEDUP_STATS_SAMTOOLS_UMITOOLS_TRANSCRIPTOME.out.bam
)
// Only run prepare_for_rsem.py on paired-end BAM files
SAMTOOLS_SORT
.out
.bam
.branch {
meta, bam ->
single_end: meta.single_end
return [ meta, bam ]
paired_end: !meta.single_end
return [ meta, bam ]
}
.set { ch_umitools_dedup_bam }
// Fix paired-end reads in name sorted BAM file
// See: https://github.com/nf-core/rnaseq/issues/828
UMITOOLS_PREPAREFORSALMON (
ch_umitools_dedup_bam.paired_end
)
ch_versions = ch_versions.mix(UMITOOLS_PREPAREFORSALMON.out.versions.first())
ch_umitools_dedup_bam
.single_end
.mix(UMITOOLS_PREPAREFORSALMON.out.bam)
.set { ch_transcriptome_bam }
}
//
// SUBWORKFLOW: Count reads from BAM alignments using Salmon
//
QUANTIFY_STAR_SALMON (
ch_transcriptome_bam,
ch_dummy_file,
PREPARE_GENOME.out.transcript_fasta,
PREPARE_GENOME.out.gtf,
true,
params.salmon_quant_libtype ?: ''
)
ch_versions = ch_versions.mix(QUANTIFY_STAR_SALMON.out.versions)
if (!params.skip_qc & !params.skip_deseq2_qc) {
DESEQ2_QC_STAR_SALMON (
QUANTIFY_STAR_SALMON.out.counts_gene_length_scaled,
ch_pca_header_multiqc,
ch_clustering_header_multiqc
)
ch_aligner_pca_multiqc = DESEQ2_QC_STAR_SALMON.out.pca_multiqc
ch_aligner_clustering_multiqc = DESEQ2_QC_STAR_SALMON.out.dists_multiqc
ch_versions = ch_versions.mix(DESEQ2_QC_STAR_SALMON.out.versions)
}
}
//
// SUBWORKFLOW: Alignment with STAR and gene/transcript quantification with RSEM
//
ch_rsem_multiqc = Channel.empty()
if (!params.skip_alignment && params.aligner == 'star_rsem') {
QUANTIFY_RSEM (
ch_filtered_reads,
PREPARE_GENOME.out.rsem_index
)
ch_genome_bam = QUANTIFY_RSEM.out.bam
ch_genome_bam_index = QUANTIFY_RSEM.out.bai
ch_samtools_stats = QUANTIFY_RSEM.out.stats
ch_samtools_flagstat = QUANTIFY_RSEM.out.flagstat
ch_samtools_idxstats = QUANTIFY_RSEM.out.idxstats
ch_star_multiqc = QUANTIFY_RSEM.out.logs
ch_rsem_multiqc = QUANTIFY_RSEM.out.stat
if (params.bam_csi_index) {
ch_genome_bam_index = QUANTIFY_RSEM.out.csi
}
ch_versions = ch_versions.mix(QUANTIFY_RSEM.out.versions)
if (!params.skip_qc & !params.skip_deseq2_qc) {
DESEQ2_QC_RSEM (
QUANTIFY_RSEM.out.merged_counts_gene,
ch_pca_header_multiqc,
ch_clustering_header_multiqc
)
ch_aligner_pca_multiqc = DESEQ2_QC_RSEM.out.pca_multiqc
ch_aligner_clustering_multiqc = DESEQ2_QC_RSEM.out.dists_multiqc
ch_versions = ch_versions.mix(DESEQ2_QC_RSEM.out.versions)
}
}
//
// SUBWORKFLOW: Alignment with HISAT2
//
ch_hisat2_multiqc = Channel.empty()
if (!params.skip_alignment && params.aligner == 'hisat2') {
FASTQ_ALIGN_HISAT2 (
ch_filtered_reads,
PREPARE_GENOME.out.hisat2_index,
PREPARE_GENOME.out.splicesites,
PREPARE_GENOME.out.fasta
)
ch_genome_bam = FASTQ_ALIGN_HISAT2.out.bam
ch_genome_bam_index = FASTQ_ALIGN_HISAT2.out.bai
ch_samtools_stats = FASTQ_ALIGN_HISAT2.out.stats
ch_samtools_flagstat = FASTQ_ALIGN_HISAT2.out.flagstat
ch_samtools_idxstats = FASTQ_ALIGN_HISAT2.out.idxstats
ch_hisat2_multiqc = FASTQ_ALIGN_HISAT2.out.summary
if (params.bam_csi_index) {
ch_genome_bam_index = FASTQ_ALIGN_HISAT2.out.csi
}
ch_versions = ch_versions.mix(FASTQ_ALIGN_HISAT2.out.versions)
//
// SUBWORKFLOW: Remove duplicate reads from BAM file based on UMIs
//
if (params.with_umi) {
BAM_DEDUP_STATS_SAMTOOLS_UMI_UMITOOLS_GENOME (
ch_genome_bam.join(ch_genome_bam_index, by: [0]),
params.umitools_dedup_stats
)
ch_genome_bam = BAM_DEDUP_STATS_SAMTOOLS_UMI_UMITOOLS_GENOME.out.bam
ch_genome_bam_index = BAM_DEDUP_STATS_SAMTOOLS_UMI_UMITOOLS_GENOME.out.bai
ch_samtools_stats = BAM_DEDUP_STATS_SAMTOOLS_UMI_UMITOOLS_GENOME.out.stats
ch_samtools_flagstat = BAM_DEDUP_STATS_SAMTOOLS_UMI_UMITOOLS_GENOME.out.flagstat
ch_samtools_idxstats = BAM_DEDUP_STATS_SAMTOOLS_UMI_UMITOOLS_GENOME.out.idxstats
if (params.bam_csi_index) {
ch_genome_bam_index = BAM_DEDUP_STATS_SAMTOOLS_UMI_UMITOOLS_GENOME.out.csi
}
ch_versions = ch_versions.mix(BAM_DEDUP_STATS_SAMTOOLS_UMI_UMITOOLS_GENOME.out.versions)
}
}
//
// Filter channels to get samples that passed STAR minimum mapping percentage
//
ch_fail_mapping_multiqc = Channel.empty()
if (!params.skip_alignment && params.aligner.contains('star')) {
ch_star_multiqc
.map { meta, align_log -> [ meta ] + WorkflowRnaseq.getStarPercentMapped(params, align_log) }
.set { ch_percent_mapped }
ch_genome_bam
.join(ch_percent_mapped, by: [0])
.map { meta, ofile, mapped, pass -> if (pass) [ meta, ofile ] }
.set { ch_genome_bam }
ch_genome_bam_index
.join(ch_percent_mapped, by: [0])
.map { meta, ofile, mapped, pass -> if (pass) [ meta, ofile ] }
.set { ch_genome_bam_index }
ch_percent_mapped
.branch { meta, mapped, pass ->
pass: pass
pass_percent_mapped[meta.id] = mapped
return [ "$meta.id\t$mapped" ]
fail: !pass
fail_percent_mapped[meta.id] = mapped
return [ "$meta.id\t$mapped" ]
}
.set { ch_pass_fail_mapped }
def header = [
"Sample",
"STAR uniquely mapped reads (%)"
]
MULTIQC_TSV_FAIL_MAPPED (
ch_pass_fail_mapped.fail.collect(),
header,
'fail_mapped_samples'
)
.set { ch_fail_mapping_multiqc }
}
//
// MODULE: Run Preseq
//
ch_preseq_multiqc = Channel.empty()
if (!params.skip_alignment && !params.skip_qc && !params.skip_preseq) {
PRESEQ_LCEXTRAP (
ch_genome_bam
)
ch_preseq_multiqc = PRESEQ_LCEXTRAP.out.lc_extrap
ch_versions = ch_versions.mix(PRESEQ_LCEXTRAP.out.versions.first())
}
//
// SUBWORKFLOW: Mark duplicate reads
//
ch_markduplicates_multiqc = Channel.empty()
if (!params.skip_alignment && !params.skip_markduplicates && !params.with_umi) {
BAM_MARKDUPLICATES_PICARD (
ch_genome_bam,
PREPARE_GENOME.out.fasta,
PREPARE_GENOME.out.fai
)
ch_genome_bam = BAM_MARKDUPLICATES_PICARD.out.bam
ch_genome_bam_index = BAM_MARKDUPLICATES_PICARD.out.bai
ch_samtools_stats = BAM_MARKDUPLICATES_PICARD.out.stats
ch_samtools_flagstat = BAM_MARKDUPLICATES_PICARD.out.flagstat
ch_samtools_idxstats = BAM_MARKDUPLICATES_PICARD.out.idxstats
ch_markduplicates_multiqc = BAM_MARKDUPLICATES_PICARD.out.metrics
if (params.bam_csi_index) {
ch_genome_bam_index = BAM_MARKDUPLICATES_PICARD.out.csi
}
ch_versions = ch_versions.mix(BAM_MARKDUPLICATES_PICARD.out.versions)
}
//
// MODULE: STRINGTIE
//
if (!params.skip_alignment && !params.skip_stringtie) {
STRINGTIE_STRINGTIE (
ch_genome_bam,
PREPARE_GENOME.out.gtf
)
ch_versions = ch_versions.mix(STRINGTIE_STRINGTIE.out.versions.first())
}
//
// MODULE: Feature biotype QC using featureCounts
//
ch_featurecounts_multiqc = Channel.empty()
if (!params.skip_alignment && !params.skip_qc && !params.skip_biotype_qc && biotype) {
PREPARE_GENOME
.out
.gtf
.map { WorkflowRnaseq.biotypeInGtf(it, biotype, log) }
.set { biotype_in_gtf }
// Prevent any samples from running if GTF file doesn't have a valid biotype
ch_genome_bam
.combine(PREPARE_GENOME.out.gtf)
.combine(biotype_in_gtf)
.filter { it[-1] }
.map { it[0..<it.size()-1] }
.set { ch_featurecounts }
SUBREAD_FEATURECOUNTS (
ch_featurecounts
)
ch_versions = ch_versions.mix(SUBREAD_FEATURECOUNTS.out.versions.first())
MULTIQC_CUSTOM_BIOTYPE (
SUBREAD_FEATURECOUNTS.out.counts,
ch_biotypes_header_multiqc
)
ch_featurecounts_multiqc = MULTIQC_CUSTOM_BIOTYPE.out.tsv
ch_versions = ch_versions.mix(MULTIQC_CUSTOM_BIOTYPE.out.versions.first())
}
//
// MODULE: Genome-wide coverage with BEDTools
//
if (!params.skip_alignment && !params.skip_bigwig) {
BEDTOOLS_GENOMECOV (
ch_genome_bam
)
ch_versions = ch_versions.mix(BEDTOOLS_GENOMECOV.out.versions.first())
//
// SUBWORKFLOW: Convert bedGraph to bigWig
//
BEDGRAPH_BEDCLIP_BEDGRAPHTOBIGWIG_FORWARD (
BEDTOOLS_GENOMECOV.out.bedgraph_forward,
PREPARE_GENOME.out.chrom_sizes
)
ch_versions = ch_versions.mix(BEDGRAPH_BEDCLIP_BEDGRAPHTOBIGWIG_FORWARD.out.versions)
BEDGRAPH_BEDCLIP_BEDGRAPHTOBIGWIG_REVERSE (
BEDTOOLS_GENOMECOV.out.bedgraph_reverse,
PREPARE_GENOME.out.chrom_sizes
)
}
//
// MODULE: Downstream QC steps
//
ch_qualimap_multiqc = Channel.empty()
ch_dupradar_multiqc = Channel.empty()
ch_bamstat_multiqc = Channel.empty()
ch_inferexperiment_multiqc = Channel.empty()
ch_innerdistance_multiqc = Channel.empty()
ch_junctionannotation_multiqc = Channel.empty()
ch_junctionsaturation_multiqc = Channel.empty()
ch_readdistribution_multiqc = Channel.empty()
ch_readduplication_multiqc = Channel.empty()
ch_fail_strand_multiqc = Channel.empty()
ch_tin_multiqc = Channel.empty()
if (!params.skip_alignment && !params.skip_qc) {
if (!params.skip_qualimap) {
QUALIMAP_RNASEQ (
ch_genome_bam,
PREPARE_GENOME.out.gtf
)
ch_qualimap_multiqc = QUALIMAP_RNASEQ.out.results
ch_versions = ch_versions.mix(QUALIMAP_RNASEQ.out.versions.first())
}
if (!params.skip_dupradar) {
DUPRADAR (
ch_genome_bam,
PREPARE_GENOME.out.gtf
)
ch_dupradar_multiqc = DUPRADAR.out.multiqc
ch_versions = ch_versions.mix(DUPRADAR.out.versions.first())
}
if (!params.skip_rseqc && rseqc_modules.size() > 0) {
BAM_RSEQC (
ch_genome_bam.join(ch_genome_bam_index, by: [0]),
PREPARE_GENOME.out.gene_bed,
rseqc_modules
)
ch_bamstat_multiqc = BAM_RSEQC.out.bamstat_txt
ch_inferexperiment_multiqc = BAM_RSEQC.out.inferexperiment_txt
ch_innerdistance_multiqc = BAM_RSEQC.out.innerdistance_freq
ch_junctionannotation_multiqc = BAM_RSEQC.out.junctionannotation_log
ch_junctionsaturation_multiqc = BAM_RSEQC.out.junctionsaturation_rscript
ch_readdistribution_multiqc = BAM_RSEQC.out.readdistribution_txt
ch_readduplication_multiqc = BAM_RSEQC.out.readduplication_pos_xls
ch_tin_multiqc = BAM_RSEQC.out.tin_txt
ch_versions = ch_versions.mix(BAM_RSEQC.out.versions)
ch_inferexperiment_multiqc
.map { meta, strand_log -> [ meta ] + WorkflowRnaseq.getInferexperimentStrandedness(strand_log, 30) }
.filter { it[0].strandedness != it[1] }
.map { meta, strandedness, sense, antisense, undetermined ->
[ "$meta.id\t$meta.strandedness\t$strandedness\t$sense\t$antisense\t$undetermined" ]
}
.set { ch_fail_strand }
def header = [
"Sample",
"Provided strandedness",
"Inferred strandedness",
"Sense (%)",
"Antisense (%)",
"Undetermined (%)"
]
MULTIQC_TSV_STRAND_CHECK (
ch_fail_strand.collect(),
header,
'fail_strand_check'
)
.set { ch_fail_strand_multiqc }
}
}
//
// SUBWORKFLOW: Pseudo-alignment and quantification with Salmon
//
ch_salmon_multiqc = Channel.empty()
ch_pseudoaligner_pca_multiqc = Channel.empty()
ch_pseudoaligner_clustering_multiqc = Channel.empty()
if (params.pseudo_aligner == 'salmon') {
QUANTIFY_SALMON (
ch_filtered_reads,
PREPARE_GENOME.out.salmon_index,
ch_dummy_file,
PREPARE_GENOME.out.gtf,
false,
params.salmon_quant_libtype ?: ''
)
ch_salmon_multiqc = QUANTIFY_SALMON.out.results
ch_versions = ch_versions.mix(QUANTIFY_SALMON.out.versions)
if (!params.skip_qc & !params.skip_deseq2_qc) {
DESEQ2_QC_SALMON (
QUANTIFY_SALMON.out.counts_gene_length_scaled,
ch_pca_header_multiqc,
ch_clustering_header_multiqc
)
ch_pseudoaligner_pca_multiqc = DESEQ2_QC_SALMON.out.pca_multiqc
ch_pseudoaligner_clustering_multiqc = DESEQ2_QC_SALMON.out.dists_multiqc
ch_versions = ch_versions.mix(DESEQ2_QC_SALMON.out.versions)
}
}
//
// MODULE: Pipeline reporting
//
CUSTOM_DUMPSOFTWAREVERSIONS (
ch_versions.unique().collectFile(name: 'collated_versions.yml')
)
//
// MODULE: MultiQC
//
if (!params.skip_multiqc) {
workflow_summary = WorkflowRnaseq.paramsSummaryMultiqc(workflow, summary_params)
ch_workflow_summary = Channel.value(workflow_summary)
methods_description = WorkflowRnaseq.methodsDescriptionText(workflow, ch_multiqc_custom_methods_description)
ch_methods_description = Channel.value(methods_description)
MULTIQC (
ch_multiqc_config,
ch_multiqc_custom_config.collect().ifEmpty([]),
CUSTOM_DUMPSOFTWAREVERSIONS.out.mqc_yml.collect(),
ch_workflow_summary.collectFile(name: 'workflow_summary_mqc.yaml'),
ch_methods_description.collectFile(name: 'methods_description_mqc.yaml'),
ch_multiqc_logo.collect().ifEmpty([]),
ch_fail_trimming_multiqc.ifEmpty([]),
ch_fail_mapping_multiqc.ifEmpty([]),
ch_fail_strand_multiqc.ifEmpty([]),
FASTQ_FASTQC_UMITOOLS_TRIMGALORE.out.fastqc_zip.collect{it[1]}.ifEmpty([]),
FASTQ_FASTQC_UMITOOLS_TRIMGALORE.out.trim_zip.collect{it[1]}.ifEmpty([]),
FASTQ_FASTQC_UMITOOLS_TRIMGALORE.out.trim_log.collect{it[1]}.ifEmpty([]),
ch_sortmerna_multiqc.collect{it[1]}.ifEmpty([]),
ch_star_multiqc.collect{it[1]}.ifEmpty([]),
ch_hisat2_multiqc.collect{it[1]}.ifEmpty([]),
ch_rsem_multiqc.collect{it[1]}.ifEmpty([]),
ch_salmon_multiqc.collect{it[1]}.ifEmpty([]),
ch_samtools_stats.collect{it[1]}.ifEmpty([]),
ch_samtools_flagstat.collect{it[1]}.ifEmpty([]),
ch_samtools_idxstats.collect{it[1]}.ifEmpty([]),
ch_markduplicates_multiqc.collect{it[1]}.ifEmpty([]),
ch_featurecounts_multiqc.collect{it[1]}.ifEmpty([]),
ch_aligner_pca_multiqc.collect().ifEmpty([]),
ch_aligner_clustering_multiqc.collect().ifEmpty([]),
ch_pseudoaligner_pca_multiqc.collect().ifEmpty([]),
ch_pseudoaligner_clustering_multiqc.collect().ifEmpty([]),
ch_preseq_multiqc.collect{it[1]}.ifEmpty([]),
ch_qualimap_multiqc.collect{it[1]}.ifEmpty([]),
ch_dupradar_multiqc.collect{it[1]}.ifEmpty([]),
ch_bamstat_multiqc.collect{it[1]}.ifEmpty([]),
ch_inferexperiment_multiqc.collect{it[1]}.ifEmpty([]),
ch_innerdistance_multiqc.collect{it[1]}.ifEmpty([]),
ch_junctionannotation_multiqc.collect{it[1]}.ifEmpty([]),
ch_junctionsaturation_multiqc.collect{it[1]}.ifEmpty([]),
ch_readdistribution_multiqc.collect{it[1]}.ifEmpty([]),
ch_readduplication_multiqc.collect{it[1]}.ifEmpty([]),
ch_tin_multiqc.collect{it[1]}.ifEmpty([])
)
multiqc_report = MULTIQC.out.report.toList()
}
}
/*
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
COMPLETION EMAIL AND SUMMARY
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
*/
workflow.onComplete {
if (params.email || params.email_on_fail) {
NfcoreTemplate.email(workflow, params, summary_params, projectDir, log, multiqc_report, fail_percent_mapped)
}
if (params.hook_url) {
NfcoreTemplate.IM_notification(workflow, params, summary_params, projectDir, log)
}
NfcoreTemplate.summary(workflow, params, log, fail_percent_mapped, pass_percent_mapped)
}
/*
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
THE END
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
*/