-
Notifications
You must be signed in to change notification settings - Fork 22
/
main.R
2014 lines (1706 loc) · 69.3 KB
/
main.R
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
#' @import logger
#' @import dplyr
#' @import Matrix
#' @importFrom data.table fread fwrite as.data.table
#' @import stringr
#' @import glue
#' @importFrom parallel mclapply
#' @import tidygraph
#' @import ggplot2
#' @import ggraph
#' @importFrom scistreer perform_nni get_mut_graph score_tree annotate_tree mut_to_tree ladderize to_phylo
#' @importFrom ggtree %<+%
#' @importFrom methods is as
#' @importFrom igraph vcount ecount E V V<- E<-
#' @import patchwork
#' @importFrom grDevices colorRampPalette
#' @importFrom stats as.dendrogram as.dist cor cutree dbinom dnbinom dnorm dpois end hclust integrate model.matrix na.omit optim p.adjust pnorm reorder rnorm setNames start t.test as.ts complete.cases is.leaf na.contiguous
#' @importFrom hahmmr logSumExp dpoilog dbbinom l_lnpois l_bbinom fit_lnpois_cpp likelihood_allele forward_back_allele run_joint_hmm_s15 run_allele_hmm_s5
#' @import tibble
#' @importFrom utils combn
#' @useDynLib numbat, .registration=TRUE
NULL
#' Run workflow to decompose tumor subclones
#'
#' @param count_mat dgCMatrix Raw count matrices where rownames are genes and column names are cells
#' @param lambdas_ref matrix Either a named vector with gene names as names and normalized expression as values, or a matrix where rownames are genes and columns are pseudobulk names
#' @param df_allele dataframe Allele counts per cell, produced by preprocess_allele
#' @param genome character Genome version (hg38, hg19, or mm10)
#' @param out_dir string Output directory
#' @param gamma numeric Dispersion parameter for the Beta-Binomial allele model
#' @param t numeric Transition probability
#' @param init_k integer Number of clusters in the initial clustering
#' @param min_cells integer Minimum number of cells to run HMM on
#' @param min_genes integer Minimum number of genes to call a segment
#' @param max_cost numeric Likelihood threshold to collapse internal branches
#' @param n_cut integer Number of cuts on the phylogeny to define subclones
#' @param tau numeric Factor to determine max_cost as a function of the number of cells (0-1)
#' @param nu numeric Phase switch rate
#' @param alpha numeric P value cutoff for diploid finding
#' @param min_overlap numeric Minimum CNV overlap threshold
#' @param max_iter integer Maximum number of iterations to run the phyologeny optimization
#' @param max_nni integer Maximum number of iterations to run NNI in the ML phylogeny inference
#' @param min_depth integer Minimum allele depth
#' @param common_diploid logical Whether to find common diploid regions in a group of peusdobulks
#' @param ncores integer Number of threads to use
#' @param ncores_nni integer Number of threads to use for NNI
#' @param max_entropy numeric Entropy threshold to filter CNVs
#' @param min_LLR numeric Minimum LLR to filter CNVs
#' @param eps numeric Convergence threshold for ML tree search
#' @param multi_allelic logical Whether to call multi-allelic CNVs
#' @param p_multi numeric P value cutoff for calling multi-allelic CNVs
#' @param use_loh logical Whether to include LOH regions in the expression baseline
#' @param skip_nj logical Whether to skip NJ tree construction and only use UPGMA
#' @param diploid_chroms vector Known diploid chromosomes
#' @param segs_loh dataframe Segments of clonal LOH to be excluded
#' @param call_clonal_loh logical Whether to call segments with clonal LOH
#' @param segs_consensus_fix dataframe Pre-determined segmentation of consensus CNVs
#' @param check_convergence logical Whether to terminate iterations based on consensus CNV convergence
#' @param random_init logical Whether to initiate phylogney using a random tree (internal use only)
#' @param exclude_neu logical Whether to exclude neutral segments from CNV retesting (internal use only)
#' @param plot logical Whether to plot results
#' @param verbose logical Verbosity
#' @return a status code
#' @export
run_numbat = function(
count_mat, lambdas_ref, df_allele, genome = 'hg38',
out_dir = tempdir(), max_iter = 2, max_nni = 100, t = 1e-5, gamma = 20, min_LLR = 5,
alpha = 1e-4, eps = 1e-5, max_entropy = 0.5, init_k = 3, min_cells = 50, tau = 0.3, nu = 1,
max_cost = ncol(count_mat) * tau, n_cut = 0, min_depth = 0, common_diploid = TRUE, min_overlap = 0.45,
ncores = 1, ncores_nni = ncores, random_init = FALSE, segs_loh = NULL, call_clonal_loh = FALSE,
verbose = TRUE, diploid_chroms = NULL, segs_consensus_fix = NULL, use_loh = NULL, min_genes = 10,
skip_nj = FALSE, multi_allelic = TRUE, p_multi = 1-alpha,
plot = TRUE, check_convergence = FALSE, exclude_neu = TRUE
) {
######### Setup output folder #########
dir.create(out_dir, showWarnings = FALSE, recursive = TRUE)
logfile = glue('{out_dir}/log.txt')
if (file.exists(logfile)) {file.remove(logfile)}
log_appender(appender_file(logfile))
######### Basic checks #########
if (genome == 'hg38') {
gtf = gtf_hg38
} else if (genome == 'hg19') {
gtf = gtf_hg19
} else if (genome == 'mm10') {
gtf = gtf_mm10
} else {
stop('Genome version must be hg38, hg19, or mm10')
}
count_mat = check_matrix(count_mat)
df_allele = annotate_genes(df_allele, gtf)
df_allele = check_allele_df(df_allele)
lambdas_ref = check_exp_ref(lambdas_ref)
# filter for annotated genes
genes_annotated = unique(gtf$gene) %>%
intersect(rownames(count_mat)) %>%
intersect(rownames(lambdas_ref))
count_mat = count_mat[genes_annotated,,drop=FALSE]
lambdas_ref = lambdas_ref[genes_annotated,,drop=FALSE]
zero_cov = names(which(colSums(count_mat) == 0))
if (length(zero_cov) > 0) {
log_message(glue('Filtering out {length(zero_cov)} cells with 0 coverage'))
count_mat = count_mat[,!colnames(count_mat) %in% zero_cov]
df_allele = df_allele %>% filter(!cell %in% zero_cov)
}
# only keep cells that have a transcriptome
df_allele = df_allele %>% filter(cell %in% colnames(count_mat))
if (nrow(df_allele) == 0){
stop('No matching cell names between count_mat and df_allele')
}
if ((!is.null(segs_loh)) & (!is.null(segs_consensus_fix))) {
stop('Cannot specify both segs_loh and segs_consensus_fix')
}
# check provided consensus CNVs
segs_consensus_fix = check_segs_fix(segs_consensus_fix)
# check clonal LOH
if (!is.null(segs_loh)) {
if (call_clonal_loh) {
stop('Cannot specify both segs_loh and call_clonal_loh')
}
segs_loh = check_segs_loh(segs_loh)
}
######### Log parameters #########
log_message(paste('\n',
glue('numbat version: ', as.character(utils::packageVersion("numbat"))),
glue('scistreer version: ', as.character(utils::packageVersion("scistreer"))),
glue('hahmmr version: ', as.character(utils::packageVersion("hahmmr"))),
'Running under parameters:',
glue('t = {t}'),
glue('alpha = {alpha}'),
glue('gamma = {gamma}'),
glue('min_cells = {min_cells}'),
glue('init_k = {init_k}'),
glue('max_cost = {max_cost}'),
glue('n_cut = {n_cut}'),
glue('max_iter = {max_iter}'),
glue('max_nni = {max_nni}'),
glue('min_depth = {min_depth}'),
glue('use_loh = {ifelse(is.null(use_loh), "auto", use_loh)}'),
glue('segs_loh = {ifelse(is.null(segs_loh), "None", "Given")}'),
glue('call_clonal_loh = {call_clonal_loh}'),
glue('segs_consensus_fix = {ifelse(is.null(segs_consensus_fix), "None", "Given")}'),
glue('multi_allelic = {multi_allelic}'),
glue('min_LLR = {min_LLR}'),
glue('min_overlap = {min_overlap}'),
glue('max_entropy = {max_entropy}'),
glue('skip_nj = {skip_nj}'),
glue('diploid_chroms = {ifelse(is.null(diploid_chroms), "None", "Given")}'),
glue('ncores = {ncores}'),
glue('ncores_nni = {ncores_nni}'),
glue('common_diploid = {common_diploid}'),
glue('tau = {tau}'),
glue('check_convergence = {check_convergence}'),
glue('plot = {plot}'),
glue('genome = {genome}'),
'Input metrics:',
glue('{ncol(count_mat)} cells'),
sep = "\n"
), verbose = verbose)
log_mem()
RhpcBLASctl::blas_set_num_threads(1)
RhpcBLASctl::omp_set_num_threads(1)
data.table::setDTthreads(1)
######## Initialization ########
if (call_clonal_loh) {
log_message('Calling segments with clonal LOH')
bulk = get_bulk(
count_mat = count_mat,
lambdas_ref = lambdas_ref,
df_allele = df_allele,
gtf = gtf,
min_depth = min_depth,
nu = nu
)
segs_loh = bulk %>% detect_clonal_loh(t = t, min_depth = min_depth)
if (!is.null(segs_loh)) {
fwrite(segs_loh, glue('{out_dir}/segs_loh.tsv'), sep = '\t')
} else {
log_message('No segments with clonal LOH detected')
}
}
sc_refs = choose_ref_cor(count_mat, lambdas_ref, gtf)
saveRDS(sc_refs, glue('{out_dir}/sc_refs.rds'))
if (random_init) {
log_message('Initializing with random tree...')
n_cells = ncol(count_mat)
dist_mat = matrix(rnorm(n_cells^2),nrow=n_cells)
colnames(dist_mat) = colnames(count_mat)
hc = hclust(as.dist(dist_mat), method = "ward.D2")
saveRDS(hc, glue('{out_dir}/hc.rds'))
# extract cell groupings
subtrees = get_nodes_celltree(hc, cutree(hc, k = init_k))
} else if (init_k == 1) {
log_message('Initializing with all-cell pseudobulk ..', verbose = verbose)
log_mem()
subtrees = list(list('cells' = colnames(count_mat), 'size' = ncol(count_mat), 'sample' = 1))
} else {
log_message('Approximating initial clusters using smoothed expression ..', verbose = verbose)
log_mem()
clust = exp_hclust(
count_mat = count_mat,
lambdas_ref = lambdas_ref,
gtf = gtf,
sc_refs = sc_refs,
ncores = ncores
)
hc = clust$hc
fwrite(
as.data.frame(clust$gexp_roll_wide) %>% tibble::rownames_to_column('cell'),
glue('{out_dir}/gexp_roll_wide.tsv.gz'),
sep = '\t',
nThread = min(4, ncores)
)
saveRDS(hc, glue('{out_dir}/hc.rds'))
# extract cell groupings
subtrees = get_nodes_celltree(hc, cutree(hc, k = init_k))
if (plot) {
p = plot_exp_roll(
gexp_roll_wide = clust$gexp_roll_wide,
hc = hc,
k = init_k,
gtf = gtf,
n_sample = 1e4
)
ggsave(glue('{out_dir}/exp_roll_clust.png'), p, width = 8, height = 4, dpi = 200)
}
}
clones = purrr::keep(subtrees, function(x) x$sample %in% 1:init_k)
normal_cells = c()
segs_consensus_old = data.frame()
######## Begin iterations ########
for (i in 1:max_iter) {
log_message(glue('Iteration {i}'), verbose = verbose)
log_mem()
subtrees = purrr::keep(subtrees, function(x) x$size > min_cells)
bulk_subtrees = make_group_bulks(
groups = subtrees,
count_mat = count_mat,
df_allele = df_allele,
lambdas_ref = lambdas_ref,
gtf = gtf,
min_depth = min_depth,
nu = nu,
segs_loh = segs_loh,
ncores = ncores)
# diagnostics
if (i == 1) {
bulk_subtrees %>% filter(sample == 0) %>% check_contam()
bulk_subtrees %>% filter(sample == 0) %>% check_exp_noise()
}
if (is.null(segs_consensus_fix)) {
bulk_subtrees = bulk_subtrees %>%
run_group_hmms(
t = t,
gamma = gamma,
alpha = alpha,
nu = nu,
min_genes = min_genes,
common_diploid = common_diploid,
diploid_chroms = diploid_chroms,
ncores = ncores,
verbose = verbose)
fwrite(bulk_subtrees, glue('{out_dir}/bulk_subtrees_{i}.tsv.gz'), sep = '\t')
if (plot) {
p = plot_bulks(bulk_subtrees, min_LLR = min_LLR, use_pos = TRUE, genome = genome)
ggsave(
glue('{out_dir}/bulk_subtrees_{i}.png'), p,
width = 13, height = 2*length(unique(bulk_subtrees$sample)), dpi = 250
)
}
# define consensus CNVs
segs_consensus = bulk_subtrees %>%
get_segs_consensus(min_LLR = min_LLR, min_overlap = min_overlap, retest = TRUE)
# check termination
if (all(segs_consensus$cnv_state_post == 'neu')) {
msg = 'No CNV remains after filtering by LLR in pseudobulks. Consider reducing min_LLR.'
log_message(msg)
return(msg)
}
# retest all segments on subtrees
bulk_subtrees = retest_bulks(
bulk_subtrees,
segs_consensus,
diploid_chroms = diploid_chroms,
gamma = gamma,
min_LLR = min_LLR,
ncores = ncores
)
fwrite(bulk_subtrees, glue('{out_dir}/bulk_subtrees_retest_{i}.tsv.gz'), sep = '\t')
# find consensus CNVs again
segs_consensus = bulk_subtrees %>%
get_segs_consensus(min_LLR = min_LLR, min_overlap = min_overlap, retest = FALSE)
# check termination again
if (all(segs_consensus$cnv_state_post == 'neu')) {
msg = 'No CNV remains after filtering by LLR in pseudobulks. Consider reducing min_LLR.'
log_message(msg)
return(msg)
}
} else {
log_message('Using fixed consensus CNVs')
segs_consensus = segs_consensus_fix
bulk_subtrees = bulk_subtrees %>%
annot_consensus(segs_consensus) %>%
annot_theta_mle() %>%
classify_alleles()
}
# retest on clones
clones = purrr::keep(clones, function(x) x$size > min_cells)
if (length(clones) == 0) {
msg = 'No clones remain after filtering by size. Consider reducing min_cells.'
log_message(msg)
return(msg)
}
bulk_clones = make_group_bulks(
groups = clones,
count_mat = count_mat,
df_allele = df_allele,
lambdas_ref = lambdas_ref,
gtf = gtf,
min_depth = min_depth,
nu = nu,
segs_loh = segs_loh,
ncores = ncores)
bulk_clones = bulk_clones %>%
run_group_hmms(
t = t,
gamma = gamma,
alpha = alpha,
nu = nu,
min_genes = min_genes,
common_diploid = common_diploid,
diploid_chroms = diploid_chroms,
ncores = ncores,
verbose = verbose,
retest = FALSE)
bulk_clones = retest_bulks(
bulk_clones,
segs_consensus,
gamma = gamma,
use_loh = use_loh,
min_LLR = min_LLR,
diploid_chroms = diploid_chroms,
ncores = ncores)
fwrite(bulk_clones, glue('{out_dir}/bulk_clones_{i}.tsv.gz'), sep = '\t')
if (plot) {
p = plot_bulks(bulk_clones, min_LLR = min_LLR, use_pos = TRUE, genome = genome)
ggsave(
glue('{out_dir}/bulk_clones_{i}.png'), p,
width = 13, height = 2*length(unique(bulk_clones$sample)), dpi = 250
)
}
# test for multi-allelic CNVs
if (multi_allelic) {
segs_consensus = test_multi_allelic(bulk_clones, segs_consensus, min_LLR = min_LLR, p_min = p_multi)
}
fwrite(segs_consensus, glue('{out_dir}/segs_consensus_{i}.tsv'), sep = '\t')
######## Evaluate CNV per cell ########
log_message('Evaluating CNV per cell ..', verbose = verbose)
log_mem()
exp_post = get_exp_post(
segs_consensus %>% mutate(cnv_state = ifelse(cnv_state == 'neu', cnv_state, cnv_state_post)),
count_mat,
lambdas_ref,
use_loh = use_loh,
segs_loh = segs_loh,
gtf = gtf,
sc_refs = sc_refs,
ncores = ncores)
haplotype_post = get_haplotype_post(
bulk_subtrees,
segs_consensus %>% mutate(cnv_state = ifelse(cnv_state == 'neu', cnv_state, cnv_state_post))
)
allele_post = get_allele_post(
df_allele,
haplotype_post,
segs_consensus %>% mutate(cnv_state = ifelse(cnv_state == 'neu', cnv_state, cnv_state_post))
)
joint_post = get_joint_post(
exp_post,
allele_post,
segs_consensus)
joint_post = joint_post %>%
group_by(seg) %>%
mutate(
avg_entropy = mean(binary_entropy(p_cnv), na.rm = TRUE)
) %>%
ungroup()
if (multi_allelic) {
log_message('Expanding allelic states..', verbose = verbose)
exp_post = expand_states(exp_post, segs_consensus)
allele_post = expand_states(allele_post, segs_consensus)
joint_post = expand_states(joint_post, segs_consensus)
}
fwrite(exp_post, glue('{out_dir}/exp_post_{i}.tsv'), sep = '\t')
fwrite(allele_post, glue('{out_dir}/allele_post_{i}.tsv'), sep = '\t')
fwrite(joint_post, glue('{out_dir}/joint_post_{i}.tsv'), sep = '\t')
######## Build phylogeny ########
log_message('Building phylogeny ..', verbose = verbose)
log_mem()
# filter CNVs
joint_post_filtered = joint_post %>%
filter(cnv_state != 'neu') %>%
filter(avg_entropy < max_entropy & LLR > min_LLR)
if (nrow(joint_post_filtered) == 0) {
msg = 'No CNV remains after filtering by entropy in single cells. Consider increasing max_entropy.'
log_message(msg)
return(msg)
} else {
n_cnv = length(unique(joint_post_filtered$seg))
log_message(glue('Using {n_cnv} CNVs to construct phylogeny'), verbose = verbose)
}
# construct genotype probability matrix
p_min = 1e-10
P = joint_post_filtered %>%
mutate(p_cnv = pmax(pmin(p_cnv, 1-p_min), p_min)) %>%
as.data.table %>%
data.table::dcast(cell ~ seg, value.var = 'p_cnv', fill = 0.5) %>%
tibble::column_to_rownames('cell') %>%
as.matrix
fwrite(
as.data.frame(P) %>% tibble::rownames_to_column('cell'),
glue('{out_dir}/geno_{i}.tsv'),
sep = '\t'
)
# contruct initial tree
dist_mat = parallelDist::parDist(rbind(P, 'outgroup' = 0), threads = ncores)
treeUPGMA = upgma(dist_mat) %>%
ape::root(outgroup = 'outgroup') %>%
ape::drop.tip('outgroup') %>%
reorder(order = 'postorder')
saveRDS(treeUPGMA, glue('{out_dir}/treeUPGMA_{i}.rds'))
UPGMA_score = score_tree(treeUPGMA, as.matrix(P))$l_tree
tree_init = treeUPGMA
# note that dist_mat gets modified by NJ
if(!skip_nj){
treeNJ = ape::nj(dist_mat) %>%
ape::root(outgroup = 'outgroup') %>%
ape::drop.tip('outgroup') %>%
reorder(order = 'postorder')
NJ_score = score_tree(treeNJ, as.matrix(P))$l_tree
if (UPGMA_score > NJ_score) {
log_message('Using UPGMA tree as seed..', verbose = verbose)
} else {
tree_init = treeNJ
log_message('Using NJ tree as seed..', verbose = verbose)
}
saveRDS(treeNJ, glue('{out_dir}/treeNJ_{i}.rds'))
} else {
log_message('Only computing UPGMA..', verbose = verbose)
log_message('Using UPGMA tree as seed..', verbose = verbose)
}
log_mem()
# maximum likelihood tree search with NNI
tree_list = perform_nni(tree_init, P, ncores = ncores_nni, eps = eps, max_iter = max_nni)
saveRDS(tree_list, glue('{out_dir}/tree_list_{i}.rds'))
treeML = tree_list[[length(tree_list)]]
saveRDS(treeML, glue('{out_dir}/treeML_{i}.rds'))
gtree = get_gtree(treeML, P, n_cut = n_cut, max_cost = max_cost)
saveRDS(gtree, glue('{out_dir}/tree_final_{i}.rds'))
G_m = get_mut_graph(gtree) %>% label_genotype()
saveRDS(G_m, glue('{out_dir}/mut_graph_{i}.rds'))
clone_post = get_clone_post(gtree, exp_post, allele_post)
fwrite(clone_post, glue('{out_dir}/clone_post_{i}.tsv'), sep = '\t')
normal_cells = clone_post %>% filter(p_cnv < 0.5) %>% pull(cell)
log_message(glue('Found {length(normal_cells)} normal cells..'), verbose = verbose)
if (plot) {
panel = plot_phylo_heatmap(
gtree,
joint_post,
segs_consensus,
clone_post,
tip_length = 0.2,
branch_width = 0.2,
line_width = 0.1,
clone_bar = TRUE
)
ggsave(glue('{out_dir}/panel_{i}.png'), panel, width = 7.5, height = 3.75, dpi = 250)
}
# form cell groupings using the obtained phylogeny
clone_to_node = setNames(V(G_m)$id, V(G_m)$clone)
subtrees = lapply(1:vcount(G_m), function(c) {
nodes = na.omit(igraph::dfs(G_m, root = c, unreachable = F)$order)
G_m %>%
igraph::as_data_frame('vertices') %>%
filter(id %in% nodes) %>%
inner_join(clone_post, by = c('GT' = 'GT_opt')) %>%
{list(sample = c, members = unique(.$GT), clones = unique(.$clone), cells = .$cell, size = length(.$cell))}
})
saveRDS(subtrees, glue('{out_dir}/subtrees_{i}.rds'))
clones = clone_post %>% split(.$clone_opt) %>%
purrr::map(function(c){list(sample = unique(c$clone_opt), members = unique(c$GT_opt), cells = c$cell, size = length(c$cell))})
saveRDS(clones, glue('{out_dir}/clones_{i}.rds'))
#### check convergence ####
if (check_convergence) {
converge = segs_equal(segs_consensus_old, segs_consensus)
if (converge) {
log_message('Convergence reached')
break
} else {
segs_consensus_old = segs_consensus
}
}
}
# Output final subclone bulk profiles - logic can be simplified
bulk_clones = make_group_bulks(
groups = clones,
count_mat = count_mat,
df_allele = df_allele,
lambdas_ref = lambdas_ref,
gtf = gtf,
min_depth = min_depth,
nu = nu,
segs_loh = segs_loh,
ncores = ncores)
bulk_clones = bulk_clones %>%
run_group_hmms(
t = t,
gamma = gamma,
alpha = alpha,
nu = nu,
min_genes = min_genes,
common_diploid = FALSE,
diploid_chroms = diploid_chroms,
ncores = ncores,
verbose = verbose,
retest = FALSE)
bulk_clones = retest_bulks(
bulk_clones,
segs_consensus,
gamma = gamma,
use_loh = use_loh,
min_LLR = min_LLR,
diploid_chroms = diploid_chroms,
ncores = ncores)
fwrite(bulk_clones, glue('{out_dir}/bulk_clones_final.tsv.gz'), sep = '\t')
if (plot) {
p = plot_bulks(bulk_clones, min_LLR = min_LLR, use_pos = TRUE, genome = genome)
ggsave(
glue('{out_dir}/bulk_clones_final.png'), p,
width = 13, height = 2*length(unique(bulk_clones$sample)), dpi = 250
)
}
log_message('All done!')
return('Success')
}
segs_equal = function(segs_1, segs_2) {
cols = c('CHROM', 'seg', 'seg_start', 'seg_end', 'cnv_state_post')
equal = isTRUE(all.equal(
segs_1 %>% select(any_of(cols)),
segs_2 %>% select(any_of(cols))
))
return(equal)
}
subtrees_equal = function(subtrees_1, subtrees_2) {
isTRUE(all.equal(subtrees_1, subtrees_2))
}
#' Run smoothed expression-based hclust
#' @param count_mat dgCMatrix Gene counts
#' @param lambdas_ref matrix Reference expression profiles
#' @param gtf dataframe Transcript GTF
#' @param sc_refs named list Reference choices for single cells
#' @param window integer Sliding window size
#' @param ncores integer Number of cores
#' @param verbose logical Verbosity
#' @keywords internal
exp_hclust = function(count_mat, lambdas_ref, gtf, sc_refs = NULL, window = 101, ncores = 1, verbose = TRUE) {
count_mat = check_matrix(count_mat)
if (is.null(sc_refs)) {
sc_refs = choose_ref_cor(count_mat, lambdas_ref, gtf)
}
lambdas_bar = get_lambdas_bar(lambdas_ref, sc_refs, verbose = FALSE)
gexp_roll_wide = smooth_expression(
count_mat,
lambdas_bar,
gtf,
window = window,
verbose = verbose
) %>% t
dist_mat = parallelDist::parDist(gexp_roll_wide, threads = ncores)
if (sum(is.na(dist_mat)) > 0) {
log_warn('NAs in distance matrix, filling with 0s. Consider filtering out cells with low coverage.')
dist_mat[is.na(dist_mat)] = 0
}
log_message('running hclust...')
hc = hclust(dist_mat, method = "ward.D2")
return(list('gexp_roll_wide' = gexp_roll_wide, 'hc' = hc))
}
#' Make a group of pseudobulks
#' @param groups list Contains fields named "sample", "cells", "size", "members"
#' @param count_mat dgCMatrix Gene counts
#' @param df_allele dataframe Alelle counts
#' @param lambdas_ref matrix Reference expression profiles
#' @param gtf dataframe Transcript GTF
#' @param min_depth integer Minimum allele depth to include
#' @param segs_loh dataframe Segments with clonal LOH to be excluded
#' @param ncores integer Number of cores
#' @return dataframe Pseudobulk profiles
#' @keywords internal
make_group_bulks = function(groups, count_mat, df_allele, lambdas_ref, gtf, min_depth = 0, nu = 1, segs_loh = NULL, ncores = NULL) {
if (length(groups) == 0) {
return(data.frame())
}
ncores = ifelse(is.null(ncores), length(groups), ncores)
results = mclapply(
groups,
mc.cores = ncores,
function(g) {
get_bulk(
count_mat = count_mat,
df_allele = df_allele,
subset = g$cells,
lambdas_ref = lambdas_ref,
gtf = gtf,
min_depth = min_depth,
nu = nu,
segs_loh = segs_loh
) %>%
mutate(
n_cells = g$size,
members = paste0(g$members, collapse = ';'),
sample = g$sample
)
})
bad = sapply(results, inherits, what = "try-error")
if (any(bad)) {
log_error(glue('job {paste(which(bad), collapse = ",")} failed'))
log_error(results[bad][[1]])
}
# unify marker index
bulks = results %>%
bind_rows() %>%
arrange(CHROM, POS) %>%
mutate(snp_id = factor(snp_id, unique(snp_id))) %>%
mutate(snp_index = as.integer(snp_id)) %>%
arrange(sample)
return(bulks)
}
#' Run multiple HMMs
#' @param bulks dataframe Pseudobulk profiles
#' @param gamma numeric Dispersion parameter for the Beta-Binomial allele model
#' @param t numeric Transition probability
#' @param alpha numeric P value cut-off to determine segment clusters in find_diploid
#' @param common_diploid logical Whether to find common diploid regions between pseudobulks
#' @param diploid_chroms character vector Known diploid chromosomes to use as baseline
#' @param retest logcial Whether to retest CNVs
#' @param run_hmm logical Whether to run HMM segments or just retest
#' @param ncores integer Number of cores
#' @param allele_only logical Whether only use allele data to run HMM
#' @keywords internal
run_group_hmms = function(
bulks, t = 1e-4, gamma = 20, alpha = 1e-4, min_genes = 10, nu = 1,
common_diploid = TRUE, diploid_chroms = NULL, allele_only = FALSE, retest = TRUE, run_hmm = TRUE,
exclude_neu = TRUE, ncores = 1, verbose = FALSE, debug = FALSE
) {
# drop samples with no allele data
bulks = bulks %>% group_by(sample) %>% filter(sum(!is.na(DP)) > 0) %>% ungroup()
if (nrow(bulks) == 0) {
return(data.frame())
}
n_groups = length(unique(bulks$sample))
if (verbose) {
log_message(glue('Running HMMs on {n_groups} cell groups..'))
}
# find common diploid region
if (!run_hmm) {
find_diploid = FALSE
} else if (common_diploid & is.null(diploid_chroms)) {
bulks = find_common_diploid(bulks, gamma = gamma, alpha = alpha, ncores = ncores)
find_diploid = FALSE
} else {
find_diploid = TRUE
}
results = mclapply(
bulks %>% split(.$sample),
mc.cores = ncores,
function(bulk) {
bulk %>% analyze_bulk(
t = t,
gamma = gamma,
nu = nu,
find_diploid = find_diploid,
run_hmm = run_hmm,
allele_only = allele_only,
diploid_chroms = diploid_chroms,
min_genes = min_genes,
retest = retest,
verbose = verbose,
exclude_neu = exclude_neu
)
})
bad = sapply(results, inherits, what = "try-error")
if (any(bad)) {
log_error(results[bad][[1]])
stop(results[bad][[1]])
}
bulks = results %>% bind_rows() %>%
group_by(seg, sample) %>%
mutate(
seg_start_index = min(snp_index),
seg_end_index = max(snp_index)
) %>%
ungroup()
return(bulks)
}
#' Extract consensus CNV segments
#' @param bulks dataframe Pseudobulks
#' @param min_LLR numeric LLR threshold to filter CNVs
#' @param min_overlap numeric Minimum overlap fraction to determine count two events as as overlapping
#' @return dataframe Consensus segments
#' @keywords internal
get_segs_consensus = function(bulks, min_LLR = 5, min_overlap = 0.45, retest = TRUE) {
if (!'sample' %in% colnames(bulks)) {
bulks$sample = 1
}
info_cols = c('sample', 'CHROM', 'seg', 'cnv_state', 'cnv_state_post',
'seg_start', 'seg_end', 'seg_start_index', 'seg_end_index',
'theta_mle', 'theta_sigma', 'phi_mle', 'phi_sigma',
'p_loh', 'p_del', 'p_amp', 'p_bamp', 'p_bdel',
'LLR', 'LLR_y', 'LLR_x', 'n_genes', 'n_snps')
# all possible aberrant segments
segs_all = bulks %>%
group_by(sample, seg, CHROM) %>%
mutate(seg_start = min(POS), seg_end = max(POS)) %>%
filter(seg_start != seg_end) %>%
ungroup() %>%
select(any_of(info_cols)) %>%
distinct() %>%
mutate(cnv_state = ifelse(LLR < min_LLR | is.na(LLR), 'neu', cnv_state))
# confident aberrant segments
segs_star = segs_all %>%
filter(cnv_state != 'neu') %>%
resolve_cnvs(
min_overlap = min_overlap
)
if (retest) {
# union of all aberrant segs
segs_cnv = segs_all %>%
filter(cnv_state != 'neu') %>%
arrange(CHROM) %>%
{GenomicRanges::GRanges(
seqnames = .$CHROM,
IRanges::IRanges(start = .$seg_start,
end = .$seg_end)
)} %>%
GenomicRanges::reduce() %>%
as.data.frame() %>%
select(CHROM = seqnames, seg_start = start, seg_end = end)
# segs to be retested
segs_retest = GenomicRanges::setdiff(
segs_cnv %>%
{GenomicRanges::GRanges(
seqnames = .$CHROM,
IRanges::IRanges(start = .$seg_start,
end = .$seg_end)
)},
segs_star %>%
{GenomicRanges::GRanges(
seqnames = .$CHROM,
IRanges::IRanges(start = .$seg_start,
end = .$seg_end)
)},
) %>%
suppressWarnings() %>%
as.data.frame() %>%
select(CHROM = seqnames, seg_start = start, seg_end = end) %>%
filter(seg_end - seg_start > 0) %>%
mutate(cnv_state = 'retest', cnv_state_post = 'retest')
} else {
segs_retest = data.frame()
}
# union of neutral segments
segs_neu = segs_all %>%
filter(cnv_state == 'neu') %>%
arrange(CHROM) %>%
{GenomicRanges::GRanges(
seqnames = .$CHROM,
IRanges::IRanges(start = .$seg_start,
end = .$seg_end)
)} %>%
GenomicRanges::reduce() %>%
as.data.frame() %>%
select(CHROM = seqnames, seg_start = start, seg_end = end) %>%
mutate(seg_length = seg_end - seg_start)
if (all(segs_all$cnv_state == 'neu')){
segs_neu = segs_neu %>%
arrange(CHROM) %>%
group_by(CHROM) %>%
mutate(
seg = paste0(CHROM, generate_postfix(1:n())),
cnv_state = 'neu',
cnv_state_post = 'neu'
) %>%
ungroup()
return(segs_neu)
}
segs_consensus = bind_rows(segs_star, segs_retest) %>%
fill_neu_segs(segs_neu) %>%
mutate(cnv_state_post = ifelse(cnv_state == 'neu', cnv_state, cnv_state_post))
return(segs_consensus)
}
#' Fill neutral regions into consensus segments
#' @param segs_consensus dataframe CNV segments from multiple samples
#' @param segs_neu dataframe Neutral segments
#' @return dataframe Collections of neutral and aberrant segments with no gaps
#' @keywords internal
fill_neu_segs = function(segs_consensus, segs_neu) {
# take complement of consensus aberrant segs
gaps = GenomicRanges::setdiff(
segs_neu %>%
{GenomicRanges::GRanges(
seqnames = .$CHROM,
IRanges::IRanges(start = .$seg_start,
end = .$seg_end)
)},
segs_consensus %>%
{GenomicRanges::GRanges(
seqnames = .$CHROM,
IRanges::IRanges(start = .$seg_start,
end = .$seg_end)
)},
) %>%
suppressWarnings() %>%
as.data.frame() %>%
select(CHROM = seqnames, seg_start = start, seg_end = end) %>%
mutate(seg_length = seg_end - seg_start) %>%
filter(seg_length > 0)
segs_consensus = segs_consensus %>%
mutate(seg_length = seg_end - seg_start) %>%
bind_rows(gaps) %>%
mutate(cnv_state = tidyr::replace_na(cnv_state, 'neu')) %>%
arrange(CHROM, seg_start) %>%
group_by(CHROM) %>%