7 static int g_min_mapQ = 0, g_flag_on = 0, g_flag_off = 0;
9 #define __g_skip_aln(b) (((b)->core.qual < g_min_mapQ) || ((b->core.flag & g_flag_on) != g_flag_on) \
10 || (b->core.flag & g_flag_off))
12 // callback function for bam_fetch()
13 static int view_func(const bam1_t *b, void *data)
15 if (!__g_skip_aln(b)) samwrite((samfile_t*)data, b);
19 static int usage(void);
21 int main_samview(int argc, char *argv[])
23 int c, is_header = 0, is_header_only = 0, is_bamin = 1, ret = 0;
24 samfile_t *in = 0, *out = 0;
25 char in_mode[4], out_mode[4], *fn_out = 0, *fn_list = 0;
27 /* parse command-line options */
28 strcpy(in_mode, "r"); strcpy(out_mode, "w");
29 while ((c = getopt(argc, argv, "Sbt:hHo:q:f:F:")) >= 0) {
31 case 'S': is_bamin = 0; break;
32 case 'b': strcat(out_mode, "b"); break;
33 case 't': fn_list = strdup(optarg); is_bamin = 0; break;
34 case 'h': is_header = 1; break;
35 case 'H': is_header_only = 1; break;
36 case 'o': fn_out = strdup(optarg); break;
37 case 'f': g_flag_on = strtol(optarg, 0, 0); break;
38 case 'F': g_flag_off = strtol(optarg, 0, 0); break;
39 case 'q': g_min_mapQ = atoi(optarg); break;
40 default: return usage();
43 if (is_header_only) is_header = 1;
44 if (is_bamin) strcat(in_mode, "b");
45 if (is_header) strcat(out_mode, "h");
46 if (argc == optind) return usage();
49 if ((in = samopen(argv[optind], in_mode, fn_list)) == 0) {
50 fprintf(stderr, "[main_samview] fail to open file for reading.\n");
53 if ((out = samopen(fn_out? fn_out : "-", out_mode, in->header)) == 0) {
54 fprintf(stderr, "[main_samview] fail to open file for writing.\n");
57 if (is_header_only) goto view_end; // no need to print alignments
59 if (argc == optind + 1) { // convert/print the entire file
60 bam1_t *b = bam_init1();
62 while ((r = samread(in, b)) >= 0) // read one alignment from `in'
64 samwrite(out, b); // write the alignment to `out'
65 if (r < -1) fprintf(stderr, "[main_samview] truncated file.\n");
67 } else { // retrieve alignments in specified regions
70 if (is_bamin) idx = bam_index_load(argv[optind]); // load BAM index
71 if (idx == 0) { // index is unavailable
72 fprintf(stderr, "[main_samview] random alignment retrieval only works for indexed BAM files.\n");
76 for (i = optind + 1; i < argc; ++i) {
78 bam_parse_region(in->header, argv[i], &tid, &beg, &end); // parse a region in the format like `chr2:100-200'
79 if (tid < 0) { // reference name is not found
80 fprintf(stderr, "[main_samview] fail to get the reference name. Continue anyway.\n");
83 bam_fetch(in->x.bam, idx, tid, beg, end, out, view_func); // fetch alignments
85 bam_index_destroy(idx); // destroy the BAM index
89 // close files, free and return
90 free(fn_list); free(fn_out);
98 fprintf(stderr, "\n");
99 fprintf(stderr, "Usage: samtools view [options] <in.bam>|<in.sam> [region1 [...]]\n\n");
100 fprintf(stderr, "Options: -b output BAM\n");
101 fprintf(stderr, " -h print header for the SAM output\n");
102 fprintf(stderr, " -H print header only (no alignments)\n");
103 fprintf(stderr, " -S input is SAM\n");
104 fprintf(stderr, " -t FILE list of reference names and lengths (force -S) [null]\n");
105 fprintf(stderr, " -o FILE output file name [stdout]\n");
106 fprintf(stderr, " -f INT required flag, 0 for unset [0]\n");
107 fprintf(stderr, " -F INT filtering flag, 0 for unset [0]\n");
108 fprintf(stderr, " -q INT minimum mapping quality [0]\n");
112 1. By default, this command assumes the file on the command line is in\n\
113 the BAM format and it prints the alignments in SAM. If `-t' is\n\
114 applied, the input file is assumed to be in the SAM format. The\n\
115 file supplied with `-t' is SPACE/TAB delimited with the first two\n\
116 fields of each line consisting of the reference name and the\n\
117 corresponding sequence length. The `.fai' file generated by `faidx'\n\
118 can be used here. This file may be empty if reads are unaligned.\n\
120 2. SAM->BAM conversion: `samtools view -bt ref.fa.fai in.sam.gz'.\n\
122 3. BAM->SAM conversion: `samtools view in.bam'.\n\
124 4. A region should be presented in one of the following formats:\n\
125 `chr1', `chr2:1,000' and `chr3:1000-2,000'. When a region is\n\
126 specified, the input alignment file must be an indexed BAM file.\n\
131 int main_import(int argc, char *argv[])
136 fprintf(stderr, "Usage: bamtk import <in.ref_list> <in.sam> <out.bam>\n");
140 argv2 = calloc(6, sizeof(char*));
141 argv2[0] = "import", argv2[1] = "-o", argv2[2] = argv[3], argv2[3] = "-bt", argv2[4] = argv[1], argv2[5] = argv[2];
142 ret = main_samview(argc2, argv2);