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 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187 1188 1189 1190 1191 1192 1193 1194 1195 1196 1197 1198 1199 1200 1201 1202 1203 1204 1205 1206 1207 1208 1209 1210 1211 1212 1213 1214 1215 1216 1217 1218 1219 1220 1221 1222 1223 1224 1225 1226 1227 1228 1229 1230 1231 1232 1233 1234 1235 1236 1237 1238 1239 1240 1241 1242 1243 1244 1245 1246 1247 1248 1249 1250 1251 1252 1253 1254 1255 1256 1257 1258 1259 1260 1261 1262 1263 1264 1265 1266 1267 1268 1269 1270 1271 1272 1273 1274 1275 1276 1277 1278 1279 1280 1281 1282 1283 1284 1285 1286 1287 1288 1289 1290 1291 1292 1293 1294 1295 1296 1297 1298 1299 1300 1301 1302 1303 1304 1305 1306 1307 1308 1309 1310 1311 1312 1313 1314 1315 1316 1317 1318 1319 1320 1321 1322 1323 1324 1325 1326 1327 1328 1329 1330 1331 1332 1333 1334 1335 1336 1337 1338 1339 1340 1341 1342 1343 1344 1345 1346 1347 1348 1349 1350 1351 1352 1353 1354 1355 1356 1357 1358 1359 1360 1361 1362 1363 1364 1365 1366 1367 1368 1369 1370 1371 1372 1373 1374 1375 1376 1377 1378 1379 1380 1381 1382 1383 1384 1385 1386 1387 1388 1389 1390 1391 1392 1393 1394 1395 1396 1397 1398 1399 1400 1401 1402 1403 1404 1405 1406 1407 1408 1409 1410 1411 1412 1413 1414 1415 1416 1417 1418 1419 1420 1421 1422 1423 1424 1425 1426 1427 1428 1429 1430 1431 1432 1433 1434 1435 1436 1437 1438 1439 1440 1441 1442 1443 1444 1445 1446 1447 1448 1449 1450 1451 1452 1453 1454 1455 1456 1457 1458 1459 1460 1461 1462 1463 1464 1465 1466 1467 1468 1469 1470 1471 1472 1473 1474 1475 1476 1477 1478 1479 1480 1481 1482 1483 1484 1485 1486 1487 1488 1489 1490 1491 1492 1493 1494 1495 1496 1497 1498 1499 1500 1501 1502 1503 1504 1505 1506 1507 1508 1509 1510 1511 1512 1513 1514 1515 1516 1517 1518 1519 1520 1521 1522 1523 1524 1525 1526 1527 1528 1529 1530 1531 1532 1533 1534 1535 1536 1537 1538 1539 1540 1541 1542 1543 1544 1545 1546 1547 1548 1549 1550 1551 1552 1553 1554 1555 1556 1557 1558 1559 1560 1561 1562 1563 1564 1565 1566 1567 1568 1569 1570 1571 1572 1573 1574 1575 1576 1577 1578 1579 1580 1581 1582 1583 1584 1585 1586 1587 1588 1589 1590 1591 1592 1593 1594 1595 1596 1597 1598 1599 1600 1601 1602 1603 1604 1605 1606 1607 1608 1609 1610 1611 1612 1613 1614 1615 1616 1617 1618 1619 1620 1621 1622 1623 1624 1625 1626 1627 1628 1629 1630 1631 1632 1633 1634 1635 1636 1637 1638 1639 1640 1641 1642 1643 1644 1645 1646 1647 1648 1649 1650 1651 1652 1653 1654 1655 1656 1657 1658 1659 1660 1661 1662 1663 1664 1665 1666 1667 1668 1669 1670 1671 1672 1673 1674 1675 1676 1677 1678 1679 1680 1681 1682 1683 1684 1685 1686 1687 1688 1689 1690 1691 1692 1693 1694 1695 1696 1697 1698 1699 1700 1701 1702 1703 1704 1705 1706 1707 1708 1709 1710 1711 1712 1713 1714 1715 1716 1717 1718 1719 1720 1721 1722 1723 1724 1725 1726 1727 1728 1729 1730 1731 1732 1733 1734 1735 1736 1737 1738 1739 1740 1741 1742 1743 1744
|
/*
Copyright (c) 2012 Erik Aronesty
Permission is hereby granted, free of charge, to any person obtaining a copy
of this software and associated documentation files (the "Software"), to deal
in the Software without restriction, including without limitation the rights
to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
copies of the Software, and to permit persons to whom the Software is
furnished to do so, subject to the following conditions:
The above copyright notice and this permission notice shall be included in
all copies or substantial portions of the Software.
THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN
THE SOFTWARE.
*/
#include <math.h>
#include <ctype.h>
#include <stdio.h>
#include <stdlib.h>
#include <unistd.h>
#include <getopt.h>
#include <string.h>
#include <errno.h>
#include <stdarg.h>
#include <gsl/gsl_randist.h>
#include <sys/stat.h>
#include <string>
#include <queue>
#include <list>
#include <google/sparse_hash_map> // or sparse_hash_set, dense_hash_map, ...
#include <google/dense_hash_map> // or sparse_hash_set, dense_hash_map, ...
#include "tidx/tidx.h"
#include "fastq-lib.h"
#define SVNREV atoi(strchr("$Revision: 632 $", ':')+1)
const char * VERSION = "0.9";
#define MIN_READ_LEN 20
#define DEFAULT_LOCII 1000000
using namespace std;
using namespace google;
void usage(FILE *f);
// #define DEBUG 1
#define meminit(l) (memset(&l,0,sizeof(l)))
#ifdef DEBUG
#define debug(s,...) fprintf(stderr,s,##__VA_ARGS__)
#else
#define debug(s,...)
#endif
#undef warn
#define warn(s,...) ++errs; fprintf(stderr,s,##__VA_ARGS__)
#define die(s,...) (fprintf(stderr,s,##__VA_ARGS__), exit(1))
#define stat_out(s,...) fprintf(stat_fout,s,##__VA_ARGS__)
#define stdev(cnt, sum, ssq) sqrt((((double)cnt)*ssq-pow((double)sum,2)) / ((double)cnt*((double)cnt-1)))
#define log10(x) (log(x)/log(10))
double quantile(const std::vector<int> &vec, double p);
double quantile(const std::vector<double> &vec, double p);
double pnorm(double x);
double qnorm(double x);
int rand_round(double x);
// basic utils
std::vector<char *> split(char* str, const char* delim);
std::string string_format(const std::string &fmt, ...);
void to_upper(const std::string str);
void rename_tmp(std::string f);
int errs=0;
extern int optind;
class Noise {
public:
Noise() {noise=0;depth=0;};
Noise(int d, double n, double q, double mq) {depth=d; noise=n;qnoise=q;mnqual=mq;};
double noise;
double qnoise;
int depth;
double mnqual;
};
double quantile_depth(const vector<Noise> &vec, double p);
bool noisebydepth (const Noise &a, const Noise &b) { return (a.depth>b.depth);}
class PileupEnt {
public:
bool is_rev;
bool is_start;
bool f;
const char *b;
int q;
int m;
int p;
};
class vcall {
public:
vcall() {base='\0'; mn_qual=mq0=fwd=rev=qual=is_ref=qual_ssq=mq_sum=mq_ssq=tail_rev=tail_fwd=0;}
char base;
bool is_ref;
int qual, fwd, rev, mq0, mn_qual, qual_ssq, mq_sum, mq_ssq, tail_rev, tail_fwd;
vector <string> seqs;
int depth() const {return fwd+rev;}
int mq_rms() const {return sqrt(mq_ssq/depth());}
int qual_rms() const {return sqrt(qual_ssq/depth());}
};
class vfinal {
public:
vfinal(vcall &c) {max_idl_cnt=0; padj=1; pcall = &c;};
vfinal & operator=(vfinal const&x) {max_idl_seq=x.max_idl_seq; max_idl_cnt=x.max_idl_cnt; padj=x.padj; pcall=x.pcall;}
vcall *pcall;
string max_idl_seq;
int max_idl_cnt;
double padj;
bool is_indel() {return max_idl_cnt > 0;};
};
bool hitolocall (const vcall &i,const vcall &j) {return ((i.depth())>(j.depth()));}
bool sortreffirst (const vfinal &i,const vfinal &j) {return (i.pcall->is_ref&&!j.pcall->is_ref)||((i.pcall->is_ref==j.pcall->is_ref) && ((i.pcall->depth())>(j.pcall->depth())));}
class Read {
public:
int MapQ;
string Seq;
Read() {MapQ=0;};
};
class PileupReads {
public:
double MeanReadLen() {return ReadBin.size() ? TotReadLen/ReadBin.size() : MIN_READ_LEN;}
int TotReadLen;
deque<Read> ReadBin;
list<Read> ReadList;
PileupReads() {TotReadLen=0;}
};
class PileupSummary {
public:
string Chr;
int Pos;
char Base;
int Depth;
int TotQual;
int NumReads;
vector<vcall> Calls;
int SkipN;
int SkipDupReads;
int SkipMinMapq;
int SkipMinQual;
int MaxDepthByPos;
int RepeatCount;
char RepeatBase;
PileupSummary(char *line, PileupReads &reads);
PileupSummary() { Base = '\0'; Pos=-1; };
};
class PileupVisitor {
public:
char InputType;
string AnnotFile; // path to file
tidx AnnotDex; // start/stop index file
char AnnotType; // b (bed) or g (gtf - preferred)
PileupReads Reads;
PileupVisitor() {InputType ='\0';}
PileupVisitor(const char *a) {InputType ='\0'; LoadIndex(a);}
void Parse(char *dat) {PileupSummary p(dat, Reads); Visit(p);};
void LoadIndex(const char *a);
virtual void Visit(PileupSummary &dat)=0;
virtual void Finish()=0;
};
class VarStatVisitor : public PileupVisitor {
public:
VarStatVisitor() : PileupVisitor() {tot_locii=0; tot_depth=0; num_reads=0;};
void Visit(PileupSummary &dat);
void Finish() {};
public:
double tot_depth;
int tot_locii;
int num_reads;
vector<Noise> stats;
vector<Noise> ins_stats;
vector<Noise> del_stats;
};
class VarCallVisitor : public PileupVisitor {
deque<PileupSummary> Win;
void VisitX(PileupSummary &dat);
public:
int WinMax;
VarCallVisitor() : PileupVisitor() {SkippedDepth=0;WinMax=0;Hets=0;Homs=0;Locii=0;};
void Visit(PileupSummary &dat);
void Finish();
int SkippedDepth;
int Locii;
int Hets;
int Homs;
};
bool hasdata(const string &file) {
struct stat st;
if (stat(file.c_str(), &st)) {
return false;
}
return st.st_size > 0;
}
int minsampdepth=20;
double pct_depth=0;
double pct_qdepth=0;
double global_error_rate=0;
double max_phred;
int total_locii=-1;
double pct_balance=0;
char *debug_xchr=NULL;
int debug_xpos=0;
int min_depth=1;
int min_mapq=0;
int min_qual=3;
int repeat_filter=7;
double artifact_filter=1;
int min_adepth=2;
int read_tail_pct=.6;
int read_tail_len=4;
int min_idepth=3;
int no_baq=0;
double zygosity=.5; // set to .1 for 1 10% admixture, or even .05 for het/admix
void parse_bams(PileupVisitor &v, int in_n, char **in, const char *ref);
FILE *noise_f=NULL, *var_f = NULL, *varsum_f = NULL, *tgt_f = NULL, *tgtsum_f = NULL, *vcf_f = NULL, *eav_f=NULL;
double alpha=.05;
int phred=33;
double phi(double x);
FILE *openordie(const char *path, const char *mode) {
FILE *f=fopen(path, mode);
if (!f) {
warn("Can't open-%s %s: %s\n", mode, path, strerror(errno));
exit(1);
}
return f;
}
int main(int argc, char **argv) {
char c;
const char *noiseout=NULL;
const char *ref=NULL;
optind = 0;
int umindepth=0;
int uminadepth=0;
int uminidepth=0;
double upctqdepth=0;
int do_stats=0;
int do_varcall=0;
char *out_prefix = NULL;
char *target_annot = NULL;
char *read_stats = NULL;
while ( (c = getopt_long(argc, argv, "?svVBhe:m:N:x:f:p:a:g:q:Q:i:o:D:R:b:L:S:",NULL,NULL)) != -1) {
switch (c) {
case 'h': usage(stdout); return 0;
case 'm': umindepth=atoi(optarg); break;
case 'q': min_qual=atoi(optarg); break;
case 'o': out_prefix=optarg; break;
case 'Q': min_mapq=atoi(optarg); break;
case 'V': printf("Version: %s.%d\n", VERSION, SVNREV); exit(0); break;
case 'R': repeat_filter=atoi(optarg); break;
case 'A': target_annot=optarg; break;
case 'a': uminadepth=atoi(optarg);break;
case 'D': artifact_filter=atof(optarg);break;
case 'i': uminidepth=atoi(optarg);break;
case 'x': {
debug_xchr=optarg;
char *p=strrchr(debug_xchr, ':');
if (!p) die("Invalid param for -x");
*p='\0';
debug_xpos=atoi(++p);
if (!p) die("Invalid param for -x, need pos");
break;
}
case 'b': pct_balance=atof(optarg)/100.0; break;
case 'B': no_baq=1; break;
case 'p': upctqdepth=atof(optarg); break;
case 'e': alpha=atof(optarg); break;
case 'g': global_error_rate=atof(optarg); break;
case 'L': total_locii=atoi(optarg); break;
case 'f': ref=optarg; break;
case 'N': noiseout=optarg; break;
case 's': do_stats=1; break;
case 'S': read_stats=optarg; break;
case 'v': do_varcall=1; break;
case '?':
if (!optopt) {
usage(stdout); return 0;
} else if (optopt && strchr("ox", optopt))
fprintf (stderr, "Option -%c requires an argument.\n", optopt);
else if (isprint(optopt))
fprintf (stderr, "Unknown option `-%c'.\n", optopt);
else
fprintf (stderr, "Unknown option character `\\x%x'.\n", optopt);
usage(stderr);
return 1;
}
}
if (!do_stats && !do_varcall || do_stats && do_varcall) {
warn("Specify -s for stats only, or -v to do variant calling\n\n");
usage(stderr);
return 1;
}
if (out_prefix) {
if (!do_varcall) {
warn("Specify -o with -v only\n\n");
usage(stderr);
return 1;
}
var_f = openordie(string_format("%s.var.tmp", out_prefix).c_str(), "w");
vcf_f = openordie(string_format("%s.vcf.tmp", out_prefix).c_str(), "w");
eav_f = openordie(string_format("%s.eav.tmp", out_prefix).c_str(), "w");
noise_f = openordie(string_format("%s.noise.tmp", out_prefix).c_str(), "w");
varsum_f = openordie(string_format("%s.varsum.tmp", out_prefix).c_str(), "w");
if (target_annot) {
tgt_f = openordie(string_format("%s.tgt.tmp", out_prefix).c_str(), "w");
tgtsum_f = openordie(string_format("%s.tgtsum.tmp", out_prefix).c_str(), "w");
}
} else {
var_f = stdout;
varsum_f = stderr;
}
if (umindepth > minsampdepth) {
minsampdepth=umindepth;
}
if (noiseout) {
noise_f = fopen(noiseout, "w");
if (!noise_f) {
warn("Can't write %s: %s\n", noiseout, strerror(errno));
exit(1);
}
}
// set argv to '-' if stdin
const char *stdv[3] = {argv[0],"-",NULL};
if (!argv[optind]) {
argc=2;
argv = (char **) stdv;
optind=1;
}
char **in=&argv[optind];
int in_n = argc-optind;
// not really random
srand(1);
max_phred = -log10(global_error_rate)*10;
if (do_stats) {
FILE *stat_fout=stdout; // stats to stdout
if (do_varcall) // unless varcalling at the same time
stat_fout=stderr;
VarStatVisitor vstat;
parse_bams(vstat, in_n, in, ref);
stat_out("version\tvarcall-%s.%d\n", VERSION, SVNREV);
stat_out("min depth\t%d\n", minsampdepth);
stat_out("alpha\t%f\n", alpha);
if (vstat.stats.size()) {
// sort by depth descending
sort(vstat.stats.begin(), vstat.stats.end(), noisebydepth);
// flip 3 and 1 because sorted in descending order for sampling (above)
double depth_q3=quantile_depth(vstat.stats, .25);
double depth_q2=quantile_depth(vstat.stats, .50);
double depth_q1=quantile_depth(vstat.stats, .75);
double depth_qx=quantile_depth(vstat.stats, .95);
// number of locii to compute error rate
int ncnt=min(100000,vstat.stats.size());
int i;
double nsum=0, nssq=0, dsum=0, dmin=vstat.stats[0].depth, qnsum=0, qnssq=0, qualsum=0;
double ins_nsum=0, ins_nssq=0, del_nsum=0, del_nssq=0;
for (i=0;i<ncnt;++i) {
if (vstat.stats[i].depth < depth_q1) {
continue;
}
nsum+=vstat.stats[i].noise;
nssq+=vstat.stats[i].noise*vstat.stats[i].noise;
dsum+=vstat.stats[i].depth;
qnsum+=vstat.stats[i].qnoise;
qnssq+=vstat.stats[i].qnoise*vstat.stats[i].qnoise;
qualsum+=vstat.stats[i].mnqual;
if (vstat.stats[i].depth < dmin) dmin = vstat.stats[i].depth;
ins_nsum+=vstat.ins_stats[i].noise;
ins_nssq+=vstat.ins_stats[i].noise*vstat.ins_stats[i].noise;
del_nsum+=vstat.del_stats[i].noise;
del_nssq+=vstat.del_stats[i].noise*vstat.del_stats[i].noise;
}
double noise_mean =nsum/ncnt;
double noise_dev = stdev(ncnt, nsum, nssq);
double qnoise_mean =qnsum/ncnt;
double qnoise_dev = stdev(ncnt, qnsum, qnssq);
double qual_mean = qualsum/ncnt;
double ins_noise_mean =ins_nsum/ncnt;
double ins_noise_dev = stdev(ncnt, ins_nsum, ins_nssq);
double del_noise_mean =del_nsum/ncnt;
double del_noise_dev = stdev(ncnt, del_nsum, del_nssq);
stat_out("qual mean\t%.4f\n", qual_mean);
stat_out("noise mean\t%.6f\n", noise_mean);
stat_out("noise dev\t%.6f\n", noise_dev);
stat_out("qnoise mean\t%.6f\n", qnoise_mean);
stat_out("qnoise dev\t%.6f\n", qnoise_dev);
stat_out("ins freq\t%.6f\n", ins_noise_mean);
stat_out("ins freq dev\t%.6f\n", ins_noise_dev);
stat_out("del freq\t%.6f\n", del_noise_mean);
stat_out("del freq dev\t%.6f\n", del_noise_dev);
if (qnoise_mean >= noise_mean ) {
stat_out("error\tpoor quality estimates\n");
}
stat_out("noise depth mean\t%.4f\n", dsum/ncnt);
stat_out("noise depth min\t%.4f\n", dmin);
stat_out("noise cnt\t%d\n", ncnt);
stat_out("depth q1\t%.4f\n", depth_q1);
stat_out("depth median\t%.4f\n", depth_q2);
stat_out("depth q3\t%.4f\n", depth_q3);
dsum=0;
for (i=0;i<vstat.stats.size();++i) {
dsum+=vstat.stats[i].depth;
}
int locii_gtmin=0;
for (i=0;i<vstat.stats.size();++i) {
if (vstat.stats[i].depth >= min_depth) {
++locii_gtmin;
}
}
stat_out("locii >= min depth\t%d\n", locii_gtmin);
stat_out("locii\t%d\n", vstat.tot_locii);
double stdevfrommean=-qnorm((alpha/locii_gtmin)/2);
stat_out("qnorm adj\t%f\n", stdevfrommean);
pct_qdepth=qnoise_mean+qnoise_dev*stdevfrommean;
stat_out("min pct qual\t%.4f\n", 100*pct_qdepth);
}
}
if (read_stats){
FILE * f = fopen(read_stats, "r");
if (!f) {
warn("File %s does not exist, quitting\n", read_stats);
exit(1);
}
line l; meminit(l);
char *val;
while(read_line(f, l)>0) {
if (val=strchr(l.s, '\t')) {
*val='\0'; ++val;
if (!strcasecmp(l.s, "min depth")) {
if (umindepth && umindepth > atoi(val)) {
fprintf(varsum_f,"warning\tsampling depth was less than variation depth\n");
}
if (!umindepth) umindepth=atoi(val);
} else if (!strcasecmp(l.s, "min pct qual")) {
if (upctqdepth<=0) upctqdepth=atof(val);
} else if (!strcasecmp(l.s, "noise mean")) {
if (global_error_rate<=0) global_error_rate=atof(val);
} else if (!strcasecmp(l.s, "locii >= min depth")) {
if (total_locii<0) total_locii=atoi(val);
} else if (!strcasecmp(l.s, "alpha")) {
if (alpha<=0) alpha=atof(val);
}
}
}
}
if (total_locii<0) total_locii=DEFAULT_LOCII;
if (total_locii==0) total_locii=1; // no adjustment
if (eav_f) {
fprintf(eav_f,"chr\tpos\tref\tdepth\tnum_states\ttop_consensus\ttop_freq\tvar_base\tvar_depth\tvar_qual\tvar_strands\tforward_strands\treverse_strands\t%cval\n",total_locii>1?'e':'p');
}
if (do_varcall) {
if (umindepth) min_depth=umindepth;
if (upctqdepth > 0) pct_qdepth=(double)upctqdepth/100;
if (uminadepth) min_adepth=uminadepth;
if (uminidepth) min_idepth=uminidepth;
if (!min_depth || (!pct_depth && !pct_qdepth)) {
fprintf(varsum_f,"warning\toutputting all variations, no minimum depths specified\n");
}
fprintf(varsum_f,"version\tvarcall-%s.%d\n", VERSION, SVNREV);
fprintf(varsum_f,"min depth\t%d\n", min_depth);
fprintf(varsum_f,"min call depth\t%d\n", min_adepth);
fprintf(varsum_f,"alpha\t%f\n", alpha);
fprintf(varsum_f,"min pct qual\t%d\n", (int)(100*pct_qdepth));
fprintf(varsum_f,"min balance\t%d\n", (int)(100*pct_balance));
fprintf(varsum_f,"artifact filter\t%f\n", artifact_filter);
fprintf(varsum_f,"min qual\t%d\n", min_qual);
fprintf(varsum_f,"min map qual\t%d\n", min_mapq);
fprintf(varsum_f,"error rate\t%f\n", global_error_rate);
fprintf(varsum_f,"locii used for adjustment\t%d\n", total_locii);
VarCallVisitor vcall;
if (repeat_filter > 0) {
fprintf(varsum_f,"homopolymer filter\t%d\n", repeat_filter);
vcall.WinMax=repeat_filter+repeat_filter+3;
} else {
vcall.WinMax=5;
}
if (vcf_f) {
// print VCF header
fprintf(vcf_f, "%s\n", "##fileformat=VCFv4.1");
}
parse_bams(vcall, in_n, in, ref);
if (vcall.InputType == 'B') {
fprintf(varsum_f,"baq correct\t%s\n", (no_baq?"no":"yes"));
}
fprintf(varsum_f,"locii\t%d\n", vcall.Locii);
fprintf(varsum_f,"hom calls\t%d\n", vcall.Homs);
fprintf(varsum_f,"het calls\t%d\n", vcall.Hets);
fprintf(varsum_f,"locii below depth\t%d\n", vcall.SkippedDepth);
if (out_prefix) {
fclose(var_f);
fclose(vcf_f);
fclose(eav_f);
fclose(noise_f);
fclose(varsum_f);
if (target_annot) {
fclose(tgt_f);
fclose(tgtsum_f);
}
rename_tmp(string_format("%s.var.tmp", out_prefix));
rename_tmp(string_format("%s.vcf.tmp", out_prefix));
rename_tmp(string_format("%s.eav.tmp", out_prefix));
rename_tmp(string_format("%s.noise.tmp", out_prefix));
rename_tmp(string_format("%s.varsum.tmp", out_prefix));
if (target_annot) {
rename_tmp(string_format("%s.tgt.tmp", out_prefix));
rename_tmp(string_format("%s.tgtsum.tmp", out_prefix));
}
}
}
}
void rename_tmp(std::string f) {
std::string notmp = f;
size_t pos = notmp.find(".tmp");
if (pos >= 0) {
notmp.replace(notmp.find(".tmp"),4,"");
rename(f.c_str(),notmp.c_str());
}
}
// normal distribution
double qnorm(double q) {
if(q == .5)
return 0;
q = 1.0 - q;
double p = (q > 0.0 && q < 0.5) ? q : (1.0 - q);
double t = sqrt(log(1.0 / pow(p, 2.0)));
double c0 = 2.515517;
double c1 = 0.802853;
double c2 = 0.010328;
double d1 = 1.432788;
double d2 = 0.189269;
double d3 = 0.001308;
double x = t - (c0 + c1 * t + c2 * pow(t, 2.0)) /
(1.0 + d1 * t + d2 * pow(t, 2.0) + d3 * pow(t, 3.0));
if(q > .5)
x *= -1.0;
return x;
}
double pnorm(double x)
{
// constants
double a1 = 0.254829592;
double a2 = -0.284496736;
double a3 = 1.421413741;
double a4 = -1.453152027;
double a5 = 1.061405429;
double p = 0.3275911;
// Save the sign of x
int sign = 1;
if (x < 0)
sign = -1;
x = fabs(x)/sqrt(2.0);
// A&S formula 7.1.26
double t = 1.0/(1.0 + p*x);
double y = 1.0 - (((((a5*t + a4)*t) + a3)*t + a2)*t + a1)*t*exp(-x*x);
return 0.5*(1.0 + sign*y);
}
void parse_bams(PileupVisitor &v, int in_n, char **in, const char *ref) {
if (!in_n) {
warn("No input files, quitting\n");
exit(1);
}
int i, bam_n=0;
for (i=0;i<in_n;++i) {
if (!strcmp(fext(in[i]), ".bam")) {
++bam_n;
}
}
if (bam_n != in_n) {
if (bam_n > 0) {
warn("Can't mix bams and other input files\n");
exit(1);
} else {
if (in_n > 1) {
warn("Can't handle multiple pileups... TODO\n");
exit(1);
} else {
warn("input\t%d pileup\n", in_n);
v.InputType='P';
}
}
} else {
warn("input\t%d bam\n", bam_n);
v.InputType='B';
}
int is_popen = 0;
FILE *fin;
if (bam_n) {
if (!ref) {
warn("Need a reference file (-f) parameter, try -h for help\n");
exit(1);
}
if (!hasdata(string(ref)+".fai")) {
int ret=system(string_format("samtools faidx '%s'", ref).c_str());
if (ret) {
warn("Need a %s.fai file, run samtools faidx\n", ref);
exit(1);
}
}
const char *nobaq = no_baq ? "-B" : "";
string mpil_cmd = string_format("samtools mpileup -Q 0 -d 100000 %s -f '%s'", nobaq, ref);
int i;
for (i=0;i<in_n;++i) {
mpil_cmd += " '";
mpil_cmd += in[i];
mpil_cmd += "' ";
}
warn("command\t%s\n", mpil_cmd.c_str());
fin = popen(mpil_cmd.c_str(), "r");
if (!fin)
exit(1);
is_popen = 1;
} else {
if (!strcmp(in[0], "-")) {
fin=stdin;
} else {
if (!strcmp(fext(in[0]), ".gz")) {
string gunz = string_format("gunzip -c '%s'", in[0]);
fin = popen(gunz.c_str(), "r");
is_popen = 1;
} else {
fin = fopen(in[0], "r");
}
if (!fin) {
warn("%s: %s", in[0], strerror(errno));
exit(1);
}
}
}
line l; meminit(l);
int cnt=0;
if (fin) {
while(read_line(fin, l)>0) {
// chr 2 G 6 ^9,^+.^*,^2,^&.^&, &.'&*- 9+*2&& 166,552,643,201,299,321
v.Parse(l.s);
++cnt;
}
v.Finish();
if (is_popen) pclose(fin); else fclose(fin);
}
if (cnt == 0) {
warn("No data in pileup, quitting\n");
exit(1);
}
}
#define T_A 0
#define T_C 1
#define T_G 2
#define T_T 3
#define T_SDEL 4
#define T_NDEL 5
#define T_INS 6
#define T_N 7
#define b2i(c) ((c)=='A'?0:(c)=='a'?0:(c)=='C'?1:(c)=='c'?1:(c)=='G'?2:(c)=='g'?2:(c)=='T'?3:(c)=='t'?3:(c)=='*'?4:(c)=='-'?5:(c)=='+'?6:7)
#define i2b(i) (i==0?'A':i==1?'C':i==2?'G':i==3?'T':i==4?'*':i==5?'-':i==6?'+':'?')
bool hitoloint (int i,int j) { return (i>j);}
int track_readlen[10000];
PileupSummary::PileupSummary(char *line, PileupReads &rds) {
vector<char *> d=split(line, "\t");
if (d.size() < 6) {
warn("Can't read pileup : %d fields, need 6 columns\n", (int) d.size());
exit(1);
}
const char * p_qual=d[5];
Chr=d[0];
Pos=atoi(d[1]);
Base=*(d[2]);
Depth = atoi(d[3]);
SkipDupReads = 0;
SkipN = 0;
SkipMinQual = 0;
SkipMinMapq = 0;
MaxDepthByPos = 0;
RepeatCount = 0;
RepeatBase = '\0';
NumReads = 0;
int i;
vector<int> depthbypos;
const char *cur_p = d[4];
list<Read>::iterator read_i = rds.ReadList.begin();
int eor=0;
for (i=0;i<Depth;++i,++read_i) {
bool sor=0;
if (*cur_p == '^') {
sor=1;
++cur_p;
Read x;
x.MapQ = *cur_p-phred;
++cur_p;
if (read_i != rds.ReadList.end()) {
++read_i;
}
read_i=rds.ReadList.insert(read_i,x);
}
if (read_i == rds.ReadList.end()) {
warn("warning\tread start without '^', partial pileup: '%s'\n", cur_p);
Read x;
x.MapQ = -1;
read_i=rds.ReadList.insert(read_i,x);
}
int pia = read_i->Seq.length()+1;
if (pia >= depthbypos.size()) {
depthbypos.resize(pia+1);
}
depthbypos[pia]++;
if (sor)
++NumReads;
char q = p_qual[i]-phred; // qual char
char mq = read_i->MapQ;
char o = *cur_p; // orig call
char c = toupper(o); // uppercase/ref
bool is_ref = 0;
if (o == '.' || o == ',') {
c = Base; // ref instead
is_ref = 1;
}
if (o == '>' || o == '<') {
c = 'N'; // no call
is_ref = 1;
}
bool skip = 0;
// probably should not be adding anything here... but the old code added 1 and floored... new code adds .5 and rounds... which is comparable
// really.. should just be adding zero, the reason the old code had it was because of a lack of max()
if (c == 'N') {
++SkipN;
skip=1;
} else if (artifact_filter > 0 && (depthbypos[pia] > max(1,rand_round(0.5+artifact_filter * (Depth/rds.MeanReadLen()))))) {
++SkipDupReads;
skip=1;
} else if (mq < min_mapq) {
++SkipMinMapq;
skip=1;
} else if (q < min_qual) {
++SkipMinQual;
skip=1;
} else {
int j = b2i(c);
if (j >= Calls.size()) {
int was = Calls.size();
Calls.resize(j+1);
int t; for (t=was;t<=j;++t) {
Calls[t].base=i2b(t);
}
}
if (is_ref)
Calls[j].is_ref = 1;
if ( o == ',' || o == 'a' || o == 'c' || o == 't' || o == 'g' ) {
++Calls[j].rev;
} else if ( c != 'N' ) {
++Calls[j].fwd;
}
Calls[j].qual+=q;
Calls[j].mn_qual+=min(mq,q);
Calls[j].mq_ssq+=mq*mq;
Calls[j].mq_sum+=mq;
Calls[j].qual_ssq+=q*q;
/*
if (pia <= read_tail_len || (rds.MeanReadLen()-pia) <= read_tail_len) {
if ( o == ',' || o == 'a' || o == 'c' || o == 't' || o == 'g' ) {
++Calls[j].tail_rev;
} else {
++Calls[j].tail_fwd;
}
}
*/
if (vcf_f) {
if (mq == 0)
Calls[j].mq0++;
}
}
if (c == '-' || c == '+') {
warn("invalid pileup, at '%s', indel not attached to read?\n", cur_p);
} else {
if (c != '*')
read_i->Seq += c;
++cur_p;
}
if (*cur_p == '+' || *cur_p == '-') {
c = *cur_p;
char *end_p;
int len = strtol(++cur_p, &end_p, 10);
string ins_seq(end_p, len);
to_upper(ins_seq);
read_i->Seq += ins_seq;
if (!skip) {
int j = b2i(c);
if (j >= Calls.size()) {
int was = Calls.size();
Calls.resize(j+1);
int t; for (t=was;t<=j;++t) {
Calls[t].base=i2b(t);
}
}
if ( o == ',' || o == 'a' || o == 'c' || o == 't' || o == 'g' ) {
++Calls[j].rev;
} else {
++Calls[j].fwd;
}
Calls[j].qual+=q;
Calls[j].mn_qual+=min(q, mq);
Calls[j].qual_ssq+=q*q;
Calls[j].mq_ssq+=mq*mq;
Calls[j].mq_sum+=mq;
Calls[j].seqs.push_back(ins_seq);
}
cur_p=end_p+len;
}
if (*cur_p == '$') {
if (read_i->MapQ > -1) {
rds.TotReadLen+=read_i->Seq.size();
rds.ReadBin.push_back(*read_i);
if (rds.ReadBin.size() > min(1000,Depth*2)) {
rds.ReadBin.pop_front();
rds.TotReadLen-=rds.ReadBin.front().Seq.size();
}
}
// printf("%d\t%s\n", read_i->MapQ, read_i->Seq.c_str());
read_i=rds.ReadList.erase(read_i);
--read_i;
++cur_p;
++eor;
}
}
if ((Depth-eor) != rds.ReadList.size()) {
warn("warning\tdepth is %d, but read list is: %d\n", Depth, (int) rds.ReadList.size());
}
if (*cur_p == '-' || *cur_p == '+') {
char *end_p;
int len = strtol(++cur_p, &end_p, 10);
// keep this
string idl(end_p, len);
cur_p=end_p+len;
}
if (*cur_p) {
warn("Failed to parse pileup %s\n", d[4]);
exit(1);
}
for (i=0;i<depthbypos.size();++i) {
if (depthbypos[i] > MaxDepthByPos) {
MaxDepthByPos = depthbypos[i];
}
}
Depth=0;
for (i=0;i<5 && i < Calls.size();++i) { // total depth (exclude inserts for tot depth, otherwise they are double-counted)
Depth+=Calls[i].depth();
}
TotQual=0;
for (i=0;i<5 && i < Calls.size();++i) { // total depth (exclude inserts for tot depth, otherwise they are double-counted)
TotQual+=Calls[i].qual;
}
}
PileupSummary JunkSummary;
void VarCallVisitor::Visit(PileupSummary &p) {
if (WinMax < 3) {
// no real window ... just go straight
VisitX(p);
return;
}
if (p.Base != '-' && p.Base != '@') {
if (Win.size() && (Win.back().Pos != (p.Pos - 1) )) {
if (Win.back().Pos < p.Pos && ((p.Pos - Win.back().Pos) <= (WinMax/2))) {
while (Win.back().Pos < (p.Pos - 1)) {
// visit/pop, add a placeholder
JunkSummary.Base = '-';
JunkSummary.Pos = Win.back().Pos + 1;
Visit(JunkSummary);
}
} else {
while (Win.size() && Win[WinMax/2].Base != '@') {
// visit/pop, but don't add anything, until it's empty
JunkSummary.Base = '@';
JunkSummary.Pos = 0;
Visit(JunkSummary);
}
}
}
}
// initialize the window with nothing, if it's not full
while (Win.size() < WinMax) {
JunkSummary.Base = '@';
JunkSummary.Pos = 0;
Win.push_back(JunkSummary);
}
Win.push_back(p);
//debug("Visit: %d\n", p.Pos);
if (Win.size() > WinMax) // queue too big? pop
Win.pop_front();
int i;
int lrc=0,rrc=0; // left repeat count, right repeat count
char lrb, rrb; // left repeat base...
int vx;
if (Win.size() < WinMax) { // small window? look at leading edge only
return;
} else {
vx = WinMax/2; // larger window? look at midpoint
}
if (Win[vx].Base == '-' || Win[vx].Base == '@')
return;
if (vx > 1) { // look left
lrb = Win[vx-1].Base;
for (i=vx-2; i >= 0; --i) { // increment repeat count
if (Win[i].Base == lrb)
++lrc;
else
break;
}
}
if (vx < (Win.size()-2)) {
rrb = Win[vx+1].Base;
for (i=vx+2; i < Win.size(); ++i) {
if (Win[i].Base == rrb)
++rrc;
else
break;
}
}
// repeat counts are now 1-based, not 0-based
++lrc;
++rrc;
// maximum repeat count and associated base
if (lrb == rrb ) {
Win[vx].RepeatCount = lrc+rrc;
Win[vx].RepeatBase = lrb;
} else if (lrc > rrc) {
Win[vx].RepeatCount = lrc;
Win[vx].RepeatBase = lrb;
} else {
Win[vx].RepeatCount = rrc;
Win[vx].RepeatBase = rrb;
}
if (debug_xpos) {
if (Win[vx].Pos == debug_xpos && !strcmp(debug_xchr,Win[vx].Chr.data())) {
fprintf(stderr,"xpos-window\t");
for (i=0;i<Win.size();++i) {
fprintf(stderr,"%c", Win[i].Base);
}
fprintf(stderr,"\n");
}
}
double drms = 0;
if (vx < Win.size()-1) {
int i;
int dminus = b2i('-');
int dstar = b2i('*');
if (Win[vx].Calls.size() > dminus && Win[vx].Calls[dminus].depth() > 0) {
if (Win[vx+1].Calls.size() > dstar && Win[vx+1].Calls[dstar].depth() > 0) {
// baq adjustment works at the 'star' not at the 'indel', so adjust qual using the next locus
double adj=Win[vx+1].Calls[dstar].qual_rms()/(double)Win[vx].Calls[dminus].qual_rms();
if (debug_xpos) {
if (Win[vx].Pos == debug_xpos && !strcmp(debug_xchr,Win[vx].Chr.data())) {
fprintf(stderr,"xpos-adj-qual\t%d to %d (%f)\n", Win[vx].Calls[dminus].qual_rms(),Win[vx+1].Calls[dstar].qual_rms(), adj);
}
}
Win[vx].Calls[dminus].qual *= adj;
Win[vx].Calls[dminus].qual_ssq *= adj;
} else {
vcall none;
if (debug_xpos) {
if (Win[vx].Pos == debug_xpos && !strcmp(debug_xchr,Win[vx].Chr.data())) {
fprintf(stderr,"xpos-skip-del-qual\t%d\n", Win[vx].Calls[dminus].depth());
}
}
Win[vx].Calls[dminus] = none;
}
}
}
VisitX(Win[vx]);
}
void VarCallVisitor::Finish() {
// finish out the rest of the pileup, with the existing window
int vx = WinMax/2+1;
while (vx < Win.size()) {
///debug("Finish: %d\n", Win[vx].Pos);
VisitX(Win[vx++]);
}
}
void VarCallVisitor::VisitX(PileupSummary &p) {
//debug("VisitX: %d\n", p.Pos);
if (debug_xpos) {
if (p.Pos != debug_xpos)
return;
if (strcmp(debug_xchr,p.Chr.data()))
return;
}
if (p.Depth < min_depth) {
if (debug_xpos) {
fprintf(stderr,"xpos-skip-depth\t%d < %d\n",p.Depth, min_depth);
fprintf(stderr,"xpos-skip-dup\t%d\n",p.SkipDupReads);
fprintf(stderr,"xpos-skip-n\t%d\n",p.SkipN);
fprintf(stderr,"xpos-skip-mapq\t%d\n",p.SkipMinMapq);
fprintf(stderr,"xpos-skip-qual\t%d\n",p.SkipMinQual);
}
++SkippedDepth;
return;
}
int ins_fwd = p.Calls.size() > 6 ? p.Calls[6].fwd : 0;
int ins_rev = p.Calls.size() > 6 ? p.Calls[6].rev : 0;
int i;
if (p.Calls.size() > 6)
p.Calls.resize(7); // toss N's before sort
sort(p.Calls.begin(), p.Calls.end(), hitolocall);
int need_out = -1;
int skipped_balance=0;
int skipped_alpha=0;
int skipped_indel=0;
int skipped_tail_hom=0;
int skipped_depth=0;
int skipped_repeat=0;
vector<vfinal> final_calls;
for (i=0;i<p.Calls.size();++i) { // all calls
// printf("CALL TOP: depth:%d base: %c, pd: %d, calls: %d\n", (int) p.Calls[i].depth(), p.Calls[i].base, p.Depth, (int) p.Calls.size());
double pct = (double) p.Calls[i].depth()/p.Depth;
double qpct = (double) p.Calls[i].qual/p.TotQual;
if (!p.Calls[i].base)
continue;
if (!p.Calls[i].depth())
continue;
double bpct = (double) min(p.Calls[i].fwd,p.Calls[i].rev)/p.Calls[i].depth();
if (pct > pct_depth && qpct >= pct_qdepth && (p.Calls[i].depth() >= min_adepth)) {
if (bpct < pct_balance) {
int fwd_adj=0, rev_adj=0;
// f=b*(f+r); r=f/b-f; adj=r-(f/b-f)
if (p.Calls[i].fwd < p.Calls[i].rev) {
rev_adj = (int) p.Calls[i].rev - ( p.Calls[i].fwd/pct_balance - p.Calls[i].fwd );
} else {
fwd_adj = (int) p.Calls[i].fwd - ( p.Calls[i].rev/pct_balance - p.Calls[i].rev );
}
if (fwd_adj + rev_adj > 1 && bpct > 0) {
// adjust call down
p.Calls[i].qual -= (rev_adj+fwd_adj)*(p.Calls[i].qual/p.Calls[i].depth());
p.Calls[i].mq_sum -= (rev_adj+fwd_adj)*(p.Calls[i].mq_sum/p.Calls[i].depth());
p.Calls[i].qual_ssq -= (rev_adj+fwd_adj)*(p.Calls[i].qual_ssq/p.Calls[i].depth());
p.Calls[i].mq_ssq -= (rev_adj+fwd_adj)*(p.Calls[i].mq_ssq/p.Calls[i].depth());
p.Calls[i].rev -= rev_adj;
p.Calls[i].fwd -= fwd_adj;
skipped_balance+=rev_adj+fwd_adj;
// fixed bpct
bpct = (double) min(p.Calls[i].fwd,p.Calls[i].rev)/p.Calls[i].depth();
} else {
// it's junk anyway
}
// fix depths after adjustment!
pct = (double) p.Calls[i].depth()/p.Depth;
qpct = (double) p.Calls[i].qual/p.TotQual;
}
}
if (pct > pct_depth && qpct >= pct_qdepth && (p.Calls[i].depth() >= min_adepth)) {
// balance is meaningless at low depths
if ((bpct >= pct_balance) || (p.Calls[i].depth()<4)) {
if (p.Calls[i].base == '+' || p.Calls[i].base == '-') {
// yuk ... time to think about a possible indel call
if (p.Calls[i].depth() >= min_idepth) {
// should really pick more than 1
// but need to allow "similar" indels to pile up
// should group into distinct bins, using some homology thing
sort(p.Calls[i].seqs.begin(), p.Calls[i].seqs.end());
string prev, maxs;
int pcnt=0, maxc=0, j;
for (j=0;j<p.Calls[i].seqs.size();++j) {
if (prev == p.Calls[i].seqs[j]) {
++pcnt;
} else {
if (pcnt > maxc) {
maxs=prev;
maxc=pcnt;
}
prev=p.Calls[i].seqs[j];
pcnt=1;
}
}
if (pcnt > maxc) {
maxs=prev;
maxc=pcnt;
}
if (maxc >= min_idepth && maxc >= min_adepth) {
// only calls 1 indel at a given position
if ((repeat_filter == 0) || (p.RepeatCount < repeat_filter)) {
// maybe use rms here... see if it helps
double mean_qual = p.Calls[i].qual/(double)p.Calls[i].depth();
double err_rate = mean_qual < max_phred ? pow(10,-mean_qual/10.0) : global_error_rate;
// expected number of non-reference = error_rate*depth
double pval=(p.Depth*err_rate==0)?0:gsl_ran_poisson_pdf(p.Calls[i].depth(), p.Depth*err_rate);
double padj=total_locii ? pval*total_locii : pval; // multiple-testing adjustment
if (padj <= alpha) {
vfinal final(p.Calls[i]);
double mq_padj=max(total_locii*pow(10,-p.Calls[i].mq_sum/10.0),padj); // never report pval as better than the total mapping quality
if (debug_xpos) fprintf(stderr,"xpos-debug-pval\tbase:%c, err:%g, pval:%g, padj:%g, mq_padj:%g, mq_sum:%d\n", p.Calls[i].base, err_rate, pval, padj, mq_padj, p.Calls[i].mq_sum);
if (mq_padj > 1) mq_padj=1;
if (need_out == -1)
need_out = i;
// printf("FINAL: depth:%d base: %s\n", (int) maxc, maxs.c_str());
final.padj=mq_padj;
final.max_idl_cnt=maxc;
final.max_idl_seq=maxs;
final_calls.push_back(final);
} else {
skipped_alpha+=p.Calls[i].depth();
}
// implicitly skip all the ohter indel calls at the same locus
skipped_indel+=p.Calls[i].depth()-maxc;
} else {
skipped_repeat+=p.Calls[i].depth();
}
} else {
skipped_indel+=p.Calls[i].depth();
}
} else {
skipped_indel+=p.Calls[i].depth();
}
} else {
if (p.Calls[i].base == '*' && (
((repeat_filter > 0) && (p.RepeatCount >= repeat_filter)) ||
(p.Calls[i].depth() < min_idepth)
)) {
skipped_indel+=p.Calls[i].depth();
} else {
// subtract inserts from reference .. perhaps > 0 is correct here....
if (p.Calls[i].is_ref && (ins_rev+ins_fwd) > max(min_idepth,min_adepth)) {
p.Calls[i].fwd-=ins_fwd;
p.Calls[i].rev-=ins_rev;
}
double mean_qual = p.Calls[i].qual/(double)p.Calls[i].depth();
/*
if ( (repeat_filter > 0) && (p.RepeatCount >= repeat_filter) ) {
p.Calls[i].fwd-=p.Calls[i].tail_fwd;
p.Calls[i].rev-=p.Calls[i].tail_rev;
skipped_tail_hom+=p.Calls[i].tail_fwd+p.Calls[i].tail_rev;
}
*/
if (p.Calls[i].depth() >= min_adepth && p.Calls[i].depth() > 0) {
double err_rate = mean_qual < max_phred ? pow(10,-mean_qual/10.0) : global_error_rate;
// expected number of non-reference bases at this position is error_rate*depth
double pval=(p.Depth*err_rate==0)?0:gsl_ran_poisson_pdf(p.Calls[i].depth(), p.Depth*err_rate);
double padj=total_locii ? pval*total_locii : pval; // multiple-testing adjustment
if (padj <= alpha) {
double mq_padj=max(total_locii*pow(10,-p.Calls[i].mq_sum/10.0),padj); // never report as better than the mapping quality
if (mq_padj > 1) mq_padj=1;
if (debug_xpos) fprintf(stderr,"xpos-debug-pval\tbase:%c, err:%g, pval:%g, padj:%g, mq_padj:%g, mq_sum:%d\n", p.Calls[i].base, err_rate, pval, padj, mq_padj, p.Calls[i].mq_sum);
if (!p.Calls[i].is_ref || debug_xpos) {
if (need_out == -1)
need_out = i;
}
vfinal final(p.Calls[i]);
final.padj=mq_padj;
final_calls.push_back(final);
} else {
skipped_alpha+=p.Calls[i].depth();
}
}
}
}
} else {
skipped_balance+=p.Calls[i].depth();
}
} else {
// depth is too low now.... technically you can just add all the rest of the calls to skipped_depth without checking
skipped_depth+=p.Calls[i].depth();
}
}
++Locii;
if (need_out>=0||debug_xpos) {
if (final_calls.size() > 1){
// printf("HERE1 %c/%c\n", final_calls[0].pcall->base, final_calls[1].pcall->base);
if(final_calls[1].pcall->is_ref) {
vfinal tmp=final_calls[1];
final_calls[1]=final_calls[0];
final_calls[0]=tmp;
// printf("HERE2 %c/%c\n", final_calls[0].pcall->base, final_calls[1].pcall->base);
}
}
// printf("allele_count: %d\n", (int) final_calls.size());
int total_call_depth=0;
int i;
for (i=0;i<final_calls.size();++i) {
total_call_depth+=final_calls[i].pcall->depth();
}
double pct_allele = 0;
if (need_out >=0) {
// more than 1 call at this position = Het
if (final_calls.size() > 1) {
if (final_calls[0].pcall->is_ref) {
pct_allele = 100.0 * final_calls[1].pcall->depth() / (double) total_call_depth;
} else {
// no reference seen... but still het?
pct_allele = 100.0 * final_calls[0].pcall->depth() / (double) total_call_depth;
}
++Hets;
} else {
pct_allele = 100.0 * final_calls[0].pcall->depth() / (double) total_call_depth;
++Homs;
}
}
if (var_f) {
int i;
string pil;
for (i=0;i<final_calls.size();++i) {
vfinal &f=final_calls[i];
if (f.is_indel()) {
pil += string_format("\t%c%s:%d,%d,%.1e",f.pcall->base,f.max_idl_seq.c_str(),f.max_idl_cnt,f.pcall->qual/f.pcall->depth(),f.padj);
} else {
pil += string_format("\t%c:%d,%d,%.1e",f.pcall->base,f.pcall->depth(),f.pcall->qual/f.pcall->depth(),f.padj);
}
}
fprintf(var_f,"%s\t%d\t%c\t%d\t%d\t%2.2f%s\n",p.Chr.c_str(), p.Pos, p.Base, p.Depth, skipped_alpha+skipped_depth+skipped_balance+p.SkipN+p.SkipDupReads+p.SkipMinMapq+p.SkipMinQual, pct_allele, pil.c_str());
}
if (vcf_f) {
for (i=0;i<final_calls.size();++i) {
vfinal &f=final_calls[i];
int qual = f.padj>0?min(40,10*(-log10(f.padj))):40;
if (f.is_indel()) {
string base;
string alt;
if (f.pcall->base =='-') {
base = p.Base + f.max_idl_seq;
alt = p.Base;
} else {
base = p.Base;
alt = p.Base + f.max_idl_seq;
}
double freq_allele = f.max_idl_cnt / (double) p.Depth;
fprintf(vcf_f,"%s\t%d\t.\t%s\t%s\t%2d\tPASS\tMQ=%d;BQ=%d;DP=%d;AF=%2.2f\n",
p.Chr.c_str(), p.Pos, base.c_str(), alt.c_str(), qual,
(int) f.pcall->mq_rms(),
(int) f.pcall->qual_rms(),
total_call_depth,
freq_allele);
} else {
char alt = f.pcall->base;
if (f.pcall->is_ref)
alt = '.';
double freq_allele = f.pcall->depth() / (double) p.Depth;
fprintf(vcf_f,"%s\t%d\t.\t%c\t%c\t%d\tPASS\tMQ=%d;BQ=%d;DP=%d;AF=%2.2f\n",
p.Chr.c_str(), p.Pos, p.Base, alt, qual,
(int) f.pcall->mq_rms(),
(int) f.pcall->qual_rms(),
total_call_depth,
freq_allele);
}
}
}
if (eav_f) {
// printf(eav_f,"chr\tpos\tref\tdepth\tnum_states\ttop_consensus\ttop_freq\tvar_base\tvar_depth\tvar_qual\tvar_strands\tforward_strands\treverse_strands\n");
string top_cons, var_base, var_depth, var_qual, var_strands, forward, reverse;
float padj=final_calls[0].padj;
if (final_calls[0].pcall->is_ref && final_calls.size() > 1) {
padj=final_calls[1].padj;
}
for (i=0;i<final_calls.size();++i) {
vfinal &f=final_calls[i];
if (i < 2) {
if (i > 0) top_cons += "/";
top_cons += f.pcall->base;
}
if (i > 0) var_base += "/";
var_base += f.pcall->base;
if (f.is_indel()) {
if (i < 2) {
top_cons += f.max_idl_seq;
}
var_base += f.max_idl_seq;
}
if (i > 0) var_depth+= ";";
var_depth+= string_format("%d",f.pcall->depth());
if (i > 0) var_qual+= ";";
var_qual+= string_format("%d",f.pcall->qual_rms());
if (i > 0) var_strands+= ";";
var_strands+= string_format("%d",(f.pcall->fwd>0)+(f.pcall->rev>0));
if (i > 0) forward += ";";
forward+= string_format("%d",f.pcall->fwd);
if (i > 0) reverse += ";";
reverse+= string_format("%d",f.pcall->rev);
}
fprintf(eav_f,"%s\t%d\t%c\t%d\t%d\t%s\t%2.2f\t%s\t%s\t%s\t%s\t%s\t%s\t%.1e\n",p.Chr.c_str(), p.Pos, p.Base, p.Depth, (int) final_calls.size(),top_cons.c_str(), pct_allele, var_base.c_str(), var_depth.c_str(), var_qual.c_str(), var_strands.c_str(), forward.c_str(), reverse.c_str(), padj);
}
if (debug_xpos) {
fprintf(stderr,"xpos-skip-dup\t%d\n",p.SkipDupReads);
fprintf(stderr,"xpos-skip-mapq\t%d\n",p.SkipMinMapq);
fprintf(stderr,"xpos-skip-qual\t%d\n",p.SkipMinQual);
fprintf(stderr,"xpos-skip-bal\t%d\n",skipped_balance);
fprintf(stderr,"xpos-skip-depth\t%d\n",skipped_depth);
fprintf(stderr,"xpos-skip-indel\t%d\n",skipped_indel);
// fprintf(stderr,"xpos-skip-tail-imbalance\t%d\n",skipped_tail_hom);
fprintf(stderr,"xpos-skip-repeat\t%d\n",skipped_repeat);
fprintf(stderr,"xpos-skip-alpha\t%d\n",skipped_alpha);
if (repeat_filter > 0) {
fprintf(stderr,"repeat-count\t%d\n",p.RepeatCount);
fprintf(stderr,"repeat-filter\t%d\n",repeat_filter);
fprintf(stderr,"repeat-base\t%c\n",p.RepeatBase);
}
exit(0);
}
}
}
void PileupVisitor::LoadIndex(const char *path) {
FILE *f = fopen(path,"r");
if (!f) {
warn("Can't open %s : %s\n", path, strerror(errno));
exit(1);
}
AnnotType = '\0';
line l; meminit(l);
int cnt=0;
while(read_line(f, l)>0) {
vector<char *> d=split(l.s, "\t");
if (d.size() < 9) {
warn("File must be a GTF or a BED: '%s'\n", path);
exit(1);
}
AnnotType = (*d[5]=='+' || *d[5] == '-') ? 'b' : '\0';
AnnotType = (*d[6]=='+' || *d[5] == '-') ? 'g' : AnnotType;
break;
}
if (!AnnotType) {
warn("File must be a GTF or a BED: '%s'\n", path);
exit(1);
}
if (!AnnotDex.read(path)) {
// void build(const char *path, const char *sep, int nchr, int nbeg, int nend, int skip_i, char skip_c);
AnnotDex.build(path, "\t", 0, 1, 2, 0, '#', 1);
}
fclose(f);
}
void VarStatVisitor::Visit(PileupSummary &p) {
tot_locii += 1;
if (p.Depth < minsampdepth)
return;
// insert and deletions have their own, separate noise levels
int ins_depth = p.Calls.size() > 6 ? p.Calls[6].depth() : 0;
int ins_qual = p.Calls.size() > 6 ? p.Calls[6].qual : 0;
double ins_noise = 0;
double ins_qnoise = 0;
if (p.Calls.size() > 1 && p.Calls[1].depth() > ins_depth && ins_depth > 0) {
ins_noise = (double) ins_depth/p.Depth;
ins_qnoise = (double) ins_qual/p.TotQual;
}
int del_depth = p.Calls.size() > 5 ? p.Calls[5].depth() : 0;
int del_qual = p.Calls.size() > 5 ? p.Calls[5].qual : 0;
double del_noise = 0;
double del_qnoise = 0;
if (p.Calls.size() > 1 && p.Calls[1].depth() > del_depth && del_depth > 0) {
del_noise = (double) del_depth/p.Depth;
del_qnoise = (double) del_qual/p.TotQual;
}
// snp's are "noise" if there are 3 alleles at a given position
int i;
if (p.Calls.size() > 5)
p.Calls.resize(5); // toss N's and inserts before sort
sort(p.Calls.begin(), p.Calls.end(), hitolocall);
double noise = p.Calls.size() > 2 ? (double) p.Calls[2].depth()/p.Depth : 0;
double qnoise = p.Calls.size() > 2 ? (double) p.Calls[2].qual/p.TotQual : 0;
double mnqual = (double)p.TotQual/p.Depth;
char pbase = p.Calls.size() > 2 ? p.Calls[2].base : '.';
if (noise_f) {
fprintf(noise_f,"%d\t%c\t%f\t%f\n", p.Depth, pbase, noise, qnoise, mnqual);
/*
if (ins_noise > 0) {
fprintf(noise_f,"%d\t%c\t%f\t%f\n", p.Depth, '+', ins_noise, ins_qnoise, mnqual);
}
if (del_noise > 0) {
fprintf(noise_f,"%d\t%c\t%f\t%f\n", p.Depth, '-', del_noise, del_qnoise, mnqual);
}
*/
}
tot_depth += p.Depth;
num_reads += p.NumReads;
stats.push_back(Noise(p.Depth, noise, qnoise, mnqual));
ins_stats.push_back(Noise(p.Depth, ins_noise, ins_qnoise, mnqual));
del_stats.push_back(Noise(p.Depth, del_noise, del_qnoise, mnqual));
}
void usage(FILE *f) {
fprintf(f,
"Usage: varcall <-s|-v> <-f REF> [options] bam1 [bam2...]\n"
"Version: %s.%d (BETA)\n"
"\n"
"Either outputs summry stats for the list of files, or performs variant calling\n"
"\n"
"Options (later options override earlier):\n"
"\n"
"-s Calculate statistics\n"
"-v Calculate variants bases on supplied parameters (see -S)\n"
"-f Reference fasta (required if using bams, ignored otherwise)\n"
"-m Min locii depth (0)\n"
"-a Min allele depth (0)\n"
"-p Min allele pct by quality (0)\n"
"-q Min qual (3)\n"
"-Q Min mapping quality (0)\n"
"-b Min pct balance (strand/total) (0)\n"
"-D FLOAT Max duplicate read fraction (depth/length per position) (1)\n"
"-B Turn off BAQ correction (false)\n"
"-R Homopolymer repeat indel filtering (8)\n"
"-e FLOAT Alpha filter to use, requires -l or -S (.05)\n"
"-g FLOAT Global minimum error rate (default: assume phred is ok)\n"
"-l INT Number of locii in total pileup used for bonferroni (1 mil)\n"
"-x CHR:POS Output this pos only, then quit\n"
"-N FIL Output noise stats to FIL\n"
"-S FIL Read in statistics and params from a previous run with -s (do this!)\n"
"-A ANNOT Calculate in-target stats using the annotation file (requires -o)\n"
"-o PREFIX Output prefix (note: overlaps with -N)\n"
"\n"
"Input files\n"
"\n"
"Files must be sorted bam files with bai index files available. Alternatively,\n"
"a single pileup file can be supplied.\n"
"\n"
"Output files\n"
"\n"
"Varcalls go to stdout. Stats go to stdout, or stderr if varcalling too\n"
"\n"
"If an output prefix is used, files are created as follows:\n"
" PREFIX.var Variant calls in tab delimited 'varcall' format\n"
" PREFIX.eav Variant calls in tab delimited 'ea-var' format\n"
" PREFIX.vcf Variant calls, in vcf format\n"
" PREFIX.varsum Summary of variant calls\n"
" PREFIX.tgt On-target stats detail\n"
" PREFIX.tgtsum Summary of on-target stats\n"
" PREFIX.noise Noise stats detail\n"
"\n"
"Stats Output:\n"
"\n"
"Contains mean, median, quartile information for depth, base quality, read len,\n"
"mapping quality, indel levels. Also estimates parameters suitable for\n"
"variant calls, and can be passed directly to this program for variant calls\n"
"\n"
"Filtering Details:\n"
"\n"
,VERSION, SVNREV);
}
std::string string_format(const std::string &fmt, ...) {
int n, size=100;
std::string str;
va_list ap;
while (1) {
str.resize(size);
va_start(ap, fmt);
int n = vsnprintf((char *)str.c_str(), size, fmt.c_str(), ap);
va_end(ap);
if (n > -1 && n < size) {
str.resize(n);
return str;
}
if (n > -1)
size=n+1;
else
size*=2;
}
}
void to_upper(const std::string str) {
std::string::iterator it;
int i;
for ( i=0;i<str.size();++i ) {
((char *)(void *)str.data())[i]=toupper(((char *)str.data())[i]);
}
}
// returns quantile depth
double quantile_depth(const std::vector<Noise> &vec, double p) {
int l = vec.size();
assert(l > 0);
double t = ((double)l-1)*p;
int it = (int) t;
int v=vec[it].depth;
if (t > (double)it) {
return (v + (t-it) * (vec[it+1].depth - v));
} else {
return v;
}
}
double quantile(const std::vector<int> &vec, double p) {
int l = vec.size();
double t = ((double)l-1)*p;
int it = (int) t;
int v=vec[it];
if (t > (double)it) {
return (v + (t-it) * (vec[it+1] - v));
} else {
return v;
}
}
double quantile(const std::vector<double> &vec, double p) {
int l = vec.size();
double t = ((double)l-1)*p;
int it = (int) t;
double v=vec[it];
if (t > (double)it) {
return (v + p * (vec[it+1] - v));
} else {
return v;
}
}
std::vector<char *> split(char* str,const char* delim)
{
char* token = strtok(str,delim);
std::vector<char *> result;
while(token != NULL)
{
result.push_back(token);
token = strtok(NULL,delim);
}
return result;
}
int rand_round(double x) {
return floor(x)+((rand()>(x-int(x))) ? 1 : 0);
//warn("rr:%f=%d\n",x);
}
|