]> git.donarmstrong.com Git - rsem.git/blob - SamParser.h
Modified the acknowledgement section of README.md
[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 lose mate info, discard. is it necessary?
182         if (!((b->core.flag & 0x0040) && (b2->core.flag & 0x0080)) && !((b->core.flag & 0x0080) && (b2->core.flag & 0x0040))) return 4;
183         // If only one mate is mapped, discard
184         if (((b->core.flag & 0x0004) && !(b2->core.flag & 0x0004)) || (!(b->core.flag & 0x0004) && (b2->core.flag & 0x0004))) return 4;
185
186         if (b->core.flag & 0x0040) {
187                 mp1 = b; mp2 = b2;
188         }
189         else  {
190                 mp1 = b2; mp2 = b;
191         }
192
193         int readType = getReadType(mp1, mp2);
194         std::string name = getName(mp1);
195
196         if (readType != 1 || (readType == 1 && read.getName().compare(name) != 0)) {
197                 val = readType;
198                 SingleRead mate1(getName(mp1), getReadSeq(mp1));
199                 SingleRead mate2(getName(mp2), getReadSeq(mp2));
200                 read = PairedEndRead(mate1, mate2);
201         }
202         else val = 5;
203
204         if (readType == 1) {
205                 if (!check(mp1) || !check(mp2)) { fprintf(stderr, "RSEM does not support gapped alignments, sorry!\n"); exit(-1); }
206
207                 if (mp1->core.tid != mp2->core.tid) {
208                         fprintf(stderr, "The two reads do not come from the same pair!");
209                         exit(-1);
210                 }
211                 //assert(mp1->core.tid == mp2->core.tid);
212                 if (getDir(mp1) > 0) {
213                         hit = PairedEndHit(transcripts.getInternalSid(mp1->core.tid + 1), mp1->core.pos, mp2->core.pos + mp2->core.l_qseq - mp1->core.pos);
214                 }
215                 else {
216                         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);
217                 }
218         }
219
220         return val;
221 }
222
223 int SamParser::parseNext(PairedEndReadQ& read, PairedEndHit& hit) {
224         int val;
225         bool canR = ((samread(sam_in, b) >= 0) && (samread(sam_in, b2) >= 0));
226         if (!canR) return -1;
227
228         general_assert((b->core.flag & 0x0001) && (b2->core.flag & 0x0001), \
229                         "One of the mate is not paired-end! (RSEM assumes the two mates of a paired-end read should be adjacent)");
230
231         bam1_t *mp1 = NULL, *mp2 = NULL;
232
233         // If lose mate info, discard. is it necessary?
234         if (!((b->core.flag & 0x0040) && (b2->core.flag & 0x0080)) && !((b->core.flag & 0x0080) && (b2->core.flag & 0x0040))) return 4;
235         // If only one mate is mapped, discard
236         if (((b->core.flag & 0x0004) && !(b2->core.flag & 0x0004)) || (!(b->core.flag & 0x0004) && (b2->core.flag & 0x0004))) return 4;
237
238         if (b->core.flag & 0x0040) {
239                 mp1 = b; mp2 = b2;
240         }
241         else  {
242                 mp1 = b2; mp2 = b;
243         }
244
245         int readType = getReadType(mp1, mp2);
246         std::string name = getName(mp1);
247
248         if (readType != 1 || (readType == 1 && read.getName().compare(name) != 0)) {
249                 val = readType;
250                 SingleReadQ mate1(getName(mp1), getReadSeq(mp1), getQScore(mp1));
251                 SingleReadQ mate2(getName(mp2), getReadSeq(mp2), getQScore(mp2));
252                 read = PairedEndReadQ(mate1, mate2);
253         }
254         else val = 5;
255
256         if (readType == 1) {
257                 if (!check(mp1) || !check(mp2)) { fprintf(stderr, "RSEM does not support gapped alignments, sorry!\n"); exit(-1); }
258
259                 if (mp1->core.tid != mp2->core.tid) {
260                         fprintf(stderr, "The two reads do not come from the same pair!");
261                         exit(-1);
262                 }
263                 //assert(mp1->core.tid == mp2->core.tid);
264                 if (getDir(mp1) > 0) {
265                         hit = PairedEndHit(transcripts.getInternalSid(mp1->core.tid + 1), mp1->core.pos, mp2->core.pos + mp2->core.l_qseq - mp1->core.pos);
266                 }
267                 else {
268                         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);
269                 }
270         }
271
272         return val;
273 }
274
275 inline std::string SamParser::getReadSeq(const bam1_t* b) {
276         uint8_t *p = bam1_seq(b);
277         std::string readseq = "";
278         char base = 0;
279
280         if (getDir(b) < 0) {
281                 for (int i = b->core.l_qseq - 1; i >= 0; i--) {
282                         switch(bam1_seqi(p, i)) {
283                         //case 0 : base = '='; break;
284                         case 1 : base = 'T'; break;
285                         case 2 : base = 'G'; break;
286                         case 4 : base = 'C'; break;
287                         case 8 : base = 'A'; break;
288                         case 15 : base = 'N'; break;
289                         default : assert(false);
290                         }
291                         readseq.append(1, base);
292                 }
293         }
294         else {
295                 for (int i = 0; i < b->core.l_qseq; i++) {
296                         switch(bam1_seqi(p, i)) {
297                         //case 0 : base = '='; break;
298                         case 1 : base = 'A'; break;
299                         case 2 : base = 'C'; break;
300                         case 4 : base = 'G'; break;
301                         case 8 : base = 'T'; break;
302                         case 15 : base = 'N'; break;
303                         default : assert(false);
304                         }
305                         readseq.append(1, base);
306                 }
307         }
308
309         return readseq;
310 }
311
312 inline std::string SamParser::getQScore(const bam1_t* b) {
313         uint8_t *p = bam1_qual(b);
314         std::string qscore = "";
315
316         if (getDir(b) > 0) {
317                 for (int i = 0; i < b->core.l_qseq; i++) {
318                         qscore.append(1, (char)(*p + 33));
319                         ++p;
320                 }
321         }
322         else {
323                 p = p + b->core.l_qseq - 1;
324                 for (int i = 0; i < b->core.l_qseq; i++) {
325                         qscore.append(1, (char)(*p + 33));
326                         --p;
327                 }
328         }
329
330         return qscore;
331 }
332
333 //0 ~ N0 , 1 ~ N1, 2 ~ N2
334 inline int SamParser::getReadType(const bam1_t* b) {
335         if (!(b->core.flag & 0x0004)) return 1;
336
337         if (!strcmp(rtTag, "")) return 0;
338
339         uint8_t *p = bam_aux_get(b, rtTag);
340         if (p == NULL) return 0;
341         return (bam_aux2i(p) > 0 ? 2 : 0);
342 }
343
344 //For paired-end reads, do not print out type 2 reads
345 inline int SamParser::getReadType(const bam1_t* b, const bam1_t* b2) {
346         if (!(b->core.flag & 0x0004) && !(b2->core.flag & 0x0004)) return 1;
347
348         if (!strcmp(rtTag, "")) return 0;
349
350         uint8_t *p = bam_aux_get(b, rtTag);
351         if (p != NULL && bam_aux2i(p) > 0) return 2;
352
353         p = bam_aux_get(b2, rtTag);
354         if (p != NULL && bam_aux2i(p) > 0) return 2;
355
356         return 0;
357 }
358
359 #endif /* SAMPARSER_H_ */