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 1745 1746 1747 1748 1749 1750 1751 1752 1753 1754 1755 1756 1757 1758 1759 1760 1761 1762 1763 1764 1765 1766 1767 1768 1769 1770 1771 1772 1773 1774 1775 1776 1777 1778 1779 1780 1781 1782 1783 1784 1785 1786 1787 1788 1789 1790 1791 1792 1793 1794 1795 1796 1797 1798 1799 1800 1801 1802 1803 1804 1805 1806 1807 1808 1809 1810 1811 1812 1813 1814 1815 1816 1817 1818 1819 1820 1821 1822 1823 1824 1825 1826 1827 1828 1829 1830 1831 1832 1833 1834 1835 1836 1837 1838 1839 1840 1841 1842 1843 1844 1845 1846 1847 1848 1849 1850 1851 1852 1853 1854 1855 1856 1857 1858 1859 1860 1861 1862 1863 1864 1865 1866 1867 1868 1869 1870 1871 1872 1873 1874 1875 1876 1877 1878 1879 1880 1881 1882 1883 1884 1885 1886 1887 1888 1889 1890 1891 1892 1893 1894 1895 1896 1897 1898 1899 1900 1901 1902 1903 1904 1905 1906 1907 1908 1909 1910 1911 1912 1913 1914 1915 1916 1917 1918 1919 1920 1921 1922 1923 1924 1925 1926 1927 1928 1929 1930 1931 1932 1933 1934 1935 1936 1937 1938 1939 1940 1941 1942 1943 1944 1945 1946 1947 1948 1949 1950 1951 1952 1953 1954 1955 1956 1957 1958 1959 1960 1961 1962 1963 1964 1965 1966 1967 1968 1969 1970 1971 1972 1973 1974 1975 1976 1977 1978 1979 1980 1981 1982 1983 1984 1985 1986 1987 1988 1989 1990 1991 1992 1993 1994 1995 1996 1997 1998 1999 2000 2001 2002 2003 2004 2005 2006 2007 2008 2009 2010 2011 2012 2013 2014 2015 2016 2017 2018 2019 2020 2021 2022 2023 2024 2025 2026 2027 2028 2029 2030 2031 2032 2033 2034 2035 2036 2037 2038 2039 2040 2041 2042 2043 2044 2045 2046 2047 2048 2049 2050 2051 2052 2053 2054 2055 2056 2057 2058 2059 2060 2061 2062 2063 2064 2065 2066 2067 2068 2069 2070 2071 2072 2073 2074 2075 2076 2077 2078 2079 2080 2081 2082 2083 2084 2085 2086 2087 2088 2089 2090 2091 2092 2093 2094 2095 2096 2097 2098 2099 2100 2101 2102 2103 2104 2105 2106 2107 2108 2109 2110 2111 2112 2113 2114 2115 2116 2117 2118 2119 2120 2121 2122 2123 2124 2125 2126 2127 2128 2129 2130 2131 2132 2133 2134 2135 2136 2137 2138 2139 2140 2141 2142 2143 2144 2145 2146 2147 2148 2149 2150 2151 2152 2153 2154 2155 2156 2157 2158 2159 2160 2161 2162 2163 2164 2165 2166 2167 2168 2169 2170 2171 2172 2173 2174 2175 2176 2177 2178 2179 2180 2181 2182 2183 2184 2185 2186 2187 2188 2189 2190 2191 2192 2193 2194 2195 2196 2197 2198 2199 2200 2201 2202 2203 2204 2205 2206 2207 2208 2209 2210 2211 2212 2213 2214 2215 2216 2217 2218 2219 2220 2221 2222 2223 2224 2225 2226 2227 2228 2229 2230 2231 2232 2233 2234 2235 2236 2237 2238 2239 2240 2241 2242 2243 2244 2245 2246 2247 2248 2249 2250 2251 2252 2253 2254 2255 2256 2257 2258 2259 2260 2261 2262 2263 2264 2265 2266 2267 2268 2269 2270 2271 2272 2273 2274 2275 2276 2277 2278 2279 2280 2281 2282
|
#include "samtools.pysam.h"
/* bam_markdup.c -- Mark duplicates from a coord sorted file that has gone
through fixmates with the mate scoring option on.
Copyright (C) 2017-2022 Genome Research Ltd.
Author: Andrew Whitwham <aw7@sanger.ac.uk>
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
Estimate library size derived from Picard DuplicationMetrics.java
Copyright (c) 2009,2018 The Broad Institute. MIT license.
*/
#include <config.h>
#include <stdlib.h>
#include <string.h>
#include <stdio.h>
#include <unistd.h>
#include <ctype.h>
#include <time.h>
#include <sys/stat.h>
#include <math.h>
#include <regex.h>
#include "htslib/thread_pool.h"
#include "htslib/sam.h"
#include "sam_opts.h"
#include "samtools.h"
#include "htslib/khash.h"
#include "htslib/klist.h"
#include "htslib/kstring.h"
#include "tmp_file.h"
#include "bam.h"
typedef struct {
samFile *in;
samFile *out;
char *prefix;
int remove_dups;
int32_t max_length;
int do_stats;
int supp;
int tag;
int opt_dist;
int no_pg;
int clear;
int mode;
int write_index;
int include_fails;
int check_chain;
char *stats_file;
char *arg_list;
char *out_fn;
regex_t *rgx;
int rgx_x;
int rgx_y;
int rgx_t;
char *barcode;
regex_t *bc_rgx;
} md_param_t;
typedef struct {
hts_pos_t this_coord;
hts_pos_t other_coord;
int32_t this_ref;
int32_t other_ref;
int32_t barcode;
int8_t single;
int8_t leftmost;
int8_t orientation;
} key_data_t;
typedef struct read_queue_s {
key_data_t pair_key;
key_data_t single_key;
bam1_t *b;
struct read_queue_s *duplicate;
hts_pos_t pos;
int dup_checked;
} read_queue_t;
typedef struct {
read_queue_t *p;
} in_hash_t;
typedef struct {
char *name;
char type;
} dup_map_t;
typedef struct {
bam1_t *b;
int64_t score;
int64_t mate_score;
long x;
long y;
int opt;
int beg;
int end;
} check_t;
typedef struct {
check_t *c;
size_t size;
size_t length;
} check_list_t;
static khint32_t do_hash(unsigned char *key, khint32_t len);
static khint_t hash_key(key_data_t key) {
int i = 0;
khint_t hash;
if (key.single) {
unsigned char sig[17];
memcpy(sig + i, &key.this_ref, 4); i += 4;
memcpy(sig + i, &key.this_coord, 8); i += 8;
memcpy(sig + i, &key.orientation, 1); i += 1;
memcpy(sig + i, &key.barcode, 4); i += 4;
hash = do_hash(sig, i);
} else {
unsigned char sig[30];
memcpy(sig + i, &key.this_ref, 4); i += 4;
memcpy(sig + i, &key.this_coord, 8); i += 8;
memcpy(sig + i, &key.other_ref, 4); i += 4;
memcpy(sig + i, &key.other_coord, 8); i += 8;
memcpy(sig + i, &key.leftmost, 1); i += 1;
memcpy(sig + i, &key.orientation, 1); i += 1;
memcpy(sig + i, &key.barcode, 4); i += 4;
hash = do_hash(sig, i);
}
return hash;
}
static int key_equal(key_data_t a, key_data_t b) {
int match = 1;
if (a.this_coord != b.this_coord)
match = 0;
else if (a.orientation != b.orientation)
match = 0;
else if (a.this_ref != b.this_ref)
match = 0;
else if (a.single != b.single)
match = 0;
else if (a.barcode != b.barcode)
match = 0;
if (!a.single) {
if (a.other_coord != b.other_coord)
match = 0;
else if (a.leftmost != b.leftmost)
match = 0;
else if (a.other_ref != b.other_ref)
match = 0;
}
return match;
}
#define __free_queue_element(p)
// Orientations (prime numbers to feed to hashing algorithm)
#define O_FF 2
#define O_RR 3
#define O_FR 5
#define O_RF 7
// Left or rightmost
#define R_LE 11
#define R_RI 13
#define BMD_WARNING_MAX 10
#define MD_MIN_QUALITY 15
// Duplicate finding mode
#define MD_MODE_TEMPLATE 0
#define MD_MODE_SEQUENCE 1
KHASH_INIT(reads, key_data_t, in_hash_t, 1, hash_key, key_equal) // read map hash
KLIST_INIT(read_queue, read_queue_t, __free_queue_element) // the reads buffer
KHASH_MAP_INIT_STR(duplicates, dup_map_t) // map of duplicates for supplementary dup id
/* The Bob Jenkins one_at_a_time hash to reduce the key to a 32 bit value. */
static khint32_t do_hash(unsigned char *key, khint32_t len) {
khint32_t hash, i;
for (hash = 0, i = 0; i < len; ++i) {
hash += key[i];
hash += (hash << 10);
hash ^= (hash >> 6);
}
hash += (hash << 3);
hash ^= (hash >> 11);
hash += (hash << 15);
return hash;
}
/* Get mate score from tag. */
static int64_t get_mate_score(bam1_t *b) {
uint8_t *data;
int64_t score;
if ((data = bam_aux_get(b, "ms"))) {
score = bam_aux2i(data);
} else {
fprintf(samtools_stderr, "[markdup] error: no ms score tag. Please run samtools fixmate on file first.\n");
return -1;
}
return score;
}
/* Calc current score from quality. */
static int64_t calc_score(bam1_t *b)
{
int64_t score = 0;
uint8_t *qual = bam_get_qual(b);
int i;
for (i = 0; i < b->core.l_qseq; i++) {
if (qual[i] >= MD_MIN_QUALITY) score += qual[i];
}
return score;
}
/* Create a signature hash of the current read and its pair.
Uses the unclipped start (or end depending on orientation),
the reference id, orientation and whether the current
read is leftmost of the pair. */
static int make_pair_key(md_param_t *param, key_data_t *key, bam1_t *bam, long *warnings) {
hts_pos_t this_coord, this_end, other_coord, other_end, leftmost;
int32_t this_ref, other_ref, barcode = 0;
int8_t orientation, left_read;
uint8_t *data;
char *cig, *bar;
long incoming_warnings = *warnings;
this_ref = bam->core.tid + 1; // avoid a 0 being put into the hash
other_ref = bam->core.mtid + 1;
this_coord = unclipped_start(bam);
this_end = unclipped_end(bam);
if ((data = bam_aux_get(bam, "MC"))) {
if (!(cig = bam_aux2Z(data))) {
fprintf(samtools_stderr, "[markdup] error: MC tag wrong type. Please use the MC tag provided by samtools fixmate.\n");
return 1;
}
other_end = unclipped_other_end(bam->core.mpos, cig);
other_coord = unclipped_other_start(bam->core.mpos, cig);
} else {
fprintf(samtools_stderr, "[markdup] error: no MC tag. Please run samtools fixmate on file first.\n");
return 1;
}
// work out orientations
if (param->mode == MD_MODE_TEMPLATE) {
if (this_ref != other_ref) {
leftmost = this_ref < other_ref;
} else {
if (bam_is_rev(bam) == bam_is_mrev(bam)) {
if (!bam_is_rev(bam)) {
leftmost = this_coord <= other_coord;
} else {
leftmost = this_end <= other_end;
}
} else {
if (bam_is_rev(bam)) {
leftmost = this_end <= other_coord;
} else {
leftmost = this_coord <= other_end;
}
}
}
// pair orientation
if (leftmost) {
if (bam_is_rev(bam) == bam_is_mrev(bam)) {
other_coord = other_end;
if (!bam_is_rev(bam)) {
if (bam->core.flag & BAM_FREAD1) {
orientation = O_FF;
} else {
orientation = O_RR;
}
} else {
if (bam->core.flag & BAM_FREAD1) {
orientation = O_RR;
} else {
orientation = O_FF;
}
}
} else {
if (!bam_is_rev(bam)) {
orientation = O_FR;
other_coord = other_end;
} else {
orientation = O_RF;
this_coord = this_end;
}
}
} else {
if (bam_is_rev(bam) == bam_is_mrev(bam)) {
this_coord = this_end;
if (!bam_is_rev(bam)) {
if (bam->core.flag & BAM_FREAD1) {
orientation = O_RR;
} else {
orientation = O_FF;
}
} else {
if (bam->core.flag & BAM_FREAD1) {
orientation = O_FF;
} else {
orientation = O_RR;
}
}
} else {
if (!bam_is_rev(bam)) {
orientation = O_RF;
other_coord = other_end;
} else {
orientation = O_FR;
this_coord = this_end;
}
}
}
} else { // MD_MODE_SEQUENCE
if (this_ref != other_ref) {
leftmost = this_ref - other_ref;
} else {
if (bam_is_rev(bam) == bam_is_mrev(bam)) {
if (!bam_is_rev(bam)) {
leftmost = this_coord - other_coord;
} else {
leftmost = this_end - other_end;
}
} else {
if (bam_is_rev(bam)) {
leftmost = this_end - other_coord;
} else {
leftmost = this_coord - other_end;
}
}
}
if (leftmost < 0) {
leftmost = 1;
} else if (leftmost > 0) {
leftmost = 0;
} else {
// tie breaks
if (bam->core.pos == bam->core.mpos) {
if (bam->core.flag & BAM_FREAD1) {
leftmost = 1;
} else {
leftmost = 0;
}
} else if (bam->core.pos < bam->core.mpos) {
leftmost = 1;
} else {
leftmost = 0;
}
}
// pair orientation
if (leftmost) {
if (bam_is_rev(bam) == bam_is_mrev(bam)) {
if (!bam_is_rev(bam)) {
orientation = O_FF;
} else {
orientation = O_RR;
}
} else {
if (!bam_is_rev(bam)) {
orientation = O_FR;
} else {
orientation = O_RF;
}
}
} else {
if (bam_is_rev(bam) == bam_is_mrev(bam)) {
if (!bam_is_rev(bam)) {
orientation = O_RR;
} else {
orientation = O_FF;
}
} else {
if (!bam_is_rev(bam)) {
orientation = O_RF;
} else {
orientation = O_FR;
}
}
}
if (!bam_is_rev(bam)) {
this_coord = unclipped_start(bam);
} else {
this_coord = unclipped_end(bam);
}
if (!bam_is_mrev(bam)) {
other_coord = unclipped_other_start(bam->core.mpos, cig);
} else {
other_coord = unclipped_other_end(bam->core.mpos, cig);
}
}
if (!leftmost)
left_read = R_RI;
else
left_read = R_LE;
if (param->barcode) {
if ((data = bam_aux_get(bam, param->barcode))) {
if (!(bar = bam_aux2Z(data))) {
(*warnings)++;
if (*warnings <= BMD_WARNING_MAX) {
fprintf(samtools_stderr, "[markdup] warning: %s tag wrong type. Aux tag needs to be a string type.\n", param->barcode);
}
} else {
barcode = do_hash((unsigned char *)bar, strlen(bar));
}
}
} else if (param->bc_rgx) {
int result;
regmatch_t matches[3];
size_t max_matches = 2;
char *qname = bam_get_qname(bam);
if ((result = regexec(param->bc_rgx, qname, max_matches, matches, 0)) == 0) {
int bc_start, bc_end;
bc_start = matches[1].rm_so;
bc_end = matches[1].rm_eo;
if (bc_start != -1) {
barcode = do_hash((unsigned char *)qname + bc_start, bc_end - bc_start);
} else {
(*warnings)++;
if (*warnings <= BMD_WARNING_MAX) {
fprintf(samtools_stderr, "[markdup] warning: barcode regex unable to match substring on %s.\n", qname);
}
}
} else {
(*warnings)++;
if (*warnings <= BMD_WARNING_MAX) {
char warn_msg[256];
regerror(result, param->bc_rgx, warn_msg, 256);
fprintf(samtools_stderr, "[markdup] warning: barcode regex match error \"%s\" on %s.\n", warn_msg, qname);
}
}
}
if ((*warnings == BMD_WARNING_MAX) && (incoming_warnings != *warnings)) {
fprintf(samtools_stderr, "[markdup] warning: %ld barcode read warnings. New warnings will not be reported.\n",
*warnings);
}
key->single = 0;
key->this_ref = this_ref;
key->this_coord = this_coord;
key->other_ref = other_ref;
key->other_coord = other_coord;
key->leftmost = left_read;
key->orientation = orientation;
key->barcode = barcode;
return 0;
}
/* Create a signature hash of single read (or read with an unmatched pair).
Uses unclipped start (or end depending on orientation), reference id,
and orientation. */
static void make_single_key(md_param_t *param, key_data_t *key, bam1_t *bam, long *warnings) {
hts_pos_t this_coord;
int32_t this_ref, barcode = 0;
int8_t orientation;
uint8_t *data;
char *bar;
long incoming_warnings = *warnings;
this_ref = bam->core.tid + 1; // avoid a 0 being put into the hash
if (bam_is_rev(bam)) {
this_coord = unclipped_end(bam);
orientation = O_RR;
} else {
this_coord = unclipped_start(bam);
orientation = O_FF;
}
if (param->barcode) {
if ((data = bam_aux_get(bam, param->barcode))) {
if (!(bar = bam_aux2Z(data))) {
(*warnings)++;
if (*warnings <= BMD_WARNING_MAX) {
fprintf(samtools_stderr, "[markdup] warning: %s tag wrong type. Aux tag needs to be a string type.\n", param->barcode);
}
} else {
barcode = do_hash((unsigned char *)bar, strlen(bar));
}
}
} else if (param->bc_rgx) {
int result;
regmatch_t matches[3];
size_t max_matches = 2;
char *qname = bam_get_qname(bam);
if ((result = regexec(param->bc_rgx, qname, max_matches, matches, 0)) == 0) {
int bc_start, bc_end;
bc_start = matches[1].rm_so;
bc_end = matches[1].rm_eo;
if (bc_start != -1) {
barcode = do_hash((unsigned char *)qname + bc_start, bc_end - bc_start);
} else {
(*warnings)++;
if (*warnings <= BMD_WARNING_MAX) {
fprintf(samtools_stderr, "[markdup] warning: barcode regex unable to match substring on %s.\n", qname);
}
}
} else {
(*warnings)++;
if (*warnings <= BMD_WARNING_MAX) {
char warn_msg[256];
regerror(result, param->bc_rgx, warn_msg, 256);
fprintf(samtools_stderr, "[markdup] warning: barcode regex match error \"%s\" on %s.\n", warn_msg, qname);
}
}
}
if ((*warnings == BMD_WARNING_MAX) && (incoming_warnings != *warnings)) {
fprintf(samtools_stderr, "[markdup] warning: %ld barcode read warnings. New warnings will not be reported.\n",
*warnings);
}
key->single = 1;
key->this_ref = this_ref;
key->this_coord = this_coord;
key->orientation = orientation;
key->barcode = barcode;
}
/* Add the duplicate name to a hash if it does not exist. */
static int add_duplicate(khash_t(duplicates) *d_hash, bam1_t *dupe, char *orig_name, char type) {
khiter_t d;
int ret;
d = kh_get(duplicates, d_hash, bam_get_qname(dupe));
if (d == kh_end(d_hash)) {
char *name = strdup(bam_get_qname(dupe));
if (name) {
d = kh_put(duplicates, d_hash, name, &ret);
} else {
ret = -1;
}
if (ret >= 0) {
if (orig_name) {
if (ret == 0) {
// replace old name
free(kh_value(d_hash, d).name);
free(name);
}
kh_value(d_hash, d).name = strdup(orig_name);
if (kh_value(d_hash, d).name == NULL) {
fprintf(samtools_stderr, "[markdup] error: unable to allocate memory for duplicate original name.\n");
return 1;
}
} else {
kh_value(d_hash, d).name = NULL;
}
kh_value(d_hash, d).type = type;
} else {
fprintf(samtools_stderr, "[markdup] error: unable to store supplementary duplicates.\n");
free(name);
return 1;
}
}
return 0;
}
/* Get coordinates from the standard Illumina style read names.
Returned values are of the x and y coordinates and a section of
the read name to test (t) for string equality e.g. lane and tile part. */
static int get_coordinates_colons(md_param_t *param, const char *qname, int *t_beg, int *t_end, long *x_coord, long *y_coord, long *warnings) {
int sep = 0;
int pos = 0;
int xpos = 0, ypos = 0;
char *end;
while (qname[pos]) {
if (qname[pos] == ':') {
sep++;
if (sep == 2) {
xpos = pos + 1;
} else if (sep == 3) {
ypos = pos + 1;
} else if (sep == 4) { // HiSeq style names
xpos = ypos;
ypos = pos + 1;
} else if (sep == 5) { // Newer Illumina format
xpos = pos + 1;
} else if (sep == 6) {
ypos = pos + 1;
}
}
pos++;
}
/* The most current Illumina read format at time of writing is:
@machine:run:flowcell:lane:tile:x:y:UMI or
@machine:run:flowcell:lane:tile:x:y
Counting the separating colons gives us a quick format check.
Older name formats have fewer elements.
*/
if (!(sep == 3 || sep == 4 || sep == 6 || sep == 7)) {
(*warnings)++;
if (*warnings <= BMD_WARNING_MAX) {
fprintf(samtools_stderr, "[markdup] warning: cannot decipher read name %s for optical duplicate marking.\n", qname);
}
return 1;
} else {
*x_coord = strtol(qname + xpos, &end, 10);
if ((qname + xpos) == end) {
(*warnings)++;
if (*warnings <= BMD_WARNING_MAX) {
fprintf(samtools_stderr, "[markdup] warning: cannot decipher x coordinate in %s .\n", qname);
}
return 1;
}
*y_coord = strtol(qname + ypos, &end, 10);
if ((qname + ypos) == end) {
(*warnings)++;
if (*warnings <= BMD_WARNING_MAX) {
fprintf(samtools_stderr, "[markdup] warning: cannot decipher y coordinate in %s .\n", qname);
}
return 1;
}
*t_beg = 0;
*t_end = xpos;
}
return 0;
}
/* Get the coordinates from the read name.
Returned values are of the x and y coordinates and an optional section of
the read name to test (t) for string equality e.g. lane and tile part. */
static inline int get_coordinates_regex(md_param_t *param, const char *qname, int *t_beg, int *t_end, long *x_coord, long *y_coord, long *warnings) {
regmatch_t matches[5];
size_t max_matches = 5;
int xpos, ypos, xend, yend, xlen, ylen;
char coord[255];
char *end;
if (!param->rgx_t)
max_matches = 4;
if (regexec(param->rgx, qname, max_matches, matches, 0))
return -1;
xpos = matches[param->rgx_x].rm_so;
ypos = matches[param->rgx_y].rm_so;
if (param->rgx_t) {
*t_beg = matches[param->rgx_t].rm_so;
*t_end = matches[param->rgx_t].rm_eo;
} else {
*t_beg = *t_end = 0;
}
if (xpos == -1 || ypos == -1 || *t_beg == -1)
return -1;
xend = matches[param->rgx_x].rm_eo;
yend = matches[param->rgx_y].rm_eo;
if ((xlen = xend - xpos) > 254) {
(*warnings)++;
if (*warnings <= BMD_WARNING_MAX) {
fprintf(samtools_stderr, "[markdup] warning: x coordinate string longer than allowed qname length in %s (%d long).\n", qname, xlen);
}
return 1;
}
strncpy(coord, qname + xpos, xlen);
coord[xlen] = '\0';
*x_coord = strtol(coord, &end, 10);
if (coord == end) {
(*warnings)++;
if (*warnings <= BMD_WARNING_MAX) {
fprintf(samtools_stderr, "[markdup] warning: cannot decipher x coordinate in %s (%s).\n", qname, coord);
}
return 1;
}
if ((ylen = yend - ypos) > 254) {
(*warnings)++;
if (*warnings <= BMD_WARNING_MAX) {
fprintf(samtools_stderr, "[markdup] warning: y coordinate string longer than allowed qname length in %s (%d long).\n", qname, ylen);
}
return 1;
}
strncpy(coord, qname + ypos, ylen);
coord[ylen] = '\0';
*y_coord = strtol(coord, &end, 10);
if (coord == end) {
(*warnings)++;
if (*warnings <= BMD_WARNING_MAX) {
fprintf(samtools_stderr, "[markdup] warning: cannot decipher y coordinate in %s (%s).\n", qname, coord);
}
return 1;
}
return 0;
}
static int get_coordinates(md_param_t *param, const char *name, int *t_beg, int *t_end, long *x_coord, long *y_coord, long *warnings) {
int ret = 1;
if (param->rgx == NULL) {
ret = get_coordinates_colons(param, name, t_beg, t_end, x_coord, y_coord, warnings);
} else {
ret = get_coordinates_regex(param, name, t_beg, t_end, x_coord, y_coord, warnings);
}
return ret;
}
/* Using the coordinates from the read name, see whether the duplicated read is
close enough (set by max_dist) to the original to be counted as optical.*/
static int is_optical_duplicate(md_param_t *param, bam1_t *ori, bam1_t *dup, long max_dist, long *warnings) {
int ret = 0;
char *original, *duplicate;
long ox, oy, dx, dy;
int o_beg = 0, o_end = 0, d_beg = 0, d_end = 0;
original = bam_get_qname(ori);
duplicate = bam_get_qname(dup);
if (get_coordinates(param, original, &o_beg, &o_end, &ox, &oy, warnings)) {
return ret;
}
if (get_coordinates(param, duplicate, &d_beg, &d_end, &dx, &dy, warnings)) {
return ret;
}
if (strncmp(original + o_beg, duplicate + d_beg, o_end - o_beg) == 0) {
long xdiff, ydiff;
if (ox > dx) {
xdiff = ox - dx;
} else {
xdiff = dx - ox;
}
if (xdiff <= max_dist) {
// still might be optical
if (oy > dy) {
ydiff = oy - dy;
} else {
ydiff = dy - oy;
}
if (ydiff <= max_dist) ret = 1;
}
}
return ret;
}
/* Using the coordinates from the Illumina read name, see whether the duplicated read is
close enough (set by max_dist) to the original to be counted as optical.
This function needs the values from the first read to be already calculated. */
static int optical_duplicate_partial(md_param_t *param, const char *name, const int o_beg, const int o_end, const long ox, const long oy, bam1_t *dup, check_t *c, long max_dist, long *warnings) {
int ret = 0;
char *duplicate;
int d_beg = 0, d_end = 0;
long dx, dy;
duplicate = bam_get_qname(dup);
if (get_coordinates(param, duplicate, &d_beg, &d_end, &dx, &dy, warnings)) {
return ret;
}
if (strncmp(name + o_beg, duplicate + d_beg, o_end - o_beg) == 0) {
// the initial parts match, look at the numbers
long xdiff, ydiff;
if (ox > dx) {
xdiff = ox - dx;
} else {
xdiff = dx - ox;
}
if (xdiff <= max_dist) {
// still might be optical
if (oy > dy) {
ydiff = oy - dy;
} else {
ydiff = dy - oy;
}
if (ydiff <= max_dist) ret = 1;
}
}
c->x = dx;
c->y = dy;
c->beg = d_beg;
c->end = d_end;
return ret;
}
/* Mark the read as a duplicate and update the duplicate hash (if needed) */
static int mark_duplicates(md_param_t *param, khash_t(duplicates) *dup_hash, bam1_t *ori, bam1_t *dup,
long *optical, long *warn) {
char dup_type = 0;
long incoming_warnings = *warn;
dup->core.flag |= BAM_FDUP;
if (param->tag) {
if (bam_aux_update_str(dup, "do", strlen(bam_get_qname(ori)) + 1, bam_get_qname(ori))) {
fprintf(samtools_stderr, "[markdup] error: unable to append 'do' tag.\n");
return -1;
}
}
if (param->opt_dist) { // mark optical duplicates
if (is_optical_duplicate(param, ori, dup, param->opt_dist, warn)) {
bam_aux_update_str(dup, "dt", 3, "SQ");
dup_type = 'O';
(*optical)++;
} else {
// not an optical duplicate
bam_aux_update_str(dup, "dt", 3, "LB");
}
}
if ((*warn == BMD_WARNING_MAX) && (incoming_warnings != *warn)) {
fprintf(samtools_stderr, "[markdup] warning: %ld decipher read name warnings. New warnings will not be reported.\n",
*warn);
}
if (param->supp) {
if (bam_aux_get(dup, "SA") || (dup->core.flag & BAM_FMUNMAP) || bam_aux_get(dup, "XA")) {
char *original = NULL;
if (param->tag) {
original = bam_get_qname(ori);
}
if (add_duplicate(dup_hash, dup, original, dup_type))
return -1;
}
}
return 0;
}
/* If the duplicate type has changed to optical then retag and duplicate hash. */
static inline int optical_retag(md_param_t *param, khash_t(duplicates) *dup_hash, bam1_t *b, int paired, long *optical_single, long *optical_pair) {
int ret = 0;
if (bam_aux_update_str(b, "dt", 3, "SQ")) {
fprintf(samtools_stderr, "[markdup] error: unable to update 'dt' tag.\n");
ret = -1;
}
if (paired) {
(*optical_pair)++;
} else {
(*optical_single)++;
}
if (param->supp) {
// Change the duplicate type
if (bam_aux_get(b, "SA") || (b->core.flag & BAM_FMUNMAP)
|| bam_aux_get(b, "XA")) {
khiter_t d;
d = kh_get(duplicates, dup_hash, bam_get_qname(b));
if (d == kh_end(dup_hash)) {
// error, name should already be in dup hash
fprintf(samtools_stderr, "[markdup] error: duplicate name %s not found in hash.\n",
bam_get_qname(b));
ret = -1;
} else {
kh_value(dup_hash, d).type = 'O';
}
}
}
return ret;
}
/* Check all duplicates of the highest quality read (the "original") for consistancy. Also
pre-calculate any values for use in check_duplicate_chain later.
Returns 0 on success, >0 on coordinate reading error (program can continue) or
<0 on an error (program should not continue. */
static int check_chain_against_original(md_param_t *param, khash_t(duplicates) *dup_hash, read_queue_t *ori,
check_list_t *list, long *warn, long *optical_single, long *optical_pair) {
int ret = 0, coord_fail = 0;
char *ori_name = bam_get_qname(ori->b);
read_queue_t *current = ori->duplicate;
int t_beg = 0, t_end = 0;
long x, y;
if (param->opt_dist) {
coord_fail = get_coordinates(param, ori_name, &t_beg, &t_end, &x, &y, warn);
}
list->length = 0;
while (current) {
check_t *c;
if (list->length >= list->size) {
check_t *tmp;
list->size *= 2;
if (!(tmp = realloc(list->c, list->size * sizeof(check_t)))) {
fprintf(samtools_stderr, "[markdup] error: Unable to expand opt check list.\n");
return -1;
}
list->c = tmp;
}
c = &list->c[list->length];
c->b = current->b;
c->x = -1;
c->y = -1;
c->opt = 0;
c->score = 0;
c->mate_score = 0;
current->dup_checked = 1;
if (param->tag) {
uint8_t *data;
// at this stage all duplicates should have a do tag
if ((data = bam_aux_get(current->b, "do")) != NULL) {
// see if we need to change the tag
char *old_name = bam_aux2Z(data);
if (old_name) {
if (strcmp(old_name, ori_name) != 0) {
if (bam_aux_update_str(current->b, "do", strlen(ori_name) + 1, (const char *)ori_name)) {
fprintf(samtools_stderr, "[markdup] error: unable to update 'do' tag.\n");
ret = -1;
break;
}
}
} else {
fprintf(samtools_stderr, "[markdup] error: 'do' tag has wrong type for read %s.\n", bam_get_qname(current->b));
ret = -1;
break;
}
}
}
if (param->opt_dist && !coord_fail) {
uint8_t *data;
char *dup_type;
int is_opt = 0;
int current_paired = (current->b->core.flag & BAM_FPAIRED) && !(current->b->core.flag & BAM_FMUNMAP);
if ((data = bam_aux_get(current->b, "dt"))) {
if ((dup_type = bam_aux2Z(data))) {
if (strcmp(dup_type, "SQ") == 0) {
c->opt = 1;
}
}
}
// need to run this to get the duplicates x and y scores
is_opt = optical_duplicate_partial(param, ori_name, t_beg, t_end, x, y, current->b, c, param->opt_dist, warn);
if (!c->opt && is_opt) {
if (optical_retag(param, dup_hash, current->b, current_paired, optical_single, optical_pair)) {
ret = -1;
break;
}
c->opt = 1;
}
c->score = calc_score(current->b);
if (current_paired) {
if ((c->mate_score = get_mate_score(current->b)) == -1) {
fprintf(samtools_stderr, "[markdup] error: no ms score tag. Please run samtools fixmate on file first.\n");
ret = -1;
break;
}
}
}
current = current->duplicate;
list->length++;
}
if (!ret && coord_fail)
ret = coord_fail;
return ret;
}
static int xcoord_sort(const void *a, const void *b) {
check_t *ac = (check_t *) a;
check_t *bc = (check_t *) b;
return (ac->x - bc->x);
}
/* Check all the duplicates against each other to see if they are optical duplicates. */
static int check_duplicate_chain(md_param_t *param, khash_t(duplicates) *dup_hash, check_list_t *list,
long *warn, long *optical_single, long *optical_pair) {
int ret = 0;
size_t curr = 0;
qsort(list->c, list->length, sizeof(list->c[0]), xcoord_sort);
while (curr < list->length - 1) {
check_t *current = &list->c[curr];
size_t count = curr;
char *cur_name = bam_get_qname(current->b);
int current_paired = (current->b->core.flag & BAM_FPAIRED) && !(current->b->core.flag & BAM_FMUNMAP);
while (++count < list->length && (list->c[count].x - current->x <= param->opt_dist)) {
// while close enough along the x coordinate
check_t *chk = &list->c[count];
if (current->opt && chk->opt)
continue;
// if both are already optical duplicates there is no need to check again, otherwise...
long ydiff;
if (current->y > chk->y) {
ydiff = current->y - chk->y;
} else {
ydiff = chk->y - current->y;
}
if (ydiff > param->opt_dist)
continue;
// the number are right, check the names
if (strncmp(cur_name + current->beg, bam_get_qname(chk->b) + chk->beg, current->end - current->beg) != 0)
continue;
// optical duplicates
int chk_dup = 0;
int chk_paired = (chk->b->core.flag & BAM_FPAIRED) && !(chk->b->core.flag & BAM_FMUNMAP);
if (current_paired != chk_paired) {
if (!chk_paired) {
// chk is single vs pair, this is a dup.
chk_dup = 1;
}
} else {
// do it by scores
int64_t cur_score, chk_score;
if ((current->b->core.flag & BAM_FQCFAIL) != (chk->b->core.flag & BAM_FQCFAIL)) {
if (current->b->core.flag & BAM_FQCFAIL) {
cur_score = 0;
chk_score = 1;
} else {
cur_score = 1;
chk_score = 0;
}
} else {
cur_score = current->score;
chk_score = chk->score;
if (current_paired) {
// they are pairs so add mate scores.
chk_score += chk->mate_score;
cur_score += current->mate_score;
}
}
if (cur_score == chk_score) {
if (strcmp(bam_get_qname(chk->b), cur_name) < 0) {
chk_score++;
} else {
chk_score--;
}
}
if (cur_score > chk_score) {
chk_dup = 1;
}
}
if (chk_dup) {
// the duplicate is the optical duplicate
if (!chk->opt) { // only change if not already an optical duplicate
if (optical_retag(param, dup_hash, chk->b, chk_paired, optical_single, optical_pair)) {
ret = -1;
goto fail;
}
chk->opt = 1;
}
} else {
if (!current->opt) {
if (optical_retag(param, dup_hash, current->b, current_paired, optical_single, optical_pair)) {
ret = -1;
goto fail;
}
current->opt = 1;
}
}
}
curr++;
}
fail:
return ret;
}
/* Where there is more than one duplicate go down the list and check for optical duplicates and change
do tags (where used) to point to original (non-duplicate) read. */
static int find_duplicate_chains(md_param_t *param, klist_t(read_queue) *read_buffer, khash_t(duplicates) *dup_hash, check_list_t *dup_list,
const hts_pos_t prev_coord, const int32_t prev_tid, long *warn, long *optical_single,
long *optical_pair, const int check_range) {
int ret = 0;
kliter_t(read_queue) *rq;
rq = kl_begin(read_buffer);
while (rq != kl_end(read_buffer)) {
read_queue_t *in_read = &kl_val(rq);
if (check_range) {
/* Just check against the moving window of reads based on coordinates and max read length. */
if (in_read->pos + param->max_length > prev_coord && in_read->b->core.tid == prev_tid && (prev_tid != -1 || prev_coord != -1)) {
break;
}
} else {
// this is the last set of results and the end entry will be blank
if (!bam_get_qname(in_read->b)) {
break;
}
}
if (!(in_read->b->core.flag & BAM_FDUP) && in_read->duplicate) { // is the head of a duplicate chain
// check against the original for tagging and optical duplication
if ((ret = check_chain_against_original(param, dup_hash, in_read, dup_list, warn, optical_single, optical_pair))) {
if (ret < 0) { // real error
ret = -1;
break;
} else { // coordinate decoding error
ret = 0;
in_read->duplicate = NULL;
continue;
}
}
// check the rest of the duplicates against each other for optical duplication
if (param->opt_dist && check_duplicate_chain(param, dup_hash, dup_list, warn, optical_single, optical_pair)) {
ret = -1;
break;
}
in_read->duplicate = NULL;
}
rq = kl_next(rq);
}
return ret;
}
/*
Function to use when estimating library size.
This is based on an approximate formula for the coverage of a set
obtained after sampling it a given number of times with replacement.
x = number of items in the set (the number of unique fragments in the library)
c = number of unique items (unique read pairs observed)
n = number of items samples (total number of read pairs)
c and n are known; x is unknown.
As n -> infinity, the coverage (c/x) can be given as:
c / x = 1 - exp(-n / x) (see https://math.stackexchange.com/questions/32800)
This needs to be solved for x, so it is rearranged to put both terms on the
left side and estimate_library_size() finds a value of x which gives a
result of zero (or as close as it can get).
*/
static inline double coverage_equation(double x, double c, double n) {
return c / x - 1 + exp(-n / x);
}
/* estimate the library size, based on the Picard code in DuplicationMetrics.java*/
static unsigned long estimate_library_size(unsigned long paired_reads, unsigned long paired_duplicate_reads, unsigned long optical) {
unsigned long estimated_size = 0;
unsigned long non_optical_pairs = (paired_reads - optical) / 2;
unsigned long unique_pairs = (paired_reads - paired_duplicate_reads) / 2;
unsigned long duplicate_pairs = (paired_duplicate_reads - optical) / 2;
if ((non_optical_pairs && duplicate_pairs && unique_pairs) && (non_optical_pairs > duplicate_pairs)) {
double m = 1;
double M = 100;
int i;
if (coverage_equation(m * (double)unique_pairs, (double)unique_pairs, (double)non_optical_pairs) < 0) {
fprintf(samtools_stderr, "[markdup] warning: unable to calculate estimated library size.\n");
return estimated_size;
}
while (coverage_equation(M * (double)unique_pairs, (double)unique_pairs, (double)non_optical_pairs) > 0) {
M *= 10;
}
for (i = 0; i < 40; i++) {
double r = (m + M) / 2;
double u = coverage_equation(r * (double)unique_pairs, (double)unique_pairs, (double)non_optical_pairs);
if (u > 0) {
m = r;
} else if (u < 0) {
M = r;
} else {
break;
}
}
estimated_size = (unsigned long)(unique_pairs * (m + M) / 2);
} else {
fprintf(samtools_stderr, "[markdup] warning: unable to calculate estimated library size."
" Read pairs %ld should be greater than duplicate pairs %ld,"
" which should both be non zero.\n",
non_optical_pairs, duplicate_pairs);
}
return estimated_size;
}
/* Compare the reads near each other (coordinate sorted) and try to spot the duplicates.
Generally the highest quality scoring is chosen as the original and all others the duplicates.
The score is based on the sum of the quality values (<= 15) of the read and its mate (if any).
While single reads are compared to only one read of a pair, the pair will chosen as the original.
The comparison is done on position and orientation, see above for details.
Marking the supplementary reads of a duplicate as also duplicates takes an extra file read/write
step. This is because the duplicate can occur before the primary read.*/
static int bam_mark_duplicates(md_param_t *param) {
bam_hdr_t *header = NULL;
khiter_t k;
khash_t(reads) *pair_hash = kh_init(reads);
khash_t(reads) *single_hash = kh_init(reads);
klist_t(read_queue) *read_buffer = kl_init(read_queue);
kliter_t(read_queue) *rq;
khash_t(duplicates) *dup_hash = kh_init(duplicates);
int32_t prev_tid;
hts_pos_t prev_coord;
read_queue_t *in_read;
int ret;
long reading, writing, excluded, duplicate, single, pair, single_dup, examined, optical, single_optical;
long np_duplicate, np_opt_duplicate;
long opt_warnings = 0, bc_warnings = 0;
tmp_file_t temp;
char *idx_fn = NULL;
int exclude = 0;
check_list_t dup_list = {NULL, 0, 0};
if (!pair_hash || !single_hash || !read_buffer || !dup_hash) {
fprintf(samtools_stderr, "[markdup] out of memory\n");
goto fail;
}
if ((header = sam_hdr_read(param->in)) == NULL) {
fprintf(samtools_stderr, "[markdup] error reading header\n");
goto fail;
}
// accept unknown, unsorted or coordinate sort order, but error on queryname sorted.
// only really works on coordinate sorted files.
kstring_t str = KS_INITIALIZE;
if (!sam_hdr_find_tag_hd(header, "SO", &str) && str.s && !strcmp(str.s, "queryname")) {
fprintf(samtools_stderr, "[markdup] error: queryname sorted, must be sorted by coordinate.\n");
ks_free(&str);
goto fail;
}
ks_free(&str);
if (!param->no_pg && sam_hdr_add_pg(header, "samtools", "VN", samtools_version(),
param->arg_list ? "CL" : NULL,
param->arg_list ? param->arg_list : NULL,
NULL) != 0) {
fprintf(samtools_stderr, "[markdup] warning: unable to add @PG line to header.\n");
}
if (sam_hdr_write(param->out, header) < 0) {
fprintf(samtools_stderr, "[markdup] error writing header.\n");
goto fail;
}
if (param->write_index) {
if (!(idx_fn = auto_index(param->out, param->out_fn, header)))
goto fail;
}
// used for coordinate order checks
prev_tid = prev_coord = 0;
// get the buffer going
in_read = kl_pushp(read_queue, read_buffer);
if (!in_read) {
fprintf(samtools_stderr, "[markdup] out of memory\n");
goto fail;
}
// handling supplementary reads needs a temporary file
if (param->supp) {
if (tmp_file_open_write(&temp, param->prefix, 1)) {
fprintf(samtools_stderr, "[markdup] error: unable to open tmp file %s.\n", param->prefix);
goto fail;
}
}
if ((in_read->b = bam_init1()) == NULL) {
fprintf(samtools_stderr, "[markdup] error: unable to allocate memory for alignment.\n");
goto fail;
}
if (param->check_chain && !(param->tag || param->opt_dist))
param->check_chain = 0;
if (param->check_chain) {
dup_list.size = 128;
dup_list.c = NULL;
if ((dup_list.c = malloc(dup_list.size * sizeof(check_t))) == NULL) {
fprintf(samtools_stderr, "[markdup] error: unable to allocate memory for dup_list.\n");
goto fail;
}
}
reading = writing = excluded = single_dup = duplicate = examined = pair = single = optical = single_optical = 0;
np_duplicate = np_opt_duplicate = 0;
while ((ret = sam_read1(param->in, header, in_read->b)) >= 0) {
int dup_checked = 0;
// do some basic coordinate order checks
if (in_read->b->core.tid >= 0) { // -1 for unmapped reads
if (in_read->b->core.tid < prev_tid ||
((in_read->b->core.tid == prev_tid) && (in_read->b->core.pos < prev_coord))) {
fprintf(samtools_stderr, "[markdup] error: not in coordinate sorted order.\n");
goto fail;
}
}
prev_coord = in_read->pos = in_read->b->core.pos;
prev_tid = in_read->b->core.tid;
in_read->pair_key.single = 1;
in_read->single_key.single = 0;
in_read->duplicate = NULL;
in_read->dup_checked = 0;
reading++;
if (param->clear && (in_read->b->core.flag & BAM_FDUP)) {
uint8_t *data;
in_read->b->core.flag ^= BAM_FDUP;
if ((data = bam_aux_get(in_read->b, "dt")) != NULL) {
bam_aux_del(in_read->b, data);
}
if ((data = bam_aux_get(in_read->b, "do")) != NULL) {
bam_aux_del(in_read->b, data);
}
}
if (param->include_fails) {
exclude |= (BAM_FSECONDARY | BAM_FSUPPLEMENTARY | BAM_FUNMAP);
} else {
exclude |= (BAM_FSECONDARY | BAM_FSUPPLEMENTARY | BAM_FUNMAP | BAM_FQCFAIL);
}
// read must not be secondary, supplementary, unmapped or (possibly) failed QC
if (!(in_read->b->core.flag & exclude)) {
examined++;
// look at the pairs first
if ((in_read->b->core.flag & BAM_FPAIRED) && !(in_read->b->core.flag & BAM_FMUNMAP)) {
int ret, mate_tmp;
key_data_t pair_key;
key_data_t single_key;
in_hash_t *bp;
if (make_pair_key(param, &pair_key, in_read->b, &bc_warnings)) {
fprintf(samtools_stderr, "[markdup] error: unable to assign pair hash key.\n");
goto fail;
}
make_single_key(param, &single_key, in_read->b, &bc_warnings);
pair++;
in_read->pos = single_key.this_coord; // cigar/orientation modified pos
// put in singles hash for checking against non paired reads
k = kh_put(reads, single_hash, single_key, &ret);
if (ret > 0) { // new
// add to single duplicate hash
bp = &kh_val(single_hash, k);
bp->p = in_read;
in_read->single_key = single_key;
} else if (ret == 0) { // exists
// look at singles only for duplication marking
bp = &kh_val(single_hash, k);
if (!(bp->p->b->core.flag & BAM_FPAIRED) || (bp->p->b->core.flag & BAM_FMUNMAP)) {
// singleton will always be marked duplicate even if
// scores more than one read of the pair
bam1_t *dup = bp->p->b;
if (param->check_chain)
in_read->duplicate = bp->p;
bp->p = in_read;
if (mark_duplicates(param, dup_hash, bp->p->b, dup, &single_optical, &opt_warnings))
goto fail;
single_dup++;
}
} else {
fprintf(samtools_stderr, "[markdup] error: single hashing failure.\n");
goto fail;
}
// now do the pair
k = kh_put(reads, pair_hash, pair_key, &ret);
if (ret > 0) { // new
// add to the pair hash
bp = &kh_val(pair_hash, k);
bp->p = in_read;
in_read->pair_key = pair_key;
} else if (ret == 0) {
int64_t old_score, new_score, tie_add = 0;
bam1_t *dup = NULL;
bp = &kh_val(pair_hash, k);
if ((bp->p->b->core.flag & BAM_FQCFAIL) != (in_read->b->core.flag & BAM_FQCFAIL)) {
if (bp->p->b->core.flag & BAM_FQCFAIL) {
old_score = 0;
new_score = 1;
} else {
old_score = 1;
new_score = 0;
}
} else {
if ((mate_tmp = get_mate_score(bp->p->b)) == -1) {
fprintf(samtools_stderr, "[markdup] error: no ms score tag. Please run samtools fixmate on file first.\n");
goto fail;
} else {
old_score = calc_score(bp->p->b) + mate_tmp;
}
if ((mate_tmp = get_mate_score(in_read->b)) == -1) {
fprintf(samtools_stderr, "[markdup] error: no ms score tag. Please run samtools fixmate on file first.\n");
goto fail;
} else {
new_score = calc_score(in_read->b) + mate_tmp;
}
}
// choose the highest score as the original
// and add it to the pair hash, mark the other as duplicate
if (new_score == old_score) {
if (strcmp(bam_get_qname(in_read->b), bam_get_qname(bp->p->b)) < 0) {
tie_add = 1;
} else {
tie_add = -1;
}
}
if (new_score + tie_add > old_score) { // swap reads
dup = bp->p->b;
if (param->check_chain) {
if (in_read->duplicate) {
read_queue_t *current = in_read->duplicate;
while (current->duplicate) {
current = current->duplicate;
}
current->duplicate = bp->p;
} else {
in_read->duplicate = bp->p;
}
}
bp->p = in_read;
} else {
if (param->check_chain) {
if (bp->p->duplicate) {
if (in_read->duplicate) {
read_queue_t *current = bp->p->duplicate;
while (current->duplicate) {
current = current->duplicate;
}
current->duplicate = in_read->duplicate;
}
in_read->duplicate = bp->p->duplicate;
}
bp->p->duplicate = in_read;
}
dup = in_read->b;
}
if (mark_duplicates(param, dup_hash, bp->p->b, dup, &optical, &opt_warnings))
goto fail;
duplicate++;
} else {
fprintf(samtools_stderr, "[markdup] error: pair hashing failure.\n");
goto fail;
}
} else { // do the single (or effectively single) reads
int ret;
key_data_t single_key;
in_hash_t *bp;
make_single_key(param, &single_key, in_read->b, &bc_warnings);
single++;
in_read->pos = single_key.this_coord; // cigar/orientation modified pos
k = kh_put(reads, single_hash, single_key, &ret);
if (ret > 0) { // new
bp = &kh_val(single_hash, k);
bp->p = in_read;
in_read->single_key = single_key;
} else if (ret == 0) { // exists
bp = &kh_val(single_hash, k);
if ((bp->p->b->core.flag & BAM_FPAIRED) && !(bp->p->b->core.flag & BAM_FMUNMAP)) {
// if matched against one of a pair just mark as duplicate
if (param->check_chain) {
if (bp->p->duplicate) {
in_read->duplicate = bp->p->duplicate;
}
bp->p->duplicate = in_read;
}
if (mark_duplicates(param, dup_hash, bp->p->b, in_read->b, &single_optical, &opt_warnings))
goto fail;
} else {
int64_t old_score, new_score;
bam1_t *dup = NULL;
old_score = calc_score(bp->p->b);
new_score = calc_score(in_read->b);
// choose the highest score as the original, add it
// to the single hash and mark the other as duplicate
if (new_score > old_score) { // swap reads
dup = bp->p->b;
if (param->check_chain)
in_read->duplicate = bp->p;
bp->p = in_read;
} else {
if (param->check_chain) {
if (bp->p->duplicate) {
in_read->duplicate = bp->p->duplicate;
}
bp->p->duplicate = in_read;
}
dup = in_read->b;
}
if (mark_duplicates(param, dup_hash, bp->p->b, dup, &single_optical, &opt_warnings))
goto fail;
}
single_dup++;
} else {
fprintf(samtools_stderr, "[markdup] error: single hashing failure.\n");
goto fail;
}
}
} else {
excluded++;
}
// loop through the stored reads and write out those we
// no longer need
rq = kl_begin(read_buffer);
while (rq != kl_end(read_buffer)) {
in_read = &kl_val(rq);
/* keep a moving window of reads based on coordinates and max read length. Any unaligned reads
should just be written as they cannot be matched as duplicates. */
if (in_read->pos + param->max_length > prev_coord && in_read->b->core.tid == prev_tid && (prev_tid != -1 || prev_coord != -1)) {
break;
}
if (!dup_checked && param->check_chain) {
// check for multiple optical duplicates of the same original read
if (find_duplicate_chains(param, read_buffer, dup_hash, &dup_list, prev_coord, prev_tid, &opt_warnings, &single_optical, &optical, 1)) {
fprintf(samtools_stderr, "[markdup] error: duplicate checking failed.\n");
goto fail;
}
dup_checked = 1;
}
if (param->check_chain && (in_read->b->core.flag & BAM_FDUP) && !in_read->dup_checked && !(in_read->b->core.flag & exclude)) {
break;
}
if (!param->remove_dups || !(in_read->b->core.flag & BAM_FDUP)) {
if (param->supp) {
if (tmp_file_write(&temp, in_read->b)) {
fprintf(samtools_stderr, "[markdup] error: writing temp output failed.\n");
goto fail;
}
} else {
if (sam_write1(param->out, header, in_read->b) < 0) {
fprintf(samtools_stderr, "[markdup] error: writing output failed.\n");
goto fail;
}
}
writing++;
}
// remove from hash
if (in_read->pair_key.single == 0) {
k = kh_get(reads, pair_hash, in_read->pair_key);
kh_del(reads, pair_hash, k);
}
if (in_read->single_key.single == 1) {
k = kh_get(reads, single_hash, in_read->single_key);
kh_del(reads, single_hash, k);
}
kl_shift(read_queue, read_buffer, NULL);
bam_destroy1(in_read->b);
rq = kl_begin(read_buffer);
}
// set the next one up for reading
in_read = kl_pushp(read_queue, read_buffer);
if (!in_read) {
fprintf(samtools_stderr, "[markdup] out of memory\n");
goto fail;
}
if ((in_read->b = bam_init1()) == NULL) {
fprintf(samtools_stderr, "[markdup] error: unable to allocate memory for alignment.\n");
goto fail;
}
}
if (ret < -1) {
fprintf(samtools_stderr, "[markdup] error: truncated input file.\n");
goto fail;
}
// one last check
if (param->tag || param->opt_dist) {
if (find_duplicate_chains(param, read_buffer, dup_hash, &dup_list, prev_coord, prev_tid, &opt_warnings, &single_optical, &optical, 0)) {
fprintf(samtools_stderr, "[markdup] error: duplicate checking failed.\n");
goto fail;
}
}
// write out the end of the list
rq = kl_begin(read_buffer);
while (rq != kl_end(read_buffer)) {
in_read = &kl_val(rq);
if (bam_get_qname(in_read->b)) { // last entry will be blank
if (!param->remove_dups || !(in_read->b->core.flag & BAM_FDUP)) {
if (param->supp) {
if (tmp_file_write(&temp, in_read->b)) {
fprintf(samtools_stderr, "[markdup] error: writing temp output failed.\n");
goto fail;
}
} else {
if (sam_write1(param->out, header, in_read->b) < 0) {
fprintf(samtools_stderr, "[markdup] error: writing output failed.\n");
goto fail;
}
}
writing++;
}
}
kl_shift(read_queue, read_buffer, NULL);
bam_destroy1(in_read->b);
rq = kl_begin(read_buffer);
}
if (param->supp) {
bam1_t *b;
if (tmp_file_end_write(&temp)) {
fprintf(samtools_stderr, "[markdup] error: unable to end tmp writing.\n");
goto fail;
}
// read data from temp file and mark duplicate supplementary alignments
if (tmp_file_begin_read(&temp)) {
goto fail;
}
b = bam_init1();
while ((ret = tmp_file_read(&temp, b)) > 0) {
if ((b->core.flag & BAM_FSUPPLEMENTARY) || (b->core.flag & BAM_FUNMAP) || (b->core.flag & BAM_FSECONDARY)) {
k = kh_get(duplicates, dup_hash, bam_get_qname(b));
if (k != kh_end(dup_hash)) {
b->core.flag |= BAM_FDUP;
np_duplicate++;
if (param->tag && kh_val(dup_hash, k).name) {
if (bam_aux_update_str(b, "do", strlen(kh_val(dup_hash, k).name) + 1, (char*)kh_val(dup_hash, k).name)) {
fprintf(samtools_stderr, "[markdup] error: unable to append supplementary 'do' tag.\n");
goto fail;
}
}
if (param->opt_dist) {
if (kh_val(dup_hash, k).type) {
bam_aux_update_str(b, "dt", 3, "SQ");
np_opt_duplicate++;
} else {
bam_aux_update_str(b, "dt", 3, "LB");
}
}
}
}
if (!param->remove_dups || !(b->core.flag & BAM_FDUP)) {
if (sam_write1(param->out, header, b) < 0) {
fprintf(samtools_stderr, "[markdup] error: writing final output failed.\n");
goto fail;
}
}
}
if (ret == -1) {
fprintf(samtools_stderr, "[markdup] error: failed to read tmp file.\n");
goto fail;
}
for (k = kh_begin(dup_hash); k != kh_end(dup_hash); ++k) {
if (kh_exist(dup_hash, k)) {
free(kh_val(dup_hash, k).name);
free((char *)kh_key(dup_hash, k));
kh_key(dup_hash, k) = NULL;
}
}
tmp_file_destroy(&temp);
bam_destroy1(b);
}
if (opt_warnings) {
fprintf(samtools_stderr, "[markdup] warning: number of failed attempts to get coordinates from read names = %ld\n",
opt_warnings);
}
if (bc_warnings) {
fprintf(samtools_stderr, "[markdup] warning: number of failed attempts to get barcodes = %ld\n", bc_warnings);
}
if (param->do_stats) {
FILE *fp;
int file_open = 0;
unsigned long els;
if (param->stats_file) {
if (NULL == (fp = fopen(param->stats_file, "w"))) {
fprintf(samtools_stderr, "[markdup] warning: cannot write stats to %s.\n", param->stats_file);
fp = samtools_stderr;
} else {
file_open = 1;
}
} else {
fp = samtools_stderr;
}
els = estimate_library_size(pair, duplicate, optical);
fprintf(fp,
"COMMAND: %s\n"
"READ: %ld\n"
"WRITTEN: %ld\n"
"EXCLUDED: %ld\n"
"EXAMINED: %ld\n"
"PAIRED: %ld\n"
"SINGLE: %ld\n"
"DUPLICATE PAIR: %ld\n"
"DUPLICATE SINGLE: %ld\n"
"DUPLICATE PAIR OPTICAL: %ld\n"
"DUPLICATE SINGLE OPTICAL: %ld\n"
"DUPLICATE NON PRIMARY: %ld\n"
"DUPLICATE NON PRIMARY OPTICAL: %ld\n"
"DUPLICATE PRIMARY TOTAL: %ld\n"
"DUPLICATE TOTAL: %ld\n"
"ESTIMATED_LIBRARY_SIZE: %ld\n", param->arg_list, reading, writing, excluded, examined, pair, single,
duplicate, single_dup, optical, single_optical, np_duplicate, np_opt_duplicate,
single_dup + duplicate, single_dup + duplicate + np_duplicate, els);
if (file_open) {
fclose(fp);
}
}
if (param->write_index) {
if (sam_idx_save(param->out) < 0) {
print_error_errno("markdup", "writing index failed");
goto fail;
}
}
if (param->check_chain && (param->tag || param->opt_dist))
free(dup_list.c);
kh_destroy(reads, pair_hash);
kh_destroy(reads, single_hash);
kl_destroy(read_queue, read_buffer);
kh_destroy(duplicates, dup_hash);
sam_hdr_destroy(header);
return 0;
fail:
for (rq = kl_begin(read_buffer); rq != kl_end(read_buffer); rq = kl_next(rq))
bam_destroy1(kl_val(rq).b);
kl_destroy(read_queue, read_buffer);
for (k = kh_begin(dup_hash); k != kh_end(dup_hash); ++k) {
if (kh_exist(dup_hash, k)) {
free((char *)kh_key(dup_hash, k));
}
}
kh_destroy(duplicates, dup_hash);
if (param->check_chain && (param->tag || param->opt_dist))
free(dup_list.c);
kh_destroy(reads, pair_hash);
kh_destroy(reads, single_hash);
sam_hdr_destroy(header);
return 1;
}
static int markdup_usage(void) {
fprintf(samtools_stderr, "\n");
fprintf(samtools_stderr, "Usage: samtools markdup <input.bam> <output.bam>\n\n");
fprintf(samtools_stderr, "Option: \n");
fprintf(samtools_stderr, " -r Remove duplicate reads\n");
fprintf(samtools_stderr, " -l INT Max read length (default 300 bases)\n");
fprintf(samtools_stderr, " -S Mark supplementary alignments of duplicates as duplicates (slower).\n");
fprintf(samtools_stderr, " -s Report stats.\n");
fprintf(samtools_stderr, " -f NAME Write stats to named file. Implies -s.\n");
fprintf(samtools_stderr, " -T PREFIX Write temporary files to PREFIX.samtools.nnnn.nnnn.tmp.\n");
fprintf(samtools_stderr, " -d INT Optical distance (if set, marks with dt tag)\n");
fprintf(samtools_stderr, " -c Clear previous duplicate settings and tags.\n");
fprintf(samtools_stderr, " -m --mode TYPE Duplicate decision method for paired reads.\n"
" TYPE = t measure positions based on template start/end (default).\n"
" s measure positions based on sequence start.\n");
fprintf(samtools_stderr, " -u Output uncompressed data\n");
fprintf(samtools_stderr, " --include-fails Include quality check failed reads.\n");
fprintf(samtools_stderr, " --no-PG Do not add a PG line\n");
fprintf(samtools_stderr, " --no-multi-dup Reduced duplicates of duplicates checking.\n");
fprintf(samtools_stderr, " --read-coords STR Regex for coords from read name.\n");
fprintf(samtools_stderr, " --coords-order STR Order of regex elements. txy (default). With t being a part of\n"
" the read names that must be equal and x/y being coordinates.\n");
fprintf(samtools_stderr, " --barcode-tag STR Use barcode a tag that duplicates much match.\n");
fprintf(samtools_stderr, " --barcode-name Use the UMI/barcode in the read name (eigth colon delimited part).\n");
fprintf(samtools_stderr, " --barcode-rgx STR Regex for barcode in the readname (alternative to --barcode-name).\n");
fprintf(samtools_stderr, " -t Mark primary duplicates with the name of the original in a \'do\' tag."
" Mainly for information and debugging.\n");
sam_global_opt_help(samtools_stderr, "-.O..@..");
fprintf(samtools_stderr, "\nThe input file must be coordinate sorted and must have gone"
" through fixmates with the mate scoring option on.\n");
return 1;
}
int bam_markdup(int argc, char **argv) {
int c, ret, bc_name = 0;
char wmode[4] = {'w', 'b', 0, 0};
sam_global_args ga = SAM_GLOBAL_ARGS_INIT;
htsThreadPool p = {NULL, 0};
kstring_t tmpprefix = {0, 0, NULL};
struct stat st;
unsigned int t;
char *regex = NULL, *bc_regex = NULL;
char *regex_order = "txy";
md_param_t param = {NULL, NULL, NULL, 0, 300, 0, 0, 0, 0, 0, 0, 0, 0, 0,
1, NULL, NULL, NULL, NULL, 0, 0, 0, NULL, NULL};
static const struct option lopts[] = {
SAM_OPT_GLOBAL_OPTIONS('-', 0, 'O', 0, 0, '@'),
{"include-fails", no_argument, NULL, 1001},
{"no-PG", no_argument, NULL, 1002},
{"mode", required_argument, NULL, 'm'},
{"no-multi-dup", no_argument, NULL, 1003},
{"read-coords", required_argument, NULL, 1004},
{"coords-order", required_argument, NULL, 1005},
{"barcode-tag", required_argument, NULL, 1006},
{"barcode-name", no_argument, NULL, 1007},
{"barcode-rgx", required_argument, NULL, 1008},
{NULL, 0, NULL, 0}
};
while ((c = getopt_long(argc, argv, "rsl:StT:O:@:f:d:cm:u", lopts, NULL)) >= 0) {
switch (c) {
case 'r': param.remove_dups = 1; break;
case 'l': param.max_length = atoi(optarg); break;
case 's': param.do_stats = 1; break;
case 'T': kputs(optarg, &tmpprefix); break;
case 'S': param.supp = 1; break;
case 't': param.tag = 1; break;
case 'f': param.stats_file = optarg; param.do_stats = 1; break;
case 'd': param.opt_dist = atoi(optarg); break;
case 'c': param.clear = 1; break;
case 'm':
if (strcmp(optarg, "t") == 0) {
param.mode = MD_MODE_TEMPLATE;
} else if (strcmp(optarg, "s") == 0) {
param.mode = MD_MODE_SEQUENCE;
} else {
fprintf(samtools_stderr, "[markdup] error: unknown mode '%s'.\n", optarg);
return markdup_usage();
}
break;
case 'u': wmode[2] = '0'; break;
case 1001: param.include_fails = 1; break;
case 1002: param.no_pg = 1; break;
case 1003: param.check_chain = 0; break;
case 1004: regex = optarg; break;
case 1005: regex_order = optarg; break;
case 1006: param.barcode = optarg; break;
case 1007: bc_name = 1; break;
case 1008: bc_name = 1, bc_regex = optarg; break;
default: if (parse_sam_global_opt(c, optarg, lopts, &ga) == 0) break;
/* else fall-through */
case '?': return markdup_usage();
}
}
if (optind + 2 > argc)
return markdup_usage();
if (param.barcode && bc_name) {
fprintf(samtools_stderr, "[markdup] Error: cannot specify --barcode-tag and "
"--barcode-name (or --barcode-rgx) at same time.\n");
return 1;
}
if (param.opt_dist < 0) param.opt_dist = 0;
if (param.max_length < 0) param.max_length = 300;
if (regex) {
int result;
// set the order the elements of the regex are assigned to.
// x and y being coordinates, t being any other important part of the read
// e.g. tile and lane
// x and y order does not matter as long as it is consistent
if ((strncmp(regex_order, "txy", 3) == 0) || (strncmp(regex_order, "tyx", 3) == 0)) {
param.rgx_t = 1;
param.rgx_x = 2;
param.rgx_y = 3;
} else if ((strncmp(regex_order, "xyt", 3) == 0) || (strncmp(regex_order, "yxt", 3) == 0)) {
param.rgx_x = 1;
param.rgx_y = 2;
param.rgx_t = 3;
} else if ((strncmp(regex_order, "xty", 3) == 0) || (strncmp(regex_order, "ytx", 3) == 0)) {
param.rgx_x = 1;
param.rgx_t = 2;
param.rgx_y = 3;
} else if ((strncmp(regex_order, "xy", 2) == 0) || (strncmp(regex_order, "yx", 2) == 0)) {
param.rgx_x = 1;
param.rgx_y = 2;
param.rgx_t = 0;
} else {
fprintf(samtools_stderr, "[markdup] error: could not recognise regex coordinate order \"%s\".\n", regex_order);
return 1;
}
if ((param.rgx = malloc(sizeof(regex_t))) == NULL) {
fprintf(samtools_stderr, "[markdup] error: could not allocate memory for regex.\n");
return 1;
}
if ((result = regcomp(param.rgx, regex, REG_EXTENDED))) {
char err_msg[256];
regerror(result, param.rgx, err_msg, 256);
fprintf(samtools_stderr, "[markdup] error: regex error \"%s\"\n", err_msg);
free(param.rgx);
return 1;
}
}
if (bc_name) {
int result;
/* From Illumina UMI documentation: "The UMI sequence is located in the
eighth colon-delimited field of the read name (QNAME)". */
char *rgx = "[0-9A-Za-z]+:[0-9]+:[0-9]+:[0-9]+:[0-9]+:[0-9]+:[0-9]+:([!-?A-~]+)";
if ((param.bc_rgx = malloc(sizeof(regex_t))) == NULL) {
fprintf(samtools_stderr, "[markdup] error: could not allocate memory for barcode regex.\n");
return 1;
}
if (bc_regex) {
rgx = bc_regex;
}
if ((result = regcomp(param.bc_rgx, rgx, REG_EXTENDED))) {
char err_msg[256];
regerror(result, param.bc_rgx, err_msg, 256);
fprintf(samtools_stderr, "[markdup] error: barcode regex error \"%s\"\n", err_msg);
free(param.bc_rgx);
return 1;
}
}
param.in = sam_open_format(argv[optind], "r", &ga.in);
if (!param.in) {
print_error_errno("markdup", "failed to open \"%s\" for input", argv[optind]);
return 1;
}
sam_open_mode(wmode + 1, argv[optind + 1], NULL);
param.out = sam_open_format(argv[optind + 1], wmode, &ga.out);
if (!param.out) {
print_error_errno("markdup", "failed to open \"%s\" for output", argv[optind + 1]);
return 1;
}
if (ga.nthreads > 0) {
if (!(p.pool = hts_tpool_init(ga.nthreads))) {
fprintf(samtools_stderr, "[markdup] error creating thread pool\n");
return 1;
}
hts_set_opt(param.in, HTS_OPT_THREAD_POOL, &p);
hts_set_opt(param.out, HTS_OPT_THREAD_POOL, &p);
}
// actual stuff happens here
// we need temp files so fix up the name here
if (tmpprefix.l == 0) {
if (strcmp(argv[optind + 1], "-") != 0)
ksprintf(&tmpprefix, "%s.", argv[optind + 1]);
else
kputc('.', &tmpprefix);
}
if (stat(tmpprefix.s, &st) == 0 && S_ISDIR(st.st_mode)) {
if (tmpprefix.s[tmpprefix.l-1] != '/') kputc('/', &tmpprefix);
}
t = ((unsigned) time(NULL)) ^ ((unsigned) clock());
ksprintf(&tmpprefix, "samtools.%d.%u.tmp", (int) getpid(), t % 10000);
param.prefix = tmpprefix.s;
param.arg_list = stringify_argv(argc + 1, argv - 1);
param.write_index = ga.write_index;
param.out_fn = argv[optind + 1];
ret = bam_mark_duplicates(¶m);
sam_close(param.in);
if (sam_close(param.out) < 0) {
fprintf(samtools_stderr, "[markdup] error closing output file\n");
ret = 1;
}
if (p.pool) hts_tpool_destroy(p.pool);
if (param.rgx) {
regfree(param.rgx);
free(param.rgx);
}
if (param.bc_rgx) {
regfree(param.bc_rgx);
free(param.bc_rgx);
}
free(param.arg_list);
free(tmpprefix.s);
sam_global_args_free(&ga);
return ret;
}
|