oarfish 0.10.0

A fast, accurate and versatile tool for long-read transcript quantification.
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
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
use crate::NamedDigestVec;
use crate::alignment_parser;
use crate::em;
use crate::kde_utils;
use crate::prog_opts::Args;
use crate::util::constants::EMPTY_READ_NAME;
use crate::util::oarfish_types::AlnInfo;
use crate::util::oarfish_types::DiscardTable;
use crate::util::oarfish_types::{
    AlignmentFilters, EMInfo, InMemoryAlignmentStore, InputSourceType, ReadChunkWithNames,
    ReadSource, TranscriptInfo,
};
use crate::util::read_function::read_short_quant_vec;
use crate::util::write_function::{write_infrep_file, write_out_prob, write_output};
use crate::util::projection::{mapping_to_genomic_alignment, projected_to_records};
use bramble_rs::ProjectionConfig;
use bramble_rs::ProjectionContext;
use bramble_rs::g2t::G2TTree;
use bramble_rs::project_group_with;
use crate::{logistic_prob, normalize_read_probs};
use arrow2::{array::Float64Array, chunk::Chunk, datatypes::Field};
use crossbeam::channel::Receiver;
use crossbeam::channel::Sender;
use crossbeam::channel::bounded;

use needletail::parse_fastx_file;
use noodles_bam as bam;
use num_format::{Locale, ToFormattedString};
use serde_json::json;
use std::io::BufRead;
use swapvec::{SwapVec, SwapVecConfig};
use tracing::{info, warn};

/// Produce a [serde_json::Value] that encodes the relevant arguments and
/// parameters of the run that we wish to record to file. Ultimately, this
/// will be written to the corresponding `meta_info.json` file for this run.
fn get_json_info(
    args: &Args,
    emi: &EMInfo,
    seqcol_digest: &NamedDigestVec,
    aln_time: std::time::Duration,
) -> serde_json::Value {
    let prob = if args.model_coverage {
        "logistic_coverage"
    } else {
        "no_coverage"
    };

    let source = if args.alignments.is_some() {
        "from_bam"
    } else {
        "from_raw_reads"
    };

    json!({
        "prob_model" : prob,
        "alignment_source" : source,
        "alignment_time" : {
            "comment" : "Time to parse (in alignment mode) or generate (in raw read mode) alignments, as well as apply filters, and compute conditional probabilities.",
            "human_time" : humantime::format_duration(aln_time).to_string(),
            "seconds" : aln_time.as_secs_f64()
        },
        "bin_width" : args.bin_width,
        "filter_options" : &emi.eq_map.filter_opts,
        "discard_table" : &emi.eq_map.discard_table,
        "alignments": &args.alignments,
        "output": args.output.as_ref().expect("present"),
        "verbose": &args.verbose,
        "single_cell": &args.single_cell,
        "quiet": &args.quiet,
        "em_max_iter": &args.max_em_iter,
        "em_convergence_thresh": &args.convergence_thresh,
        "threads": &args.threads,
        "filter_group": &args.filter_group,
        "write_assignment_probs": &emi.eq_map.filter_opts.write_assignment_probs_type,
        "short_quant": &args.short_quant,
        "num_bootstraps": &args.num_bootstraps,
        "digest": seqcol_digest.to_json()
    })
}

