File: AlignmentLibrary.hpp

package info (click to toggle)
salmon 1.10.3%2Bds1-3
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 35,148 kB
  • sloc: cpp: 200,707; ansic: 171,082; sh: 859; python: 792; makefile: 238
file content (683 lines) | stat: -rw-r--r-- 24,440 bytes parent folder | download | duplicates (3)
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
#ifndef ALIGNMENT_LIBRARY_HPP
#define ALIGNMENT_LIBRARY_HPP

extern "C" {
#include "io_lib/os.h"
#include "io_lib/scram.h"
#undef max
#undef min
}

// Our includes
#include "AlignmentGroup.hpp"
#include "AlignmentModel.hpp"
#include "ONTAlignmentModel.hpp"
#include "BAMQueue.hpp"
#include "BAMUtils.hpp"
#include "ClusterForest.hpp"
#include "DistributionUtils.hpp"
#include "EquivalenceClassBuilder.hpp"
#include "FASTAParser.hpp"
#include "FragmentLengthDistribution.hpp"
#include "FragmentStartPositionDistribution.hpp"
#include "GCFragModel.hpp"
#include "LibraryFormat.hpp"
#include "LibraryTypeDetector.hpp"
#include "ReadKmerDist.hpp"
#include "SBModel.hpp"
#include "SalmonOpts.hpp"
#include "SalmonUtils.hpp"
#include "SimplePosBias.hpp"
#include "SpinLock.hpp" // From pufferfish, with try_lock
#include "Transcript.hpp"
#include "concurrentqueue.h"
#include "parallel_hashmap/phmap.h"

// Boost includes
#include <boost/filesystem.hpp>

// Standard includes
#include <functional>
#include <memory>
#include <vector>
#include <stdexcept>

template <typename T> class NullFragmentFilter;

/**
 *  This class represents a library of alignments used to quantify
 *  a set of target transcripts.  The AlignmentLibrary contains info
 *  about both the alignment file and the target sequence (transcripts).
 *  It is used to group them together and track information about them
 *  during the quantification procedure.
 */
