forked from smithlabcode/preseq
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathpreseq.cpp
1755 lines (1478 loc) · 58.2 KB
/
preseq.cpp
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
/*
* preseq: to predict properties of genomic sequencing libraries
*
* Copyright (C) 2013-2015 University of Southern California and
* Andrew D. Smith and Timothy Daley
*
* Authors: Timothy Daley, Chao Deng, Victoria Helus, and Andrew Smith
*
* This program is free software: you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* This program is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with this program. If not, see <http://www.gnu.org/licenses/>.
*/
#include <fstream>
#include <numeric>
#include <vector>
#include <iomanip>
#include <queue>
#include <sys/types.h>
#include <unistd.h>
#include <cstring>
#include <tr1/unordered_map>
#include <cmath>
#include <fstream>
#include <iostream>
#include <sstream>
#include <gsl/gsl_cdf.h>
#include <gsl/gsl_randist.h>
#include <gsl/gsl_statistics_double.h>
#include <gsl/gsl_sf_gamma.h>
#include <OptionParser.hpp>
#include <smithlab_utils.hpp>
#include <GenomicRegion.hpp>
#include <RNG.hpp>
#include <smithlab_os.hpp>
#define PRESEQ_VERSION "2.0.3"
// AS: might not be good to depend on mapped read here
// TD: if we're including gc_extrap, we need the dependence
#include "continued_fraction.hpp"
#include "load_data_for_complexity.hpp"
#include "moment_sequence.hpp"
using std::string;
using std::min;
using std::vector;
using std::endl;
using std::cerr;
using std::max;
using std::ifstream;
using std::isfinite;
using std::setw;
using std::fixed;
using std::setprecision;
using std::tr1::unordered_map;
/////////////////////////////////////////////////////////
// Confidence interval stuff
/*
static inline double
alpha_log_confint_multiplier(const double estimate,
const double variance, const double alpha) {
const double inv_norm_alpha = gsl_cdf_ugaussian_Qinv(alpha/2.0);
return exp(inv_norm_alpha*
sqrt(log(1.0 + variance/pow(estimate, 2))));
}
*/
static void
median_and_ci(const vector<double> &estimates,
const double ci_level,
double &median_estimate,
double &lower_ci_estimate,
double &upper_ci_estimate){
assert(!estimates.empty());
const double alpha = 1.0 - ci_level;
const size_t n_est = estimates.size();
vector<double> sorted_estimates(estimates);
sort(sorted_estimates.begin(), sorted_estimates.end());
median_estimate =
gsl_stats_median_from_sorted_data(&sorted_estimates[0],
1, n_est);
lower_ci_estimate =
gsl_stats_quantile_from_sorted_data(&sorted_estimates[0],
1, n_est, alpha/2);
upper_ci_estimate =
gsl_stats_quantile_from_sorted_data(&sorted_estimates[0],
1, n_est, 1.0 - alpha/2);
}
static void
vector_median_and_ci(const vector<vector<double> > &bootstrap_estimates,
const double ci_level,
vector<double> &yield_estimates,
vector<double> &lower_ci_lognormal,
vector<double> &upper_ci_lognormal) {
yield_estimates.clear();
lower_ci_lognormal.clear();
upper_ci_lognormal.clear();
assert(!bootstrap_estimates.empty());
const size_t n_est = bootstrap_estimates.size();
vector<double> estimates_row(bootstrap_estimates.size(), 0.0);
for (size_t i = 0; i < bootstrap_estimates[0].size(); i++) {
// estimates is in wrong order, work locally on const val
for (size_t k = 0; k < n_est; ++k)
estimates_row[k] = bootstrap_estimates[k][i];
double median_estimate, lower_ci_estimate, upper_ci_estimate;
median_and_ci(estimates_row, ci_level, median_estimate,
lower_ci_estimate, upper_ci_estimate);
sort(estimates_row.begin(), estimates_row.end());
yield_estimates.push_back(median_estimate);
lower_ci_lognormal.push_back(lower_ci_estimate);
upper_ci_lognormal.push_back(upper_ci_estimate);
}
}
void
log_mean(const bool VERBOSE,
const vector<double> &estimates,
const double c_level,
double &log_mean,
double &log_lower_ci,
double &log_upper_ci){
vector<double> log_estimates(estimates);
for(size_t i = 0; i < log_estimates.size(); i++)
log_estimates[i] = log(log_estimates[i]);
log_mean = exp(gsl_stats_mean(&log_estimates[0], 1,
log_estimates.size()) );
double log_std_dev = std::sqrt(gsl_stats_variance(&log_estimates[0], 1,
log_estimates.size()) );
const double inv_norm_alpha = gsl_cdf_ugaussian_Qinv((1.0 - c_level)/2.0);
log_lower_ci = exp(log(log_mean) - inv_norm_alpha*log_std_dev);
log_upper_ci = exp(log(log_mean) + inv_norm_alpha*log_std_dev);
}
void
mean_and_ci(const vector<double> &estimates,
const double ci_level,
double &mean_estimate,
double &lower_ci_estimate,
double &upper_ci_estimate){
assert(!estimates.empty());
const double alpha = 1.0 - ci_level;
const size_t n_est = estimates.size();
vector<double> sorted_estimates(estimates);
sort(sorted_estimates.begin(), sorted_estimates.end());
mean_estimate =
gsl_stats_mean(&sorted_estimates[0], 1, n_est);
lower_ci_estimate =
gsl_stats_quantile_from_sorted_data(&sorted_estimates[0],
1, n_est, alpha/2);
upper_ci_estimate =
gsl_stats_quantile_from_sorted_data(&sorted_estimates[0],
1, n_est, 1.0 - alpha/2);
}
////////////////////////////////////////////////////////////////////////
////////////////////////////////////////////////////////////////////////
////////////////////////////////////////////////////////////////////////
////////////////////////////////////////////////////////////////////////
/////
///// EXTRAP MODE BELOW HERE
/////
// vals_hist[j] = n_{j} = # (counts = j)
// vals_hist_distinct_counts[k] = kth index j s.t. vals_hist[j] > 0
// stores kth index of vals_hist that is positive
// distinct_counts_hist[k] = vals_hist[vals_hist_distinct_counts[k]]
// stores the kth positive value of vals_hist
void
resample_hist(const gsl_rng *rng, const vector<size_t> &vals_hist_distinct_counts,
const vector<double> &distinct_counts_hist,
vector<double> &out_hist) {
vector<unsigned int> sample_distinct_counts_hist(distinct_counts_hist.size(), 0);
const unsigned int distinct =
static_cast<unsigned int>(accumulate(distinct_counts_hist.begin(),
distinct_counts_hist.end(), 0.0));
gsl_ran_multinomial(rng, distinct_counts_hist.size(), distinct,
&distinct_counts_hist.front(),
&sample_distinct_counts_hist.front());
out_hist.clear();
out_hist.resize(vals_hist_distinct_counts.back() + 1, 0.0);
for(size_t i = 0; i < sample_distinct_counts_hist.size(); i++)
out_hist[vals_hist_distinct_counts[i]] =
static_cast<double>(sample_distinct_counts_hist[i]);
}
// interpolate by explicit calculating the expectation
// for sampling without replacement;
// see K.L Heck 1975
// N total sample size; S the total number of distincts
// n sub sample size
static double
interpolate_distinct(vector<double> &hist, size_t N,
size_t S, const size_t n) {
double denom = gsl_sf_lngamma(N + 1) - gsl_sf_lngamma(n + 1) - gsl_sf_lngamma(N - n + 1);
vector<double> numer(hist.size(), 0);
for (size_t i = 1; i < hist.size(); i++) {
// N - i -n + 1 should be greater than 0
if (N < i + n) {
numer[i] = 0;
} else {
numer[i] = gsl_sf_lngamma(N - i + 1) - gsl_sf_lngamma(n + 1) - gsl_sf_lngamma(N - i - n + 1);
numer[i] = exp(numer[i] - denom) * hist[i];
}
}
return S - accumulate(numer.begin(), numer.end(), 0);
}
// check if estimates are finite, increasing, and concave
static bool
check_yield_estimates(const vector<double> &estimates) {
if (estimates.empty())
return false;
// make sure that the estimate is increasing in the time_step and is
// below the initial distinct per step_size
if (!isfinite(accumulate(estimates.begin(), estimates.end(), 0.0)))
return false;
for (size_t i = 1; i < estimates.size(); ++i)
if ((estimates[i] < estimates[i - 1]) ||
(i >= 2 && (estimates[i] - estimates[i - 1] >
estimates[i - 1] - estimates[i - 2])) ||
(estimates[i] < 0.0))
return false;
return true;
}
void
extrap_bootstrap(const bool VERBOSE, const bool DEFECTS,
const unsigned long int seed,
const vector<double> &orig_hist,
const size_t bootstraps, const size_t orig_max_terms,
const int diagonal, const double bin_step_size,
const double max_extrapolation, const size_t max_iter,
vector< vector<double> > &bootstrap_estimates) {
// clear returning vectors
bootstrap_estimates.clear();
//setup rng
srand(time(0) + getpid());
gsl_rng_env_setup();
gsl_rng *rng = gsl_rng_alloc(gsl_rng_default);
gsl_rng_set(rng, seed);
double vals_sum = 0.0;
for(size_t i = 0; i < orig_hist.size(); i++)
vals_sum += orig_hist[i]*i;
const double initial_distinct
= accumulate(orig_hist.begin(), orig_hist.end(), 0.0);
vector<size_t> orig_hist_distinct_counts;
vector<double> distinct_orig_hist;
for (size_t i = 0; i < orig_hist.size(); i++){
if (orig_hist[i] > 0) {
orig_hist_distinct_counts.push_back(i);
distinct_orig_hist.push_back(orig_hist[i]);
}
}
for (size_t iter = 0;
(iter < max_iter && bootstrap_estimates.size() < bootstraps);
++iter) {
vector<double> yield_vector;
vector<double> hist;
resample_hist(rng, orig_hist_distinct_counts, distinct_orig_hist, hist);
double sample_vals_sum = 0.0;
for(size_t i = 0; i < hist.size(); i++)
sample_vals_sum += i*hist[i];
//resize boot_hist to remove excess zeros
while (hist.back() == 0)
hist.pop_back();
// compute complexity curve by random sampling w/out replacement
const size_t upper_limit = static_cast<size_t>(sample_vals_sum);
const size_t distinct = static_cast<size_t>(accumulate(hist.begin(), hist.end(), 0.0));
const size_t step = static_cast<size_t>(bin_step_size);
size_t sample = step;
while(sample < upper_limit){
yield_vector.push_back(interpolate_distinct(hist, upper_limit, distinct, sample));
sample += step;
}
// ENSURE THAT THE MAX TERMS ARE ACCEPTABLE
size_t counts_before_first_zero = 1;
while (counts_before_first_zero < hist.size() &&
hist[counts_before_first_zero] > 0)
++counts_before_first_zero;
size_t max_terms = std::min(orig_max_terms, counts_before_first_zero - 1);
// refit curve for lower bound (degree of approx is 1 less than
// max_terms)
max_terms = max_terms - (max_terms % 2 == 1);
// defect mode, simple extrapolation
if(DEFECTS){
vector<double> ps_coeffs;
for (size_t j = 1; j <= max_terms; j++)
ps_coeffs.push_back(hist[j]*std::pow((double)(-1), (int)(j + 1)) );
const ContinuedFraction
defect_cf(ps_coeffs, diagonal, max_terms);
double sample_size = static_cast<double>(sample);
while(sample_size < max_extrapolation){
double t = (sample_size - sample_vals_sum)/sample_vals_sum;
assert(t >= 0.0);
yield_vector.push_back(initial_distinct + t*defect_cf(t));
sample_size += bin_step_size;
}
// no checking of curve in defect mode
bootstrap_estimates.push_back(yield_vector);
if (VERBOSE) cerr << '.';
}
else{
//refit curve for lower bound
const ContinuedFractionApproximation
lower_cfa(diagonal, max_terms);
const ContinuedFraction
lower_cf(lower_cfa.optimal_cont_frac_distinct(hist));
//extrapolate the curve start
if (lower_cf.is_valid()){
double sample_size = static_cast<double>(sample);
while(sample_size < max_extrapolation){
double t = (sample_size - sample_vals_sum)/sample_vals_sum;
assert(t >= 0.0);
yield_vector.push_back(initial_distinct + t*lower_cf(t));
sample_size += bin_step_size;
}
// SANITY CHECK
if (check_yield_estimates(yield_vector)) {
bootstrap_estimates.push_back(yield_vector);
if (VERBOSE) cerr << '.';
}
else if (VERBOSE){
cerr << "_";
}
}
else if (VERBOSE){
cerr << "_";
}
}
}
if (VERBOSE)
cerr << endl;
if (bootstrap_estimates.size() < bootstraps)
throw SMITHLABException("too many defects in the approximation, consider running in defect mode");
}
static bool
extrap_single_estimate(const bool VERBOSE, const bool DEFECTS,
vector<double> &hist,
size_t max_terms, const int diagonal,
const double step_size,
const double max_extrapolation,
vector<double> &yield_estimate) {
yield_estimate.clear();
double vals_sum = 0.0;
for(size_t i = 0; i < hist.size(); i++)
vals_sum += i*hist[i];
const double initial_distinct
= accumulate(hist.begin(), hist.end(), 0.0);
// interpolate complexity curve by random sampling w/out replacement
size_t upper_limit = static_cast<size_t>(vals_sum);
size_t step = static_cast<size_t>(step_size);
size_t sample = step;
while (sample < upper_limit){
yield_estimate.push_back(
interpolate_distinct(hist, upper_limit,
static_cast<size_t>(initial_distinct), sample));
sample += step;
}
// ENSURE THAT THE MAX TERMS ARE ACCEPTABLE
size_t counts_before_first_zero = 1;
while (counts_before_first_zero < hist.size() &&
hist[counts_before_first_zero] > 0)
++counts_before_first_zero;
// Ensure we are not using a zero term
max_terms = std::min(max_terms, counts_before_first_zero - 1);
// refit curve for lower bound (degree of approx is 1 less than
// max_terms)
max_terms = max_terms - (max_terms % 2 == 1);
if(DEFECTS){
vector<double> ps_coeffs;
for (size_t j = 1; j <= max_terms; j++)
ps_coeffs.push_back(hist[j]*std::pow((double)(-1), (int)(j + 1)) );
const ContinuedFraction
defect_cf(ps_coeffs, diagonal, max_terms);
double sample_size = static_cast<double>(sample);
while(sample_size < max_extrapolation){
const double one_minus_fold_extrap
= (sample_size - vals_sum)/vals_sum;
assert(one_minus_fold_extrap >= 0.0);
double tmp = one_minus_fold_extrap*defect_cf(one_minus_fold_extrap);
yield_estimate.push_back(initial_distinct + tmp);
sample_size += step_size;
}
if (VERBOSE) {
if(defect_cf.offset_coeffs.size() > 0){
cerr << "CF_OFFSET_COEFF_ESTIMATES" << endl;
copy(defect_cf.offset_coeffs.begin(), defect_cf.offset_coeffs.end(),
std::ostream_iterator<double>(cerr, "\n"));
}
if(defect_cf.cf_coeffs.size() > 0){
cerr << "CF_COEFF_ESTIMATES" << endl;
copy(defect_cf.cf_coeffs.begin(), defect_cf.cf_coeffs.end(),
std::ostream_iterator<double>(cerr, "\n"));
}
}
// NO FAIL! DEFECT MODE DOESN'T CARE ABOUT FAILURE
}
else{
const ContinuedFractionApproximation
lower_cfa(diagonal, max_terms);
const ContinuedFraction
lower_cf(lower_cfa.optimal_cont_frac_distinct(hist));
// extrapolate curve
if (lower_cf.is_valid()){
double sample_size = static_cast<double>(sample);
while(sample_size < max_extrapolation){
const double one_minus_fold_extrap
= (sample_size - vals_sum)/vals_sum;
assert(one_minus_fold_extrap >= 0.0);
double tmp = one_minus_fold_extrap*lower_cf(one_minus_fold_extrap);
yield_estimate.push_back(initial_distinct + tmp);
sample_size += step_size;
}
}
else{
// FAIL!
// lower_cf unacceptable, need to bootstrap to obtain estimates
return false;
}
if (VERBOSE) {
if(lower_cf.offset_coeffs.size() > 0){
cerr << "CF_OFFSET_COEFF_ESTIMATES" << endl;
copy(lower_cf.offset_coeffs.begin(), lower_cf.offset_coeffs.end(),
std::ostream_iterator<double>(cerr, "\n"));
}
if(lower_cf.cf_coeffs.size() > 0){
cerr << "CF_COEFF_ESTIMATES" << endl;
copy(lower_cf.cf_coeffs.begin(), lower_cf.cf_coeffs.end(),
std::ostream_iterator<double>(cerr, "\n"));
}
}
}
// SUCCESS!!
return true;
}
static double
GoodToulmin2xExtrap(const vector<double> &counts_hist){
double two_fold_extrap = 0.0;
for(size_t i = 0; i < counts_hist.size(); i++)
two_fold_extrap += pow(-1.0, i + 1)*counts_hist[i];
return two_fold_extrap;
}
static void
write_predicted_complexity_curve(const string outfile,
const double c_level, const double step_size,
const vector<double> &yield_estimates,
const vector<double> &yield_lower_ci_lognormal,
const vector<double> &yield_upper_ci_lognormal) {
std::ofstream of;
if (!outfile.empty()) of.open(outfile.c_str());
std::ostream out(outfile.empty() ? std::cout.rdbuf() : of.rdbuf());
out << "TOTAL_READS\tEXPECTED_DISTINCT\t"
<< "LOWER_" << c_level << "CI\t"
<< "UPPER_" << c_level << "CI" << endl;
out.setf(std::ios_base::fixed, std::ios_base::floatfield);
out.precision(1);
out << 0 << '\t' << 0 << '\t' << 0 << '\t' << 0 << endl;
for (size_t i = 0; i < yield_estimates.size(); ++i)
out << (i + 1)*step_size << '\t'
<< yield_estimates[i] << '\t'
<< yield_lower_ci_lognormal[i] << '\t'
<< yield_upper_ci_lognormal[i] << endl;
}
static void
write_predicted_coverage_curve(const string outfile,
const double c_level,
const double base_step_size,
const size_t bin_size,
const vector<double> &coverage_estimates,
const vector<double> &coverage_lower_ci_lognormal,
const vector<double> &coverage_upper_ci_lognormal) {
std::ofstream of;
if (!outfile.empty()) of.open(outfile.c_str());
std::ostream out(outfile.empty() ? std::cout.rdbuf() : of.rdbuf());
out << "TOTAL_BASES\tEXPECTED_COVERED_BASES\t"
<< "LOWER_" << 100*c_level << "%CI\t"
<< "UPPER_" << 100*c_level << "%CI" << endl;
out.setf(std::ios_base::fixed, std::ios_base::floatfield);
out.precision(1);
out << 0 << '\t' << 0 << '\t' << 0 << '\t' << 0 << endl;
for (size_t i = 0; i < coverage_estimates.size(); ++i)
out << (i + 1)*base_step_size << '\t'
<< coverage_estimates[i]*bin_size << '\t'
<< coverage_lower_ci_lognormal[i]*bin_size << '\t'
<< coverage_upper_ci_lognormal[i]*bin_size << endl;
}
static int
lc_extrap(const int argc, const char **argv) {
try {
const size_t MIN_REQUIRED_COUNTS = 4;
/* FILES */
string outfile;
size_t orig_max_terms = 100;
double max_extrapolation = 1.0e10;
double step_size = 1e6;
size_t bootstraps = 100;
int diagonal = 0;
double c_level = 0.95;
unsigned long int seed = 0;
/* FLAGS */
bool VERBOSE = false;
bool VALS_INPUT = false;
bool PAIRED_END = false;
bool HIST_INPUT = false;
bool SINGLE_ESTIMATE = false;
bool DEFECTS = false;
#ifdef HAVE_SAMTOOLS
bool BAM_FORMAT_INPUT = false;
size_t MAX_SEGMENT_LENGTH = 5000;
#endif
/********** GET COMMAND LINE ARGUMENTS FOR LC EXTRAP ***********/
OptionParser opt_parse(strip_path(argv[1]),
"", "<sorted-bed-file>");
opt_parse.add_opt("output", 'o', "yield output file (default: stdout)",
false , outfile);
opt_parse.add_opt("extrap",'e',"maximum extrapolation "
"(default: " + toa(max_extrapolation) + ")",
false, max_extrapolation);
opt_parse.add_opt("step",'s',"step size in extrapolations "
"(default: " + toa(step_size) + ")",
false, step_size);
opt_parse.add_opt("bootstraps",'n',"number of bootstraps "
"(default: " + toa(bootstraps) + "), ",
false, bootstraps);
opt_parse.add_opt("cval", 'c', "level for confidence intervals "
"(default: " + toa(c_level) + ")", false, c_level);
opt_parse.add_opt("terms",'x',"maximum number of terms", false,
orig_max_terms);
opt_parse.add_opt("verbose", 'v', "print more information",
false, VERBOSE);
#ifdef HAVE_SAMTOOLS
opt_parse.add_opt("bam", 'B', "input is in BAM format",
false, BAM_FORMAT_INPUT);
opt_parse.add_opt("seg_len", 'l', "maximum segment length when merging "
"paired end bam reads (default: "
+ toa(MAX_SEGMENT_LENGTH) + ")",
false, MAX_SEGMENT_LENGTH);
#endif
opt_parse.add_opt("pe", 'P', "input is paired end read file",
false, PAIRED_END);
opt_parse.add_opt("vals", 'V',
"input is a text file containing only the observed counts",
false, VALS_INPUT);
opt_parse.add_opt("hist", 'H',
"input is a text file containing the observed histogram",
false, HIST_INPUT);
opt_parse.add_opt("quick",'Q',
"quick mode, estimate yield without bootstrapping for confidence intervals",
false, SINGLE_ESTIMATE);
opt_parse.add_opt("defects", 'D',
"defects mode to extrapolate without testing for defects",
false, DEFECTS);
opt_parse.add_opt("seed", 'r', "seed for random number generator",
false, seed);
vector<string> leftover_args;
opt_parse.parse(argc-1, argv+1, leftover_args);
if (argc == 2 || opt_parse.help_requested()) {
cerr << opt_parse.help_message() << endl;
return EXIT_SUCCESS;
}
if (opt_parse.about_requested()) {
cerr << opt_parse.about_message() << endl;
return EXIT_SUCCESS;
}
if (opt_parse.option_missing()) {
cerr << opt_parse.option_missing_message() << endl;
return EXIT_SUCCESS;
}
if (leftover_args.empty()) {
cerr << opt_parse.help_message() << endl;
return EXIT_SUCCESS;
}
const string input_file_name = leftover_args.front();
/******************************************************************/
// if seed is not set, make it random
if(seed == 0){
seed = rand();
}
vector<double> counts_hist;
size_t n_reads = 0;
// LOAD VALUES
if(HIST_INPUT){
if(VERBOSE)
cerr << "HIST_INPUT" << endl;
n_reads = load_histogram(input_file_name, counts_hist);
}
else if(VALS_INPUT){
if(VERBOSE)
cerr << "VALS_INPUT" << endl;
n_reads = load_counts(input_file_name, counts_hist);
}
#ifdef HAVE_SAMTOOLS
else if (BAM_FORMAT_INPUT && PAIRED_END){
if(VERBOSE)
cerr << "PAIRED_END_BAM_INPUT" << endl;
const size_t MAX_READS_TO_HOLD = 5000000;
size_t n_paired = 0;
size_t n_mates = 0;
n_reads = load_counts_BAM_pe(VERBOSE, input_file_name,
MAX_SEGMENT_LENGTH,
MAX_READS_TO_HOLD, n_paired,
n_mates, counts_hist);
if(VERBOSE){
cerr << "MERGED PAIRED END READS = " << n_paired << endl;
cerr << "MATES PROCESSED = " << n_mates << endl;
}
}
else if(BAM_FORMAT_INPUT){
if(VERBOSE)
cerr << "BAM_INPUT" << endl;
n_reads = load_counts_BAM_se(input_file_name, counts_hist);
}
#endif
else if(PAIRED_END){
if(VERBOSE)
cerr << "PAIRED_END_BED_INPUT" << endl;
n_reads = load_counts_BED_pe(input_file_name, counts_hist);
}
else{ // default is single end bed file
if(VERBOSE)
cerr << "BED_INPUT" << endl;
n_reads = load_counts_BED_se(input_file_name, counts_hist);
}
const size_t max_observed_count = counts_hist.size() - 1;
const double distinct_reads = accumulate(counts_hist.begin(),
counts_hist.end(), 0.0);
// ENSURE THAT THE MAX TERMS ARE ACCEPTABLE
size_t counts_before_first_zero = 1;
while (counts_before_first_zero < counts_hist.size() &&
counts_hist[counts_before_first_zero] > 0)
++counts_before_first_zero;
orig_max_terms = std::min(orig_max_terms, counts_before_first_zero - 1);
orig_max_terms = orig_max_terms - (orig_max_terms % 2 == 1);
const size_t distinct_counts =
static_cast<size_t>(std::count_if(counts_hist.begin(), counts_hist.end(),
bind2nd(std::greater<double>(), 0.0)));
if (VERBOSE)
cerr << "TOTAL READS = " << n_reads << endl
<< "DISTINCT READS = " << distinct_reads << endl
<< "DISTINCT COUNTS = " << distinct_counts << endl
<< "MAX COUNT = " << max_observed_count << endl
<< "COUNTS OF 1 = " << counts_hist[1] << endl
<< "MAX TERMS = " << orig_max_terms << endl;
if (VERBOSE) {
// OUTPUT THE ORIGINAL HISTOGRAM
cerr << "OBSERVED COUNTS (" << counts_hist.size() << ")" << endl;
for (size_t i = 0; i < counts_hist.size(); i++)
if (counts_hist[i] > 0)
cerr << i << '\t' << static_cast<size_t>(counts_hist[i]) << endl;
cerr << endl;
}
// check to make sure library is not overly saturated
const double two_fold_extrap = GoodToulmin2xExtrap(counts_hist);
if(two_fold_extrap < 0.0)
throw SMITHLABException("Library expected to saturate in doubling of "
"size, unable to extrapolate");
size_t total_reads = 0;
for(size_t i = 0; i < counts_hist.size(); i++){
total_reads += i*counts_hist[i];
}
//assert(total_reads == n_reads);
// catch if all reads are distinct
if (orig_max_terms < MIN_REQUIRED_COUNTS)
throw SMITHLABException("max count before zero is les than min required "
"count (4), sample not sufficiently deep or "
"duplicates removed");
/////////////////////////////////////////////////////////////////////
/////////////////////////////////////////////////////////////////////
/////////////////////////////////////////////////////////////////////
// ESTIMATE COMPLEXITY CURVE
if(VERBOSE)
cerr << "[ESTIMATING YIELD CURVE]" << endl;
vector<double> yield_estimates;
if(SINGLE_ESTIMATE){
bool SINGLE_ESTIMATE_SUCCESS =
extrap_single_estimate(VERBOSE, DEFECTS, counts_hist, orig_max_terms,
diagonal, step_size, max_extrapolation,
yield_estimates);
// IF FAILURE, EXIT
if(!SINGLE_ESTIMATE_SUCCESS)
throw SMITHLABException("SINGLE ESTIMATE FAILED, NEED TO RUN "
"FULL MODE FOR ESTIMATES");
std::ofstream of;
if (!outfile.empty()) of.open(outfile.c_str());
std::ostream out(outfile.empty() ? std::cout.rdbuf() : of.rdbuf());
out << "TOTAL_READS\tEXPECTED_DISTINCT" << endl;
out.setf(std::ios_base::fixed, std::ios_base::floatfield);
out.precision(1);
out << 0 << '\t' << 0 << endl;
for (size_t i = 0; i < yield_estimates.size(); ++i)
out << (i + 1)*step_size << '\t'
<< yield_estimates[i] << endl;
}
else{
if (VERBOSE)
cerr << "[BOOTSTRAPPING HISTOGRAM]" << endl;
const size_t max_iter = 10*bootstraps;
vector<vector <double> > bootstrap_estimates;
extrap_bootstrap(VERBOSE, DEFECTS, seed, counts_hist, bootstraps,
orig_max_terms, diagonal, step_size, max_extrapolation,
max_iter, bootstrap_estimates);
/////////////////////////////////////////////////////////////////////
if (VERBOSE)
cerr << "[COMPUTING CONFIDENCE INTERVALS]" << endl;
// yield ci
vector<double> yield_upper_ci_lognormal, yield_lower_ci_lognormal;
vector_median_and_ci(bootstrap_estimates, c_level, yield_estimates,
yield_lower_ci_lognormal, yield_upper_ci_lognormal);
/////////////////////////////////////////////////////////////////////
if (VERBOSE)
cerr << "[WRITING OUTPUT]" << endl;
write_predicted_complexity_curve(outfile, c_level, step_size,
yield_estimates, yield_lower_ci_lognormal,
yield_upper_ci_lognormal);
}
}
catch (SMITHLABException &e) {
cerr << "ERROR:\t" << e.what() << endl;
return EXIT_FAILURE;
}
catch (std::bad_alloc &ba) {
cerr << "ERROR: could not allocate memory" << endl;
return EXIT_FAILURE;
}
return EXIT_SUCCESS;
}
///////////////////////////////////////////////////////////////////////
//////////////////////////////////////////////////////////////////////
///// GC_EXTRAP: predicting genomic coverage
/////
static int
gc_extrap(const int argc, const char **argv) {
try {
const size_t MIN_REQUIRED_COUNTS = 4;
int diagonal = 0;
size_t orig_max_terms = 100;
size_t bin_size = 10;
bool VERBOSE = false;
string outfile;
double base_step_size = 1.0e8;
size_t max_width = 10000;
bool SINGLE_ESTIMATE = false;
double max_extrapolation = 1.0e12;
size_t bootstraps = 100;
unsigned long int seed = 0;
bool DEFECTS = false;
bool NO_SEQUENCE = false;
double c_level = 0.95;
// ********* GET COMMAND LINE ARGUMENTS FOR GC EXTRAP **********
OptionParser opt_parse(strip_path(argv[1]),
"", "<sorted-mapped-read-file>");
opt_parse.add_opt("output", 'o', "coverage yield output file (default: stdout)",
false , outfile);
opt_parse.add_opt("max_width", 'w', "max fragment length, "
"set equal to read length for single end reads",
false, max_width);
opt_parse.add_opt("bin_size", 'b', "bin size "
"(default: " + toa(bin_size) + ")",
false, bin_size);
opt_parse.add_opt("extrap",'e',"maximum extrapolation in base pairs"
"(default: " + toa(max_extrapolation) + ")",
false, max_extrapolation);
opt_parse.add_opt("step",'s',"step size in bases between extrapolations "
"(default: " + toa(base_step_size) + ")",
false, base_step_size);
opt_parse.add_opt("bootstraps",'n',"number of bootstraps "
"(default: " + toa(bootstraps) + "), ",
false, bootstraps);
opt_parse.add_opt("cval", 'c', "level for confidence intervals "
"(default: " + toa(c_level) + ")", false, c_level);
opt_parse.add_opt("terms",'x',"maximum number of terms",
false, orig_max_terms);
opt_parse.add_opt("verbose", 'v', "print more information",
false, VERBOSE);
opt_parse.add_opt("bed", 'D',
"input is in bed format without sequence information",
false, NO_SEQUENCE);
opt_parse.add_opt("quick",'Q',
"quick mode: run gc_extrap without "
"bootstrapping for confidence intervals",
false, SINGLE_ESTIMATE);
opt_parse.add_opt("defects", 'D',
"defects mode to extrapolate without testing for defects",
false, DEFECTS);
opt_parse.add_opt("seed", 'r', "seed for random number generator",
false, seed);
vector<string> leftover_args;
opt_parse.parse(argc-1, argv+1, leftover_args);
if (argc == 2 || opt_parse.help_requested()) {
cerr << opt_parse.help_message() << endl;
return EXIT_SUCCESS;
}
if (opt_parse.about_requested()) {
cerr << opt_parse.about_message() << endl;
return EXIT_SUCCESS;
}
if (opt_parse.option_missing()) {
cerr << opt_parse.option_missing_message() << endl;
return EXIT_SUCCESS;
}
if (leftover_args.empty()) {
cerr << opt_parse.help_message() << endl;
return EXIT_SUCCESS;
}
const string input_file_name = leftover_args.front();
// ****************************************************************
// if seed is not set, set it to random
if(seed == 0){
seed = rand();
}
vector<double> coverage_hist;
size_t n_reads = 0;
if(VERBOSE)
cerr << "LOADING READS" << endl;
if(NO_SEQUENCE){
if(VERBOSE)
cerr << "BED FORMAT" << endl;
n_reads = load_coverage_counts_GR(input_file_name, bin_size,
max_width, coverage_hist);
}
else{
if(VERBOSE)
cerr << "MAPPED READ FORMAT" << endl;
n_reads = load_coverage_counts_MR(VERBOSE, input_file_name, bin_size,
max_width, coverage_hist);
}
double total_bins = 0.0;
for(size_t i = 0; i < coverage_hist.size(); i++)
total_bins += coverage_hist[i]*i;
const double distinct_bins =
accumulate(coverage_hist.begin(), coverage_hist.end(), 0.0);
const double avg_bins_per_read = total_bins/n_reads;
double bin_step_size = base_step_size/bin_size;
const size_t max_observed_count = coverage_hist.size() - 1;
// ENSURE THAT THE MAX TERMS ARE ACCEPTABLE
size_t counts_before_first_zero = 1;
while (counts_before_first_zero < coverage_hist.size() &&
coverage_hist[counts_before_first_zero] > 0)
++counts_before_first_zero;
orig_max_terms = std::min(orig_max_terms, counts_before_first_zero - 1);
if (VERBOSE)
cerr << "TOTAL READS = " << n_reads << endl
<< "BASE STEP SIZE = " << base_step_size << endl
<< "BIN STEP SIZE = " << bin_step_size << endl
<< "TOTAL BINS = " << total_bins << endl
<< "BINS PER READ = " << avg_bins_per_read << endl
<< "DISTINCT BINS = " << distinct_bins << endl
<< "TOTAL BASES = " << total_bins*bin_size << endl
<< "TOTAL COVERED BASES = " << distinct_bins*bin_size << endl
<< "MAX COVERAGE COUNT = " << max_observed_count << endl
<< "COUNTS OF 1 = " << coverage_hist[1] << endl;