-
Notifications
You must be signed in to change notification settings - Fork 37
/
hts_utils.h
347 lines (280 loc) · 7.78 KB
/
hts_utils.h
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
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
/* The MIT License
Copyright (c) 2013 Adrian Tan <atks@umich.edu>
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.
*/
#ifndef HTS_UTILS_H
#define HTS_UTILS_H
#include <string>
#include <vector>
#include <cstdlib>
#include <cstdint>
#include <cstring>
#include <cmath>
#include <cfloat>
#include <vector>
#include <map>
#include <queue>
#include "htslib/kstring.h"
#include "htslib/khash.h"
#include "htslib/hts.h"
#include "htslib/sam.h"
#include "htslib/vcf.h"
#include "htslib/vcfutils.h"
#include "utils.h"
/**************
*BAM HDR UTILS
**************/
/**
* Copies contigs found in bam header to bcf header.
*/
void bam_hdr_transfer_contigs_to_bcf_hdr(const bam_hdr_t *sh, bcf_hdr_t *vh);
/**
* Get number of sequences.
*/
#define bam_hdr_get_n_targets(h) ((h)->n_targets)
/**
* Get list of sequence names.
*/
#define bam_hdr_get_target_name(h) ((h)->target_name)
/**
* Get list of sequence lenghts.
*/
#define bam_hdr_get_target_len(h) ((h)->target_len)
/**********
*BAM UTILS
**********/
//used when a base on a read is not aligned to the human genome reference
//in particular for soft clips and insertions
#define BAM_READ_INDEX_NA -1
/**
* Gets the start position of the first mapped base in the sequence.
*/
#define bam_get_chrom(h, s) ((h)->target_name[(s)->core.tid])
/**
* Gets the 1 based start position of the first mapped base in the read.
*/
#define bam_get_pos1(s) ((s)->core.pos + 1)
/**
* Gets the 0 based start position of the first mapped base in the read.
*/
#define bam_get_pos0(s) ((s)->core.pos)
/**
* Gets the end position of the last mapped base in the read.
*/
int32_t bam_get_end_pos1(bam1_t *srec);
/**
* Gets the template ID of the paired read.
*/
#define bam_get_mtid(s) ((s)->core.mtid)
/**
* Gets the start position of the first mapped base in the read.
*/
#define bam_get_mpos1(s) ((s)->core.mpos)
/**
* Gets the read sequence from a bam record
*/
void bam_get_seq_string(bam1_t *s, kstring_t *seq);
/**
* Gets the base qualities from a bam record
*/
void bam_get_qual_string(bam1_t *s, kstring_t *qual);
/**
* Gets the cigar sequence from a bam record
*/
#define bam_get_n_cigar_op(b) ((b)->core.n_cigar)
/**
* Gets the cigar from a BAM record
*/
void bam_get_cigar_string(bam1_t *s, kstring_t *cigar);
/**
* Gets the cigar string from a bam record
*/
void bam_get_cigar_expanded_string(bam1_t *s, kstring_t *cigar_string);
/**
* Is this sequence the first read?
*/
#define bam_is_fread1(s) (((s)->core.flag&BAM_FREAD1) != 0)
/**
* Is this sequence the second read?
*/
#define bam_is_fread2(s) (((s)->core.flag&BAM_FREAD2) != 0)
/**
* Gets the base in the read that is mapped to a genomic position.
* Returns -1 if it does not exists.
*/
void bam_get_base_and_qual(bam1_t *srec, uint32_t pos, char& base, char& qual, int32_t& rpos);
/**
* Gets the base in the read that is mapped to a genomic position.
* Returns true if successful, false if the location is a deletion or not on the read.
*/
void bam_get_base_and_qual_and_read_and_qual(bam1_t *s, uint32_t pos, char& base, char& qual, int32_t& rpos, kstring_t* readseq, kstring_t* readqual);
/**
* Gets flag
*/
#define bam_get_flag(s) ((s)->core.flag)
/**
* Get map quality
*/
#define bam_get_mapq(s) ((s)->core.qual)
/**
*Get tid - e.g. chromosome id
*/
#define bam_get_tid(s) ((s)->core.tid)
/**
* Get read length
*/
#define bam_get_l_qseq(s) ((s)->core.l_qseq)
/**************
*BCF HDR UTILS
**************/
/**
* Copies contigs found in bcf header to another bcf header.
*/
void bcf_hdr_transfer_contigs(const bcf_hdr_t *sh, bcf_hdr_t *vh);
/**
* Get samples from bcf header
*/
#define bcf_hdr_get_samples(h) ((h)->samples)
/**
* Get ith sample name from bcf header
*/
#define bcf_hdr_get_sample_name(h, i) ((h)->samples[i])
/**
* Get number of samples in bcf header
*/
int32_t bcf_hdr_get_n_sample(bcf_hdr_t *h);
/**
* Reads header of a VCF file and returns the bcf header object.
* This wraps around vcf_hdr_read from the original htslib to
* allow for an alternative header file to be read in.
*
* this searches for the alternative header saved as <filename>.hdr
*/
bcf_hdr_t *bcf_alt_hdr_read(htsFile *fp);
/**
* Subsets a record by samples.
* @imap - indices the subsetted samples
*/
int bcf_hdr_subset_samples(const bcf_hdr_t *h, bcf1_t *v, std::vector<int32_t>& imap);
/**********
*BCF UTILS
**********/
/**
* Gets a string representation of a variant.
*/
void bcf_variant2string(bcf_hdr_t *h, bcf1_t *v, kstring_t *var);
/**
* Gets a sorted string representation of a variant.
*/
void bcf_variant2string_sorted(bcf_hdr_t *h, bcf1_t *v, kstring_t *var);
/**
* Gets a sorted string representation of the alleles of a variant.
*/
void bcf_alleles2string_sorted(bcf_hdr_t *h, bcf1_t *v, kstring_t *var);
/**
* Prints a VCF record to STDERR.
*/
void bcf_print(bcf_hdr_t *h, bcf1_t *v);
/**
* Prints a VCF record in compact string representation to STDERR.
*/
void bcf_print_lite(bcf_hdr_t *h, bcf1_t *v);
/**
* Prints a VCF record in compact string representation to STDERR with alleles sorted.
*/
void bcf_print_lite_sorted(bcf_hdr_t *h, bcf1_t *v);
/**
* Prints a VCF record in compact string representation to STDERR with a new line.
*/
void bcf_print_liten(bcf_hdr_t *h, bcf1_t *v);
/**
* Get number of chromosomes
*/
#define bcf_get_n_chrom(h) ((h)->n[BCF_DT_CTG])
/**
* Get chromosome name
*/
#define bcf_get_chrom(h, v) ((h)->id[BCF_DT_CTG][(v)->rid].key)
/**
* Set chromosome name
*/
void bcf_set_chrom(bcf_hdr_t *h, bcf1_t *v, const char* chrom);
/**
* Get RID
*/
#define bcf_get_rid(v) ((v)->rid)
/**
* Set RID
*/
#define bcf_set_rid(v, c) ((v)->rid=(c))
/**
* Check if variant is passed
*/
bool bcf_is_passed(bcf_hdr_t *h, bcf1_t *v);
/**
* Get 0-based position
*/
#define bcf_get_pos0(v) ((v)->pos)
/**
* Set 0-based position
*/
#define bcf_set_pos0(v, p) ((v)->pos = (p))
/**
* Get 1-based position
*/
#define bcf_get_pos1(v) ((v)->pos+1)
/**
* Get 1-based end position
*/
#define bcf_get_end_pos1(v) ((v)->pos + strlen((v)->d.allele[0]))
/**
* Set 1-based position
*/
#define bcf_set_pos1(v, p) ((v)->pos = (p)-1);
/**
* Get allele
*/
#define bcf_get_allele(v) ((v)->d.allele)
/**
* Get n_alleles of a bcf record
*/
#define bcf_get_n_allele(v) ((v)->n_allele)
/**
* Get reference base for a SNP, assumes the record is a biallelic SNP
*/
#define bcf_get_snp_ref(v) ((v)->d.allele[0][0])
/**
* Get alternative base for a SNP, assumes the record is a biallelic SNP
*/
#define bcf_get_snp_alt(v) ((v)->d.allele[1][0])
/**
* Get reference allele
*/
#define bcf_get_ref(v) ((v)->d.allele[0])
/**
* Get alternative allele
*/
#define bcf_get_alt(v, i) ((v)->d.allele[i])
/**
* Get variant type
*/
#define bcf_get_var_type(v) ((v)->d.var_type)
/**
* Set number of samples in bcf record
*/
#define bcf_set_n_sample(v, n) ((v)->n_sample = (n));
#endif