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
use std::collections::{HashMap, HashSet};
use std::fmt::Display;
use std::hash::Hash;
use std::ops::Deref;
use std::str;
use anyhow::Result;
use bio::stats::{
bayesian::bayes_factors::evidence::KassRaftery, bayesian::bayes_factors::BayesFactor, LogProb,
PHREDProb, Prob,
};
use counter::Counter;
use half::f16;
use itertools::join;
use itertools::Itertools;
use ordered_float::NotNan;
use rust_htslib::bcf::Read;
use rust_htslib::{bam, bam::record::Cigar, bcf};
use crate::variants::model;
use crate::Event;
pub(crate) mod adaptive_integration;
pub(crate) mod anonymize;
pub(crate) mod collect_variants;
pub(crate) mod comparison;
pub(crate) mod homopolymers;
pub(crate) mod log2_fold_change;
pub(crate) mod variant_buffer;
pub(crate) use collect_variants::collect_variants;
pub(crate) const NUMERICAL_EPSILON: f64 = 1e-3;
lazy_static! {
pub(crate) static ref PROB_05: LogProb = LogProb::from(Prob(0.5f64));
pub(crate) static ref PROB_033: LogProb = LogProb::from(Prob(1.0 / 3.0));
pub(crate) static ref PROB_025: LogProb = LogProb::from(Prob(0.25));
pub(crate) static ref PROB_095: LogProb = LogProb::from(Prob(0.95));
pub(crate) static ref PROB_09: LogProb = LogProb::from(Prob(0.9));
}
pub(crate) fn aux_tag_strand_info(record: &bam::Record) -> Option<&[u8]> {
if let Ok(bam::record::Aux::String(strand_info)) = record.aux(b"SI") {
Some(strand_info.as_bytes())
} else {
None
}
}
pub(crate) fn is_sv_bcf(reader: &bcf::Reader) -> bool {
for rec in reader.header().header_records() {
if let bcf::header::HeaderRecord::Info { values, .. } = rec {
if values.get("ID").map_or(false, |id| id == "SVTYPE") {
return true;
}
}
}
false
}
pub fn is_bnd(record: &mut bcf::Record) -> Result<bool> {
Ok(record
.info(b"SVTYPE")
.string()?
.map_or(false, |entries| entries[0] == b"BND"))
}
pub(crate) fn info_tag_svtype(record: &mut bcf::Record) -> Result<Option<Vec<u8>>> {
Ok(record.info(b"SVTYPE").string()?.map(|v| v[0].to_owned()))
}
pub(crate) fn info_tag_event(record: &mut bcf::Record) -> Result<Option<Vec<u8>>> {
Ok(record.info(b"EVENT").string()?.map(|v| v[0].to_owned()))
}
pub(crate) fn info_tag_mateid(record: &mut bcf::Record) -> Result<Option<Vec<u8>>> {
Ok(record.info(b"MATEID").string()?.map(|v| v[0].to_owned()))
}
pub(crate) fn is_reverse_strand(record: &bam::Record) -> bool {
record.flags() & 0x10 != 0
}
pub(crate) fn contains_indel_op(record: &bam::Record) -> bool {
record
.cigar_cached()
.expect("bug: cigar accessed before caching")
.iter()
.any(|op| matches!(op, Cigar::Ins(_) | Cigar::Del(_)))
}
#[derive(new, Getters, CopyGetters, Debug)]
pub(crate) struct GenomicLocus {
#[getset(get = "pub")]
chrom: Vec<u8>,
#[getset(get_copy = "pub")]
pos: u32,
}
pub(crate) fn generalized_cigar<T: Hash + Eq + Clone + Display, F, K>(
items: impl Iterator<Item = T>,
keep_order: bool,
aux_sort: F,
) -> String
where
F: FnMut(&(T, usize)) -> K,
K: Ord,
{
if keep_order {
join(
items
.map(|item| (item, 1))
.coalesce(|(a, n), (b, m)| {
if a == b {
Ok((a, n + m))
} else {
Err(((a, n), (b, m)))
}
})
.map(|(item, count)| format!("{}{}", count, item)),
"",
)
} else {
let items: Counter<T> = items.collect();
join(
items
.most_common()
.into_iter()
.sorted_by_key(aux_sort)
.map(|(item, count)| format!("{}{}", count, item)),
"",
)
}
}
pub(crate) fn bayes_factor_to_letter(bayes_factor: BayesFactor) -> char {
match bayes_factor.evidence_kass_raftery() {
KassRaftery::Barely => 'B',
KassRaftery::None if relative_eq!(*bayes_factor, 1.0) => 'E',
KassRaftery::None => 'N',
KassRaftery::Positive => 'P',
KassRaftery::Strong => 'S',
KassRaftery::VeryStrong => 'V',
}
}
pub(crate) fn tags_prob_sum(
record: &mut bcf::Record,
tags: &[String],
vartype: Option<&model::VariantType>,
) -> Result<Vec<Option<LogProb>>> {
let mut skips = SimpleCounter::default();
let variants = collect_variants(record, false, Some(&mut skips))?;
let mut tags_probs_out = vec![Vec::new(); variants.len()];
for tag in tags {
if let Some(tags_probs_in) = (record.info(tag.as_bytes()).float())? {
for (i, (variant, tag_prob)) in variants.iter().zip(tags_probs_in.iter()).enumerate() {
if (vartype.is_some() && !variant.is_type(vartype.unwrap())) || tag_prob.is_nan() {
continue;
}
tags_probs_out[i].push(LogProb::from(PHREDProb(*tag_prob as f64)));
}
}
}
Ok(tags_probs_out
.into_iter()
.map(|probs| {
if !probs.is_empty() {
Some(LogProb::ln_sum_exp(&probs).cap_numerical_overshoot(NUMERICAL_EPSILON))
} else {
None
}
})
.collect_vec())
}
pub(crate) fn events_to_tags<E>(events: &[E]) -> Vec<String>
where
E: Event,
{
events.iter().map(|e| e.tag_name("PROB")).collect_vec()
}
pub(crate) fn collect_prob_dist<E>(
calls: &mut bcf::Reader,
events: &[E],
vartype: Option<&model::VariantType>,
) -> Result<Vec<NotNan<f64>>>
where
E: Event,
{
let mut visited_breakend_events = HashSet::new();
let mut record = calls.empty_record();
let mut prob_dist = Vec::new();
let tags = events_to_tags(events);
loop {
match calls.read(&mut record) {
None => break,
Some(res) => res?,
}
if let Ok(Some(event)) = info_tag_event(&mut record) {
if visited_breakend_events.contains(&event) {
continue;
} else {
visited_breakend_events.insert(event.to_owned());
}
}
for p in (tags_prob_sum(&mut record, &tags, vartype)?)
.into_iter()
.flatten()
{
prob_dist.push(NotNan::new(*p)?);
}
}
prob_dist.sort();
Ok(prob_dist)
}
pub(crate) fn filter_by_threshold<E: Event>(
calls: &mut bcf::Reader,
threshold: Option<LogProb>,
out: &mut bcf::Writer,
events: &[E],
vartype: Option<&model::VariantType>,
) -> Result<()> {
let mut breakend_event_decisions = HashMap::new();
let tags = events.iter().map(|e| e.tag_name("PROB")).collect_vec();
let filter = |record: &mut bcf::Record| -> Result<Vec<bool>> {
let bnd_event = info_tag_event(record).ok().flatten();
let keep = if let Some(event) = bnd_event.as_ref() {
breakend_event_decisions.get(event).cloned()
} else {
None
};
let probs = tags_prob_sum(record, &tags, vartype)?;
assert!(
bnd_event.is_none() || probs.len() == 1,
"breakend events may only occur in single variant records"
);
Ok(probs
.into_iter()
.map(|p| {
if let Some(keep) = keep {
keep
} else {
let keep = match (p, threshold) {
(Some(p), Some(threshold))
if p > threshold || relative_eq!(*p, *threshold) =>
{
true
}
(Some(_), None) => true,
_ => false,
};
if let Some(event) = bnd_event.as_ref() {
breakend_event_decisions.insert(event.to_owned(), keep);
}
keep
}
})
.collect())
};
filter_calls(calls, out, filter)
}
pub(crate) fn filter_calls<F, I, II>(
calls: &mut bcf::Reader,
out: &mut bcf::Writer,
mut filter: F,
) -> Result<()>
where
F: FnMut(&mut bcf::Record) -> Result<II>,
I: Iterator<Item = bool>,
II: IntoIterator<Item = bool, IntoIter = I>,
{
let mut record = calls.empty_record();
let mut i = 1;
loop {
match calls.read(&mut record) {
None => return Ok(()),
Some(res) => res?,
}
let mut remove = vec![false];
remove.extend(filter(&mut record)?.into_iter().map(|keep| !keep));
assert_eq!(
remove.len(),
record.allele_count() as usize,
"bug: filter passed to filter_calls has to return a bool for each alt allele at record {}.",
i,
);
if !remove[1..].iter().all(|r| *r) {
record.remove_alleles(&remove)?;
out.write(&record)?;
}
i += 1;
}
}
pub(crate) fn is_phred_scaled(inbcf: &bcf::Reader) -> bool {
get_event_tags(inbcf)
.iter()
.all(|(_, d)| d.ends_with("(PHRED)") || !d.ends_with(')'))
}
pub(crate) fn get_event_tags(inbcf: &bcf::Reader) -> Vec<(String, String)> {
inbcf
.header()
.header_records()
.into_iter()
.filter_map(|rec| {
if let bcf::header::HeaderRecord::Info { values, .. } = rec {
let id = values["ID"].clone();
if id.starts_with("PROB_") {
let description = values["Description"].clone();
let description = description.trim_matches('"').into();
return Some((id, description));
}
}
None
})
.collect_vec()
}
#[derive(Debug, Clone, Copy, Serialize, Deserialize)]
pub(crate) enum MiniLogProb {
F16(f16),
F32(f32),
}
impl MiniLogProb {
pub(crate) fn new(prob: LogProb) -> Self {
let half = f16::from_f64(*prob);
let proj = half.to_f64();
if *prob < -10.0 && proj.floor() as i64 == prob.floor() as i64 {
MiniLogProb::F16(half)
} else {
MiniLogProb::F32(*prob as f32)
}
}
pub(crate) fn to_logprob(&self) -> LogProb {
LogProb(match self {
MiniLogProb::F16(p) => p.to_f64(),
MiniLogProb::F32(p) => *p as f64,
})
}
}
#[cfg(test)]
mod tests {
use super::*;
use crate::variants::model::VariantType;
use crate::ComplementEvent;
use crate::SimpleEvent;
use bio::stats::{LogProb, Prob};
use rust_htslib::bcf::{self, Read};
#[test]
fn test_tags_prob_sum() {
let test_file = "tests/resources/test_tags_prob_sum/overshoot.vcf";
let mut overshoot_calls = bcf::Reader::from_path(test_file).unwrap();
let mut record = overshoot_calls.empty_record();
match overshoot_calls.read(&mut record) {
None => panic!("BCF reading error: seems empty"),
Some(res) => res.unwrap(),
}
let alt_tags = [
String::from("PROB_ADO_TO_REF"),
String::from("PROB_ADO_TO_ALT"),
String::from("PROB_HOM_ALT"),
String::from("PROB_HET"),
String::from("PROB_ERR_REF"),
];
let snv = VariantType::Snv;
if let Ok(prob_sum) = tags_prob_sum(&mut record, &alt_tags, Some(&snv)) {
assert_eq!(LogProb::ln_one(), prob_sum[0].unwrap());
} else {
panic!("tags_prob_sum(&overshoot_calls, &alt_events, &snv) returned Error")
}
}
#[test]
fn test_collect_prob_dist() {
let events = vec![
SimpleEvent {
name: "germline".to_owned(),
},
SimpleEvent {
name: "somatic".to_owned(),
},
];
let absent_event = vec![ComplementEvent {
name: "absent".to_owned(),
}];
let test_file = "tests/resources/test_collect_prob_dist/min.calls.vcf";
let del = VariantType::Deletion(None);
let mut del_calls_1 = bcf::Reader::from_path(test_file).unwrap();
let prob_del = collect_prob_dist(&mut del_calls_1, &events, Some(&del)).unwrap();
println!("prob_del[0]: {:?}", prob_del[0].into_inner());
assert_eq!(prob_del.len(), 1);
assert_relative_eq!(prob_del[0].into_inner(), Prob(0.8).ln(), epsilon = 0.000005);
let mut del_calls_2 = bcf::Reader::from_path(test_file).unwrap();
let prob_del_abs = collect_prob_dist(&mut del_calls_2, &absent_event, Some(&del)).unwrap();
assert_eq!(prob_del_abs.len(), 1);
assert_relative_eq!(
prob_del_abs[0].into_inner(),
Prob(0.2).ln(),
epsilon = 0.000005
);
let ins = VariantType::Insertion(None);
let mut ins_calls_1 = bcf::Reader::from_path(test_file).unwrap();
let prob_ins = collect_prob_dist(&mut ins_calls_1, &events, Some(&ins)).unwrap();
assert_eq!(prob_ins.len(), 1);
assert_relative_eq!(prob_ins[0].into_inner(), Prob(0.2).ln(), epsilon = 0.000005);
let mut ins_calls_2 = bcf::Reader::from_path(test_file).unwrap();
let prob_ins_abs = collect_prob_dist(&mut ins_calls_2, &absent_event, Some(&ins)).unwrap();
assert_eq!(prob_ins_abs.len(), 1);
assert_relative_eq!(
prob_ins_abs[0].into_inner(),
Prob(0.8).ln(),
epsilon = 0.000005
);
}
#[test]
fn test_filter_by_threshold() {
}
}
#[derive(CopyGetters)]
pub(crate) struct SimpleCounter<T>
where
T: Eq + Hash,
{
inner: HashMap<T, usize>,
#[getset(get_copy = "pub(crate)")]
total_count: usize,
}
impl<T> SimpleCounter<T>
where
T: Eq + Hash,
{
pub(crate) fn incr(&mut self, event: T) {
self.total_count += 1;
*self.inner.entry(event).or_insert(0) += 1;
}
}
impl<T> Default for SimpleCounter<T>
where
T: Eq + Hash,
{
fn default() -> Self {
SimpleCounter {
inner: HashMap::new(),
total_count: 0,
}
}
}
impl<T> Deref for SimpleCounter<T>
where
T: Eq + Hash,
{
type Target = HashMap<T, usize>;
fn deref(&self) -> &Self::Target {
&self.inner
}
}