#[allow(clippy::too_many_arguments)]
fn perform_inference_and_write_output(
    header: &noodles_sam::header::Header,
    store: &mut InMemoryAlignmentStore,
    name_vec: Option<SwapVec<String>>,
    txps: &mut [TranscriptInfo],
    txps_name: &[String],
    seqcol_digest: NamedDigestVec,
    aln_time: std::time::Duration,
    args: &Args,
) -> anyhow::Result<()> {
    // print discard table information in which the user might be interested.
    info!("\ndiscard_table: \n{}\n", store.discard_table.to_table());

    // if we are using the KDE, create that here.
    let kde_opt: Option<kders::kde::KDEModel> = if args.use_kde {
        Some(kde_utils::get_kde_model(txps, store)?)
    } else {
        None
    };

    if store.filter_opts.model_coverage {
        //obtaining the Cumulative Distribution Function (CDF) for each transcript
        logistic_prob(txps, args.growth_rate, &args.bin_width, args.threads);
        //Normalize the probabilities for the records of each read
        normalize_read_probs(store, txps, &args.bin_width);
    }

    info!(
        "Total number of alignment records : {}",
        store.total_len().to_formatted_string(&Locale::en)
    );
    info!(
        "number of aligned reads : {}",
        store.num_aligned_reads().to_formatted_string(&Locale::en)
    );
    info!(
        "number of unique alignments : {}",
        store.unique_alignments().to_formatted_string(&Locale::en)
    );

    // if we are seeding the quantification estimates with short read
    // abundances, then read those in here.
    let init_abundances = args.short_quant.as_ref().map(|sr_path| {
        read_short_quant_vec(sr_path, txps_name).unwrap_or_else(|e| panic!("{}", e))
    });

    // wrap up all of the relevant information we need for estimation
    // in an EMInfo struct and then call the EM algorithm.
    let emi = EMInfo {
        eq_map: store,
        txp_info: txps,
        max_iter: args.max_em_iter,
        convergence_thresh: args.convergence_thresh,
        init_abundances,
        kde_model: kde_opt,
    };

    if args.use_kde {
        /*
        // run EM for model train iterations
        let orig_iter = emi.max_iter;
        emi.max_iter = 10;
        let counts = em::em(&emi, args.threads);
        // relearn the kde
        let new_model =
        kde_utils::refresh_kde_model(&txps, &store, &emi.kde_model.unwrap(), &counts);
        info!("refreshed KDE model");
        emi.kde_model = Some(new_model?);
        emi.max_iter = orig_iter;
        */
    }

    let counts = if args.threads > 4 {
        em::em_par(&emi, args.threads)
    } else {
        em::em(&emi, args.threads)
    };

    let aux_txp_counts = crate::util::aux_counts::get_aux_counts(store, txps)?;

    // prepare the JSON object we'll write
    // to meta_info.json
    let json_info = get_json_info(args, &emi, &seqcol_digest, aln_time);

    // write the output
    write_output(
        args.output.as_ref().expect("present"),
        json_info,
        header,
        &counts,
        &aux_txp_counts,
    )?;

    // if the user requested bootstrap replicates,
    // compute and write those out now.
    if args.num_bootstraps > 0 {
        let breps = em::bootstrap(&emi, args.num_bootstraps, args.threads);

        let mut new_arrays = vec![];
        let mut bs_fields = vec![];
        for (i, b) in breps.into_iter().enumerate() {
            let bs_array = Float64Array::from_vec(b);
            bs_fields.push(Field::new(
                format!("bootstrap.{}", i),
                bs_array.data_type().clone(),
                false,
            ));
            new_arrays.push(bs_array.boxed());
        }
        let chunk = Chunk::new(new_arrays);
        write_infrep_file(args.output.as_ref().expect("present"), bs_fields, chunk)?;
    }

    if args.write_assignment_probs.is_some() {
        let name_vec = name_vec
            .expect("cannot write assignment probabilities without valid vector of read names");
        write_out_prob(
            args.output.as_ref().expect("present"),
            &emi,
            &counts,
            name_vec,
            txps_name,
            args.display_thresh,
        )?;
    }

    Ok(())
}

pub fn quantify_bulk_alignments_from_bam<R: BufRead>(
    header: &noodles_sam::Header,
    filter_opts: AlignmentFilters,
    reader: &mut bam::io::Reader<R>,
    txps: &mut [TranscriptInfo],
    txps_name: &[String],
    args: &Args,
    seqcol_digest: NamedDigestVec,
) -> anyhow::Result<()> {
    let mut name_vec = if filter_opts.write_assignment_probs {
        Some(SwapVec::<String>::with_config(SwapVecConfig {
            swap_after: Default::default(),
            batch_size: Default::default(),
            compression: Some(swapvec::Compression::Lz4),
        }))
    } else {
        None
    };
    // now parse the actual alignments for the reads and store the results
    // in our in-memory stor
    let mut store = InMemoryAlignmentStore::new(filter_opts, header);
    let read_aln_start = std::time::SystemTime::now();
    alignment_parser::parse_alignments(
        &mut store,
        &mut name_vec,
        header,
        reader,
        txps,
        args.sort_check_num,
        args.quiet,
    )?;
    let read_aln_time = read_aln_start.elapsed()?;
    info!(
        "Parsing of alignments from input took: {}",
        humantime::format_duration(read_aln_time).to_string()
    );

    perform_inference_and_write_output(
        header,
        &mut store,
        name_vec,
        txps,
        txps_name,
        seqcol_digest,
        read_aln_time,
        args,
    )
}

