]> git.donarmstrong.com Git - rsem.git/blobdiff - SamParser.h
Added error detection for cases such as a read's two mates having different names...
[rsem.git] / SamParser.h
index 335984e7394f2645ef4ee178032c57799c69d312..4747ef4281af593f714450d68ccfc65090662891 100644 (file)
@@ -12,9 +12,7 @@
 #include "sam/sam.h"
 
 #include "utils.h"
-
-#include "Transcript.h"
-#include "Transcripts.h"
+#include "my_assert.h"
 
 #include "SingleRead.h"
 #include "SingleReadQ.h"
 #include "SingleHit.h"
 #include "PairedEndHit.h"
 
+#include "Transcripts.h"
+
 class SamParser {
 public:
-       SamParser(char, const char*, Transcripts&, const char* = 0);
+       SamParser(char, const char*, const char*, Transcripts&, const char*);
        ~SamParser();
 
        /**
@@ -51,6 +51,7 @@ private:
        bam_header_t *header;
        bam1_t *b, *b2;
 
+       Transcripts& transcripts;
 
        //tag used by aligner
        static char rtTag[STRLEN];
@@ -61,7 +62,7 @@ private:
        }
 
        std::string getName(const bam1_t* b) {
-               return std::string((char*)bam1_qname(b));
+               return std::string(bam1_qname(b));
        }
 
        std::string getReadSeq(const bam1_t*);
@@ -79,33 +80,23 @@ private:
 char SamParser::rtTag[STRLEN] = ""; // default : no tag, thus no Type 2 reads
 
 // aux, if not 0, points to the file name of fn_list
-SamParser::SamParser(char inpType, const char* inpF, Transcripts& transcripts, const char* aux) {
+SamParser::SamParser(char inpType, const char* inpF, const char* aux, Transcripts& transcripts, const char* imdName)
+       : transcripts(transcripts)
+{
        switch(inpType) {
        case 'b': sam_in = samopen(inpF, "rb", aux); break;
        case 's': sam_in = samopen(inpF, "r", aux); break;
        default: assert(false);
        }
 
-       if (sam_in == 0) { fprintf(stderr, "Cannot open %s! It may not exist.\n", inpF); exit(-1); }
-    header = sam_in->header;
-    if (header == 0) { fprintf(stderr, "Fail to parse sam header!\n"); exit(-1); }
-
-    // Check if the reference used for aligner is the transcript set RSEM generated
-    if (transcripts.getM() != header->n_targets) {
-       fprintf(stderr, "Number of transcripts does not match! Please align reads against the transcript set and use RSEM generated reference for your aligner!\n");
-       exit(-1);
-    }
-    for (int i = 0; i < header->n_targets; i++) {
-       const Transcript& transcript = transcripts.getTranscriptAt(i + 1);
-       // If update int to long, chance the (int) conversion
-       if (transcript.getTranscriptID().compare(header->target_name[i]) != 0 || transcript.getLength() != (int)header->target_len[i]) {
-               fprintf(stderr, "Transcript information does not match! Please align reads against the transcript set and use RSEM generated reference for your aligner!\n");
-               exit(-1);
-       }
-    }
-
-    b = bam_init1();
-    b2 = bam_init1();
+       general_assert(sam_in != 0, "Cannot open " + cstrtos(inpF) + "! It may not exist.");
+       header = sam_in->header;
+       general_assert(header != 0, "Fail to parse sam header!");
+
+       transcripts.buildMappings(header->n_targets, header->target_name, imdName);
+
+       b = bam_init1();
+       b2 = bam_init1();
 }
 
 SamParser::~SamParser() {
@@ -121,8 +112,7 @@ int SamParser::parseNext(SingleRead& read, SingleHit& hit) {
        bool canR = (samread(sam_in, b) >= 0);
        if (!canR) return -1;
 
-       if (b->core.flag & 0x0001) { fprintf(stderr, "Find a paired end read in the file!\n"); exit(-1); }
-       //(b->core.flag & 0x0100) &&  && !(b->core.flag & 0x0004)
+       general_assert(!(b->core.flag & 0x0001), "Find a paired end read in the file!");
 
        int readType = getReadType(b);
        std::string name = getName(b);
@@ -137,10 +127,10 @@ int SamParser::parseNext(SingleRead& read, SingleHit& hit) {
                if (!check(b)) { fprintf(stderr, "RSEM does not support gapped alignments, sorry!\n"); exit(-1); }
 
                if (getDir(b) > 0) {
-                       hit = SingleHit(b->core.tid + 1, b->core.pos);
+                       hit = SingleHit(transcripts.getInternalSid(b->core.tid + 1), b->core.pos);
                }
                else {
-                       hit = SingleHit(-(b->core.tid + 1), header->target_len[b->core.tid] - b->core.pos - b->core.l_qseq);
+                       hit = SingleHit(-transcripts.getInternalSid(b->core.tid + 1), header->target_len[b->core.tid] - b->core.pos - b->core.l_qseq);
                }
        }
 
@@ -152,8 +142,7 @@ int SamParser::parseNext(SingleReadQ& read, SingleHit& hit) {
        bool canR = (samread(sam_in, b) >= 0);
        if (!canR) return -1;
 
-       if (b->core.flag & 0x0001) { fprintf(stderr, "Find a paired end read in the file!\n"); exit(-1); }
-       //assert(!(b->core.flag & 0x0001)); //(b->core.flag & 0x0100) &&  && !(b->core.flag & 0x0004)
+       general_assert(!(b->core.flag & 0x0001), "Find a paired end read in the file!");
 
        int readType = getReadType(b);
        std::string name = getName(b);
@@ -168,10 +157,10 @@ int SamParser::parseNext(SingleReadQ& read, SingleHit& hit) {
                if (!check(b)) { fprintf(stderr, "RSEM does not support gapped alignments, sorry!\n"); exit(-1); }
 
                if (getDir(b) > 0) {
-                       hit = SingleHit(b->core.tid + 1, b->core.pos);
+                       hit = SingleHit(transcripts.getInternalSid(b->core.tid + 1), b->core.pos);
                }
                else {
-                       hit = SingleHit(-(b->core.tid + 1), header->target_len[b->core.tid] - b->core.pos - b->core.l_qseq);
+                       hit = SingleHit(-transcripts.getInternalSid(b->core.tid + 1), header->target_len[b->core.tid] - b->core.pos - b->core.l_qseq);
                }
        }
 
@@ -184,21 +173,24 @@ int SamParser::parseNext(PairedEndRead& read, PairedEndHit& hit) {
        bool canR = ((samread(sam_in, b) >= 0) && (samread(sam_in, b2) >= 0));
        if (!canR) return -1;
 
-       if (!((b->core.flag & 0x0001) && (b2->core.flag & 0x0001))) {
-               fprintf(stderr, "One of the mate is not paired-end! (RSEM assumes the two mates of a paired-end read should be adjacent)\n");
-               exit(-1);
-       }
-       //assert((b->core.flag & 0x0001) && (b2->core.flag & 0x0001));
+       general_assert((b->core.flag & 0x0001) && (b2->core.flag & 0x0001), \
+                       "One of the mate is not paired-end! (RSEM assumes the two mates of a paired-end read should be adjacent)");
 
        bam1_t *mp1 = NULL, *mp2 = NULL;
 
-       if ((b->core.flag & 0x0040) && (b2->core.flag & 0x0080)) {
+       // If lose mate info, discard. is it necessary?
+       if (!((b->core.flag & 0x0040) && (b2->core.flag & 0x0080)) && !((b->core.flag & 0x0080) && (b2->core.flag & 0x0040))) return 4;
+       // If only one mate is mapped, discard
+       if (((b->core.flag & 0x0004) && !(b2->core.flag & 0x0004)) || (!(b->core.flag & 0x0004) && (b2->core.flag & 0x0004))) return 4;
+
+       if (b->core.flag & 0x0040) {
                mp1 = b; mp2 = b2;
        }
-       else if ((b->core.flag & 0x0080) && (b2->core.flag & 0x0040)) {
+       else  {
                mp1 = b2; mp2 = b;
        }
-       else return 4; // If lose mate info, discard. is it necessary?
+
+       general_assert(!strcmp(bam1_qname(mp1), bam1_qname(mp2)), "Detected a read pair whose two mates have different names: " + getName(mp1) + " , " + getName(mp2) + " !");
 
        int readType = getReadType(mp1, mp2);
        std::string name = getName(mp1);
@@ -220,10 +212,10 @@ int SamParser::parseNext(PairedEndRead& read, PairedEndHit& hit) {
                }
                //assert(mp1->core.tid == mp2->core.tid);
                if (getDir(mp1) > 0) {
-                       hit = PairedEndHit(mp1->core.tid + 1, mp1->core.pos, mp2->core.pos + mp2->core.l_qseq - mp1->core.pos);
+                       hit = PairedEndHit(transcripts.getInternalSid(mp1->core.tid + 1), mp1->core.pos, mp2->core.pos + mp2->core.l_qseq - mp1->core.pos);
                }
                else {
-                       hit = PairedEndHit(-(mp1->core.tid + 1), header->target_len[mp1->core.tid] - mp1->core.pos - mp1->core.l_qseq, mp1->core.pos + mp1->core.l_qseq - mp2->core.pos);
+                       hit = PairedEndHit(-transcripts.getInternalSid(mp1->core.tid + 1), header->target_len[mp1->core.tid] - mp1->core.pos - mp1->core.l_qseq, mp1->core.pos + mp1->core.l_qseq - mp2->core.pos);
                }
        }
 
@@ -235,21 +227,24 @@ int SamParser::parseNext(PairedEndReadQ& read, PairedEndHit& hit) {
        bool canR = ((samread(sam_in, b) >= 0) && (samread(sam_in, b2) >= 0));
        if (!canR) return -1;
 
-       if (!((b->core.flag & 0x0001) && (b2->core.flag & 0x0001))) {
-               fprintf(stderr, "One of the mate is not paired-end! (RSEM assumes the two mates of a paired-end read should be adjacent)\n");
-               exit(-1);
-       }
-       //assert((b->core.flag & 0x0001) && (b2->core.flag & 0x0001));
+       general_assert((b->core.flag & 0x0001) && (b2->core.flag & 0x0001), \
+                       "One of the mate is not paired-end! (RSEM assumes the two mates of a paired-end read should be adjacent)");
 
        bam1_t *mp1 = NULL, *mp2 = NULL;
 
-       if ((b->core.flag & 0x0040) && (b2->core.flag & 0x0080)) {
+       // If lose mate info, discard. is it necessary?
+       if (!((b->core.flag & 0x0040) && (b2->core.flag & 0x0080)) && !((b->core.flag & 0x0080) && (b2->core.flag & 0x0040))) return 4;
+       // If only one mate is mapped, discard
+       if (((b->core.flag & 0x0004) && !(b2->core.flag & 0x0004)) || (!(b->core.flag & 0x0004) && (b2->core.flag & 0x0004))) return 4;
+
+       if (b->core.flag & 0x0040) {
                mp1 = b; mp2 = b2;
        }
-       else if ((b->core.flag & 0x0080) && (b2->core.flag & 0x0040)) {
+       else  {
                mp1 = b2; mp2 = b;
        }
-       else return 4;
+
+       general_assert(!strcmp(bam1_qname(mp1), bam1_qname(mp2)), "Detected a read pair whose two mates have different names: " + getName(mp1) + " , " + getName(mp2) + " !");
 
        int readType = getReadType(mp1, mp2);
        std::string name = getName(mp1);
@@ -271,10 +266,10 @@ int SamParser::parseNext(PairedEndReadQ& read, PairedEndHit& hit) {
                }
                //assert(mp1->core.tid == mp2->core.tid);
                if (getDir(mp1) > 0) {
-                       hit = PairedEndHit(mp1->core.tid + 1, mp1->core.pos, mp2->core.pos + mp2->core.l_qseq - mp1->core.pos);
+                       hit = PairedEndHit(transcripts.getInternalSid(mp1->core.tid + 1), mp1->core.pos, mp2->core.pos + mp2->core.l_qseq - mp1->core.pos);
                }
                else {
-                       hit = PairedEndHit(-(mp1->core.tid + 1), header->target_len[mp1->core.tid] - mp1->core.pos - mp1->core.l_qseq, mp1->core.pos + mp1->core.l_qseq - mp2->core.pos);
+                       hit = PairedEndHit(-transcripts.getInternalSid(mp1->core.tid + 1), header->target_len[mp1->core.tid] - mp1->core.pos - mp1->core.l_qseq, mp1->core.pos + mp1->core.l_qseq - mp2->core.pos);
                }
        }
 
@@ -350,10 +345,17 @@ inline int SamParser::getReadType(const bam1_t* b) {
        return (bam_aux2i(p) > 0 ? 2 : 0);
 }
 
-
 //For paired-end reads, do not print out type 2 reads
 inline int SamParser::getReadType(const bam1_t* b, const bam1_t* b2) {
-       if ((b->core.flag & 0x0002) && (b2->core.flag & 0x0002)) return 1;
+       if (!(b->core.flag & 0x0004) && !(b2->core.flag & 0x0004)) return 1;
+
+       if (!strcmp(rtTag, "")) return 0;
+
+       uint8_t *p = bam_aux_get(b, rtTag);
+       if (p != NULL && bam_aux2i(p) > 0) return 2;
+
+       p = bam_aux_get(b2, rtTag);
+       if (p != NULL && bam_aux2i(p) > 0) return 2;
 
        return 0;
 }