]> git.donarmstrong.com Git - rsem.git/blob - BamWriter.h
Allowed > 2^31 hits
[rsem.git] / BamWriter.h
1 #ifndef BAMWRITER_H_
2 #define BAMWRITER_H_
3
4 #include<cmath>
5 #include<cstdio>
6 #include<cstring>
7 #include<cassert>
8 #include<string>
9 #include<sstream>
10 #include<iostream>
11
12 #include <stdint.h>
13 #include "sam/bam.h"
14 #include "sam/sam.h"
15 #include "sam_rsem_aux.h"
16 #include "sam_rsem_cvt.h"
17
18 #include "utils.h"
19
20 #include "SingleHit.h"
21 #include "PairedEndHit.h"
22
23 #include "HitWrapper.h"
24 #include "Transcript.h"
25 #include "Transcripts.h"
26
27 class BamWriter {
28 public:
29         BamWriter(char, const char*, const char*, const char*, Transcripts&);
30         ~BamWriter();
31
32         void work(HitWrapper<SingleHit>);
33         void work(HitWrapper<PairedEndHit>);
34 private:
35         samfile_t *in, *out;
36         Transcripts& transcripts;
37
38         //convert bam1_t
39         void convert(bam1_t*, double);
40 };
41
42 //fn_list can be NULL
43 BamWriter::BamWriter(char inpType, const char* inpF, const char* fn_list, const char* outF, Transcripts& transcripts)
44         : transcripts(transcripts)
45 {
46         switch(inpType) {
47         case 's': in = samopen(inpF, "r", fn_list); break;
48         case 'b': in = samopen(inpF, "rb", fn_list); break;
49         default: assert(false);
50         }
51         assert(in != 0);
52
53         //build mappings from external sid to internal sid
54         transcripts.buildMappings(in->header->n_targets, in->header->target_name);
55
56         //generate output's header
57         bam_header_t *out_header = bam_header_dwt(in->header);
58         for (int i = 0; i < out_header->n_targets; i++) {
59                 out_header->target_len[i] = transcripts.getTranscriptViaEid(i + 1).getLength();  // transcript length without poly(A) tail
60         }
61
62         std::ostringstream strout;
63         strout<<"@HD\tVN:1.4\tSO:unknown\n@PG\tID:RSEM\n";
64         std::string content = strout.str();
65         append_header_text(out_header, content.c_str(), content.length());
66
67         out = samopen(outF, "wb", out_header);
68         assert(out != 0);
69
70         bam_header_destroy(out_header);
71 }
72
73 BamWriter::~BamWriter() {
74         samclose(in);
75         samclose(out);
76 }
77
78 void BamWriter::work(HitWrapper<SingleHit> wrapper) {
79         bam1_t *b;
80         SingleHit *hit;
81
82         HIT_INT_TYPE cnt = 0;
83
84         b = bam_init1();
85
86         while (samread(in, b) >= 0) {
87                 ++cnt;
88                 if (verbose && cnt % 1000000 == 0) { std::cout<< cnt<< "alignment lines are loaded!"<< std::endl; }
89
90                 if (b->core.flag & 0x0004) continue;
91
92                 hit = wrapper.getNextHit();
93                 assert(hit != NULL);
94
95                 assert(transcripts.getInternalSid(b->core.tid + 1) == hit->getSid());
96                 convert(b, hit->getConPrb());
97                 if (b->core.qual > 0) samwrite(out, b); // output only when MAPQ > 0
98         }
99
100         assert(wrapper.getNextHit() == NULL);
101
102         bam_destroy1(b);
103         if (verbose) { std::cout<< "Bam output file is generated!"<< std::endl; }
104 }
105
106 void BamWriter::work(HitWrapper<PairedEndHit> wrapper) {
107         bam1_t *b, *b2;
108         PairedEndHit *hit;
109
110         int cnt = 0;
111
112         b = bam_init1();
113         b2 = bam_init1();
114
115         while (samread(in, b) >= 0 && samread(in, b2) >= 0) {
116                 cnt += 2;
117                 if (verbose && cnt % 1000000 == 0) { std::cout<< cnt<< "alignment lines are loaded!"<< std::endl; }
118                 //mate info is not complete, skip
119                 if (!(((b->core.flag & 0x0040) && (b2->core.flag & 0x0080)) || ((b->core.flag & 0x0080) && (b2->core.flag & 0x0040)))) continue;
120                 //unalignable reads, skip
121                 if ((b->core.flag & 0x0004) || (b2->core.flag & 0x0004)) continue;
122
123                 //swap if b is mate 2
124                 if (b->core.flag & 0x0080) {
125                         assert(b2->core.flag & 0x0040);
126                         bam1_t *tmp = b;
127                         b = b2; b2 = tmp;
128                 }
129
130                 hit = wrapper.getNextHit();
131                 assert(hit != NULL);
132
133                 assert(transcripts.getInternalSid(b->core.tid + 1) == hit->getSid());
134                 assert(transcripts.getInternalSid(b2->core.tid + 1) == hit->getSid());
135
136                 convert(b, hit->getConPrb());
137                 convert(b2, hit->getConPrb());
138
139                 b->core.mpos = b2->core.pos;
140                 b2->core.mpos = b->core.pos;
141
142                 if (b->core.qual > 0) {
143                         samwrite(out, b);
144                         samwrite(out, b2);
145                 }
146         }
147
148         assert(wrapper.getNextHit() == NULL);
149
150         bam_destroy1(b);
151         bam_destroy1(b2);
152
153         if (verbose) { std::cout<< "Bam output file is generated!"<< std::endl; }
154 }
155
156 void BamWriter::convert(bam1_t *b, double prb) {
157         const Transcript& transcript = transcripts.getTranscriptViaEid(b->core.tid + 1);
158
159         int pos = b->core.pos;
160         int readlen = b->core.l_qseq;
161
162         std::vector<uint32_t> data;
163         data.clear();
164
165         int core_pos, core_n_cigar;
166         std::vector<Interval> vec;
167         vec.assign(1, Interval(1, transcript.getLength()));
168         // make an artificial chromosome coordinates for the transcript to get new CIGAR strings
169         tr2chr(Transcript("", "", "", '+', vec, ""), pos + 1, pos + readlen, core_pos, core_n_cigar, data);
170         assert(core_pos >= 0);
171
172         int rest_len = b->data_len - b->core.l_qname - b->core.n_cigar * 4;
173         b->data_len = b->core.l_qname + core_n_cigar * 4 + rest_len;
174         expand_data_size(b);
175         uint8_t* pt = b->data + b->core.l_qname;
176         memmove(pt + core_n_cigar * 4, pt + b->core.n_cigar * 4, rest_len);
177         for (int i = 0; i < core_n_cigar; i++) { memmove(pt, &data[i], 4); pt += 4; }
178
179         b->core.pos = core_pos;
180         b->core.n_cigar = core_n_cigar;
181         b->core.qual = getMAPQ(prb);
182         b->core.bin = bam_reg2bin(b->core.pos, bam_calend(&(b->core), bam1_cigar(b)));
183
184         float val = (float)prb;
185         bam_aux_append(b, "ZW", 'f', bam_aux_type2size('f'), (uint8_t*)&val);
186 }
187
188 #endif /* BAMWRITER_H_ */