/// Genome-BAM mode: parse a name-collated, spliced genome BAM, project each
/// read's alignments onto the transcriptome with bramble, and quantify.
///
/// `txp_header` / `txps` / `txps_name` describe the *transcriptome* (built from
/// the annotation via the projection bridge); `genome_header` is the input BAM's
/// (chromosome) header used to read records and resolve reference ids.
#[allow(clippy::too_many_arguments)]
pub fn quantify_genome_alignments_from_bam<R: BufRead>(
    genome_header: &noodles_sam::Header,
    txp_header: &noodles_sam::Header,
    g2t: &G2TTree,
    proj_config: &ProjectionConfig,
    filter_opts: AlignmentFilters,
    reader: &mut bam::io::Reader<R>,
    txps: &mut [TranscriptInfo],
    txps_name: &[String],
    args: &Args,
    seqcol_digest: NamedDigestVec,
) -> anyhow::Result<()> {
    let mut name_vec = if filter_opts.write_assignment_probs {
        Some(SwapVec::<String>::with_config(SwapVecConfig {
            swap_after: Default::default(),
            batch_size: Default::default(),
            compression: Some(swapvec::Compression::Lz4),
        }))
    } else {
        None
    };

    // the store keys alignments by transcript, so it lives over the
    // transcriptome header.
    let mut store = InMemoryAlignmentStore::new(filter_opts, txp_header);
    let read_aln_start = std::time::SystemTime::now();
    alignment_parser::parse_genome_alignments(
        &mut store,
        &mut name_vec,
        genome_header,
        g2t,
        proj_config,
        args.projected_prob_beta,
        args.projected_prob_source,
        reader,
        txps,
        args.sort_check_num,
        args.quiet,
    )?;
    let read_aln_time = read_aln_start.elapsed()?;
    info!(
        "Parsing and projection of genome alignments took: {}",
        humantime::format_duration(read_aln_time).to_string()
    );

    perform_inference_and_write_output(
        txp_header,
        &mut store,
        name_vec,
        txps,
        txps_name,
        seqcol_digest,
        read_aln_time,
        args,
    )
}