template <typename FragT, typename EQBuilderT, typename AlignModelT> class AlignmentLibrary {

public:
  AlignmentLibrary(std::vector<boost::filesystem::path>& alnFiles,
                   const boost::filesystem::path& transcriptFile,
                   LibraryFormat libFmt, SalmonOpts& salmonOpts)
      : alignmentFiles_(alnFiles), transcriptFile_(transcriptFile),
        libFmt_(libFmt), transcripts_(std::vector<Transcript>()),
        fragStartDists_(5), posBiasFW_(5), posBiasRC_(5), posBiasExpectFW_(5),
        posBiasExpectRC_(5), /*seqBiasModel_(1.0),*/
        eqBuilder_(salmonOpts.jointLog, salmonOpts.maxHashResizeThreads), quantificationPasses_(0),
        expectedBias_(constExprPow(4, readBias_[0].getK()), 1.0),
        expectedGC_(salmonOpts.numConditionalGCBins, salmonOpts.numFragGCBins,
                    distribution_utils::DistributionSpace::LOG),
        observedGC_(salmonOpts.numConditionalGCBins, salmonOpts.numFragGCBins,
                    distribution_utils::DistributionSpace::LOG) {
    namespace bfs = boost::filesystem;

    // Make sure the alignment file exists.
    for (auto& alignmentFile : alignmentFiles_) {
      if (!bfs::exists(alignmentFile)) {
        std::stringstream ss;
        ss << "The provided alignment file: " << alignmentFile
           << " does not exist!\n";
        throw std::invalid_argument(ss.str());
      }
    }

    // Make sure the transcript file exists.
    if (!bfs::exists(transcriptFile_)) {
      std::stringstream ss;
      ss << "The provided transcript file: " << transcriptFile_
         << " does not exist!\n";
      throw std::invalid_argument(ss.str());
    }

    // The alignment file existed, so create the alignment queue
    size_t numParseThreads = salmonOpts.numParseThreads;
    std::cerr << "parseThreads = " << numParseThreads << "\n";
    bq = std::unique_ptr<BAMQueue<FragT>>(
        new BAMQueue<FragT>(alnFiles, libFmt_, numParseThreads,
                            salmonOpts.mappingCacheMemoryLimit));

    std::cerr << "Checking that provided alignment files have consistent "
                 "headers . . . ";
    if (!salmon::utils::headersAreConsistent(bq->headers())) {
      std::stringstream ss;
      ss << "\nThe multiple alignment files provided had inconsistent "
            "headers.\n";
      ss << "Currently, we require that if multiple SAM/BAM files are "
            "provided,\n";
      ss << "they must have identical @SQ records.\n";
      throw std::invalid_argument(ss.str());
    }
    std::cerr << "done\n";

    SAM_hdr* header = bq->header();

    // Figure out aligner information from the header if we can
    aligner_ = salmon::bam_utils::inferAlignerFromHeader(header);

    // in this case check for decoys and make a list of their names
    phmap::flat_hash_set<std::string> decoys;
    if (aligner_ == salmon::bam_utils::AlignerDetails::PUFFERFISH) {
     // for each reference
     for (decltype(header->nref) i = 0; i < header->nref; ++i) {
       // for each tag 
       SAM_hdr_tag *tag;
	     for (tag = header->ref[i].tag; tag; tag = tag->next) {
         // if this tag marks it as a decoy
         if ((tag->len == 4) and (std::strncmp(tag->str, "DS:D", 4) == 0)) {
           decoys.insert(header->ref[i].name);
           break;
         } // end if decoy tag

       } // end for each tag
      } // end for each referecne
    }
    
    if (!decoys.empty()) {
        bq->forceEndParsing();
        bq.reset();
        salmonOpts.jointLog->error(
        "Salmon is being run in alignment-mode with a SAM/BAM file that contains decoy\n"
        "sequences (marked as such during salmon indexing). This SAM/BAM file had {}\n"
        "such sequences tagged in the header. Since alignments to decoys are not\n"
        "intended for decoy-level quantification, this functionality is not currently\n"
        "supported.  If you wish to run salmon with this SAM/BAM file, please \n"
        "filter out / remove decoy transcripts (those tagged with `DS:D`) from the \n"
        "header, and all SAM/BAM records that represent alignments to decoys \n"
        "(those tagged with `XT:A:D`). If you believe you are receiving this message\n"
        "in error, please report this issue on GitHub.", decoys.size());
        salmonOpts.jointLog->flush();
        std::stringstream ss;
        ss << "\nCannot quantify from SAM/BAM file containing decoy transcripts or alignment records!\n";
        throw std::runtime_error(ss.str());
    }

    // The transcript file existed, so load up the transcripts
    double alpha = 0.005;
    // we know how many we will have, so reserve the space for 
    // them.
    transcripts_.reserve(header->nref);
    for (decltype(header->nref) i = 0; i < header->nref; ++i) {
      transcripts_.emplace_back(i, header->ref[i].name, header->ref[i].len,
                                alpha);
    }

    FASTAParser fp(transcriptFile.string());

    fmt::print(stderr, "Populating targets from aln = {}, fasta = {} . . .",
               alnFiles.front(), transcriptFile_);
    fp.populateTargets(transcripts_, salmonOpts);
    /*
for (auto& txp : transcripts_) {
    // Length classes taken from
    // ======
    // Roberts, Adam, et al.
    // "Improving RNA-Seq expression estimates by correcting for fragment bias."
    // Genome Biol 12.3 (2011): R22.
    // ======
    // perhaps, define these in a more data-driven way
    if (txp.RefLength <= 1334) {
        txp.lengthClassIndex(0);
    } else if (txp.RefLength <= 2104) {
        txp.lengthClassIndex(0);
    } else if (txp.RefLength <= 2988) {
        txp.lengthClassIndex(0);
    } else if (txp.RefLength <= 4389) {
        txp.lengthClassIndex(0);
    } else {
        txp.lengthClassIndex(0);
    }
}
*/

    std::vector<uint32_t> lengths;
    lengths.reserve(transcripts_.size());
    for (auto& txp : transcripts_) {
      lengths.push_back(txp.RefLength);
    }
    setTranscriptLengthClasses_(lengths, posBiasFW_.size());

    fmt::print(stderr, "done\n");

    // Create the cluster forest for this set of transcripts
    clusters_.reset(new ClusterForest(transcripts_.size(), transcripts_));

    // Initialize the fragment length distribution
    size_t maxFragLen = salmonOpts.fragLenDistMax;
    double meanFragLen = salmonOpts.fragLenDistPriorMean;
    double fragLenStd = salmonOpts.fragLenDistPriorSD;
    size_t fragLenKernelN = 4;
    double fragLenKernelP = 0.5;
    flDist_.reset(new FragmentLengthDistribution(1.0, maxFragLen, meanFragLen,
                                                 fragLenStd, fragLenKernelN,
                                                 fragLenKernelP, 1));

    alnMod_.reset(new AlignModelT(1.0, salmonOpts.numErrorBins));
    alnMod_->setLogger(salmonOpts.jointLog);

    if (libFmt.type == ReadType::SINGLE_END) {
      // Convert the PMF to non-log scale
      std::vector<double> logPMF;
      size_t minVal;
      size_t maxVal;
      flDist_->dumpPMF(logPMF, minVal, maxVal);
      double sum = salmon::math::LOG_0;
      for (auto v : logPMF) {
        sum = salmon::math::logAdd(sum, v);
      }
      for (auto& v : logPMF) {
        v -= sum;
      }

      // Create the non-logged distribution.
      // Here, we multiply by 100 to discourage small
      // numbers in the correctionFactorsfromCounts call
      // below.
      std::vector<double> pmf(maxVal + 1, 0.0);
      for (size_t i = minVal; i < maxVal; ++i) {
        pmf[i] = 100.0 * std::exp(logPMF[i - minVal]);
      }

      using distribution_utils::DistributionSpace;
      // We compute the factors in linear space (since we've de-logged the pmf)
      conditionalMeans_ = distribution_utils::correctionFactorsFromMass(
          pmf, DistributionSpace::LINEAR);
    }

    salmon::utils::markAuxiliaryTargets(salmonOpts, transcripts_);

    // Start parsing the alignments
    NullFragmentFilter<FragT>* nff = nullptr;
    bool onlyProcessAmbiguousAlignments = false;
    bq->start(nff, onlyProcessAmbiguousAlignments);
  }

  AlignmentLibrary(std::vector<boost::filesystem::path>& alnFiles,
                   LibraryFormat libFmt, SalmonOpts& salmonOpts,
                   bool /*eqClassMode_*/, std::vector<std::string>& tnames,
                   std::vector<double>& tefflens)
      : alignmentFiles_(alnFiles),
        libFmt_(libFmt), transcripts_(std::vector<Transcript>()),
        fragStartDists_(5), posBiasFW_(5), posBiasRC_(5), posBiasExpectFW_(5),
        posBiasExpectRC_(5), /*seqBiasModel_(1.0),*/
        eqBuilder_(salmonOpts.jointLog, salmonOpts.maxHashResizeThreads), quantificationPasses_(0),
        expectedBias_(constExprPow(4, readBias_[0].getK()), 1.0),
        expectedGC_(salmonOpts.numConditionalGCBins, salmonOpts.numFragGCBins,
                    distribution_utils::DistributionSpace::LOG),
        observedGC_(salmonOpts.numConditionalGCBins, salmonOpts.numFragGCBins,
                    distribution_utils::DistributionSpace::LOG) {
    namespace bfs = boost::filesystem;

    // Make sure the alignment file exists.
    for (auto& alignmentFile : alignmentFiles_) {
      if (!bfs::exists(alignmentFile)) {
        std::stringstream ss;
        ss << "The provided eqClass file: " << alignmentFile
           << " does not exist!\n";
        throw std::invalid_argument(ss.str());
      }
    }

    // The transcript file existed, so load up the transcripts
    double alpha = 0.005;
    for (size_t i = 0; i < tnames.size(); ++i) {
      transcripts_.emplace_back(i, tnames[i].c_str(), tefflens[i], true, alpha);
    }

    // Initialize the fragment length distribution
    size_t maxFragLen = salmonOpts.fragLenDistMax;
    double meanFragLen = salmonOpts.fragLenDistPriorMean;
    double fragLenStd = salmonOpts.fragLenDistPriorSD;
    size_t fragLenKernelN = 4;
    double fragLenKernelP = 0.5;
    flDist_.reset(new FragmentLengthDistribution(1.0, maxFragLen, meanFragLen,
                                                 fragLenStd, fragLenKernelN,
                                                 fragLenKernelP, 1));

    alnMod_.reset(new AlignModelT(1.0, salmonOpts.numErrorBins));
    alnMod_->setLogger(salmonOpts.jointLog);
    salmon::utils::markAuxiliaryTargets(salmonOpts, transcripts_);
  }

  EQBuilderT& equivalenceClassBuilder() { return eqBuilder_; }

  std::string getIndexSeqHash256() const { return ""; }
  std::string getIndexNameHash256() const { return ""; }
  std::string getIndexSeqHash512() const { return ""; }
  std::string getIndexNameHash512() const { return ""; }
  std::string getIndexDecoySeqHash256() const { return ""; }
  std::string getIndexDecoyNameHash256() const { return ""; }

  /**
   * Return true if this read library is for paired-end reads and false
   * otherwise.
   */
  bool isPairedEnd() { return (libFmt_.type == ReadType::PAIRED_END); }

  salmon::bam_utils::AlignerDetails getAlignerType() const { return aligner_; }

  // TODO: Make same as mapping-based
  void updateTranscriptLengthsAtomic(std::atomic<bool>& done) {
    if (sl_.try_lock()) {
      if (!done) {

        auto fld = flDist_.get();
        // Convert the PMF to non-log scale
        std::vector<double> logPMF;
        size_t minVal;
        size_t maxVal;
        fld->dumpPMF(logPMF, minVal, maxVal);
        double sum = salmon::math::LOG_0;
        for (auto v : logPMF) {
          sum = salmon::math::logAdd(sum, v);
        }
        for (auto& v : logPMF) {
          v -= sum;
        }

        // Create the non-logged distribution.
        // Here, we multiply by 100 to discourage small
        // numbers in the correctionFactorsfromCounts call
        // below.
        std::vector<double> pmf(maxVal + 1, 0.0);
        for (size_t i = minVal; i < maxVal; ++i) {
          pmf[i] = 100.0 * std::exp(logPMF[i - minVal]);
        }

        using distribution_utils::DistributionSpace;
        // We compute the factors in linear space (since we've de-logged the
        // pmf)
        auto correctionFactors = distribution_utils::correctionFactorsFromMass(
            pmf, DistributionSpace::LINEAR);
        // Since we'll continue treating effective lengths in log space,
        // populate them as such
        distribution_utils::computeSmoothedEffectiveLengths(
            pmf.size(), transcripts_, correctionFactors,
            DistributionSpace::LOG);

        /*
                // Update the effective length of *every* transcript
                for( auto& t : transcripts_ ) {
                    t.updateEffectiveLength(logPMF, logFLDMean, minVal, maxVal);
                }
        */
        // then declare that we are done
        done = true;
      }
      sl_.unlock();
    }
  }

  const std::vector<double>& condMeans() const { return conditionalMeans_; }

  std::vector<Transcript>& transcripts() { return transcripts_; }
  const std::vector<Transcript>& transcripts() const { return transcripts_; }

  inline bool getAlignmentGroup(AlignmentGroup<FragT>*& ag) {
    return bq->getAlignmentGroup(ag);
  }

  // inline t_pool* threadPool() { return threadPool_.get(); }

  inline SAM_hdr* header() { return bq->header(); }

  inline std::vector<FragmentStartPositionDistribution>&
  fragmentStartPositionDistributions() {
    return fragStartDists_;
  }

  inline FragmentLengthDistribution* fragmentLengthDistribution() const {
    return flDist_.get();
  }

  inline AlignModelT& alignmentModel() { return *alnMod_.get(); }

  // SequenceBiasModel& sequenceBiasModel() { return seqBiasModel_; }

  //    inline oneapi::tbb::concurrent_queue<FragT*>& fragmentQueue() {
  inline oneapi::tbb::concurrent_queue<FragT*>& fragmentQueue() {
    return bq->getFragmentQueue();
  }

  //    inline oneapi::tbb::concurrent_bounded_queue<AlignmentGroup<FragT*>*>&
  //    alignmentGroupQueue() {
  inline moodycamel::ConcurrentQueue<AlignmentGroup<FragT*>*>&
  alignmentGroupQueue() {
    return bq->getAlignmentGroupQueue();
  }

  inline BAMQueue<FragT>& getAlignmentGroupQueue() { return *bq.get(); }

  inline size_t upperBoundHits() { return bq->numMappedFragments(); }
  inline size_t numObservedFragments() const {
    return bq->numObservedFragments();
  }
  inline size_t numMappedFragments() const { return bq->numMappedFragments(); }
  inline size_t numUniquelyMappedFragments() {
    return bq->numUniquelyMappedFragments();
  }
  inline double effectiveMappingRate() const {
    return static_cast<double>(numMappedFragments()) / numObservedFragments();
  }

  // const boost::filesystem::path& alignmentFile() { return alignmentFile_; }

  ClusterForest& clusterForest() { return *clusters_.get(); }

  template <typename FilterT>
  bool reset(bool incPasses = true, FilterT filter = nullptr,
             bool onlyProcessAmbiguousAlignments = false) {
    namespace bfs = boost::filesystem;

    for (auto& alignmentFile : alignmentFiles_) {
      if (!bfs::is_regular_file(alignmentFile)) {
        return false;
      }
    }

    bq->reset();
    bq->start(filter, onlyProcessAmbiguousAlignments);
    if (incPasses) {
      quantificationPasses_++;
      fmt::print(stderr, "Current iteration = {}\n", quantificationPasses_);
    }
    return true;
  }

  inline LibraryFormat format() { return libFmt_; }
  inline const LibraryFormat format() const { return libFmt_; }

  /**
   * If this is set, attempt to automatically detect this library's type
   */
  void enableAutodetect() {
    // if auto detection is not already enabled, and we're enabling it
    if (!detector_) {
      detector_.reset(new LibraryTypeDetector(libFmt_.type));
    }
  }

  bool autoDetect() const { return (detector_.get() != nullptr); }

  LibraryTypeDetector* getDetector() { return detector_.get(); }

  LibraryFormat& getFormat() { return libFmt_; }
  const LibraryFormat& getFormat() const { return libFmt_; }

  void setGCFracForward(double fracForward) { gcFracFwd_ = fracForward; }

  double gcFracFwd() const { return gcFracFwd_; }
  double gcFracRC() const { return 1.0 - gcFracFwd_; }

  std::vector<double>& expectedSeqBias() { return expectedBias_; }

  const std::vector<double>& expectedSeqBias() const { return expectedBias_; }

  void setExpectedGCBias(const GCFragModel& expectedBiasIn) {
    expectedGC_ = expectedBiasIn;
  }

  GCFragModel& expectedGCBias() { return expectedGC_; }

  const GCFragModel& expectedGCBias() const { return expectedGC_; }

  const GCFragModel& observedGC() const { return observedGC_; }

  GCFragModel& observedGC() { return observedGC_; }

  std::vector<SimplePosBias>& posBias(salmon::utils::Direction dir) {
    return (dir == salmon::utils::Direction::FORWARD) ? posBiasFW_ : posBiasRC_;
  }
  const std::vector<SimplePosBias>&
  posBias(salmon::utils::Direction dir) const {
    return (dir == salmon::utils::Direction::FORWARD) ? posBiasFW_ : posBiasRC_;
  }

  std::vector<SimplePosBias>& posBiasExpected(salmon::utils::Direction dir) {
    return (dir == salmon::utils::Direction::FORWARD) ? posBiasExpectFW_
                                                      : posBiasExpectRC_;
  }

  const std::vector<SimplePosBias>&
  posBiasExpected(salmon::utils::Direction dir) const {
    return (dir == salmon::utils::Direction::FORWARD) ? posBiasExpectFW_
                                                      : posBiasExpectRC_;
  }

  ReadKmerDist<6, std::atomic<uint32_t>>&
  readBias(salmon::utils::Direction dir) {
    return (dir == salmon::utils::Direction::FORWARD) ? readBias_[0]
                                                      : readBias_[1];
  }
  const ReadKmerDist<6, std::atomic<uint32_t>>&
  readBias(salmon::utils::Direction dir) const {
    return (dir == salmon::utils::Direction::FORWARD) ? readBias_[0]
                                                      : readBias_[1];
  }

  SBModel& readBiasModelObserved(salmon::utils::Direction dir) {
    return (dir == salmon::utils::Direction::FORWARD)
               ? readBiasModelObserved_[0]
               : readBiasModelObserved_[1];
  }
  const SBModel& readBiasModelObserved(salmon::utils::Direction dir) const {
    return (dir == salmon::utils::Direction::FORWARD)
               ? readBiasModelObserved_[0]
               : readBiasModelObserved_[1];
  }

  SBModel& readBiasModelExpected(salmon::utils::Direction dir) {
    return (dir == salmon::utils::Direction::FORWARD)
               ? readBiasModelExpected_[0]
               : readBiasModelExpected_[1];
  }
  const SBModel& readBiasModelExpected(salmon::utils::Direction dir) const {
    return (dir == salmon::utils::Direction::FORWARD)
               ? readBiasModelExpected_[0]
               : readBiasModelExpected_[1];
  }

  void setReadBiasModelExpected(SBModel&& model, salmon::utils::Direction dir) {
    size_t idx = (dir == salmon::utils::Direction::FORWARD) ? 0 : 1;
    readBiasModelExpected_[idx] = std::move(model);
  }

  const std::vector<uint32_t>& getLengthQuantiles() const {
    return lengthQuantiles_;
  }

  uint64_t getNumDecoys() const {
    return numDecoys_;
  }

  salmon::utils::DuplicateTargetStatus index_retains_duplicates() const { 
    return salmon::utils::DuplicateTargetStatus::UNKNOWN; 
  }

private:

  void setTranscriptLengthClasses_(std::vector<uint32_t>& lengths,
                                   size_t nbins) {
    auto n = lengths.size();
    if (n > nbins) {
      lengthQuantiles_.clear();
      lengthQuantiles_.reserve(nbins);

      size_t step = lengths.size() / nbins;
      size_t cumStep = 0;
      for (size_t i = 0; i < nbins; ++i) {
        cumStep += step;
        size_t ind = std::min(cumStep, n - 1);
        std::nth_element(lengths.begin(), lengths.begin() + ind, lengths.end());
        // Find the proper quantile
        lengthQuantiles_.push_back(*(lengths.begin() + ind));
      }
    } else {
      lengthQuantiles_.clear();
      lengthQuantiles_.reserve(n);
      std::sort(lengths.begin(), lengths.end());
      for (auto l : lengths) {
        lengthQuantiles_.push_back(l);
      }
      posBiasFW_.resize(n);
      posBiasRC_.resize(n);
      posBiasExpectFW_.resize(n);
      posBiasExpectRC_.resize(n);
    }

    auto qb = lengthQuantiles_.begin();
    auto qe = lengthQuantiles_.end();
    auto maxQuant = std::distance(qb, qe) - 1;
    for (auto& t : transcripts_) {
      auto ind = std::min(
          maxQuant, std::distance(qb, std::upper_bound(qb, qe, t.RefLength)));
      // the index is the smallest quantile longer than this length
      t.lengthClassIndex(ind);
    }
  }

  /**
   * The file from which the alignments will be read.
   * This can be a SAM or BAM file, and can be a regular
   * file or a fifo.
   */
  std::vector<boost::filesystem::path> alignmentFiles_;
  /**
   * The file from which the transcripts are read.
   * This is expected to be a FASTA format file.
   */
  boost::filesystem::path transcriptFile_;
  /**
   * Describes the expected format of the sequencing
   * fragment library.
   */
  LibraryFormat libFmt_;
  /**
   * The targets (transcripts) to be quantified.
   */
  std::vector<Transcript> transcripts_;
  /**
   * A pointer to the queue from which the fragments
   * will be read.
   */
  // std::unique_ptr<t_pool, std::function<void(t_pool*)>> threadPool_;
  std::unique_ptr<BAMQueue<FragT>> bq;

  // SequenceBiasModel seqBiasModel_;

  /**
   * The cluster forest maintains the dynamic relationship
   * defined by transcripts and reads --- if two transcripts
   * share an ambiguously mapped read, then they are placed
   * in the same cluster.
   */
  std::unique_ptr<ClusterForest> clusters_;

  /**
   * The emperical fragment start position distribution
   */
  std::vector<FragmentStartPositionDistribution> fragStartDists_;

  /**
   * The emperical fragment length distribution.
   *
   */
  std::unique_ptr<FragmentLengthDistribution> flDist_;
  /**
   * The emperical error model
   */
  std::unique_ptr<AlignModelT> alnMod_;

  /** Keeps track of the number of passes that have been
   *  made through the alignment file.
   */
  size_t quantificationPasses_;
  SpinLock sl_;
  EQBuilderT eqBuilder_;

  /** Positional bias things**/
  std::vector<uint32_t> lengthQuantiles_;
  std::vector<SimplePosBias> posBiasFW_;
  std::vector<SimplePosBias> posBiasRC_;
  std::vector<SimplePosBias> posBiasExpectFW_;
  std::vector<SimplePosBias> posBiasExpectRC_;

  /** GC-fragment bias things **/
  // One bin for each percentage GC content
  double gcFracFwd_;
  GCFragModel observedGC_;
  GCFragModel expectedGC_;

  // Since multiple threads can touch this dist, we
  // need atomic counters.
  std::array<ReadKmerDist<6, std::atomic<uint32_t>>, 2> readBias_;
  std::array<SBModel, 2> readBiasModelObserved_;
  std::array<SBModel, 2> readBiasModelExpected_;

  // ReadKmerDist<6, std::atomic<uint32_t>> readBias_;
  std::vector<double> expectedBias_;
  std::unique_ptr<LibraryTypeDetector> detector_{nullptr};
  std::vector<double> conditionalMeans_;

  salmon::bam_utils::AlignerDetails aligner_{salmon::bam_utils::AlignerDetails::UNKNOWN};
  uint64_t numDecoys_{0};
};

#endif // ALIGNMENT_LIBRARY_HPP