15 #include "sam_rsem_aux.h"
16 #include "sam_rsem_cvt.h"
20 #include "SingleHit.h"
21 #include "PairedEndHit.h"
23 #include "HitWrapper.h"
24 #include "Transcript.h"
25 #include "Transcripts.h"
29 BamWriter(char, const char*, const char*, const char*, Transcripts&);
32 void work(HitWrapper<SingleHit>);
33 void work(HitWrapper<PairedEndHit>);
36 Transcripts& transcripts;
39 void convert(bam1_t*, double);
43 BamWriter::BamWriter(char inpType, const char* inpF, const char* fn_list, const char* outF, Transcripts& transcripts)
44 : transcripts(transcripts)
47 case 's': in = samopen(inpF, "r", fn_list); break;
48 case 'b': in = samopen(inpF, "rb", fn_list); break;
49 default: assert(false);
53 //build mappings from external sid to internal sid
54 transcripts.buildMappings(in->header->n_targets, in->header->target_name);
56 //generate output's header
57 bam_header_t *out_header = bam_header_dwt(in->header);
59 std::ostringstream strout;
60 strout<<"@HD\tVN:1.4\tSO:unknown\n@PG\tID:RSEM\n";
61 std::string content = strout.str();
62 append_header_text(out_header, content.c_str(), content.length());
64 out = samopen(outF, "wb", out_header);
67 bam_header_destroy(out_header);
70 BamWriter::~BamWriter() {
75 void BamWriter::work(HitWrapper<SingleHit> wrapper) {
83 while (samread(in, b) >= 0) {
85 if (verbose && cnt % 1000000 == 0) { std::cout<< cnt<< " alignment lines are loaded!"<< std::endl; }
87 if (!(b->core.flag & 0x0004)) {
88 hit = wrapper.getNextHit();
91 assert(transcripts.getInternalSid(b->core.tid + 1) == hit->getSid());
92 convert(b, hit->getConPrb());
95 //if (b->core.qual > 0) samwrite(out, b); // output only when MAPQ > 0
99 assert(wrapper.getNextHit() == NULL);
102 if (verbose) { std::cout<< "Bam output file is generated!"<< std::endl; }
105 void BamWriter::work(HitWrapper<PairedEndHit> wrapper) {
109 HIT_INT_TYPE cnt = 0;
114 while (samread(in, b) >= 0 && samread(in, b2) >= 0) {
116 if (verbose && cnt % 1000000 == 0) { std::cout<< cnt<< " alignment lines are loaded!"<< std::endl; }
118 assert((b->core.flag & 0x0001) && (b2->core.flag & 0x0001));
119 assert(((b->core.flag & 0x0040) && (b2->core.flag & 0x0080)) || ((b->core.flag & 0x0080) && (b2->core.flag & 0x0040)));
121 //unalignable reads, skip
122 bool notgood = (b->core.flag & 0x0004) || (b2->core.flag & 0x0004);
125 //swap if b is mate 2
126 if (b->core.flag & 0x0080) {
127 assert(b2->core.flag & 0x0040);
132 hit = wrapper.getNextHit();
135 assert(transcripts.getInternalSid(b->core.tid + 1) == hit->getSid());
136 assert(transcripts.getInternalSid(b2->core.tid + 1) == hit->getSid());
138 convert(b, hit->getConPrb());
139 convert(b2, hit->getConPrb());
141 // omit "b->core.mpos = b2->core.pos; b2->core.mpos = b->core.pos" for the reason that it is possible that only one mate is alignable
148 assert(wrapper.getNextHit() == NULL);
153 if (verbose) { std::cout<< "Bam output file is generated!"<< std::endl; }
156 void BamWriter::convert(bam1_t *b, double prb) {
157 b->core.qual = getMAPQ(prb);
158 float val = (float)prb;
159 bam_aux_append(b, "ZW", 'f', bam_aux_type2size('f'), (uint8_t*)&val);
162 #endif /* BAMWRITER_H_ */