/// Genome-read mode: spliced-align raw reads to the genome with the rammap
/// aligner, project each read's mappings onto the transcriptome with bramble,
/// and quantify.
///
/// This mirrors [`quantify_bulk_alignments_raw_reads`] (same producer → mapper →
/// consumer pipeline), but the mapper projects each read's genomic mappings onto
/// the transcripts (via bramble) rather than mapping directly to a transcriptome.
/// `txp_header` / `txps` / `txps_name` describe the transcriptome built from the
/// annotation; `aligner` is a spliced genome aligner; `g2t` the genome→
/// transcriptome index over the same reference order as the aligner targets.
#[allow(clippy::too_many_arguments)]
#[allow(unused_mut)]
pub fn quantify_genome_raw_reads(
    txp_header: &noodles_sam::Header,
    mut aligner: crate::util::mapper::Mapper,
    g2t: &G2TTree,
    proj_config: &ProjectionConfig,
    filter_opts: AlignmentFilters,
    read_paths: &[std::path::PathBuf],
    txps: &mut [TranscriptInfo],
    txps_name: &[String],
    args: &Args,
    seqcol_digest: NamedDigestVec,
) -> anyhow::Result<()> {
    // shared, read-only view of the transcript info (for the projected filters).
    let mut txp_info_view: Vec<TranscriptInfo> = Vec::with_capacity(txps.len());
    for ti in txps.iter() {
        txp_info_view.push(ti.clone());
    }

    let map_threads = args.threads.saturating_sub(2).max(1);

    let beta = args.projected_prob_beta;
    let use_fasta = proj_config.use_fasta;
    let prob_source = args.projected_prob_source;

    type ReadGroup = ReadChunkWithNames;
    type AlignmentGroupInfo = (Vec<AlnInfo>, Vec<f32>, Vec<usize>, Option<Vec<String>>);

    let (read_sender, read_receiver): (Sender<ReadGroup>, Receiver<ReadGroup>) =
        bounded(args.threads * 10);

    const READ_CHUNK_SIZE: usize = 200;
    let mut rpaths = vec![];
    read_paths.clone_into(&mut rpaths);

    // Producer thread: read sequences and send them to the channel.
    let producer = std::thread::spawn(move || {
        let mut ctr = 0_usize;
        let mut chunk_size = 0_usize;
        let mut read_chunk = ReadChunkWithNames::new();

        let mark_chunk = |chunk_size: &mut usize,
                          ctr: &mut usize,
                          read_chunk: &mut ReadGroup,
                          read_sender: &Sender<ReadGroup>| {
            *chunk_size += 1;
            *ctr += 1;
            if *chunk_size >= READ_CHUNK_SIZE {
                read_sender
                    .send(read_chunk.clone())
                    .expect("Error sending sequence");
                read_chunk.clear();
                *chunk_size = 0;
            }
        };

        for read_path in rpaths {
            match get_source_type(&read_path) {
                InputSourceType::Ubam => {
                    let mut reader = std::fs::File::open(read_path)
                        .map(bam::io::Reader::new)
                        .expect("could not create BAM reader");
                    let header = reader.read_header().expect("could not read BAM header");
                    for result in reader.record_bufs(&header) {
                        let record = result.expect("Error reading ubam record");
                        record.add_to_read_group(&mut read_chunk);
                        mark_chunk(&mut chunk_size, &mut ctr, &mut read_chunk, &read_sender);
                    }
                }
                s @ (InputSourceType::Fastx | InputSourceType::Unknown) => {
                    if matches!(s, InputSourceType::Unknown) {
                        warn!(
                            "could not determine input file type for {} from suffix; assuming (possibly gzipped) fastx",
                            read_path.display()
                        );
                    }
                    let mut reader =
                        parse_fastx_file(read_path).expect("valid path/file to read sequences");
                    while let Some(result) = reader.next() {
                        let record = result.expect("Error reading record");
                        record.add_to_read_group(&mut read_chunk);
                        mark_chunk(&mut chunk_size, &mut ctr, &mut read_chunk, &read_sender);
                    }
                }
            }
        }
        if chunk_size > 0 {
            read_sender
                .send(read_chunk)
                .expect("Error sending sequence");
        }
        ctr
    });

    let (mut store, name_vec, aln_time) = std::thread::scope(
        |s| -> anyhow::Result<(
            InMemoryAlignmentStore,
            Option<SwapVec<String>>,
            std::time::Duration,
        )> {
            const ALN_GROUP_CHUNK_LIMIT: usize = 100;

            let (aln_group_sender, aln_group_receiver): (
                Sender<AlignmentGroupInfo>,
                Receiver<AlignmentGroupInfo>,
            ) = bounded(args.threads * 100);

            let write_assignment_probs: bool = args.write_assignment_probs.is_some();
            // Diagnostics: reads the aligner maps to the genome vs reads that
            // produce >=1 projected transcriptome alignment (projection loss).
            let n_genome_mapped = std::sync::Arc::new(std::sync::atomic::AtomicU64::new(0));
            let n_projected = std::sync::Arc::new(std::sync::atomic::AtomicU64::new(0));
            // Mapper threads: align each read to the genome, then project its
            // mappings onto the transcriptome and filter.
            let consumers: Vec<_> = (0..map_threads)
                .map(|_| {
                    let receiver = read_receiver.clone();
                    let filter = filter_opts.clone();
                    let loc_aligner = aligner.clone();
                    let my_txp_info_view = &txp_info_view;
                    let aln_group_sender = aln_group_sender.clone();
                    let n_genome_mapped = std::sync::Arc::clone(&n_genome_mapped);
                    let n_projected = std::sync::Arc::clone(&n_projected);

                    s.spawn(move || {
                        let mut discard_table = DiscardTable::new();

                        let mut chunk_size = 0_usize;
                        let mut aln_group_alns: Vec<AlnInfo> = Vec::new();
                        let mut aln_group_probs: Vec<f32> = Vec::new();
                        let mut aln_group_boundaries: Vec<usize> = Vec::new();
                        let mut aln_group_read_names = write_assignment_probs.then(Vec::new);
                        aln_group_boundaries.push(0);
                        // reused across all reads this worker projects (avoids
                        // per-read allocation of bramble's projection scratch).
                        let mut pctx = ProjectionContext::new();
                        // Reused across every read this worker projects, to avoid
                        // re-allocating the per-read alignment/score scratch vectors.
                        let mut galns: Vec<bramble_rs::GenomicAlignment> = Vec::new();
                        let mut src_scores: Vec<i32> = Vec::new();

                        for read_chunk in receiver {
                            for (name, seq) in read_chunk.iter() {
                                let map_res_opt =
                                    crate::util::mapper::map_read(&loc_aligner, name, seq);
                                let Ok(mappings) = map_res_opt else {
                                    warn!("Error encountered mapping read: {}", map_res_opt.unwrap_err());
                                    continue;
                                };

                                let query_name = String::from_utf8_lossy(name).into_owned();
                                // build GenomicAlignments and a parallel vector of
                                // their alignment scores (used by the
                                // score/combined probability sources).
                                galns.clear();
                                src_scores.clear();
                                for m in mappings.iter() {
                                    if let Some(ga) =
                                        mapping_to_genomic_alignment(m, &query_name, seq.len())
                                    {
                                        galns.push(ga);
                                        src_scores.push(crate::util::mapper::alignment_score(m));
                                    }
                                }
                                if galns.is_empty() {
                                    continue;
                                }
                                n_genome_mapped
                                    .fetch_add(1, std::sync::atomic::Ordering::Relaxed);
                                // attach the read sequence to the first alignment when
                                // soft-clip rescue is enabled (bramble shares one seq
                                // across the read group). bramble expects the seq in
                                // forward-reference orientation (as a BAM stores SEQ),
                                // so reverse-complement it when the canonical mapping is
                                // on the reverse strand — the read is held here in its
                                // original FASTQ orientation.
                                if use_fasta {
                                    galns[0].sequence = Some(if galns[0].is_reverse {
                                        crate::util::projection::revcomp(seq)
                                    } else {
                                        seq.to_vec()
                                    });
                                }

                                let projected =
                                    project_group_with(&galns, g2t, proj_config, &mut pctx);
                                if projected.is_empty() {
                                    continue;
                                }
                                n_projected.fetch_add(1, std::sync::atomic::Ordering::Relaxed);
                                let recs = projected_to_records(&projected, &src_scores);
                                let (ag, aprobs) = filter.filter_projected(
                                    &mut discard_table,
                                    my_txp_info_view,
                                    &recs,
                                    seq.len(),
                                    beta,
                                    prob_source,
                                );

                                if !ag.is_empty() {
                                    aln_group_alns.extend_from_slice(&ag);
                                    aln_group_probs.extend_from_slice(&aprobs);
                                    aln_group_boundaries.push(aln_group_alns.len());
                                    if let Some(ref mut names_vec) = aln_group_read_names {
                                        names_vec.push(query_name);
                                    }
                                    chunk_size += 1;
                                }

                                if chunk_size >= ALN_GROUP_CHUNK_LIMIT {
                                    aln_group_sender
                                        .send((
                                            aln_group_alns.clone(),
                                            aln_group_probs.clone(),
                                            aln_group_boundaries.clone(),
                                            aln_group_read_names,
                                        ))
                                        .expect("Error sending alignment group");
                                    aln_group_alns.clear();
                                    aln_group_probs.clear();
                                    aln_group_boundaries.clear();
                                    aln_group_boundaries.push(0);
                                    aln_group_read_names = write_assignment_probs.then(Vec::new);
                                    chunk_size = 0;
                                }
                            }
                        }
                        if chunk_size > 0 {
                            aln_group_sender
                                .send((
                                    aln_group_alns,
                                    aln_group_probs,
                                    aln_group_boundaries,
                                    aln_group_read_names,
                                ))
                                .expect("Error sending alignment group");
                        }
                        discard_table
                    })
                })
                .collect();

            #[allow(clippy::useless_asref)]
            let txps_mut = txps.as_mut();
            let filter_opts_store = filter_opts.clone();
            let aln_group_consumer = s.spawn(move || {
                let mut name_vec = if filter_opts_store.write_assignment_probs {
                    Some(SwapVec::<String>::with_config(SwapVecConfig {
                        swap_after: Default::default(),
                        batch_size: Default::default(),
                        compression: Some(swapvec::Compression::Lz4),
                    }))
                } else {
                    None
                };

                let mut store = InMemoryAlignmentStore::new(filter_opts_store, txp_header);

                let pb = if args.quiet {
                    indicatif::ProgressBar::hidden()
                } else {
                    indicatif::ProgressBar::new_spinner().with_message("Number of reads mapped")
                };
                pb.set_style(
                    indicatif::ProgressStyle::with_template(
                        "[{elapsed_precise}] {spinner:4.green/blue} {msg} {human_pos:>12}",
                    )
                    .unwrap()
                    .tick_chars("⠁⠁⠉⠙⠚⠒⠂⠂⠒⠲⠴⠤⠄⠄⠤⠠⠠⠤⠦⠖⠒⠐⠐⠒⠓⠋⠉⠈⠈"),
                );
                pb.set_draw_target(indicatif::ProgressDrawTarget::stderr_with_hz(4));

                for (ags, aprobs, aln_boundaries, read_names) in aln_group_receiver {
                    let mut reversed_read_names = if let Some(mut names_vec) = read_names {
                        names_vec.reverse();
                        Some(names_vec)
                    } else {
                        None
                    };

                    for window in aln_boundaries.windows(2) {
                        pb.inc(1);
                        let ag = &ags[window[0]..window[1]];
                        let as_probs = &aprobs[window[0]..window[1]];
                        let read_name_opt = if let Some(ref mut names_vec) = reversed_read_names {
                            names_vec.pop()
                        } else {
                            None
                        };

                        if store.add_filtered_group(ag, as_probs, txps_mut) {
                            if let Some(ref mut nvec) = name_vec {
                                let read_name =
                                    read_name_opt.unwrap_or(EMPTY_READ_NAME.to_string());
                                nvec.push(read_name)
                                    .expect("cannot push name to read name vector");
                            }
                            if ag.len() == 1 {
                                store.inc_unique_alignments();
                            }
                        }
                    }
                }
                pb.finish_with_message("Finished aligning and projecting reads.");
                (store, name_vec)
            });

            let aln_start = std::time::SystemTime::now();
            let total_reads = producer.join().expect("Producer thread panicked");

            let mut discard_tables: Vec<DiscardTable> = Vec::with_capacity(map_threads);
            for consumer in consumers {
                discard_tables.push(consumer.join().expect("Consumer thread panicked"));
            }

            drop(aln_group_sender);

            let (mut store, name_vec) = aln_group_consumer
                .join()
                .expect("Alignment group consumer panicked");

            info!(
                "Parsed {} total reads",
                total_reads.to_formatted_string(&Locale::en)
            );
            let n_gm = n_genome_mapped.load(std::sync::atomic::Ordering::Relaxed);
            let n_pr = n_projected.load(std::sync::atomic::Ordering::Relaxed);
            info!(
                "reads aligned to genome: {}; of those, projected to >=1 transcript: {} ({:.2}%); projection-dropped: {}",
                n_gm.to_formatted_string(&Locale::en),
                n_pr.to_formatted_string(&Locale::en),
                if n_gm > 0 { 100.0 * (n_pr as f64) / (n_gm as f64) } else { 0.0 },
                (n_gm - n_pr).to_formatted_string(&Locale::en),
            );
            let aln_time = aln_start.elapsed()?;
            info!(
                "Spliced genome alignment + projection of raw reads took: {}",
                humantime::format_duration(aln_time).to_string()
            );

            for dt in &discard_tables {
                store.aggregate_discard_table(dt);
            }
            Ok((store, name_vec, aln_time))
        },
    )?;

    perform_inference_and_write_output(
        txp_header,
        &mut store,
        name_vec,
        txps,
        txps_name,
        seqcol_digest,
        aln_time,
        args,
    )
}

