-
Notifications
You must be signed in to change notification settings - Fork 2
/
af.py
876 lines (604 loc) · 23 KB
/
af.py
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
# -*- coding: utf-8 -*-
"""
Allele frequency calculations.
See also the examples at:
- http://nbviewer.ipython.org/github/alimanfoo/anhima/blob/master/examples/af.ipynb
""" # noqa
from __future__ import division, print_function, absolute_import
# third party dependencies
import numpy as np
# internal dependencies
import anhima
def _check_genotypes(genotypes):
"""
Internal function to check the genotypes input argument meets
expectations.
"""
genotypes = np.asarray(genotypes)
assert genotypes.ndim >= 2
if genotypes.ndim == 2:
# assume haploid, add ploidy dimension
genotypes = genotypes[..., np.newaxis]
return genotypes
def is_variant(genotypes):
"""Find variants with at least one non-reference allele observation.
Parameters
----------
genotypes : array_like
An array of shape (n_variants, n_samples, ploidy) where each
element of the array is an integer corresponding to an allele index
(-1 = missing, 0 = reference allele, 1 = first alternate allele,
2 = second alternate allele, etc.).
Returns
-------
is_variant : ndarray, bool
An array of shape (n_variants,) where an element is True if there
are at least `min_ac` non-reference alleles found for the corresponding
variant.
Notes
-----
Applicable to polyploid genotype calls.
Applicable to multiallelic variants.
"""
# check inputs
genotypes = _check_genotypes(genotypes)
# aggregate over samples and ploidy dimensions
out = np.sum(genotypes > 0, axis=(1, 2)) >= 1
return out
def count_variant(genotypes):
"""Count variants with at least one non-reference allele observed.
Parameters
----------
genotypes : array_like
An array of shape (n_variants, n_samples, ploidy) where each
element of the array is an integer corresponding to an allele index
(-1 = missing, 0 = reference allele, 1 = first alternate allele,
2 = second alternate allele, etc.).
min_ac : int, optional
The minimum number of non-reference alleles required to consider
variant.
Returns
-------
n : int
The number of variants.
Notes
-----
Applicable to polyploid genotype calls.
Applicable to multiallelic variants.
"""
return np.count_nonzero(is_variant(genotypes))
def is_non_variant(genotypes):
"""Find variants with no non-reference alleles.
Parameters
----------
genotypes : array_like
An array of shape (n_variants, n_samples, ploidy) where each
element of the array is an integer corresponding to an allele index
(-1 = missing, 0 = reference allele, 1 = first alternate allele,
2 = second alternate allele, etc.).
Returns
-------
is_non_variant : ndarray, bool
An array of shape (n_variants,) where an element is True if there
are no non-reference alleles found for the corresponding variant.
Notes
-----
Applicable to polyploid genotype calls.
Applicable to multiallelic variants.
"""
# check inputs
genotypes = _check_genotypes(genotypes)
# aggregate over samples and ploidy dimensions
out = np.all(genotypes <= 0, axis=(1, 2))
return out
def count_non_variant(genotypes):
"""Count variants with no non-reference alleles.
Parameters
----------
genotypes : array_like
An array of shape (n_variants, n_samples, ploidy) where each
element of the array is an integer corresponding to an allele index
(-1 = missing, 0 = reference allele, 1 = first alternate allele,
2 = second alternate allele, etc.).
Returns
-------
n : int
The number of variants.
Notes
-----
Applicable to polyploid genotype calls.
Applicable to multiallelic variants.
"""
return np.count_nonzero(is_non_variant(genotypes))
def is_singleton(genotypes, allele=1):
"""Find variants with only a single instance of `allele` observed.
Parameters
----------
genotypes : array_like
An array of shape (n_variants, n_samples, ploidy) where each
element of the array is an integer corresponding to an allele index
(-1 = missing, 0 = reference allele, 1 = first alternate allele,
2 = second alternate allele, etc.).
allele : int, optional
The allele to find singletons of.
Returns
-------
is_singleton : ndarray, bool
An array of shape (n_variants,) where an element is True if there
is a single instance of `allele` observed.
Notes
-----
Applicable to polyploid genotype calls.
Applicable to multiallelic variants, but note this function checks for a
specific `allele`.
"""
# check inputs
genotypes = _check_genotypes(genotypes)
# aggregate over samples and ploidy dimensions
out = np.sum(genotypes == allele, axis=(1, 2)) == 1
return out
def count_singletons(genotypes, allele=1):
"""Count variants with only a single instance of `allele` observed.
Parameters
----------
genotypes : array_like
An array of shape (n_variants, n_samples, ploidy) where each
element of the array is an integer corresponding to an allele index
(-1 = missing, 0 = reference allele, 1 = first alternate allele,
2 = second alternate allele, etc.).
allele : int, optional
The allele to find singletons of.
Returns
-------
n : int
The number of variants.
Notes
-----
Applicable to polyploid genotype calls.
Applicable to multiallelic variants, but note this function checks for a
specific `allele`.
"""
return np.count_nonzero(is_singleton(genotypes, allele))
def is_doubleton(genotypes, allele=1):
"""Find variants with only two instances of `allele` observed.
Parameters
----------
genotypes : array_like
An array of shape (n_variants, n_samples, ploidy) where each
element of the array is an integer corresponding to an allele index
(-1 = missing, 0 = reference allele, 1 = first alternate allele,
2 = second alternate allele, etc.).
allele : int, optional
The allele to find doubletons of.
Returns
-------
is_doubleton : ndarray, bool
An array of shape (n_variants,) where an element is True if there
are exactly two instances of `allele` observed.
Notes
-----
Applicable to polyploid genotype calls.
Applicable to multiallelic variants, but note this function checks for a
specific `allele`.
"""
# check inputs
genotypes = _check_genotypes(genotypes)
# aggregate over samples and ploidy dimensions
out = np.sum(genotypes == allele, axis=(1, 2)) == 2
return out
def count_doubletons(genotypes, allele=1):
"""Count variants with only two instances of `allele` observed.
Parameters
----------
genotypes : array_like
An array of shape (n_variants, n_samples, ploidy) where each
element of the array is an integer corresponding to an allele index
(-1 = missing, 0 = reference allele, 1 = first alternate allele,
2 = second alternate allele, etc.).
allele : int, optional
The allele to find doubletons of.
Returns
-------
n : int
The number of variants.
Notes
-----
Applicable to polyploid genotype calls.
Applicable to multiallelic variants, but note this function checks for a
specific `allele`.
"""
return np.count_nonzero(is_doubleton(genotypes, allele))
def allele_number(genotypes):
"""Count the number of non-missing allele calls per variant.
Parameters
----------
genotypes : array_like
An array of shape (n_variants, n_samples, ploidy) where each
element of the array is an integer corresponding to an allele index
(-1 = missing, 0 = reference allele, 1 = first alternate allele,
2 = second alternate allele, etc.).
Returns
-------
an : ndarray, int
An array of shape (n_variants,) counting the total number of
non-missing alleles observed.
Notes
-----
Applicable to polyploid genotype calls.
Applicable to multiallelic variants.
"""
# check inputs
genotypes = _check_genotypes(genotypes)
# aggregate over samples and ploidy dimensions
an = np.sum(genotypes >= 0, axis=(1, 2))
return an
def allele_count(genotypes, allele=1):
"""Calculate number of observations of the given allele per variant.
Parameters
----------
genotypes : array_like
An array of shape (n_variants, n_samples, ploidy) where each
element of the array is an integer corresponding to an allele index
(-1 = missing, 0 = reference allele, 1 = first alternate allele,
2 = second alternate allele, etc.).
allele : int, optional
The allele to count.
Returns
-------
ac : ndarray, int
An array of shape (n_variants,) counting the number of
times the given `allele` was observed.
Notes
-----
Applicable to polyploid genotype calls.
Applicable to multiallelic variants, but note that this function
calculates the frequency of a specific `allele`.
"""
# check inputs
genotypes = _check_genotypes(genotypes)
# aggregate over samples and ploidy dimensions
ac = np.sum(genotypes == allele, axis=(1, 2))
return ac
def allele_frequency(genotypes, allele=1):
"""Calculate frequency of the given allele per variant.
Parameters
----------
genotypes : array_like
An array of shape (n_variants, n_samples, ploidy) where each
element of the array is an integer corresponding to an allele index
(-1 = missing, 0 = reference allele, 1 = first alternate allele,
2 = second alternate allele, etc.).
allele : int, optional
The allele to calculate the frequency of.
Returns
-------
an : ndarray, int
An array of shape (n_variants,) counting the total number of
non-missing alleles observed.
ac : ndarray, int
An array of shape (n_variants,) counting the number of
times the given `allele` was observed.
af : ndarray, float
An array of shape (n_variants,) containing the allele frequency.
Notes
-----
Applicable to polyploid genotype calls.
Applicable to multiallelic variants, but note that this function
calculates the frequency of a specific `allele`.
"""
# count non-missing alleles
an = allele_number(genotypes)
# count alleles
ac = allele_count(genotypes, allele=allele)
# calculate allele frequency, accounting for missingness
err = np.seterr(invalid='ignore')
af = np.where(an > 0, ac / an, 0)
np.seterr(**err)
return an, ac, af
def allele_counts(genotypes, alleles=None):
"""Calculate allele counts per variant.
Parameters
----------
genotypes : array_like
An array of shape (n_variants, n_samples, ploidy) where each
element of the array is an integer corresponding to an allele index
(-1 = missing, 0 = reference allele, 1 = first alternate allele,
2 = second alternate allele, etc.).
alleles : sequence of ints, optional
The alleles to count. If not specified, all alleles will be counted.
Returns
-------
ac : ndarray, int
An array of shape (n_variants, n_alleles) counting the number of
times the given `alleles` were observed.
Notes
-----
Applicable to polyploid genotype calls.
Applicable to multiallelic variants.
"""
# check input
genotypes = _check_genotypes(genotypes)
# determine number of variants
n_variants = genotypes.shape[0]
# if alleles not specified, count all alleles
if alleles is None:
m = np.amax(genotypes)
alleles = range(m+1)
# count alleles
ac = np.zeros((n_variants, len(alleles)), dtype='i4')
for i, allele in enumerate(alleles):
np.sum(genotypes == allele, axis=(1, 2), out=ac[:, i])
return ac
def allele_frequencies(genotypes, alleles=None):
"""Calculate allele frequencies per variant.
Parameters
----------
genotypes : array_like
An array of shape (n_variants, n_samples, ploidy) where each
element of the array is an integer corresponding to an allele index
(-1 = missing, 0 = reference allele, 1 = first alternate allele,
2 = second alternate allele, etc.).
alleles : sequence of ints, optional
The alleles to calculate the frequency of. If not specified, all
alleles will be counted.
Returns
-------
an : ndarray, int
An array of shape (n_variants,) counting the total number of
non-missing alleles observed.
ac : ndarray, int
An array of shape (n_variants, n_alleles) counting the number of
times the given `alleles` were observed.
af : ndarray, float
An array of shape (n_variants, n_alleles) containing the allele
frequencies.
Notes
-----
Applicable to polyploid genotype calls.
Applicable to multiallelic variants.
"""
# count non-missing alleles
an = allele_number(genotypes)[:, np.newaxis]
# count alleles
ac = allele_counts(genotypes, alleles=alleles)
# calculate allele frequencies, accounting for missingness
err = np.seterr(invalid='ignore')
af = np.where(an > 0, ac / an, 0)
np.seterr(**err)
return an, ac, af
def allelism(genotypes):
"""Determine the number of distinct alleles found for each variant.
Parameters
----------
genotypes : array_like
An array of shape (n_variants, n_samples, ploidy) where each
element of the array is an integer corresponding to an allele index
(-1 = missing, 0 = reference allele, 1 = first alternate allele,
2 = second alternate allele, etc.).
Returns
-------
n : ndarray, int
An array of shape (n_variants,) where an element holds the allelism
of the corresponding variant.
See Also
--------
max_allele
Notes
-----
Applicable to polyploid genotype calls.
Applicable to multiallelic variants.
"""
# check inputs
genotypes = _check_genotypes(genotypes)
# calculate allele counts
ac = allele_counts(genotypes)
# count alleles present
n = np.sum(ac > 0, axis=1)
return n
def is_non_segregating(genotypes, allele=None):
"""Find non-segregating variants (fixed for a single allele).
Parameters
----------
genotypes : array_like
An array of shape (n_variants, n_samples, ploidy) where each
element of the array is an integer corresponding to an allele index
(-1 = missing, 0 = reference allele, 1 = first alternate allele,
2 = second alternate allele, etc.).
allele : int, optional
If given, find variants fixed with respect to `allele`. Otherwise
find variants fixed for any allele.
Returns
-------
is_non_segregating : ndarray, bool
An array of shape (n_variants,) where an element is True if all
genotype calls for the corresponding variant are either missing or
equal to the same allele.
Notes
-----
Applicable to polyploid genotype calls.
Applicable to multiallelic variants.
"""
# check inputs
genotypes = _check_genotypes(genotypes)
if allele is None:
# count distinct alleles
n_alleles = allelism(genotypes)
# find fixed variants
out = n_alleles == 1
else:
# find fixed variants with respect to a specific allele
out = np.all((genotypes < 0) | (genotypes == allele), axis=(1, 2))
return out
def count_non_segregating(genotypes, allele=None):
"""Count non-segregating variants.
Parameters
----------
genotypes : array_like
An array of shape (n_variants, n_samples, ploidy) where each
element of the array is an integer corresponding to an allele index
(-1 = missing, 0 = reference allele, 1 = first alternate allele,
2 = second alternate allele, etc.).
allele : int, optional
If given, find variants fixed with respect to `allele`.
Returns
-------
n : int
The number of variants.
Notes
-----
Applicable to polyploid genotype calls.
Applicable to multiallelic variants.
"""
return np.count_nonzero(is_non_segregating(genotypes, allele=allele))
def is_segregating(genotypes):
"""Find segregating variants (where more than one allele is observed).
Parameters
----------
genotypes : array_like
An array of shape (n_variants, n_samples, ploidy) where each
element of the array is an integer corresponding to an allele index
(-1 = missing, 0 = reference allele, 1 = first alternate allele,
2 = second alternate allele, etc.).
Returns
-------
is_segregating : ndarray, bool
An array of shape (n_variants,) where an element is True if more
than one allele is found for the given variant.
Notes
-----
Applicable to polyploid genotype calls.
Applicable to multiallelic variants.
"""
# check inputs
genotypes = _check_genotypes(genotypes)
# count distinct alleles
n_alleles = allelism(genotypes)
# find segregating variants
out = n_alleles > 1
return out
def count_segregating(genotypes):
"""Count segregating variants.
Parameters
----------
genotypes : array_like
An array of shape (n_variants, n_samples, ploidy) where each
element of the array is an integer corresponding to an allele index
(-1 = missing, 0 = reference allele, 1 = first alternate allele,
2 = second alternate allele, etc.).
Returns
-------
n : int
The number of variants.
Notes
-----
Applicable to polyploid genotype calls.
Applicable to multiallelic variants.
"""
return np.count_nonzero(is_segregating(genotypes))
def maximum_likelihood_ancestry(genotypes, qa, qb, filter_size=0):
"""Given alternate allele frequencies in two populations `qa` and `qb`,
predict the ancestry for a set of `genotypes`.
Parameters
----------
genotypes : array_like
An array of diploid genotype calls of shape (n_variants, n_samples,
2) where each element of the array is an integer corresponding to an
allele index (-1 = missing, 0 = reference allele, 1 = first alternate
allele, 2 = second alternate allele, etc.).
qa : array_like, float
A 1-dimensional array of shape (n_variants, ) containing alternate
allele frequencies for population A.
qb : array_like, float
A 1-dimensional array of shape (n_variants, ) containing alternate
allele frequencies for population B.
filter_size : int, optional
Sum likelihoods in a moving window of size `filter_size`.
Returns
-------
ancestry : ndarray, int, shape (n_variants, n_samples)
An array containing the ancestry predictions, where 0 = AA (both
alleles derive from population A), 1 = AB (hybrid ancestry) and 2 =
BB (both alleles derive from population B), and -1 = ambiguous (models
are equally likely).
confidence : ndarray, float, shape (n_variants, n_samples)
The confidence in the ancestry prediction (natural logarithm of the
likelihood ratio for the two most likely models).
Notes
-----
Where allele frequencies are similar between populations A and B,
ancestry predictions will have low confidence, because different ancestry
models will have similar likelihoods. Greater confidence will be obtained by
filtering variants to select those where the difference in allele
frequencies is greater. E.g.::
>>> flt = np.abs(qa - qb) > .5
>>> genotypes_flt = genotypes[flt]
>>> qa_flt = qa[flt]
>>> qb_flt = qb[flt]
>>> ancestry, confidence = maximum_likelihood_ancestry(genotypes_flt, qa_flt, qb_flt)
""" # noqa
# check inputs
genotypes = _check_genotypes(genotypes)
# require biallelic genotypes
assert np.amax(genotypes) < 2
n_variants, n_samples, ploidy = genotypes.shape
# require diploid genotypes
assert ploidy == 2
qa = np.asarray(qa)
qb = np.asarray(qb)
assert qa.ndim == qb.ndim == 1
assert n_variants == qa.shape[0] == qb.shape[0]
# calculate reference allele frequencies, assuming biallelic variants
pa = 1 - qa
pb = 1 - qb
# work around zero frequencies which cause problems when calculating logs
pa[pa == 0] = np.exp(-250)
qa[qa == 0] = np.exp(-250)
pb[pb == 0] = np.exp(-250)
qb[qb == 0] = np.exp(-250)
# calculate likelihoods
logpa = np.log(pa)
logqa = np.log(qa)
logpb = np.log(pb)
logqb = np.log(qb)
# set up likelihoods array
n_models = 3
n_gn_states = 3
log_likelihoods = np.empty((n_variants, n_samples, n_models, n_gn_states),
dtype='f8')
# probability of genotype (e.g., 0 = hom ref) given model (e.g., 0 = aa)
log_likelihoods[:, :, 0, 0] = (2 * logpa)[:, np.newaxis]
log_likelihoods[:, :, 1, 0] = (np.log(2) + logpa + logqa)[:, np.newaxis]
log_likelihoods[:, :, 2, 0] = (2 * logqa)[:, np.newaxis]
log_likelihoods[:, :, 0, 1] = (logpa + logpb)[:, np.newaxis]
log_likelihoods[:, :, 1, 1] = (np.logaddexp(logpa + logqb,
logqa + logpb)[:, np.newaxis])
log_likelihoods[:, :, 2, 1] = (logqa + logqb)[:, np.newaxis]
log_likelihoods[:, :, 0, 2] = (2 * logpb)[:, np.newaxis]
log_likelihoods[:, :, 1, 2] = (np.log(2) + logpb + logqb)[:, np.newaxis]
log_likelihoods[:, :, 2, 2] = (2 * logqb)[:, np.newaxis]
# transform genotypes for convenience
gn = anhima.gt.as_012(genotypes)
# calculate actual model likelihoods for each genotype call
model_likelihoods = np.empty((n_variants, n_samples, n_models), dtype='f8')
model_likelihoods.fill(-250)
for model in 0, 1, 2:
for gn_state in 0, 1, 2:
model_likelihoods[:, :, model][gn == gn_state] = \
log_likelihoods[:, :, model, gn_state][gn == gn_state]
# optionally combine likelihoods in a moving window
if filter_size:
model_likelihoods = np.apply_along_axis(np.convolve,
0,
model_likelihoods,
np.ones((filter_size,)))
# remove edges
model_likelihoods = \
model_likelihoods[filter_size//2:-1*(filter_size//2), ...]
# predict ancestry as model with highest likelihood
ancestry = np.argmax(model_likelihoods, axis=2)
# calculate confidence by comparing first and second most likely models
model_likelihoods.sort(axis=2)
confidence = model_likelihoods[:, :, 2] - model_likelihoods[:, :, 1]
# recind prediction where confidence is zero (models are equally likely)
ancestry[confidence == 0] = -1
return ancestry, confidence