]> git.donarmstrong.com Git - samtools.git/blobdiff - bam_rmdup.c
* samtools-0.1.5-7 (r403)
[samtools.git] / bam_rmdup.c
index 321938f3215506dbf54ebc76666c40165d2024a6..1fa6cad5d7582b9b7a1d148c47df4e00aa5b638e 100644 (file)
@@ -1,7 +1,6 @@
 #include <stdlib.h>
 #include <string.h>
 #include <stdio.h>
-#include <assert.h>
 #include <zlib.h>
 #include "bam.h"
 
@@ -76,6 +75,11 @@ void bam_rmdup_core(bamFile in, bamFile out)
                                        fprintf(stderr, "[bam_rmdup_core] %llu unmatched pairs\n", (long long)kh_size(del_set));
                                        clear_del_set(del_set);
                                }
+                               if ((int)c->tid == -1) { // append unmapped reads
+                                       bam_write1(out, b);
+                                       while (bam_read1(in, b) >= 0) bam_write1(out, b);
+                                       break;
+                               }
                                last_tid = c->tid;
                                fprintf(stderr, "[bam_rmdup_core] processing reference %s...\n", header->target_name[c->tid]);
                        }