X-Git-Url: http://woldlab.caltech.edu/gitweb/?p=samtools.git;a=blobdiff_plain;f=sam_view.c;h=d0fdad25445508e90052970a90de92e1f0f61ae5;hp=113c6c437479237791cee16164d6fbc06cd7af8b;hb=85bb95099e58e20cc03456b7528248f7baed4db4;hpb=d363084f0412f3bcdeb0304aeb0974c9a10c7649 diff --git a/sam_view.c b/sam_view.c index 113c6c4..d0fdad2 100644 --- a/sam_view.c +++ b/sam_view.c @@ -2,26 +2,56 @@ #include #include #include +#include +#include "sam_header.h" #include "sam.h" #include "faidx.h" +#include "khash.h" +KHASH_SET_INIT_STR(rg) +typedef khash_t(rg) *rghash_t; + +rghash_t g_rghash = 0; static int g_min_mapQ = 0, g_flag_on = 0, g_flag_off = 0; static char *g_library, *g_rg; +static int g_sol2sanger_tbl[128]; + +static void sol2sanger(bam1_t *b) +{ + int l; + uint8_t *qual = bam1_qual(b); + if (g_sol2sanger_tbl[30] == 0) { + for (l = 0; l != 128; ++l) { + g_sol2sanger_tbl[l] = (int)(10.0 * log(1.0 + pow(10.0, (l - 64 + 33) / 10.0)) / log(10.0) + .499); + if (g_sol2sanger_tbl[l] >= 93) g_sol2sanger_tbl[l] = 93; + } + } + for (l = 0; l < b->core.l_qseq; ++l) { + int q = qual[l]; + if (q > 127) q = 127; + qual[l] = g_sol2sanger_tbl[q]; + } +} static inline int __g_skip_aln(const bam_header_t *h, const bam1_t *b) { if (b->core.qual < g_min_mapQ || ((b->core.flag & g_flag_on) != g_flag_on) || (b->core.flag & g_flag_off)) return 1; - if (g_library || g_rg) { + if (g_rg || g_rghash) { uint8_t *s = bam_aux_get(b, "RG"); if (s) { - if (g_rg && strcmp(g_rg, (char*)(s + 1)) == 0) return 0; - if (g_library) { - const char *p = bam_strmap_get(h->rg2lib, (char*)(s + 1)); - return (p && strcmp(p, g_library) == 0)? 0 : 1; - } return 1; - } else return 1; - } else return 0; + if (g_rg) return (strcmp(g_rg, (char*)(s + 1)) == 0)? 0 : 1; + if (g_rghash) { + khint_t k = kh_get(rg, g_rghash, (char*)(s + 1)); + return (k != kh_end(g_rghash))? 0 : 1; + } + } + } + if (g_library) { + const char *p = bam_get_library((bam_header_t*)h, b); + return (p && strcmp(p, g_library) == 0)? 0 : 1; + } + return 0; } // callback function for bam_fetch() @@ -36,15 +66,16 @@ static int usage(int is_long_help); int main_samview(int argc, char *argv[]) { - int c, is_header = 0, is_header_only = 0, is_bamin = 1, ret = 0, is_uncompressed = 0, is_bamout = 0; + int c, is_header = 0, is_header_only = 0, is_bamin = 1, ret = 0, is_uncompressed = 0, is_bamout = 0, slx2sngr = 0; int of_type = BAM_OFDEC, is_long_help = 0; samfile_t *in = 0, *out = 0; - char in_mode[5], out_mode[5], *fn_out = 0, *fn_list = 0, *fn_ref = 0; + char in_mode[5], out_mode[5], *fn_out = 0, *fn_list = 0, *fn_ref = 0, *fn_rg = 0; /* parse command-line options */ strcpy(in_mode, "r"); strcpy(out_mode, "w"); - while ((c = getopt(argc, argv, "Sbt:hHo:q:f:F:ul:r:xX?T:")) >= 0) { + while ((c = getopt(argc, argv, "Sbt:hHo:q:f:F:ul:r:xX?T:CR:")) >= 0) { switch (c) { + case 'C': slx2sngr = 1; break; case 'S': is_bamin = 0; break; case 'b': is_bamout = 1; break; case 't': fn_list = strdup(optarg); is_bamin = 0; break; @@ -57,6 +88,7 @@ int main_samview(int argc, char *argv[]) case 'u': is_uncompressed = 1; break; case 'l': g_library = strdup(optarg); break; case 'r': g_rg = strdup(optarg); break; + case 'R': fn_rg = strdup(optarg); break; case 'x': of_type = BAM_OFHEX; break; case 'X': of_type = BAM_OFSTR; break; case '?': is_long_help = 1; break; @@ -74,21 +106,36 @@ int main_samview(int argc, char *argv[]) if (is_bamin) strcat(in_mode, "b"); if (is_header) strcat(out_mode, "h"); if (is_uncompressed) strcat(out_mode, "u"); - if (argc == optind) return usage(is_long_help); + if (argc == optind) return usage(is_long_help); // potential memory leak... + + // read the list of read groups + if (fn_rg) { + FILE *fp_rg; + char buf[1024]; + int ret; + g_rghash = kh_init(rg); + fp_rg = fopen(fn_rg, "r"); + while (!feof(fp_rg) && fscanf(fp_rg, "%s", buf) > 0) // this is not a good style, but bear me... + kh_put(rg, g_rghash, strdup(buf), &ret); // we'd better check duplicates... + fclose(fp_rg); + } // generate the fn_list if necessary if (fn_list == 0 && fn_ref) fn_list = samfaipath(fn_ref); // open file handlers if ((in = samopen(argv[optind], in_mode, fn_list)) == 0) { - fprintf(stderr, "[main_samview] fail to open file for reading.\n"); + fprintf(stderr, "[main_samview] fail to open \"%s\" for reading.\n", argv[optind]); + ret = 1; goto view_end; } if (in->header == 0) { - fprintf(stderr, "[main_samview] fail to read the header.\n"); + fprintf(stderr, "[main_samview] fail to read the header from \"%s\".\n", argv[optind]); + ret = 1; goto view_end; } if ((out = samopen(fn_out? fn_out : "-", out_mode, in->header)) == 0) { - fprintf(stderr, "[main_samview] fail to open file for writing.\n"); + fprintf(stderr, "[main_samview] fail to open \"%s\" for writing.\n", fn_out? fn_out : "standard output"); + ret = 1; goto view_end; } if (is_header_only) goto view_end; // no need to print alignments @@ -96,10 +143,16 @@ int main_samview(int argc, char *argv[]) if (argc == optind + 1) { // convert/print the entire file bam1_t *b = bam_init1(); int r; - while ((r = samread(in, b)) >= 0) // read one alignment from `in' - if (!__g_skip_aln(in->header, b)) + while ((r = samread(in, b)) >= 0) { // read one alignment from `in' + if (!__g_skip_aln(in->header, b)) { + if (slx2sngr) sol2sanger(b); samwrite(out, b); // write the alignment to `out' - if (r < -1) fprintf(stderr, "[main_samview] truncated file.\n"); + } + } + if (r < -1) { + fprintf(stderr, "[main_samview] truncated file.\n"); + ret = 1; + } bam_destroy1(b); } else { // retrieve alignments in specified regions int i; @@ -114,17 +167,28 @@ int main_samview(int argc, char *argv[]) int tid, beg, end; bam_parse_region(in->header, argv[i], &tid, &beg, &end); // parse a region in the format like `chr2:100-200' if (tid < 0) { // reference name is not found - fprintf(stderr, "[main_samview] fail to get the reference name. Continue anyway.\n"); + fprintf(stderr, "[main_samview] region \"%s\" specifies an unknown reference name. Continue anyway.\n", argv[i]); continue; } - bam_fetch(in->x.bam, idx, tid, beg, end, out, view_func); // fetch alignments + // fetch alignments + if (bam_fetch(in->x.bam, idx, tid, beg, end, out, view_func) < 0) { + fprintf(stderr, "[main_samview] retrieval of region \"%s\" failed due to truncated file or corrupt BAM index file\n", argv[i]); + ret = 1; + break; + } } bam_index_destroy(idx); // destroy the BAM index } view_end: // close files, free and return - free(fn_list); free(fn_ref); free(fn_out); free(g_library); free(g_rg); + free(fn_list); free(fn_ref); free(fn_out); free(g_library); free(g_rg); free(fn_rg); + if (g_rghash) { + khint_t k; + for (k = 0; k < kh_end(g_rghash); ++k) + if (kh_exist(g_rghash, k)) free((char*)kh_key(g_rghash, k)); + kh_destroy(rg, g_rghash); + } samclose(in); samclose(out); return ret; @@ -144,6 +208,7 @@ static int usage(int is_long_help) fprintf(stderr, " -t FILE list of reference names and lengths (force -S) [null]\n"); fprintf(stderr, " -T FILE reference sequence file (force -S) [null]\n"); fprintf(stderr, " -o FILE output file name [stdout]\n"); + fprintf(stderr, " -R FILE list of read groups to be outputted [null]\n"); fprintf(stderr, " -f INT required flag, 0 for unset [0]\n"); fprintf(stderr, " -F INT filtering flag, 0 for unset [0]\n"); fprintf(stderr, " -q INT minimum mapping quality [0]\n");