-
Notifications
You must be signed in to change notification settings - Fork 23
/
mocha.c
3751 lines (3446 loc) · 163 KB
/
mocha.c
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
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
/* The MIT License
Copyright (C) 2015-2024 Giulio Genovese
Author: Giulio Genovese <[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 <stdio.h>
#include <unistd.h>
#include <getopt.h>
#include <errno.h>
#include <ctype.h>
#include <math.h>
#include <htslib/kseq.h>
#include <htslib/vcf.h>
#include <htslib/synced_bcf_reader.h>
#include <htslib/ksort.h>
#include "regidx.h" // cannot use htslib/regdix.h see http://github.com/samtools/htslib/pull/761
#include "kmin.h"
#include "mocha.h"
#include "genome_rules.h"
#include "beta_binom.h"
#include "bcftools.h"
#include "filter.h"
#include "tsv2vcf.h"
#define MOCHA_VERSION "2024-09-27"
/****************************************
* CONSTANT DEFINITIONS *
****************************************/
#define SIGN(x) (((x) > 0) - ((x) < 0))
#define BDEV_LRR_BAF_DFLT "-2.0,-4.0,-6.0,10.0,6.0,4.0"
#define BDEV_BAF_PHASE_DFLT "6.0,8.0,10.0,15.0,20.0,30.0,50.0,80.0,130.0,210.0,340.0,550.0"
#define MIN_DST_DFLT "400"
#define ADJ_BAF_LRR_DFLT "5"
#define REGRESS_BAF_LRR_DFLT "15"
#define LRR_GC_ORDER_DFLT "2"
#define MAX_ORDER 5
#define XY_MAJOR_PL_DFLT "65.0"
#define XY_MINOR_PL_DFLT "35.0"
#define AUTO_TEL_PL_DFLT "20.0"
#define CHRX_TEL_PL_DFLT "8.0"
#define CHRY_TEL_PL_DFLT "6.0"
#define ERR_PL_DFLT "15.0"
#define FLIP_PL_DFLT "20.0"
#define SHORT_ARM_CHRS_DFLT "13,14,15,21,22,chr13,chr14,chr15,chr21,chr22"
#define LRR_BIAS_DFLT "0.2"
// http://www.illumina.com/documents/products/technotes/technote_cnv_algorithms.pdf
#define LRR_HAP2DIP_DFLT "0.45"
#define FLT_INCLUDE (1 << 0)
#define FLT_EXCLUDE (1 << 1)
#define WGS_DATA (1 << 2)
#define NO_LOG (1 << 3)
#define NO_ANNOT (1 << 4)
#define USE_SHORT_ARMS (1 << 5)
#define USE_CENTROMERES (1 << 6)
#define USE_MALES_XTR (1 << 7)
#define USE_MALES_PAR2 (1 << 8)
#define USE_NO_RULES_CHRS (1 << 9)
#define LRR 0
#define BAF 1
#define AD0 0
#define AD1 1
#define LDEV 0
#define BDEV 1
#define LRR_BAF 0
#define BAF_PHASE 1
#define MOCHA_UNDET 0
#define MOCHA_LOSS 1
#define MOCHA_GAIN 2
#define MOCHA_CNLOH 3
#define MOCHA_CNP_LOSS 4
#define MOCHA_CNP_GAIN 5
#define MOCHA_CNP_CNV 6
#define MOCHA_NOT 0
#define MOCHA_CEN 1
#define MOCHA_ARM 2
#define MOCHA_TEL 3
#define GT_NC 0
#define GT_AA 1
#define GT_AB 2
#define GT_BB 3
/****************************************
* DATA STRUCTURES *
****************************************/
typedef struct {
int pos;
int allele_a;
int allele_b;
float adjust[3][3]; // shift
} locus_t;
typedef struct {
int allele_a_id, allele_b_id, gc_id, gt_id, ad_id, baf_id, lrr_id;
float xy_major_log_prb;
float xy_minor_log_prb;
float auto_tel_log_prb;
float chrX_tel_log_prb;
float chrY_tel_log_prb;
float err_log_prb;
float flip_log_prb;
float *bdev_lrr_baf, *bdev_baf_phase;
int bdev_lrr_baf_n, bdev_baf_phase_n;
int min_dst;
float lrr_cutoff;
float lrr_hap2dip;
float lrr_bias;
int adj_baf_lrr;
int regress_baf_lrr;
int lrr_gc_order;
int flags;
int filter_logic;
filter_t *filter;
genome_rules_t *genome_rules;
regidx_t *cnp_idx;
regitr_t *cnp_itr;
regidx_t *mhc_idx;
regidx_t *kir_idx;
int rid;
int n;
locus_t *locus_arr;
int m_locus;
float *gc_arr;
int m_gc;
int n_flipped;
} model_t;
typedef struct {
int sample_idx;
int computed_gender;
int rid;
int beg_pos;
int end_pos;
int length;
int8_t p_arm;
int8_t q_arm;
int n_sites;
int n_hets;
int n50_hets;
float bdev;
float bdev_se;
float ldev;
float ldev_se;
float lod_lrr_baf;
float lod_baf_phase;
int n_flips;
float baf_conc;
float lod_baf_conc;
int8_t type;
float cf;
} mocha_t;
typedef struct {
int n, m;
mocha_t *a;
} mocha_table_t;
typedef struct {
float call_rate;
float lrr_median;
float lrr_sd;
float lrr_auto;
float dispersion; // either rho(AD0, AD1) for WGS model or sd(BAF)
float baf_conc;
float baf_auto;
float lrr_gc_rel_ess;
float coeffs[MAX_ORDER + 1];
} stats_t;
typedef struct {
int idx;
int computed_gender;
float adjlrr_sd;
int n_sites;
int n_missing_gts;
int n_hets;
int x_nonpar_n_hets;
int par1_n_hets;
int xtr_n_hets;
int par2_n_hets;
float x_nonpar_dispersion; // either rho(AD0, AD1) for WGS model or sd(BAF)
float x_nonpar_lrr_median;
float y_nonpar_lrr_median;
float mt_lrr_median;
stats_t stats;
stats_t *stats_arr;
int m_stats, n_stats;
int n;
int *vcf_imap_arr;
int m_vcf_imap;
int16_t *data_arr[2];
int m_data[2];
int8_t *phase_arr;
int m_phase;
} sample_t;
/****************************************
* INLINE FUNCTIONS AND CONSTANTS *
****************************************/
// this macro from ksort.h defines the function
// void ks_introsort_int(size_t n, int a[]);
KSORT_INIT_GENERIC(int)
static inline float sqf(float x) { return x * x; }
static inline double sq(double x) { return x * x; }
// the x == y is necessary in case x == -INFINITY
static inline float log_mean_expf(float x, float y) {
return x == y ? x : (x > y ? x + logf(1.0f + expf(y - x)) : y + logf(1.0f + expf(x - y))) - (float)M_LN2;
}
static beta_binom_t *beta_binom_null, *beta_binom_alt;
static inline float beta_binom_log_lkl(const beta_binom_t *self, int16_t ad0, int16_t ad1) {
return ad0 == bcf_int16_missing || ad1 == bcf_int16_missing ? 0.0f : beta_binom_log_unsafe(self, ad0, ad1);
}
/****************************************
* CONVERT FLOAT TO INT16 AND VICEVERSA *
****************************************/
#define INT16_SCALE 1000 // BAF values from Illumina are scaled to 1000
static inline int16_t float_to_int16(float in) {
return isnan(in) ? bcf_int16_missing : (int16_t)roundf(INT16_SCALE * in);
}
static inline float int16_to_float(int16_t in) { return in == bcf_int16_missing ? NAN : ((float)in) / INT16_SCALE; }
/******************************************
* LRR AND COVERAGE POLYNOMIAL REGRESSION *
******************************************/
// http://stackoverflow.com/questions/5083465/fast-efficient-least-squares-fit-algorithm-in-c
// http://github.com/natedomin/polyfit/blob/master/polyfit.c
// http://software.intel.com/sites/products/documentation/doclib/mkl_sa/11/mkl_lapack_examples/sgels_ex.c.htm
// this function needs to use doubles internally when dealing with WGS data
static int polyfit(const float *lrr, const float *gc, int n, const int *imap, int order, float *coeffs) {
int i, j, k, m = order + 1;
if (n < m || order > MAX_ORDER) return -1;
double B[MAX_ORDER + 1] = {0.0};
double P[((MAX_ORDER + 1) * 2) + 1] = {0.0};
double A[(MAX_ORDER + 1) * 2 * (MAX_ORDER + 1)] = {0.0};
// identify the column vector
for (i = 0; i < n; i++) {
float x = imap ? gc[imap[i]] : gc[i];
float y = lrr[i];
if (isnan(x) || isnan(y)) continue;
float powx = 1.0f;
for (j = 0; j < m; j++) {
B[j] += (double)(y * powx);
powx *= x;
}
}
// initialize the PowX array
P[0] = (float)n;
// compute the sum of the powers of X
for (i = 0; i < n; i++) {
float x = imap ? gc[imap[i]] : gc[i];
if (isnan(x)) continue;
float powx = x;
for (j = 1; j < ((2 * m) + 1); j++) {
P[j] += (double)powx;
powx *= x;
}
}
// initialize the reduction matrix
for (i = 0; i < m; i++) {
for (j = 0; j < m; j++) {
A[(i * (2 * m)) + j] = P[i + j];
}
A[(i * (2 * m)) + (i + m)] = 1.0;
}
// move the identity matrix portion of the redux matrix to the
// left side (find the inverse of the left side of the redux matrix)
for (i = 0; i < m; i++) {
double x = A[(i * (2 * m)) + i];
if (x != 0) {
for (k = 0; k < (2 * m); k++) {
A[(i * (2 * m)) + k] /= x;
}
for (j = 0; j < m; j++) {
if (i != j) {
double y = A[(j * (2 * m)) + i];
for (k = 0; k < (2 * m); k++) {
A[(j * (2 * m)) + k] -= y * A[(i * (2 * m)) + k];
}
}
}
} else {
// cannot work with singular matrices
return -1;
}
}
// calculate coefficients
for (i = 0; i < m; i++) {
for (j = 0; j < m; j++) {
double x = 0.0;
for (k = 0; k < m; k++) {
x += (A[(i * (2 * m)) + (k + m)] * B[k]);
}
coeffs[i] = (float)x;
}
}
return 0;
}
static void ad_to_lrr_baf(const int16_t *ad0, const int16_t *ad1, float *lrr, float *baf, int n) {
// this function keeps a list of logarithms of integers to minimize log calls
static float *logf_arr = NULL;
static int n_logf = 0, m_logf = 0;
if (ad0 == NULL && ad1 == NULL && n == 0) {
free(logf_arr);
return;
}
int i, j;
for (i = 0; i < n; i++) {
if (ad0[i] == bcf_int16_missing && ad1[i] == bcf_int16_missing) {
lrr[i] = NAN;
baf[i] = NAN;
continue;
}
int cov = (int)(ad0[i] == bcf_int16_missing ? 0 : ad0[i]) + (int)(ad1[i] == bcf_int16_missing ? 0 : ad1[i]);
if (cov == 0) {
lrr[i] = 0;
baf[i] = NAN;
} else {
if (cov > n_logf) {
hts_expand(float, cov, m_logf, logf_arr);
for (j = n_logf; j < cov; j++) logf_arr[j] = logf(j + 1);
n_logf = cov;
}
lrr[i] = logf_arr[cov - 1];
baf[i] = (ad0[i] == bcf_int16_missing || ad1[i] == bcf_int16_missing) ? NAN : (float)ad1[i] / (float)cov;
}
}
}
static void adjust_lrr(float *lrr, const float *gc, int n, const int *imap, const float *coeffs, int order) {
int i, j;
for (i = 0; i < n; i++) {
float x = imap ? gc[imap[i]] : gc[i];
float powx = 1.0f;
for (j = 0; j <= order; j++) {
lrr[i] -= coeffs[j] * powx;
powx *= x;
}
}
}
// computes total sum of squares
// this function needs to use doubles internally when dealing with WGS data
static float get_tss(const float *v, int n) {
double mean = 0.0;
int i, j = 0;
for (i = 0; i < n; i++) {
if (!isnan(v[i])) {
mean += (double)v[i];
j++;
}
}
if (j <= 1) return NAN;
mean /= (double)j;
double tss = 0.0;
for (i = 0; i < n; i++) {
if (!isnan(v[i])) tss += sq((double)v[i] - mean);
}
return (float)tss;
}
/*********************************
* HMM AND OPTIMIZATION METHODS *
*********************************/
// compute Viterbi path from probabilities
static int8_t *retrace_viterbi(int T, int N, const float *log_prb, const int8_t *ptr) {
int i, t;
int8_t *path = (int8_t *)malloc(T * sizeof(int8_t));
// initialize last path state
path[T - 1] = 0;
for (i = 1; i < N; i++)
if (log_prb[(int)path[T - 1]] < log_prb[i]) path[T - 1] = (int8_t)i;
// compute best path by tracing back the Markov chain
for (t = T - 1; t > 0; t--) path[t - 1] = ptr[(t - 1) * N + (int)path[t]];
return path;
}
// rescale Viterbi log probabilities to avoid underflow issues
static inline void rescale_log_prb(float *log_prb, int n) {
float max = -INFINITY;
int i;
for (i = 0; i < n; i++) max = max > log_prb[i] ? max : log_prb[i];
for (i = 0; i < n; i++) log_prb[i] -= max;
}
// compute the Viterbi path from BAF
// n is the length of the hidden Markov model
// m is the number of possible BAF deviations
static int8_t *log_viterbi_run(const float *emis_log_lkl, int T, int m, float xy_major_log_prb, float xy_minor_log_prb,
float tel_log_prb, float flip_log_prb, int last_p, int first_q) {
int t, i, j, changeidx;
// determine the number of hidden states based on whether phase information is used
int N = 1 + m + (isnan(flip_log_prb) ? 0 : m);
// allocate memory necessary for running the algorithm
float *log_prb = (float *)malloc(N * sizeof(float));
float *new_log_prb = (float *)malloc(N * sizeof(float));
int8_t *ptr = (int8_t *)malloc(N * (T - 1) * sizeof(int8_t));
int8_t *path;
// initialize and rescale the first state
log_prb[0] = emis_log_lkl[0];
for (i = 1; i < N; i++) log_prb[i] = tel_log_prb + emis_log_lkl[i];
rescale_log_prb(log_prb, N);
// compute best probabilities at each position
for (t = 1; t < T; t++) {
// this causes a penalty for mosaic chromosomal calls across the centromeres
float exit_log_prb = t > last_p ? xy_major_log_prb : xy_minor_log_prb;
float enter_log_prb = t < first_q ? xy_major_log_prb : xy_minor_log_prb;
for (i = 0; i < N; i++) {
new_log_prb[i] = log_prb[i];
ptr[(t - 1) * N + i] = (int8_t)i;
}
// compute whether a state switch should be considered for null state
for (i = 1; i < N; i++) {
if (new_log_prb[0] < log_prb[i] + exit_log_prb) {
new_log_prb[0] = log_prb[i] + exit_log_prb;
ptr[(t - 1) * N] = ptr[(t - 1) * N + i];
}
if (new_log_prb[i] < log_prb[0] + enter_log_prb) {
new_log_prb[i] = log_prb[0] + enter_log_prb;
ptr[(t - 1) * N + i] = ptr[(t - 1) * N];
}
}
// compute whether a state switch should be considered for each other state
// it will run twice if and only if phasing is used
for (j = 0; j == 0 || (!isnan(flip_log_prb) && j == m); j += m) {
float change_log_prb = log_prb[0] + enter_log_prb;
changeidx = 0;
for (i = 0; i < m; i++) {
if (change_log_prb < log_prb[1 + j + i] + (xy_major_log_prb + xy_minor_log_prb) * 0.75f) {
change_log_prb = log_prb[1 + j + i] + (xy_major_log_prb + xy_minor_log_prb) * 0.75f;
changeidx = 1 + j + i;
}
}
for (i = 0; i < m; i++) {
if (new_log_prb[1 + j + i] < change_log_prb) {
new_log_prb[1 + j + i] = change_log_prb;
ptr[(t - 1) * N + 1 + j + i] = ptr[(t - 1) * N + changeidx];
}
}
}
// compute whether a phase flip should be considered for non-null states
if (!isnan(flip_log_prb)) {
for (i = 0; i < m; i++) {
if (new_log_prb[1 + i] < new_log_prb[1 + m + i] + flip_log_prb) {
new_log_prb[1 + i] = new_log_prb[1 + m + i] + flip_log_prb;
ptr[(t - 1) * N + 1 + i] = ptr[(t - 1) * N + 1 + m + i];
}
if (new_log_prb[1 + m + i] < new_log_prb[1 + i] + flip_log_prb) {
new_log_prb[1 + m + i] = new_log_prb[1 + i] + flip_log_prb;
ptr[(t - 1) * N + 1 + m + i] = ptr[(t - 1) * N + 1 + i];
}
}
}
// update and rescale the current state
new_log_prb[0] += emis_log_lkl[t * N];
for (i = 0; i < m; i++) {
new_log_prb[1 + i] += emis_log_lkl[t * N + 1 + i];
if (!isnan(flip_log_prb)) new_log_prb[1 + m + i] += emis_log_lkl[t * N + 1 + m + i];
}
for (i = 0; i < N; i++) log_prb[i] = new_log_prb[i];
rescale_log_prb(log_prb, N);
}
// add closing cost to the last state
for (i = 1; i < N; i++) log_prb[i] += tel_log_prb;
rescale_log_prb(log_prb, N);
path = retrace_viterbi(T, N, log_prb, ptr);
// free memory
free(log_prb);
free(new_log_prb);
free(ptr);
// symmetrize the path
if (!isnan(flip_log_prb))
for (i = 0; i < T; i++)
if (path[i] > m) path[i] = (int8_t)m - path[i];
return path;
}
/*********************************
* LRR AND BAF LIKELIHOODS *
*********************************/
// rescale emission probabilities to avoid problems with outliers
// TODO find a better name for this function
static void rescale_emis_log_lkl(float *log_prb, int n, float err_log_prb) {
float min_thr = -INFINITY;
int i;
for (i = 0; i < n; i++)
if (min_thr < log_prb[i]) min_thr = log_prb[i];
min_thr += err_log_prb;
for (i = 0; i < n; i++)
if (log_prb[i] < min_thr) log_prb[i] = min_thr;
}
static inline float norm_log_lkl(float x, float m, float s, float w) {
return isnan(x) ? 0.0f : -0.5f * sqf((x - m) / s) * w;
}
// lrr_bias is used in a different way from what done by Petr Danecek in bcftools/vcfcnv.c
static inline float lrr_baf_log_lkl(float lrr, float baf, float ldev, float bdev, float lrr_sd, float baf_sd,
float lrr_bias) {
return norm_log_lkl(lrr, ldev, lrr_sd, lrr_bias)
+ log_mean_expf(norm_log_lkl(baf - 0.5f, bdev, baf_sd, 1.0f), norm_log_lkl(baf - 0.5f, -bdev, baf_sd, 1.0f));
}
static inline float baf_phase_log_lkl(float baf, int8_t phase, float bdev, float baf_sd) {
return phase == 0 ? log_mean_expf(norm_log_lkl(baf - 0.5f, bdev, baf_sd, 1.0f),
norm_log_lkl(baf - 0.5f, -bdev, baf_sd, 1.0f))
: norm_log_lkl(baf - 0.5f, (float)SIGN(phase) * bdev, baf_sd, 1.0f);
}
// precomupute emission probabilities
static float *lrr_baf_emis_log_lkl(const float *lrr, const float *baf, int T, const int *imap, float err_log_prb,
float lrr_bias, float lrr_hap2dip, float lrr_sd, float baf_sd,
const float *bdev_lrr_baf_arr, int m) {
float *ldev = (float *)malloc(m * sizeof(float));
int i, t;
for (i = 0; i < m; i++) ldev[i] = -logf(1.0f - 2.0f * bdev_lrr_baf_arr[i]) / (float)M_LN2 * lrr_hap2dip;
int N = 1 + 2 * m;
float *emis_log_lkl = (float *)malloc(N * T * sizeof(float));
for (t = 0; t < T; t++) {
float x = imap ? lrr[imap[t]] : lrr[t];
float y = imap ? baf[imap[t]] : baf[t];
emis_log_lkl[t * N] = lrr_baf_log_lkl(x, y, 0.0f, 0.0f, lrr_sd, baf_sd, lrr_bias);
for (i = 0; i < m; i++) {
emis_log_lkl[t * N + 1 + i] = lrr_baf_log_lkl(x, y, ldev[i], bdev_lrr_baf_arr[i], lrr_sd, baf_sd, lrr_bias);
}
// add states to distinguish LRR waves from true mosaic gains/losses
for (i = 0; i < m; i++) {
if (bdev_lrr_baf_arr[i] < -1.0f / 6.0f || bdev_lrr_baf_arr[i] >= 1.0f / 6.0f)
emis_log_lkl[t * N + 1 + m + i] = emis_log_lkl[t * N] + err_log_prb;
else
emis_log_lkl[t * N + 1 + m + i] =
lrr_baf_log_lkl(x, y, bdev_lrr_baf_arr[i], 0.0f, lrr_sd, baf_sd, lrr_bias);
}
rescale_emis_log_lkl(&emis_log_lkl[t * N], N, err_log_prb);
}
free(ldev);
return emis_log_lkl;
}
// precomupute emission probabilities
static float *baf_phase_emis_log_lkl(const float *baf, const int8_t *gt_phase, int T, const int *imap,
float err_log_prb, float baf_sd, const float *bdev, int m) {
int i, t, N = 1 + 2 * m;
float *emis_log_lkl = (float *)malloc(N * T * sizeof(float));
for (t = 0; t < T; t++) {
float x = imap ? baf[imap[t]] : baf[t];
int8_t p = imap ? gt_phase[imap[t]] : gt_phase[t];
emis_log_lkl[t * N] = baf_phase_log_lkl(x, (int8_t)1, 0.0f, baf_sd);
for (i = 0; i < m; i++) {
emis_log_lkl[t * N + 1 + i] = baf_phase_log_lkl(x, p, bdev[i], baf_sd);
if (p == 0)
emis_log_lkl[t * N + 1 + m + i] = emis_log_lkl[t * N + 1 + i];
else
emis_log_lkl[t * N + 1 + m + i] = baf_phase_log_lkl(x, p, -bdev[i], baf_sd);
}
rescale_emis_log_lkl(&emis_log_lkl[t * N], N, err_log_prb);
}
return emis_log_lkl;
}
static int cnp_edge_is_not_cn2_lrr_baf(const float *lrr, const float *baf, int n, int a, int b, float xy_log_prb,
float err_log_prb, float lrr_bias, float lrr_hap2dip, float lrr_sd, float baf_sd,
float ldev, float bdev) {
// test left edge
float sum_log_lkl = 0.0f;
int i;
for (i = a - 1; i >= 0; i--) {
float log_lkl = lrr_baf_log_lkl(lrr[i], baf[i], ldev, bdev, lrr_sd, baf_sd, lrr_bias)
- lrr_baf_log_lkl(lrr[i], baf[i], 0.0f, 0.0f, lrr_sd, baf_sd, lrr_bias);
if (!isnan(err_log_prb)) {
if (log_lkl < err_log_prb)
log_lkl = err_log_prb;
else if (log_lkl > -err_log_prb)
log_lkl = -err_log_prb;
}
sum_log_lkl += log_lkl;
if (sum_log_lkl > -xy_log_prb) return -1;
if (sum_log_lkl < xy_log_prb) break;
}
// test right edge
sum_log_lkl = 0.0f;
for (i = b + 1; i < n; i++) {
float log_lkl = lrr_baf_log_lkl(lrr[i], baf[i], ldev, bdev, lrr_sd, baf_sd, lrr_bias)
- lrr_baf_log_lkl(lrr[i], baf[i], 0, 0, lrr_sd, baf_sd, lrr_bias);
if (!isnan(err_log_prb)) {
if (log_lkl < err_log_prb)
log_lkl = err_log_prb;
else if (log_lkl > -err_log_prb)
log_lkl = -err_log_prb;
}
sum_log_lkl += log_lkl;
if (sum_log_lkl > -xy_log_prb) return -1;
if (sum_log_lkl < xy_log_prb) break;
}
return 0;
}
// return the LOD likelihood for a segment
typedef struct {
const float *lrr_arr;
const float *baf_arr;
int n;
const int *imap;
float err_log_prb;
float lrr_bias;
float lrr_hap2dip;
float lrr_sd;
float baf_sd;
} minus_lrr_baf_lod_t;
static double minus_lrr_baf_lod(double bdev_lrr_baf, void *ap) {
minus_lrr_baf_lod_t *data = (minus_lrr_baf_lod_t *)ap;
if (data->n == 0 || bdev_lrr_baf < -0.5 || bdev_lrr_baf > 0.25) return INFINITY; // kmin_brent does not handle NAN
float ldev = -logf(1.0f - 2.0f * (float)bdev_lrr_baf) / (float)M_LN2 * data->lrr_hap2dip;
float ret = 0.0f;
int i;
for (i = 0; i < data->n; i++) {
float lrr = data->imap ? data->lrr_arr[data->imap[i]] : data->lrr_arr[i];
float baf = data->imap ? data->baf_arr[data->imap[i]] : data->baf_arr[i];
float log_lkl = lrr_baf_log_lkl(lrr, baf, ldev, (float)bdev_lrr_baf, data->lrr_sd, data->baf_sd, data->lrr_bias)
- lrr_baf_log_lkl(lrr, baf, 0.0f, 0.0f, data->lrr_sd, data->baf_sd, data->lrr_bias);
if (!isnan(data->err_log_prb)) {
if (log_lkl < data->err_log_prb)
log_lkl = data->err_log_prb;
else if (log_lkl > -data->err_log_prb)
log_lkl = -data->err_log_prb;
}
ret += log_lkl;
}
return -(double)ret * M_LOG10E;
}
// return the LOD likelihood for a segment
typedef struct {
const float *baf_arr;
const int8_t *gt_phase;
int n;
const int *imap;
const int8_t *as;
float err_log_prb;
float baf_sd;
} minus_baf_phase_lod_t;
static double minus_baf_phase_lod(double bdev, void *ap) {
minus_baf_phase_lod_t *data = (minus_baf_phase_lod_t *)ap;
if (data->n == 0 || bdev < 0.0 || bdev > 0.5) return INFINITY; // kmin_brent does not handle NAN
float ret = 0.0f;
int i;
for (i = 0; i < data->n; i++) {
float baf = data->imap ? data->baf_arr[data->imap[i]] : data->baf_arr[i];
int8_t p = data->imap ? data->gt_phase[data->imap[i]] : data->gt_phase[i];
if (data->as) p *= (int8_t)SIGN(data->as[i]); // notice as has no imap
float log_lkl =
baf_phase_log_lkl(baf, p, (float)bdev, data->baf_sd) - baf_phase_log_lkl(baf, 0, 0.0f, data->baf_sd);
if (!isnan(data->err_log_prb)) {
if (log_lkl < data->err_log_prb)
log_lkl = data->err_log_prb;
else if (log_lkl > -data->err_log_prb)
log_lkl = -data->err_log_prb;
}
ret += log_lkl;
}
return -(double)ret * M_LOG10E;
}
// TODO find a better title for this function
static float compare_models(const float *baf, const int8_t *gt_phase, int n, const int *imap, float xy_log_prb,
float err_log_prb, float flip_log_prb, float tel_log_prb, float baf_sd, const float *bdev,
int m) {
if (n == 0) return NAN;
float *emis_log_lkl = baf_phase_emis_log_lkl(baf, gt_phase, n, imap, err_log_prb, baf_sd, bdev, m);
int8_t *path = log_viterbi_run(emis_log_lkl, n, m, xy_log_prb, xy_log_prb, tel_log_prb, flip_log_prb, 0, 0);
free(emis_log_lkl);
int n_flips = 0;
int i;
for (i = 1; i < n; i++)
if (path[i - 1] && path[i] && path[i - 1] != path[i]) n_flips++;
minus_baf_phase_lod_t data = {baf, gt_phase, n, imap, path, err_log_prb, baf_sd};
double x, fx = kmin_brent(minus_baf_phase_lod, 0.1, 0.2, (void *)&data, KMIN_EPS, &x);
free(path);
return -(float)fx + (float)n_flips * flip_log_prb * (float)M_LOG10E;
}
typedef struct {
const float *baf_arr;
int n;
const int *imap;
float baf_sd;
} minus_baf_log_lkl_t;
static double minus_baf_log_lkl(double bdev, void *ap) {
minus_baf_log_lkl_t *data = (minus_baf_log_lkl_t *)ap;
if (data->n == 0 || bdev < 0.0 || bdev > 0.5) return INFINITY; // kmin_brent does not handle NAN
double ret = 0.0;
int i;
for (i = 0; i < data->n; i++) {
float baf = data->imap ? data->baf_arr[data->imap[i]] : data->baf_arr[i];
if (isnan(baf)) continue;
float log_lkl = log_mean_expf(norm_log_lkl(baf - 0.5f, (float)bdev, data->baf_sd, 1.0f),
norm_log_lkl(baf - 0.5f, -(float)bdev, data->baf_sd, 1.0f));
ret += (double)log_lkl;
}
return -ret * M_LOG10E;
}
static float get_sample_mean(const float *v, int n, const int *imap) {
float mean = 0.0;
int i, j = 0;
for (i = 0; i < n; i++) {
float tmp = imap ? v[imap[i]] : v[i];
if (!isnan(tmp)) {
mean += tmp;
j++;
}
}
if (j <= 1) return NAN;
return mean /= (float)j;
}
static float get_baf_bdev(const float *baf_arr, int n, const int *imap, float baf_sd) {
double bdev = 0.0;
int i, j = 0;
for (i = 0; i < n; i++) {
float baf = imap ? baf_arr[imap[i]] : baf_arr[i];
if (isnan(baf)) continue;
bdev += fabs((double)baf - 0.5);
j++;
}
if (j == 0) return NAN;
bdev /= j;
// simple method to compute bdev should work well for germline duplications
if ((float)bdev > 2.0f * baf_sd) return (float)bdev;
minus_baf_log_lkl_t data = {baf_arr, n, imap, baf_sd};
kmin_brent(minus_baf_log_lkl, 0.1, 0.2, (void *)&data, KMIN_EPS, &bdev);
return (float)bdev < 1e-4 ? (float)NAN : (float)bdev;
}
static void get_max_sum(const int16_t *ad0, const int16_t *ad1, int n, const int *imap, int *n1, int *n2) {
*n1 = 0;
*n2 = 0;
int i;
for (i = 0; i < n; i++) {
int a = imap ? ad0[imap[i]] : ad0[i];
int b = imap ? ad1[imap[i]] : ad1[i];
if (a != bcf_int16_missing && b != bcf_int16_missing) {
if (a > *n1) *n1 = a;
if (b > *n1) *n1 = b;
if (a + b > *n2) *n2 = a + b;
}
}
}
typedef struct {
const int16_t *ad0_arr;
const int16_t *ad1_arr;
int n;
const int *imap;
float ad_rho;
} minus_ad_log_lkl_t;
static double minus_ad_log_lkl(double bdev, void *ap) {
minus_ad_log_lkl_t *data = (minus_ad_log_lkl_t *)ap;
if (data->n == 0 || bdev < 0.0 || bdev > 0.5) return INFINITY; // kmin_brent does not handle NAN
int n1, n2;
get_max_sum(data->ad0_arr, data->ad1_arr, data->n, data->imap, &n1, &n2);
beta_binom_update(beta_binom_alt, 0.5f + (float)bdev, data->ad_rho, n1, n2);
double ret = 0.0;
int i;
for (i = 0; i < data->n; i++) {
int16_t ad0 = data->imap ? data->ad0_arr[data->imap[i]] : data->ad0_arr[i];
int16_t ad1 = data->imap ? data->ad1_arr[data->imap[i]] : data->ad1_arr[i];
float log_lkl =
log_mean_expf(beta_binom_log_lkl(beta_binom_alt, ad0, ad1), beta_binom_log_lkl(beta_binom_alt, ad1, ad0));
ret += (double)log_lkl;
}
return -(double)ret * M_LOG10E;
}
static float get_ad_bdev(const int16_t *ad0_arr, const int16_t *ad1_arr, int n, const int *imap, float ad_rho) {
double bdev = 0.0;
minus_ad_log_lkl_t data = {ad0_arr, ad1_arr, n, imap, ad_rho};
kmin_brent(minus_ad_log_lkl, 0.1, 0.2, (void *)&data, KMIN_EPS, &bdev);
return (float)bdev < 1e-4 ? (float)NAN : (float)bdev;
}
/*********************************
* WGS AD LIKELIHOODS *
*********************************/
static inline float lrr_ad_log_lkl(float lrr, int16_t ad0, int16_t ad1, float ldev, float lrr_sd, float lrr_bias,
const beta_binom_t *beta_binom) {
return norm_log_lkl(lrr, ldev, lrr_sd, lrr_bias)
+ log_mean_expf(beta_binom_log_lkl(beta_binom, ad0, ad1), beta_binom_log_lkl(beta_binom, ad1, ad0));
}
static inline float ad_phase_log_lkl(int16_t ad0, int16_t ad1, int8_t phase, const beta_binom_t *beta_binom) {
return phase == 0
? log_mean_expf(beta_binom_log_lkl(beta_binom, ad0, ad1), beta_binom_log_lkl(beta_binom, ad1, ad0))
: (phase > 0 ? beta_binom_log_lkl(beta_binom, ad0, ad1) : beta_binom_log_lkl(beta_binom, ad1, ad0));
}
static float *lrr_ad_emis_log_lkl(const float *lrr, const int16_t *ad0, const int16_t *ad1, int T, const int *imap,
float err_log_prb, float lrr_bias, float lrr_hap2dip, float lrr_sd, float ad_rho,
const float *bdev_lrr_baf_arr, int m) {
int i, t, N = 1 + 2 * m;
int n1, n2;
get_max_sum(ad0, ad1, T, NULL, &n1, &n2);
float *emis_log_lkl = (float *)malloc(N * T * sizeof(float));
for (i = 0; i < 1 + m; i++) {
float ldev = i == 0 ? 0.0f : -logf(1.0f - 2.0f * bdev_lrr_baf_arr[i - 1]) / (float)M_LN2 * lrr_hap2dip;
float bdev = i == 0 ? 0.0f : fabsf(bdev_lrr_baf_arr[i - 1]);
beta_binom_t *beta_binom = i == 0 ? beta_binom_null : beta_binom_alt;
beta_binom_update(beta_binom, 0.5f + bdev, ad_rho, n1, n2);
for (t = 0; t < T; t++) {
float x = imap ? lrr[imap[t]] : lrr[t];
int16_t a = imap ? ad0[imap[t]] : ad0[t];
int16_t b = imap ? ad1[imap[t]] : ad1[t];
emis_log_lkl[t * N + i] = lrr_ad_log_lkl(x, a, b, ldev, lrr_sd, lrr_bias, beta_binom);
}
}
// generate states that should attract LRR waves with no BAF signal
for (i = 0; i < m; i++) {
// do not make extreme states compete with LRR waves
if (bdev_lrr_baf_arr[i] < -1.0f / 6.0f || bdev_lrr_baf_arr[i] >= 1.0f / 6.0f) {
for (t = 0; t < T; t++) emis_log_lkl[t * N + 1 + m + i] = emis_log_lkl[t * N] + err_log_prb;
} else {
float ldev = -logf(1.0f - 2.0f * bdev_lrr_baf_arr[i]) / (float)M_LN2 * lrr_hap2dip;
for (t = 0; t < T; t++) {
float x = imap ? lrr[imap[t]] : lrr[t];
int16_t a = imap ? ad0[imap[t]] : ad0[t];
int16_t b = imap ? ad1[imap[t]] : ad1[t];
emis_log_lkl[t * N + 1 + m + i] = lrr_ad_log_lkl(x, a, b, ldev, lrr_sd, lrr_bias, beta_binom_null);
}
}
}
for (t = 0; t < T; t++) rescale_emis_log_lkl(&emis_log_lkl[t * N], N, err_log_prb);
return emis_log_lkl;
}
static float *ad_phase_emis_log_lkl(const int16_t *ad0, const int16_t *ad1, const int8_t *gt_phase, int T,
const int *imap, float err_log_prb, float ad_rho, const float *bdev_arr, int m) {
int i, t, N = 1 + 2 * m;
float *emis_log_lkl = (float *)malloc(N * T * sizeof(float));
for (i = 0; i < 1 + m; i++) {
float bdev = i == 0 ? 0.0f : bdev_arr[i - 1];
// TODO this function should come out of the loop
int n1, n2;
get_max_sum(ad0, ad1, T, imap, &n1, &n2);
beta_binom_t *beta_binom = i == 0 ? beta_binom_null : beta_binom_alt;
beta_binom_update(beta_binom, 0.5f + bdev, ad_rho, n1, n2);
for (t = 0; t < T; t++) {
int16_t a = imap ? ad0[imap[t]] : ad0[t];
int16_t b = imap ? ad1[imap[t]] : ad1[t];
int8_t p = imap ? gt_phase[imap[t]] : gt_phase[t];
emis_log_lkl[t * N + i] = ad_phase_log_lkl(a, b, p, beta_binom);
if (i > 0) emis_log_lkl[t * N + m + i] = ad_phase_log_lkl(b, a, p, beta_binom);
}
}
for (t = 0; t < T; t++) rescale_emis_log_lkl(&emis_log_lkl[t * N], N, err_log_prb);
return emis_log_lkl;
}
static int cnp_edge_is_not_cn2_lrr_ad(const float *lrr, int16_t *ad0, int16_t *ad1, int n, int a, int b,
float xy_log_prb, float err_log_prb, float lrr_bias, float lrr_hap2dip,
float lrr_sd, float ad_rho, float ldev, float bdev) {
int n1, n2;
get_max_sum(ad0, ad1, n, NULL, &n1, &n2);
beta_binom_update(beta_binom_null, 0.5f, ad_rho, n1, n2);
beta_binom_update(beta_binom_alt, 0.5f + bdev, ad_rho, n1, n2);
// test left edge
float sum_log_lkl = 0.0f;
int i;
for (i = a - 1; i >= 0; i--) {
float log_lkl = lrr_ad_log_lkl(lrr[i], ad0[i], ad1[i], ldev, lrr_sd, lrr_bias, beta_binom_alt)
- lrr_ad_log_lkl(lrr[i], ad0[i], ad1[i], 0.0f, lrr_sd, lrr_bias, beta_binom_null);
if (!isnan(err_log_prb)) {
if (log_lkl < err_log_prb)
log_lkl = err_log_prb;
else if (log_lkl > -err_log_prb)
log_lkl = -err_log_prb;
}
sum_log_lkl += log_lkl;
if (sum_log_lkl > -xy_log_prb) return -1;
if (sum_log_lkl < xy_log_prb) break;
}
// test right edge
sum_log_lkl = 0.0f;
for (i = b + 1; i < n; i++) {
float log_lkl = lrr_ad_log_lkl(lrr[i], ad0[i], ad1[i], ldev, lrr_sd, lrr_bias, beta_binom_alt)
- lrr_ad_log_lkl(lrr[i], ad0[i], ad1[i], 0.0f, lrr_sd, lrr_bias, beta_binom_null);
if (!isnan(err_log_prb)) {
if (log_lkl < err_log_prb)
log_lkl = err_log_prb;
else if (log_lkl > -err_log_prb)
log_lkl = -err_log_prb;
}
sum_log_lkl += log_lkl;
if (sum_log_lkl > -xy_log_prb) return -1;
if (sum_log_lkl < xy_log_prb) break;
}
return 0;
}
// return the LOD likelihood for a segment
typedef struct {
const float *lrr_arr;
const int16_t *ad0_arr;
const int16_t *ad1_arr;
int n;
const int *imap;
float err_log_prb;
float lrr_bias;
float lrr_hap2dip;
float lrr_sd;
float ad_rho;
} minus_lrr_ad_lod_t;