forked from pkrusche/vt
-
Notifications
You must be signed in to change notification settings - Fork 0
/
bcf_single_genotyping_buffered_reader.cpp
261 lines (205 loc) · 6.56 KB
/
bcf_single_genotyping_buffered_reader.cpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
/* The MIT License
Copyright (c) 2015 Adrian Tan <[email protected]>
Permission is hereby granted, free of charge, to any person obtaining a copy
of this software and associated documentation files (the "Software"), to deal
in the Software without restriction, including without limitation the rights
to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
copies of the Software, and to permit persons to whom the Software is
furnished to do so, subject to the following conditions:
The above copyright notice and this permission notice shall be included in
all copies or substantial portions of the Software.
THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN
THE SOFTWARE.
*/
#include "bcf_single_genotyping_buffered_reader.h"
/**
* Constructor.
*/
BCFSingleGenotypingBufferedReader::BCFSingleGenotypingBufferedReader(std::string input_vcf_file, std::vector<GenomeInterval>& intervals, std::string output_vcf_file)
{
/////////////////////
//io initialization//
/////////////////////
odr = new BCFOrderedReader(input_vcf_file, intervals);
//////////////////////////
//options initialization//
//////////////////////////
output_annotations = false;
////////////////////////
//stats initialization//
////////////////////////
no_snps_genotyped = 0;
no_indels_genotyped = 0;
no_vntrs_genotyped = 0;
////////////////////////
//tools initialization//
////////////////////////
vm = new VariantManip();
// fai = fai_load(ref_fasta_file.c_str());
// if (fai==NULL)
// {
// fprintf(stderr, "[%s:%d %s] Cannot load genome index: %s\n", __FILE__, __LINE__, __FUNCTION__, ref_fasta_file.c_str());
// exit(1);
// }
}
/**
* Collects sufficient statistics from read for variants to be genotyped.
*
* The VCF records in the buffer must never occur before
*/
void BCFSingleGenotypingBufferedReader::process_read(bam_hdr_t *h, bam1_t *s)
{
//wrap bam1_t in AugmentBAMRecord
as.initialize(h, s);
uint32_t tid = bam_get_tid(s);
uint32_t beg1 = as.beg1;
uint32_t end1 = as.end1;
//collect statistics for variant records that are in the buffer and overlap with the read
GenotypingRecord* g;
for (std::list<GenotypingRecord*>::iterator i=buffer.begin(); i!=buffer.end(); ++i)
{
g = *i;
// std::cerr << g->pos1 << " " << g->beg1 << " " << g->end1 << " ";
//same chromosome
if (tid==g->rid)
{
if (end1 < g->beg1)
{
//can terminate
return;
}
else if (beg1 > g->end1)
{
//this should not occur if the buffer was flushed before invoking process read
continue;
}
//else if (beg1 <= g->beg1 && g->end1 <= end1)
else if (beg1 <= g->pos1 && g->pos1 <= end1)
{
// collect_sufficient_statistics(*i, as);
}
else
{
//other types of overlap, just ignore
}
// std::cerr << "\n";
}
//prior chromosome
else if (tid<g->rid)
{
//this should not occur if the buffer was flushed before invoking process read
return;
}
//latter chromosome
else if (tid>g->rid)
{
//in case if the buffer has a VCF record later in the list which matches it
continue;
}
}
//you will only reach here if a read occurs after or overlaps the last record in the buffer
//adding new VCF records and collecting statistics if necessary
bcf1_t *v = bcf_init();
while (odr->read(v))
{
int32_t vtype = vm->classify_variant(odr->hdr, v, variant);
g = create_genotyping_record(odr->hdr, v, 2, variant);
buffer.push_back(g);
if (tid==g->rid)
{
//if (end1>=g->beg1 && pos1<=g->end1)
if (beg1 <= g->pos1 && g->pos1 <= end1)
{
// collect_sufficient_statistics(g, as);
}
}
//VCF record occurs after the read
if (tid < g->rid || end1 < g->beg1)
{
return;
}
else
{
v = bcf_init();
}
}
//this means end of file
bcf_destroy(v);
}
/**
* Flush records.
*/
void BCFSingleGenotypingBufferedReader::flush(bam_hdr_t *h, bam1_t *s, bool flush_all)
{
if (flush_all)
{
//read all the remaining from the reference genotyping file
bcf1_t *v = bcf_init();
while (odr->read(v))
{
int32_t vtype = vm->classify_variant(odr->hdr, v, variant);
GenotypingRecord* g = create_genotyping_record(odr->hdr, v, 2, variant);
buffer.push_back(g);
v = bcf_init();
}
bcf_destroy(v);
GenotypingRecord* g;
while (!buffer.empty())
{
g = buffer.front();
// genotype_and_print(g);
delete g;
buffer.pop_front();
}
}
else
{
//std::cerr << "partial flush\n";
uint32_t tid = bam_get_tid(s);
GenotypingRecord* g;
while (!buffer.empty())
{
g = buffer.front();
if (tid==g->rid)
{
if (bam_get_pos1(s) > g->end1)
{
// genotype_and_print(g);
delete g;
buffer.pop_front();
}
else
{
return;
}
}
else if (tid>g->rid)
{
// genotype_and_print(g);
delete g;
buffer.pop_front();
}
else
{
return;
}
}
}
}
/**
* Create appropriate genotyping record.
*/
GenotypingRecord* BCFSingleGenotypingBufferedReader::create_genotyping_record(bcf_hdr_t* h, bcf1_t* v, uint32_t ploidy, Variant& variant)
{
GenotypingRecord* g = NULL;
if (variant.type==VT_SNP)
{
g = new SNPGenotypingRecord(h, v, 1, 2, NULL);
}
return g;
}