]> git.donarmstrong.com Git - rsem.git/blob - SamParser.h
Improved warning messages in rsem-gen-transcript-plots
[rsem.git] / SamParser.h
1 /* ReadType here means if the read is unalignable, alignable or aligned too much. It is NOT siheaderngle read or paired-end read */
2 #ifndef SAMPARSER_H_
3 #define SAMPARSER_H_
4
5 #include<cstdio>
6 #include<cstring>
7 #include<cstdlib>
8 #include<cassert>
9 #include<string>
10
11 #include "sam/bam.h"
12 #include "sam/sam.h"
13
14 #include "utils.h"
15 #include "my_assert.h"
16
17 #include "SingleRead.h"
18 #include "SingleReadQ.h"
19 #include "PairedEndRead.h"
20 #include "PairedEndReadQ.h"
21 #include "SingleHit.h"
22 #include "PairedEndHit.h"
23
24 #include "Transcripts.h"
25
26 class SamParser {
27 public:
28         SamParser(char, const char*, Transcripts&, const char* = 0);
29         ~SamParser();
30
31         /**
32          * return value
33          * -1 : no more alignment
34          * 0 : new read , type 0
35          * 1 : new read , type 1 with alignment
36          * 2 : new read , type 2
37          * 4 : discard this read
38          * 5 : new alignment but same read
39          */
40         int parseNext(SingleRead&, SingleHit&);
41         int parseNext(SingleReadQ&, SingleHit&);
42         int parseNext(PairedEndRead&, PairedEndHit&);
43         int parseNext(PairedEndReadQ&, PairedEndHit&);
44
45         static void setReadTypeTag(const char* tag) {
46                 strcpy(rtTag, tag);
47         }
48
49 private:
50         samfile_t *sam_in;
51         bam_header_t *header;
52         bam1_t *b, *b2;
53
54         Transcripts& transcripts;
55
56         //tag used by aligner
57         static char rtTag[STRLEN];
58
59         //1 +, -1 - here
60         int getDir(const bam1_t* b) {
61                 return ((b->core.flag & 0x0010) ? -1 : 1);
62         }
63
64         std::string getName(const bam1_t* b) {
65                 return std::string((char*)bam1_qname(b));
66         }
67
68         std::string getReadSeq(const bam1_t*);
69         std::string getQScore(const bam1_t*);
70
71         //0 ~ N0 1 ~ N1 2 ~ N2
72         int getReadType(const bam1_t*);
73         int getReadType(const bam1_t*, const bam1_t*); // for paired-end reads
74
75         bool check(bam1_t *b) {
76                 return (b->core.n_cigar == 1) && ((*bam1_cigar(b) & BAM_CIGAR_MASK) == BAM_CMATCH) && (b->core.l_qseq == (int32_t)(*bam1_cigar(b) >> BAM_CIGAR_SHIFT));
77         }
78 };
79
80 char SamParser::rtTag[STRLEN] = ""; // default : no tag, thus no Type 2 reads
81
82 // aux, if not 0, points to the file name of fn_list
83 SamParser::SamParser(char inpType, const char* inpF, Transcripts& transcripts, const char* aux)
84         : transcripts(transcripts)
85 {
86         switch(inpType) {
87         case 'b': sam_in = samopen(inpF, "rb", aux); break;
88         case 's': sam_in = samopen(inpF, "r", aux); break;
89         default: assert(false);
90         }
91
92         general_assert(sam_in != 0, "Cannot open " + cstrtos(inpF) + "! It may not exist.");
93         header = sam_in->header;
94         general_assert(header != 0, "Fail to parse sam header!");
95
96         transcripts.buildMappings(header->n_targets, header->target_name);
97
98         b = bam_init1();
99         b2 = bam_init1();
100 }
101
102 SamParser::~SamParser() {
103         samclose(sam_in);
104         bam_destroy1(b);
105         bam_destroy1(b2);
106 }
107
108 // If sam_read1 returns 0 , what does it mean?
109 //Assume b.core.tid is 0-based
110 int SamParser::parseNext(SingleRead& read, SingleHit& hit) {
111         int val; // return value
112         bool canR = (samread(sam_in, b) >= 0);
113         if (!canR) return -1;
114
115         general_assert(!(b->core.flag & 0x0001), "Find a paired end read in the file!");
116
117         int readType = getReadType(b);
118         std::string name = getName(b);
119
120         if (readType != 1 || (readType == 1 && read.getName().compare(name) != 0)) {
121                 val = readType;
122                 read = SingleRead(name, getReadSeq(b));
123         }
124         else val = 5;
125
126         if (readType == 1) {
127                 if (!check(b)) { fprintf(stderr, "RSEM does not support gapped alignments, sorry!\n"); exit(-1); }
128
129                 if (getDir(b) > 0) {
130                         hit = SingleHit(transcripts.getInternalSid(b->core.tid + 1), b->core.pos);
131                 }
132                 else {
133                         hit = SingleHit(-transcripts.getInternalSid(b->core.tid + 1), header->target_len[b->core.tid] - b->core.pos - b->core.l_qseq);
134                 }
135         }
136
137         return val;
138 }
139
140 int SamParser::parseNext(SingleReadQ& read, SingleHit& hit) {
141         int val;
142         bool canR = (samread(sam_in, b) >= 0);
143         if (!canR) return -1;
144
145         general_assert(!(b->core.flag & 0x0001), "Find a paired end read in the file!");
146
147         int readType = getReadType(b);
148         std::string name = getName(b);
149
150         if (readType != 1 || (readType == 1 && read.getName().compare(name) != 0)) {
151                 val = readType;
152                 read = SingleReadQ(name, getReadSeq(b), getQScore(b));
153         }
154         else val = 5;
155
156         if (readType == 1) {
157                 if (!check(b)) { fprintf(stderr, "RSEM does not support gapped alignments, sorry!\n"); exit(-1); }
158
159                 if (getDir(b) > 0) {
160                         hit = SingleHit(transcripts.getInternalSid(b->core.tid + 1), b->core.pos);
161                 }
162                 else {
163                         hit = SingleHit(-transcripts.getInternalSid(b->core.tid + 1), header->target_len[b->core.tid] - b->core.pos - b->core.l_qseq);
164                 }
165         }
166
167         return val;
168 }
169
170 //Assume whether aligned or not , two mates of paired-end reads are always get together
171 int SamParser::parseNext(PairedEndRead& read, PairedEndHit& hit) {
172         int val;
173         bool canR = ((samread(sam_in, b) >= 0) && (samread(sam_in, b2) >= 0));
174         if (!canR) return -1;
175
176         general_assert((b->core.flag & 0x0001) && (b2->core.flag & 0x0001), \
177                         "One of the mate is not paired-end! (RSEM assumes the two mates of a paired-end read should be adjacent)");
178
179         bam1_t *mp1 = NULL, *mp2 = NULL;
180
181         if ((b->core.flag & 0x0040) && (b2->core.flag & 0x0080)) {
182                 mp1 = b; mp2 = b2;
183         }
184         else if ((b->core.flag & 0x0080) && (b2->core.flag & 0x0040)) {
185                 mp1 = b2; mp2 = b;
186         }
187         else return 4; // If lose mate info, discard. is it necessary?
188
189         int readType = getReadType(mp1, mp2);
190         std::string name = getName(mp1);
191
192         if (readType != 1 || (readType == 1 && read.getName().compare(name) != 0)) {
193                 val = readType;
194                 SingleRead mate1(getName(mp1), getReadSeq(mp1));
195                 SingleRead mate2(getName(mp2), getReadSeq(mp2));
196                 read = PairedEndRead(mate1, mate2);
197         }
198         else val = 5;
199
200         if (readType == 1) {
201                 if (!check(mp1) || !check(mp2)) { fprintf(stderr, "RSEM does not support gapped alignments, sorry!\n"); exit(-1); }
202
203                 if (mp1->core.tid != mp2->core.tid) {
204                         fprintf(stderr, "The two reads do not come from the same pair!");
205                         exit(-1);
206                 }
207                 //assert(mp1->core.tid == mp2->core.tid);
208                 if (getDir(mp1) > 0) {
209                         hit = PairedEndHit(transcripts.getInternalSid(mp1->core.tid + 1), mp1->core.pos, mp2->core.pos + mp2->core.l_qseq - mp1->core.pos);
210                 }
211                 else {
212                         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);
213                 }
214         }
215
216         return val;
217 }
218
219 int SamParser::parseNext(PairedEndReadQ& read, PairedEndHit& hit) {
220         int val;
221         bool canR = ((samread(sam_in, b) >= 0) && (samread(sam_in, b2) >= 0));
222         if (!canR) return -1;
223
224         general_assert((b->core.flag & 0x0001) && (b2->core.flag & 0x0001), \
225                         "One of the mate is not paired-end! (RSEM assumes the two mates of a paired-end read should be adjacent)");
226
227         bam1_t *mp1 = NULL, *mp2 = NULL;
228
229         if ((b->core.flag & 0x0040) && (b2->core.flag & 0x0080)) {
230                 mp1 = b; mp2 = b2;
231         }
232         else if ((b->core.flag & 0x0080) && (b2->core.flag & 0x0040)) {
233                 mp1 = b2; mp2 = b;
234         }
235         else return 4;
236
237         int readType = getReadType(mp1, mp2);
238         std::string name = getName(mp1);
239
240         if (readType != 1 || (readType == 1 && read.getName().compare(name) != 0)) {
241                 val = readType;
242                 SingleReadQ mate1(getName(mp1), getReadSeq(mp1), getQScore(mp1));
243                 SingleReadQ mate2(getName(mp2), getReadSeq(mp2), getQScore(mp2));
244                 read = PairedEndReadQ(mate1, mate2);
245         }
246         else val = 5;
247
248         if (readType == 1) {
249                 if (!check(mp1) || !check(mp2)) { fprintf(stderr, "RSEM does not support gapped alignments, sorry!\n"); exit(-1); }
250
251                 if (mp1->core.tid != mp2->core.tid) {
252                         fprintf(stderr, "The two reads do not come from the same pair!");
253                         exit(-1);
254                 }
255                 //assert(mp1->core.tid == mp2->core.tid);
256                 if (getDir(mp1) > 0) {
257                         hit = PairedEndHit(transcripts.getInternalSid(mp1->core.tid + 1), mp1->core.pos, mp2->core.pos + mp2->core.l_qseq - mp1->core.pos);
258                 }
259                 else {
260                         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);
261                 }
262         }
263
264         return val;
265 }
266
267 inline std::string SamParser::getReadSeq(const bam1_t* b) {
268         uint8_t *p = bam1_seq(b);
269         std::string readseq = "";
270         char base = 0;
271
272         if (getDir(b) < 0) {
273                 for (int i = b->core.l_qseq - 1; i >= 0; i--) {
274                         switch(bam1_seqi(p, i)) {
275                         //case 0 : base = '='; break;
276                         case 1 : base = 'T'; break;
277                         case 2 : base = 'G'; break;
278                         case 4 : base = 'C'; break;
279                         case 8 : base = 'A'; break;
280                         case 15 : base = 'N'; break;
281                         default : assert(false);
282                         }
283                         readseq.append(1, base);
284                 }
285         }
286         else {
287                 for (int i = 0; i < b->core.l_qseq; i++) {
288                         switch(bam1_seqi(p, i)) {
289                         //case 0 : base = '='; break;
290                         case 1 : base = 'A'; break;
291                         case 2 : base = 'C'; break;
292                         case 4 : base = 'G'; break;
293                         case 8 : base = 'T'; break;
294                         case 15 : base = 'N'; break;
295                         default : assert(false);
296                         }
297                         readseq.append(1, base);
298                 }
299         }
300
301         return readseq;
302 }
303
304 inline std::string SamParser::getQScore(const bam1_t* b) {
305         uint8_t *p = bam1_qual(b);
306         std::string qscore = "";
307
308         if (getDir(b) > 0) {
309                 for (int i = 0; i < b->core.l_qseq; i++) {
310                         qscore.append(1, (char)(*p + 33));
311                         ++p;
312                 }
313         }
314         else {
315                 p = p + b->core.l_qseq - 1;
316                 for (int i = 0; i < b->core.l_qseq; i++) {
317                         qscore.append(1, (char)(*p + 33));
318                         --p;
319                 }
320         }
321
322         return qscore;
323 }
324
325 //0 ~ N0 , 1 ~ N1, 2 ~ N2
326 inline int SamParser::getReadType(const bam1_t* b) {
327         if (!(b->core.flag & 0x0004)) return 1;
328
329         if (!strcmp(rtTag, "")) return 0;
330
331         uint8_t *p = bam_aux_get(b, rtTag);
332         if (p == NULL) return 0;
333         return (bam_aux2i(p) > 0 ? 2 : 0);
334 }
335
336
337 //For paired-end reads, do not print out type 2 reads
338 inline int SamParser::getReadType(const bam1_t* b, const bam1_t* b2) {
339         if ((b->core.flag & 0x0002) && (b2->core.flag & 0x0002)) return 1;
340
341         return 0;
342 }
343
344 #endif /* SAMPARSER_H_ */