-
Notifications
You must be signed in to change notification settings - Fork 23
/
Copy pathmocha.h
334 lines (304 loc) · 15.3 KB
/
mocha.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
/* The MIT License
Copyright (C) 2018-2024 Giulio Genovese
Author: Giulio Genovese <giulio.genovese@gmail.com>
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.
*/
/*
Computes beta-binomial likelihoods using recursively pre-calculated tables
*/
#ifndef __MOCHA_H__
#define __MOCHA_H__
#include <htslib/kseq.h>
#include <htslib/ksort.h>
#include <htslib/vcf.h>
#include "tsv2vcf.h"
#include "bcftools.h"
#define GENDER_UNKNOWN 0
#define GENDER_MALE 1
#define GENDER_FEMALE 2
#define GENDER_KLINEFELTER 3
/****************************************
* TSV FUNCTIONS *
****************************************/
// adapted from Petr Danecek's implementation of tsv_init() in bcftools/tsv2vcf.c
tsv_t *tsv_init_delimiter(const char *str, char delimiter) {
tsv_t *tsv = (tsv_t *)calloc(1, sizeof(tsv_t));
kstring_t tmp = {0, 0, 0};
const char *ss = str, *se = ss;
tsv->ncols = 0;
while (*ss) {
if (delimiter == '\0')
while (*se && !isspace(*se)) se++;
else
while (*se && *se != delimiter) se++;
tsv->ncols++;
tsv->cols = (tsv_col_t *)realloc(tsv->cols, sizeof(tsv_col_t) * tsv->ncols);
tsv->cols[tsv->ncols - 1].name = NULL;
tsv->cols[tsv->ncols - 1].setter = NULL;
tmp.l = 0;
kputsn(ss, se - ss, &tmp);
if (strcasecmp("-", tmp.s)) tsv->cols[tsv->ncols - 1].name = strdup(tmp.s);
if (!*se) break;
se++;
if (delimiter == '\0')
while (*se && isspace(*se)) se++;
ss = se;
}
free(tmp.s);
return tsv;
}
int tsv_parse_delimiter(tsv_t *tsv, bcf1_t *rec, char *str, char delimiter) {
int status = 0;
tsv->icol = 0;
tsv->ss = tsv->se = str;
while (*tsv->ss && tsv->icol < tsv->ncols) {
if (delimiter == '\0')
while (*tsv->se && !isspace(*tsv->se)) tsv->se++;
else
while (*tsv->se && *tsv->se != delimiter) tsv->se++;
if (tsv->cols[tsv->icol].setter) {
int ret = tsv->cols[tsv->icol].setter(tsv, rec, tsv->cols[tsv->icol].usr);
if (ret < 0) return -1;
status++;
}
if (*tsv->se) {
tsv->se++;
if (delimiter == '\0')
while (*tsv->se && isspace(*tsv->se)) tsv->se++;
}
tsv->ss = tsv->se;
tsv->icol++;
}
return status ? 0 : -1;
}
int tsv_read_sample_id(tsv_t *tsv, bcf1_t *rec, void *usr) {
bcf_hdr_t *hdr = (bcf_hdr_t *)rec;
int *idx = (int *)usr;
char tmp = *tsv->se;
*tsv->se = 0;
*idx = bcf_hdr_id2int(hdr, BCF_DT_SAMPLE, tsv->ss);
*tsv->se = tmp;
return 0;
}
int tsv_read_computed_gender(tsv_t *tsv, bcf1_t *rec, void *usr) {
int *computed_gender = (int *)usr;
if (toupper(*tsv->ss) == 'M') {
*computed_gender = GENDER_MALE;
} else if (toupper(*tsv->ss) == 'F') {
*computed_gender = GENDER_FEMALE;
} else if (toupper(*tsv->ss) == 'K') {
*computed_gender = GENDER_KLINEFELTER;
} else if (toupper(*tsv->ss) == 'U' || toupper(*tsv->ss) == 'N') {
*computed_gender = GENDER_UNKNOWN;
} else {
char *endptr;
*computed_gender = (int)strtol(tsv->ss, &endptr, 0);
if (tsv->ss == endptr)
error(
"Could not parse gender %.*s\n(Acceptable values: 1/M/m = male, 2/F/f = female, 0/U/u/N/n = missing)\n",
(int)(tsv->se - tsv->ss), tsv->ss);
}
return 0;
}
/****************************************
* BASIC STATISTICS FUNCTIONS *
****************************************/
// this macro from ksort.h defines the function
// float ks_ksmall_float(size_t n, float arr[], size_t kk);
KSORT_INIT_GENERIC(float)
// compute the median of a vector using the ksort library (with iterator)
float get_median(const float *v, int n, const int *imap) {
if (n == 0) return NAN;
float tmp, *w = (float *)malloc(n * sizeof(float));
int i, j = 0;
for (i = 0; i < n; i++) {
tmp = imap ? v[imap[i]] : v[i];
if (!isnan(tmp)) w[j++] = tmp;
}
if (j == 0) {
free(w);
return NAN;
}
float ret = ks_ksmall_float((size_t)j, w, (size_t)j / 2);
if (j % 2 == 0) ret = (ret + w[j / 2 - 1]) * 0.5f;
free(w);
return ret;
}
// rho = xyss / sqrtf(xss * yss)
// m = xyss / xss
// b = ym - m * xm;
static inline int get_cov(const float *x, const float *y, int n, const int *imap, float *xss, float *yss, float *xyss) {
if (n < 2) return -1;
float xm = 0.0f;
float ym = 0.0f;
int i;
for (i = 0; i < n; i++) {
int idx = imap ? imap[i] : i;
xm += x[idx];
ym += y[idx];
}
xm /= n;
ym /= n;
*xss = 0.0f;
*yss = 0.0f;
*xyss = 0.0f;
for (i = 0; i < n; i++) {
int idx = imap ? imap[i] : i;
float xd = x[idx] - xm;
float yd = y[idx] - ym;
*xss += xd * xd;
*yss += yd * yd;
*xyss += xd * yd;
}
return 0;
}
/****************************************
* EXTRACT DATA FROM VCF FUNCTIONS *
****************************************/
// retrieve unphased genotype alleles information from BCF record
// assumes little endian architecture
static inline int bcf_get_unphased_genotype_alleles(const bcf_fmt_t *fmt, int16_t *gt0_arr, int16_t *gt1_arr,
int nsmpl) {
if (!fmt || fmt->n != 2) return 0;
int i;
#define BRANCH(type_t, bcf_type_vector_end) \
{ \
type_t *p = (type_t *)fmt->p; \
for (i = 0; i < nsmpl; i++, p += 2) { \
if (p[0] == bcf_type_vector_end || bcf_gt_is_missing(p[0]) || p[1] == bcf_type_vector_end \
|| bcf_gt_is_missing(p[1])) { \
gt0_arr[i] = bcf_int16_missing; \
gt1_arr[i] = bcf_int16_missing; \
} else { \
int swap = bcf_gt_allele(p[0]) > bcf_gt_allele(p[1]); \
gt0_arr[i] = (int16_t)bcf_gt_allele(p[swap]); \
gt1_arr[i] = (int16_t)bcf_gt_allele(p[!swap]); \
} \
} \
}
switch (fmt->type) {
case BCF_BT_INT8:
BRANCH(int8_t, bcf_int8_vector_end);
break;
case BCF_BT_INT16:
BRANCH(int16_t, bcf_int16_vector_end);
break;
case BCF_BT_INT32:
BRANCH(int32_t, bcf_int32_vector_end);
break;
default:
error("Unexpected type %d\n", fmt->type);
}
#undef BRANCH
return 1;
}
// retrieve phase information from BCF record
// bcf_int8_missing if genotype is missing
// bcf_int8_vector_end if phase does not apply
// 0 if phase is not available
// 1 if higher number allele received from the mother
// -1 if higher number allele received from the father
// assumes little endian architecture
static inline int bcf_get_genotype_phase(const bcf_fmt_t *fmt, int8_t *gt_phase_arr, int nsmpl) {
if (!fmt || fmt->n != 2) return 0;
int i;
#define BRANCH(type_t, bcf_type_vector_end) \
{ \
type_t *p = (type_t *)fmt->p; \
for (i = 0; i < nsmpl; i++, p += 2) { \
if (p[0] == bcf_type_vector_end || bcf_gt_is_missing(p[0]) || p[1] == bcf_type_vector_end \
|| bcf_gt_is_missing(p[1])) { \
gt_phase_arr[i] = bcf_int8_missing; \
} else { \
type_t gt0 = bcf_gt_allele(p[0]); \
type_t gt1 = bcf_gt_allele(p[1]); \
if (gt0 == gt1) \
gt_phase_arr[i] = bcf_int8_vector_end; \
else if (!bcf_gt_is_phased(p[1])) \
gt_phase_arr[i] = 0; \
else if (gt1 > gt0) \
gt_phase_arr[i] = 1; \
else \
gt_phase_arr[i] = -1; \
} \
} \
}
switch (fmt->type) {
case BCF_BT_INT8:
BRANCH(int8_t, bcf_int8_vector_end);
break;
case BCF_BT_INT16:
BRANCH(int16_t, bcf_int16_vector_end);
break;
case BCF_BT_INT32:
BRANCH(int32_t, bcf_int32_vector_end);
break;
default:
error("Unexpected type %d\n", fmt->type);
}
#undef BRANCH
return 1;
}
// retrive allelic depth information from BCF record
// assumes little endian architecture
static inline int bcf_get_allelic_depth(const bcf_fmt_t *fmt, const int16_t *gt0_arr, const int16_t *gt1_arr,
int16_t *ad0_arr, int16_t *ad1_arr, int nsmpl) {
if (!fmt) return 0;
int i, nalleles = fmt->n;
#define BRANCH(type_t, bcf_type_vector_end, bcf_type_missing) \
{ \
type_t *p = (type_t *)fmt->p; \
for (i = 0; i < nsmpl; i++, p += nalleles) { \
if ((gt0_arr[i] != bcf_int16_missing && (int)gt0_arr[i] >= nalleles) \
|| (gt1_arr[i] != bcf_int16_missing && (int)gt1_arr[i] >= nalleles)) \
error("Error: found VCF record with GT alleles %d and %d and %d number of alleles\n", gt0_arr[i], \
gt1_arr[i], nalleles); \
if (gt0_arr[i] == bcf_int16_missing || gt1_arr[i] == bcf_int16_missing \
|| p[gt0_arr[i]] == bcf_type_vector_end || p[gt0_arr[i]] == bcf_type_missing \
|| p[gt1_arr[i]] == bcf_type_vector_end || p[gt1_arr[i]] == bcf_type_missing) { \
ad0_arr[i] = bcf_int16_missing; \
ad1_arr[i] = bcf_int16_missing; \
} else { \
type_t gt0 = gt0_arr[i] > 0; \
type_t gt1 = gt1_arr[i] > 0; \
if (gt0 == gt1) { \
ad0_arr[i] = (int16_t)(gt0_arr[i] == gt1_arr[i] ? p[gt0_arr[i]] : p[gt0_arr[i]] + p[gt1_arr[i]]); \
ad1_arr[i] = bcf_int16_missing; \
} else { \
ad0_arr[i] = (int16_t)p[gt0_arr[i]]; \
ad1_arr[i] = (int16_t)p[gt1_arr[i]]; \
} \
} \
} \
}
switch (fmt->type) {
case BCF_BT_INT8:
BRANCH(int8_t, bcf_int8_vector_end, bcf_int8_missing);
break;
case BCF_BT_INT16:
BRANCH(int16_t, bcf_int16_vector_end, bcf_int16_missing);
break;
case BCF_BT_INT32:
BRANCH(int32_t, bcf_int32_vector_end, bcf_int32_missing);
break;
default:
error("Unexpected type %d\n", fmt->type);
}
#undef BRANCH
return 1;
}
#endif