fn get_source_type(pb: &std::path::Path) -> InputSourceType {
    let faq_endings = vec![
        ".fasta",
        ".fastq",
        ".FASTA",
        ".FASTQ",
        ".fa",
        ".fq",
        ".FA",
        ".FQ",
        ".fasta.gz",
        ".fastq.gz",
        ".FASTA.GZ",
        ".FASTQ.GZ",
        ".fa.gz",
        ".fq.gz",
        ".FA.GZ",
        ".FQ.GZ",
    ];
    let ubam_endings = vec![".bam", ".BAM", ".ubam", ".UBAM"];
    if let Some(ps) = pb.to_str() {
        for fe in faq_endings {
            if ps.ends_with(fe) {
                return InputSourceType::Fastx;
            }
        }
        for be in ubam_endings {
            if ps.ends_with(be) {
                return InputSourceType::Ubam;
            }
        }
    }

    InputSourceType::Unknown
}

#[allow(clippy::too_many_arguments)]
#[allow(unused_mut)]
pub fn quantify_bulk_alignments_raw_reads(
    header: &noodles_sam::Header,
    mut aligner: crate::util::mapper::Mapper,
    filter_opts: AlignmentFilters,
    read_paths: &[std::path::PathBuf],
    txps: &mut [TranscriptInfo],
    txps_name: &[String],
    args: &Args,
    seqcol_digest: NamedDigestVec,
) -> anyhow::Result<()> {
    // now parse the actual alignments for the reads and store the results
    // in our in-memory stor

    // we will take only shared refs to this
    let mut txp_info_view: Vec<TranscriptInfo> = Vec::with_capacity(txps.len());
    for ti in txps.iter() {
        txp_info_view.push(ti.clone());
    }

    // at least one mapping thread, otherwise everything but the fastx parser
    // and the in memory alignment store populator
    let map_threads = args.threads.saturating_sub(2).max(1);

    type ReadGroup = ReadChunkWithNames;
    type AlignmentGroupInfo = (Vec<AlnInfo>, Vec<f32>, Vec<usize>, Option<Vec<String>>);

    let (read_sender, read_receiver): (Sender<ReadGroup>, Receiver<ReadGroup>) =
        bounded(args.threads * 10);

    const READ_CHUNK_SIZE: usize = 200;
    let mut rpaths = vec![];
    read_paths.clone_into(&mut rpaths);

    // Producer thread: reads sequences and sends them to the channel
    let producer = std::thread::spawn(move || {
        let mut ctr = 0_usize;
        let mut chunk_size = 0_usize;
        let mut read_chunk = ReadChunkWithNames::new();

        // work shared between the two different
        // source types
        let mark_chunk = |chunk_size: &mut usize,
                          ctr: &mut usize,
                          read_chunk: &mut ReadGroup,
                          read_sender: &Sender<ReadGroup>| {
            *chunk_size += 1;
            *ctr += 1;
            if *chunk_size >= READ_CHUNK_SIZE {
                read_sender
                    .send(read_chunk.clone())
                    .expect("Error sending sequence");
                // prepare for the next chunk
                read_chunk.clear();
                *chunk_size = 0;
            }
        };

        // read from either a UBAM or (possibly compressed) FASTX file
        for read_path in rpaths {
            match get_source_type(&read_path) {
                InputSourceType::Ubam => {
                    let mut reader = std::fs::File::open(read_path)
                        .map(bam::io::Reader::new)
                        .expect("could not create BAM reader");
                    let header = reader.read_header().expect("could not read BAM header");
                    for result in reader.record_bufs(&header) {
                        let record = result.expect("Error reading ubam record");
                        record.add_to_read_group(&mut read_chunk);
                        mark_chunk(&mut chunk_size, &mut ctr, &mut read_chunk, &read_sender);
                    }
                }
                s @ (InputSourceType::Fastx | InputSourceType::Unknown) => {
                    if matches!(s, InputSourceType::Unknown) {
                        warn!(
                            "could not determine input file type for {} from suffix; assuming (possibly gzipped) fastx",
                            read_path.display()
                        );
                    }
                    let mut reader =
                        parse_fastx_file(read_path).expect("valid path/file to read sequences");
                    while let Some(result) = reader.next() {
                        let record = result.expect("Error reading record");
                        record.add_to_read_group(&mut read_chunk);
                        mark_chunk(&mut chunk_size, &mut ctr, &mut read_chunk, &read_sender);
                    }
                }
            }
        }
        // if any reads remain, send them off
        if chunk_size > 0 {
            read_sender
                .send(read_chunk)
                .expect("Error sending sequence");
        }
        ctr
    });

    // we need the scope here so we can borrow the relevant non-'static data
    let (mut store, name_vec, aln_time) = std::thread::scope(
        |s| -> anyhow::Result<(
            InMemoryAlignmentStore,
            Option<SwapVec<String>>,
            std::time::Duration,
        )> {
            const ALN_GROUP_CHUNK_LIMIT: usize = 100;

            let (aln_group_sender, aln_group_receiver): (
                Sender<AlignmentGroupInfo>,
                Receiver<AlignmentGroupInfo>,
            ) = bounded(args.threads * 100);

            // Consumer threads: receive sequences and perform alignment
            let write_assignment_probs: bool = args.write_assignment_probs.is_some();
            let consumers: Vec<_> = (0..map_threads)
                .map(|_| {
                    let receiver = read_receiver.clone();
                    let mut filter = filter_opts.clone();
                    let loc_aligner = aligner.clone();

                    let my_txp_info_view = &txp_info_view;
                    let aln_group_sender = aln_group_sender.clone();
                    s.spawn(move || {
                        let mut discard_table = DiscardTable::new();

                        let mut chunk_size = 0_usize;
                        let mut aln_group_alns: Vec<AlnInfo> = Vec::new();
                        let mut aln_group_probs: Vec<f32> = Vec::new();
                        let mut aln_group_boundaries: Vec<usize> = Vec::new();
                        let mut aln_group_read_names = write_assignment_probs.then(Vec::new);
                        aln_group_boundaries.push(0);

                        // get the next chunk of reads
                        for read_chunk in receiver {
                            // iterate over every read
                            for (name, seq) in read_chunk.iter() {
                                // map the next read, with cigar string
                                let map_res_opt =
                                    crate::util::mapper::map_read(&loc_aligner, name, seq);
                                if let Ok(mut mappings) = map_res_opt {
                                    let (ag, aprobs) = filter.filter(
                                        &mut discard_table,
                                        header,
                                        my_txp_info_view,
                                        &mut mappings,
                                    );

                                    if !ag.is_empty() {
                                        aln_group_alns.extend_from_slice(&ag);
                                        aln_group_probs.extend_from_slice(&aprobs);
                                        aln_group_boundaries.push(aln_group_alns.len());
                                        // if we are storing read names
                                        if let Some(ref mut names_vec) = aln_group_read_names {
                                            let name_str =
                                                String::from_utf8_lossy(name).into_owned();
                                            names_vec.push(name_str);
                                        }
                                        chunk_size += 1;
                                    }
                                    if chunk_size >= ALN_GROUP_CHUNK_LIMIT {
                                        aln_group_sender
                                            .send((
                                                aln_group_alns.clone(),
                                                aln_group_probs.clone(),
                                                aln_group_boundaries.clone(),
                                                aln_group_read_names,
                                            ))
                                            .expect("Error sending alignment group");
                                        aln_group_alns.clear();
                                        aln_group_probs.clear();
                                        aln_group_boundaries.clear();
                                        aln_group_boundaries.push(0);
                                        aln_group_read_names =
                                            write_assignment_probs.then(Vec::new);
                                        chunk_size = 0;
                                    }
                                } else {
                                    warn!(
                                        "Error encountered mappread_ing read : {}",
                                        map_res_opt.unwrap_err()
                                    );
                                }
                            }
                        }
                        if chunk_size > 0 {
                            aln_group_sender
                                .send((
                                    aln_group_alns,
                                    aln_group_probs,
                                    aln_group_boundaries,
                                    aln_group_read_names,
                                ))
                                .expect("Error sending alignment group");
                        }
                        discard_table
                    })
                })
                .collect();

            #[allow(clippy::useless_asref)]
            let txps_mut = txps.as_mut();
            let filter_opts_store = filter_opts.clone();
            let aln_group_consumer = s.spawn(move || {
                let mut name_vec = if filter_opts_store.write_assignment_probs {
                    Some(SwapVec::<String>::with_config(SwapVecConfig {
                        swap_after: Default::default(),
                        batch_size: Default::default(),
                        compression: Some(swapvec::Compression::Lz4),
                    }))
                } else {
                    None
                };

                let mut store = InMemoryAlignmentStore::new(filter_opts_store, header);

                let pb = if args.quiet {
                    indicatif::ProgressBar::hidden()
                } else {
                    indicatif::ProgressBar::new_spinner().with_message("Number of reads mapped")
                };

                pb.set_style(
                    indicatif::ProgressStyle::with_template(
                        "[{elapsed_precise}] {spinner:4.green/blue} {msg} {human_pos:>12}",
                    )
                    .unwrap()
                    .tick_chars("⠁⠁⠉⠙⠚⠒⠂⠂⠒⠲⠴⠤⠄⠄⠤⠠⠠⠤⠦⠖⠒⠐⠐⠒⠓⠋⠉⠈⠈"),
                );
                pb.set_draw_target(indicatif::ProgressDrawTarget::stderr_with_hz(4));

                for (ags, aprobs, aln_boundaries, read_names) in aln_group_receiver {
                    // if we are getting read names out then we are going to "reverse" them
                    // here so that we can simply pop the strings off the back to get them
                    // in order. We do this since we cannot otherwise "move" a string out of a
                    // Vec.
                    let mut reversed_read_names = if let Some(mut names_vec) = read_names {
                        names_vec.reverse();
                        Some(names_vec)
                    } else {
                        None
                    };

                    for window in aln_boundaries.windows(2) {
                        pb.inc(1);
                        let group_start = window[0];
                        let group_end = window[1];
                        let ag = &ags[group_start..group_end];
                        let as_probs = &aprobs[group_start..group_end];
                        let read_name_opt = if let Some(ref mut names_vec) = reversed_read_names {
                            names_vec.pop()
                        } else {
                            None
                        };

                        if store.add_filtered_group(ag, as_probs, txps_mut) {
                            if let Some(ref mut nvec) = name_vec {
                                let read_name =
                                    read_name_opt.unwrap_or(EMPTY_READ_NAME.to_string());
                                nvec.push(read_name)
                                    .expect("cannot push name to read name vector");
                            }
                            if ag.len() == 1 {
                                store.inc_unique_alignments();
                            }
                        }
                    }
                }
                pb.finish_with_message("Finished aligning reads.");
                (store, name_vec)
            });

            let aln_start = std::time::SystemTime::now();
            // Wait for the producer to finish reading
            let total_reads = producer.join().expect("Producer thread panicked");

            let mut discard_tables: Vec<DiscardTable> = Vec::with_capacity(map_threads);
            for consumer in consumers {
                let dt = consumer.join().expect("Consumer thread panicked");
                discard_tables.push(dt);
            }

            drop(aln_group_sender);

            let (mut store, name_vec) = aln_group_consumer
                .join()
                .expect("Alignment group consumer panicked");

            info!(
                "Parsed {} total reads",
                total_reads.to_formatted_string(&Locale::en)
            );

            let aln_time = aln_start.elapsed()?;
            info!(
                "Alignment of raw reads using rammap took: {}",
                humantime::format_duration(aln_time).to_string()
            );

            for dt in &discard_tables {
                store.aggregate_discard_table(dt);
            }
            Ok((store, name_vec, aln_time))
        },
    )?;

    perform_inference_and_write_output(
        header,
        &mut store,
        name_vec,
        txps,
        txps_name,
        seqcol_digest,
        aln_time,
        args,
    )
}