-
Notifications
You must be signed in to change notification settings - Fork 11
/
nuc_dynamics.py
2022 lines (1433 loc) · 67.5 KB
/
nuc_dynamics.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
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
# Python functions for NucDynamics
import os
import sys
import numpy as np
import multiprocessing
import subprocess
import traceback
import gzip
QUIET = False # Global verbosity flag
LOG_FILE_OBJ = None # Created if needed
PROG_NAME = 'nuc_dynamics'
VERSION = '1.3.0'
DESCRIPTION = 'Single-cell Hi-C genome and chromosome structure calculation module for Nuc3D and NucTools'
FILE_BUFFER_SIZE = 2**16
N3D = 'n3d'
PDB = 'pdb'
FORMATS = [N3D, PDB]
MAX_CORES = multiprocessing.cpu_count()
DEFAULT_STRUCT_AMBIG_SIZE_MB = 0.4
DEFAULT_REMOVE_MODELS_SIZE_MB = 1
DEFAULT_STRUCT_AMBIG_SIZE = DEFAULT_STRUCT_AMBIG_SIZE_MB * int(1e6)
DEFAULT_REMOVE_MODELS_SIZE = DEFAULT_REMOVE_MODELS_SIZE_MB * int(1e6)
Restraint = np.dtype([('indices', 'int32', 2), ('dists', 'float64', 2),
('ambiguity', 'int32'), ('weight', 'float64')])
Contact = np.dtype([('pos', 'uint32', 2), ('ambiguity', 'int32'), ('number', 'uint32')])
def report(msg, prefix=''):
if prefix:
msg = '%s: %s' % (prefix, msg)
if LOG_FILE_OBJ:
LOG_FILE_OBJ.write(msg + '\n')
LOG_FILE_OBJ.flush()
if not QUIET:
print(msg)
else:
print(msg)
def info(msg, prefix='INFO'):
report(msg, prefix)
def warn(msg, prefix='WARNING'):
report(msg, prefix)
def critical(msg, prefix='ABORT'):
report(msg, prefix)
sys.exit(0)
def _parallel_func_wrapper(queue, target_func, proc_data, common_args):
for t, data_item in proc_data:
result = target_func(data_item, *common_args)
if queue:
queue.put((t, result))
elif isinstance(result, Exception):
raise(result)
def open_file(file_path, mode=None, buffer_size=FILE_BUFFER_SIZE, gzip_exts=('.gz','.gzip'), partial=False):
"""
GZIP agnostic file opening
"""
import io
if os.path.splitext(file_path)[1].lower() in gzip_exts:
if mode and 'w' in mode:
file_obj = io.BufferedWriter(gzip.open(file_path, mode), buffer_size)
else:
if partial:
file_obj = io.BufferedReader(gzip.open(file_path, mode or 'rb'), buffer_size)
else:
try:
file_obj = subprocess.Popen(['zcat', file_path], stdout=subprocess.PIPE).stdout
except OSError:
file_obj = io.BufferedReader(gzip.open(file_path, mode or 'rb'), buffer_size)
if sys.version_info.major > 2:
file_obj = io.TextIOWrapper(file_obj, encoding="utf-8")
else:
if sys.version_info.major > 2:
file_obj = open(file_path, mode or 'rU', buffer_size, encoding='utf-8')
else:
file_obj = open(file_path, mode or 'rU', buffer_size)
return file_obj
def parallel_split_job(target_func, split_data, common_args, num_cpu=MAX_CORES, collect_output=True):
num_tasks = len(split_data)
num_process = min(num_cpu, num_tasks)
if num_process == 1:
results = []
for data_item in split_data:
result = target_func(data_item, *common_args)
if isinstance(result, Exception):
raise(result)
results.append(result)
return results
processes = []
if collect_output:
queue = multiprocessing.Queue() # Queue will collect parallel process output
else:
queue = None
for p in range(num_process):
# Task IDs and data for each task
# Each process can have multiple tasks if there are more tasks than processes/cpus
proc_data = [(t, data_item) for t, data_item in enumerate(split_data) if t % num_cpu == p]
args = (queue, target_func, proc_data, common_args)
proc = multiprocessing.Process(target=_parallel_func_wrapper, args=args)
processes.append(proc)
for proc in processes:
proc.start()
if queue:
results = [None] * num_tasks
for i in range(num_tasks):
t, result = queue.get() # Asynchronous fetch output: whichever process completes a task first
if isinstance(result, Exception):
warn(result)
critical('\n* * * * C/Cython code may need to be recompiled. Try running "python setup_cython.py build_ext --inplace" * * * *\n')
results[t] = result
queue.close()
return results
else:
for proc in processes: # Asynchromous wait and no output captured
proc.join()
def load_ncc_file(file_path, active_only=True, nmax=int(1e6)):
"""Load chromosome and contact data from NCC format file, as output from NucProcess"""
from numpy import array
contact_dict = {}
ambig_group = 0
with open_file(file_path) as file_obj:
n_active = 0
for i, line in enumerate(file_obj):
if line.startswith('#'):
continue
row = line.split()
if not '.' in row[12]:
critical('The old version of the NCC format is not supported by the version of NucDynamics. Please use the latest version of NucProcess')
ambig_group = 0
file_obj.seek(0)
for i, line in enumerate(file_obj):
if line.startswith('#'):
continue
chr_a, f_start_a, f_end_a, start_a, end_a, strand_a, chr_b, f_start_b, f_end_b, start_b, end_b, strand_b, ambig_code, pair_id, swap_pair = line.split()
if int(float(ambig_code)) > 0:
ambig_group += 1
if active_only and ambig_code.endswith('.0'): # inactive
continue
pos_a = int(f_start_a if strand_a == '+' else f_end_a)
pos_b = int(f_start_b if strand_b == '+' else f_end_b)
if chr_a > chr_b:
chr_a, chr_b = chr_b, chr_a
pos_a, pos_b = pos_b, pos_a
if chr_a not in contact_dict:
contact_dict[chr_a] = {}
if chr_b not in contact_dict[chr_a]:
contact_dict[chr_a][chr_b] = []
contact_dict[chr_a][chr_b].append(((pos_a, pos_b), ambig_group, i)) # Uses original file line number
n_active += 1
if n_active > nmax:
critical('Too many contacts in ncc file (> {:,}), this code is meant for single-cell data'.format(nmax))
for chr_a in contact_dict:
for chr_b in contact_dict[chr_a]:
contact_dict[chr_a][chr_b] = array(contact_dict[chr_a][chr_b], dtype=Contact)
return contact_dict
def save_ncc_file(file_path_in, name, contact_dict, particle_size=None, offset=0):
if particle_size is not None:
name = '%s_%d' % (name, int(particle_size/1000))
file_root, file_ext = os.path.splitext(file_path_in)
if file_ext.lower() in ('.gz','.gzip'):
file_root, file_ext = os.path.splitext(file_root)
file_path_out ='%s_%s.ncc' % (file_root, name)
active_numbers = set()
for chr_a in contact_dict:
for chr_b in contact_dict[chr_a]:
active_numbers |= set(contact_dict[chr_a][chr_b]['number'])
with open_file(file_path_in) as file_obj_in:
with open_file(file_path_out, 'w') as file_obj_out:
n = offset
for i, line in enumerate(file_obj_in):
if line.startswith('#'):
file_obj_out.write(line)
else:
row = line.split()
ag = row[12]
activ = '1' if i in active_numbers else '0'
row[12] = ag[:-1] + activ
line = ' '.join(row) + '\n'
file_obj_out.write(line)
def calc_limits(contact_dict):
chromo_limits = {}
for chrs, contacts in flatten_dict(contact_dict).items():
if len(contacts) < 1:
continue
for chromo, seq_pos in zip(chrs, contacts['pos'].T):
min_pos = seq_pos.min()
max_pos = seq_pos.max()
prev_min, prev_max = chromo_limits.get(chromo, (min_pos, max_pos))
chromo_limits[chromo] = [min(prev_min, min_pos), max(prev_max, max_pos)]
return chromo_limits
def export_n3d_coords(file_path, coords_dict, seq_pos_dict):
file_obj = open(file_path, 'w')
write = file_obj.write
for chromo in seq_pos_dict:
chromo_coords = coords_dict[chromo]
chromo_seq_pos = seq_pos_dict[chromo]
num_models = len(chromo_coords)
num_coords = len(chromo_seq_pos)
line = '%s\t%d\t%d\n' % (chromo, num_coords, num_models)
write(line)
for j in range(num_coords):
data = chromo_coords[:,j].ravel().tolist()
data = '\t'.join('%.8f' % d for d in data)
line = '%d\t%s\n' % (chromo_seq_pos[j], data)
write(line)
file_obj.close()
def export_pdb_coords(file_path, coords_dict, seq_pos_dict, particle_size, scale=1.0, extended=True):
"""
Write chromosome particle coordinates as a PDB format file
"""
alc = ' '
ins = ' '
prefix = 'HETATM'
line_format = '%-80.80s\n'
if extended:
pdb_format = '%-6.6s%5.1d %4.4s%s%3.3s %s%4.1d%s %8.3f%8.3f%8.3f%6.2f%6.2f %2.2s %10d\n'
ter_format = '%-6.6s%5.1d %s %s%4.1d%s %10d\n'
else:
pdb_format = '%-6.6s%5.1d %4.4s%s%3.3s %s%4.1d%s %8.3f%8.3f%8.3f%6.2f%6.2f %2.2s \n'
ter_format = '%-6.6s%5.1d %s %s%4.1d%s \n'
file_obj = open(file_path, 'w')
write = file_obj.write
chromosomes = list(seq_pos_dict.keys())
sort_chromos = []
for chromo in chromosomes:
if chromo[:3] == 'chr':
key = chromo[3:]
else:
key = chromo
if key.isdigit():
key = '%03d' % int(key)
sort_chromos.append((key, chromo))
sort_chromos.sort()
sort_chromos = [x[1] for x in sort_chromos]
num_models = len(coords_dict[chromosomes[0]])
title = 'NucDynamics genome structure export'
write(line_format % 'TITLE %s' % title)
write(line_format % 'REMARK 210')
write(line_format % 'REMARK 210 Atom type C is used for all particles')
write(line_format % 'REMARK 210 Atom number increases every %s bases' % particle_size)
write(line_format % 'REMARK 210 Residue code indicates chromosome')
write(line_format % 'REMARK 210 Residue number represents which sequence Mb the atom is in')
write(line_format % 'REMARK 210 Chain letter is different every chromosome, where Chr1=a, Chr2=b etc.')
if extended:
file_obj.write(line_format % 'REMARK 210 Extended PDB format with particle seq. pos. in last column')
file_obj.write(line_format % 'REMARK 210')
pos_chromo = {}
for m in range(num_models):
line = 'MODEL %4d' % (m+1)
write(line_format % line)
c = 0
j = 1
seqPrev = None
for k, chromo in enumerate(sort_chromos):
chain_code = chr(ord('a')+k)
tlc = chromo
while len(tlc) < 2:
tlc = '_' + tlc
if len(tlc) == 2:
tlc = 'C' + tlc
if len(tlc) > 3:
tlc = tlc[:3]
chromo_model_coords = coords_dict[chromo][m]
if not len(chromo_model_coords):
continue
pos = seq_pos_dict[chromo]
for i, seqPos in enumerate(pos):
c += 1
seqMb = int(seqPos//1e6) + 1
if seqMb == seqPrev:
j += 1
else:
j = 1
el = 'C'
a = 'C%d' % j
aName = '%-3s' % a
x, y, z = chromo_model_coords[i] #XYZ coordinates
x *= scale
y *= scale
z *= scale
seqPrev = seqMb
pos_chromo[c] = chromo
if extended:
line = pdb_format % (prefix,c,aName,alc,tlc,chain_code,seqMb,ins,x,y,z,0.0,0.0,el,seqPos)
else:
line = pdb_format % (prefix,c,aName,alc,tlc,chain_code,seqMb,ins,x,y,z,0.0,0.0,el)
write(line)
write(line_format % 'ENDMDL')
for i in range(c-2):
if pos_chromo[i+1] == pos_chromo[i+2]:
line = 'CONECT%5.1d%5.1d' % (i+1, i+2)
write(line_format % line)
write(line_format % 'END')
file_obj.close()
def remove_ambiguous_contacts(contact_dict):
ambiguity_list = []
for contacts in flatten_dict(contact_dict).values():
ambiguity_list.extend(contacts['ambiguity'])
ambiguity_array = np.array(sorted(ambiguity_list))
_, inverse, counts = np.unique(
ambiguity_array, return_inverse=True, return_counts=True
)
ambiguity_set = set(ambiguity_array[counts[inverse] == 1])
from collections import defaultdict
new_contacts = defaultdict(dict)
for (chromoA, chromoB), contacts in flatten_dict(contact_dict).items():
contacts = [c for c in contacts if c['ambiguity'] in ambiguity_set]
if contacts:
new_contacts[chromoA][chromoB] = np.array(contacts, dtype=Contact)
return dict(new_contacts)
def homologous_pair(chromoA, chromoB):
# slightly dangerous code, relies on chromoA/B being of form 'chr11.a', etc.
nA = chromoA.rfind('.')
nB = chromoB.rfind('.')
if nA < 0 or nB < 0:
return False
return (chromoA[:nA] == chromoB[:nB]) and (chromoA[nA+1:] != chromoB[nB+1:])
def remove_homologous_pairs(contact_dict):
from collections import defaultdict
new_contacts = defaultdict(dict)
for (chromoA, chromoB), contacts in flatten_dict(contact_dict).items():
if not homologous_pair(chromoA, chromoB):
new_contacts[chromoA][chromoB] = contacts
return dict(new_contacts)
DEFAULT_DISAMBIG_MIN = 3
DEFAULT_DISAMBIG_FRAC = 0.02
def resolve_homolog_ambiguous(contact_dict, min_disambig=DEFAULT_DISAMBIG_MIN, frac_disambig=DEFAULT_DISAMBIG_FRAC):
from collections import defaultdict
group_sizes = defaultdict(int)
chromo_pair_groups = defaultdict(set)
# Consider all trans chromosome pairs
# If the pair has fewer than min_disambig unambigous contacts it can have no homologous chromosome ambigous contacts
# This filtering is only for homologous chromosome ambiguity, not positional ambiguity
# Fetch the sizes of all ambiguity groups and which group goes with which chromo pair
for chr_a in contact_dict:
for chr_b in contact_dict[chr_a]:
contacts = contact_dict[chr_a][chr_b].T
for contact in contacts:
ambig_group = contact['ambiguity']
group_sizes[ambig_group] += 1
if chr_a > chr_b:
chr_a, chr_b = chr_b, chr_a
chromo_key = (chr_a, chr_b)
chromo_pair_groups[chromo_key].add(ambig_group)
n_groups = len(group_sizes)
# Count total unambiguous groups for each chromo pair
pair_stats = {}
unambig_fracs = []
unambig_counts = []
for chromo_key in chromo_pair_groups:
pair_groups = chromo_pair_groups[chromo_key]
n_unambig = 0
for ambig_group in pair_groups:
count = group_sizes[ambig_group]
if count == 1:
n_unambig += 1
frac_unambig = n_unambig/float(len(pair_groups))
if 0.0 < frac_unambig < 1.0:
unambig_counts.append(n_unambig)
unambig_fracs.append(frac_unambig)
pair_stats[chromo_key] = (n_unambig, frac_unambig)
# Calculate threshold statistic and excluded pairs
p90 = np.percentile(unambig_counts, [90.0])[0] # Threshold to seek only chromo pairs where counts are high
median_frac = np.median([unambig_fracs[i] for i, v in enumerate(unambig_counts) if v >= p90]) # An estimate of the average unambiguous proportion
non_contact_chromos = set()
for chromo_key in pair_stats:
n_unambig, frac_unambig = pair_stats[chromo_key]
chr_a, chr_b = chromo_key
if (n_unambig < min_disambig) or (frac_unambig < (frac_disambig * median_frac)):
chr_a, chr_b = chromo_key
if chr_a != chr_b:
non_contact_chromos.add(chromo_key)
new_contact_dict = defaultdict(dict)
for chr_a in contact_dict:
for chr_b in contact_dict[chr_a]:
contacts = contact_dict[chr_a][chr_b].T
contacts_new = []
for contact in contacts:
ambig_group = contact['ambiguity']
group_size = group_sizes[ambig_group]
if chr_a > chr_b:
chr_a, chr_b = chr_b, chr_a
chromo_key = (chr_a, chr_b)
if chromo_key in non_contact_chromos:
if group_size == 1: # Keep unambiguous for excluded pairs unless homologs
not_homolog = chr_a.split('.')[0] != chr_b.split('.')[0]
if not_homolog:
contacts_new.append(contact)
else:
contacts_new.append(contact)
if contacts_new:
contacts = np.array(contacts_new, dtype=Contact)
new_contact_dict[chr_a][chr_b] = np.array(contacts, dtype= Contact)
return new_contact_dict
def resolve_3d_ambiguous(contact_dict, seq_pos_dict, coords_dict, percentile_threshold=98.0, min_bead_sep=3):
from collections import defaultdict
ambig_groups = defaultdict(list)
contact_pos = defaultdict(list)
for chr_a in contact_dict:
if chr_a not in seq_pos_dict:
continue
for chr_b in contact_dict[chr_a]:
if chr_b not in seq_pos_dict:
continue
contacts = contact_dict[chr_a][chr_b].T
for contact in contacts:
pos_a, pos_b = contact['pos']
ambig_group = contact['ambiguity']
# None vals below will be filled with structure distances and bead separations
ambig_groups[ambig_group].append([chr_a, pos_a, chr_b, pos_b, None, None])
contact_pos[chr_a].append(pos_a)
contact_pos[chr_b].append(pos_b)
ambig_groups = {g:ambig_groups[g] for g in ambig_groups if len(ambig_groups[g]) > 1}
info('Get contacts particle indices')
pos_idx = {}
for chromo in contact_pos:
if chromo in coords_dict:
pos_idx[chromo] = cps = {}
seq_pos = seq_pos_dict[chromo]
j_max = len(seq_pos)-1
idx = np.searchsorted(seq_pos, contact_pos[chromo], side='right')
for i, pos in enumerate(contact_pos[chromo]):
j = min(j_max, idx[i])
cps[pos] = j
info('Get distances')
closest_nz = []
min_dists = {}
min_seps = {}
for ambig_group in ambig_groups:
dists = []
bead_seps = []
for vals in ambig_groups[ambig_group]:
chr_a, pos_a, chr_b, pos_b, null, null2 = vals
coords_a = coords_dict[chr_a]
coords_b = coords_dict[chr_b]
i = pos_idx[chr_a][pos_a]
j = pos_idx[chr_b][pos_b]
deltas = coords_a[:,i] - coords_b[:,j]
dists2 = (deltas * deltas).sum(axis=1)
dist = np.sqrt(dists2.mean())
dists.append(dist)
if chr_a == chr_b:
bead_sep = abs(j-i)
else:
bead_sep = 99999 # Just something huge for trans
bead_seps.append(bead_sep)
vals[4] = dist
vals[5] = bead_sep
min_sep = min(bead_seps)
min_seps[ambig_group] = min_sep
min_dist = min(dists)
min_dists[ambig_group] = min_dist
if min_dist > 0.0:
closest_nz.append(min_dist)
remove = set()
if closest_nz:
dist_thresh = np.percentile(closest_nz, [percentile_threshold])[0]
info('Filter ambiguity groups. Threshold = %.2f' % dist_thresh)
for ambig_group in ambig_groups:
min_dist = min_dists[ambig_group]
min_sep = min_seps[ambig_group]
for chr_a, pos_a, chr_b, pos_b, dist, bead_sep in ambig_groups[ambig_group]:
if (min_dist == 0) and (chr_a != chr_b):
remove.add((ambig_group, chr_a, pos_a, chr_b, pos_b))
continue
if (min_sep < min_bead_sep) and (chr_a != chr_b):
remove.add((ambig_group, chr_a, pos_a, chr_b, pos_b))
continue
if dist != min_dist:
if dist > min_dist + dist_thresh/2:
remove.add((ambig_group, chr_a, pos_a, chr_b, pos_b))
continue
if dist > dist_thresh:
remove.add((ambig_group, chr_a, pos_a, chr_b, pos_b))
continue
new_contact_dict = defaultdict(dict)
for (chromoA, chromoB), contacts in flatten_dict(contact_dict).items():
contacts = [c for c in contacts if (c['ambiguity'], chromoA, c['pos'][0], chromoB, c['pos'][1]) not in remove]
if contacts:
new_contact_dict[chromoA][chromoB] = np.array(contacts, dtype=Contact)
return dict(new_contact_dict)
def between(x, a, b):
return (a < x) & (x < b)
def initial_clean_contacts(contact_dict, threshold_cis=int(2e6), threshold_trans=int(10e6), pos_error=100, ignore=()):
"""
Select only unambigous contacts that are within a given sequence separation of another
contact, for the same chromosome pair.
"""
info('Cleaning initial contact list')
from numpy import array, zeros, eye
from collections import defaultdict
new_contacts = defaultdict(dict)
ag_size = defaultdict(int)
contact_pos = defaultdict(list)
for chr_pair, contacts in flatten_dict(contact_dict).items():
for ag in contacts['ambiguity']:
ag_size[ag] += 1
for (chromoA, chromoB), contacts in flatten_dict(contact_dict).items():
# Select only unambiguous
unambig = np.array([ag_size[ag] == 1 for ag in contacts['ambiguity']])
contacts = contacts[unambig]
if chromoA == chromoB:
threshold = threshold_cis
else:
threshold = threshold_trans
# positions is N x 2 matrix, where there are N contacts (and the 2 is because of pos_a, pos_b)
positions = contacts['pos'].astype('int32')
if chromoA in ignore or chromoB in ignore:
new_contacts[chromoA][chromoB] = contacts
continue
if len(positions) == 0: # Sometimes empty e.g. for MT, Y chromos
continue
# positions[:, 0] and positions[:, 1] have shape N
# the [None, ...] converts that into shape 1 x N (but it looks like you don't need the ...)
pos_a = positions[:, 0][None, ...]
pos_b = positions[:, 1][None, ...]
nondiagonal = ~eye(len(positions), dtype='bool') # False on diagonal, True off diagonal
supported = zeros(len(positions), dtype='bool')
# do not need nondiagonal in first expression below as long as pos_error > 0
# the diagonal in the second expression so as not to include i if a_i happens to be close to b_i
# abs(pos_a - pos_b.T) gives N x N matrix with difference in positions between a_i and b_j
# we want both ends in a contact to be close to the ends in another contact
supported |= (between(abs(pos_a - pos_a.T), pos_error, threshold) &
between(abs(pos_b - pos_b.T), pos_error, threshold)).any(axis=0)
supported |= (between(abs(pos_a - pos_b.T), pos_error, threshold) &
between(abs(pos_b - pos_a.T), pos_error, threshold) &
nondiagonal).any(axis=0)
if not supported.any():
continue
new_contacts[chromoA][chromoB] = contacts[supported]
return dict(new_contacts)
def remove_violated_contacts(contact_dict, coords_dict, particle_seq_pos, threshold=5.0):
"""
Remove contacts whith structure distances that exceed a given threshold
"""
from numpy import int32, sqrt, array
from collections import defaultdict
new_contacts = defaultdict(dict)
for chr_a in contact_dict:
if chr_a not in coords_dict:
continue
for chr_b, contacts in contact_dict[chr_a].items():
if chr_b not in coords_dict:
continue
contact_pos_a = contacts['pos'][:, 0].astype(int32)
contact_pos_b = contacts['pos'][:, 1].astype(int32)
coords_a = coords_dict[chr_a]
coords_b = coords_dict[chr_b]
struc_dists = []
for m in range(len(coords_a)):
coord_data_a = get_interpolated_coords(coords_a[m], contact_pos_a, particle_seq_pos[chr_a])
coord_data_b = get_interpolated_coords(coords_b[m], contact_pos_b, particle_seq_pos[chr_b])
deltas = coord_data_a - coord_data_b
dists = sqrt((deltas*deltas).sum(axis=1))
struc_dists.append(dists)
# Average over all conformational models
struc_dists = array(struc_dists).T.mean(axis=1)
# Select contacts with distances below distance threshold
indices = (struc_dists < threshold).nonzero()[0]
new_contacts[chr_a][chr_b] = contacts[indices]
return dict(new_contacts)
def get_random_coords(shape, radius=10.0):
"""
Get random, uniformly sampled coorinate positions, restricted to
a nD-sphere of given radius
"""
from numpy import random
from numpy.linalg import norm
u = random.uniform(size=shape[:-1])
x = random.normal(size=shape)
scaling = (radius * u ** (1/shape[-1])) / norm(x, axis=-1)
return scaling[..., None] * x
def get_interpolated_coords(coords, pos, prev_pos):
from numpy import interp, apply_along_axis
from functools import partial
coords = apply_along_axis(partial(interp, pos, prev_pos), -2, coords)
return coords
def pack_chromo_coords(coords_dict, chromosomes):
"""
Place chromosome 3D coordinates stored in a dictionary keyed by
chromosome name into a single, ordered array. The chromosomes argument
is required to set the correct array storage order.
"""
chromo_num_particles = [len(coords_dict[chromo][0]) for chromo in chromosomes]
n_particles = sum(chromo_num_particles)
n_models = len(coords_dict[chromosomes[0]])
coords = np.empty((n_models, n_particles, 3), float)
j = 0
for i, chromo in enumerate(chromosomes):
span = chromo_num_particles[i]
coords[:,j:j+span] = coords_dict[chromo]
j += span
return coords
def unpack_chromo_coords(coords, chromosomes, seq_pos_dict):
"""
Exctract coords for multiple chromosomes stored in a single array into
a dictionary, keyed by chromosome name. The chromosomes argument is required
to get the correct array storage order.
"""
chromo_num_particles = [len(seq_pos_dict[chromo]) for chromo in chromosomes]
n_seq_pos = sum(chromo_num_particles)
n_models, n_particles, dims = coords.shape
if n_seq_pos != n_particles:
msg = 'Model coordinates must be an array of num models x %d' % (n_seq_pos,)
raise(Exception(msg))
coords_dict = {}
j = 0
for i, chromo in enumerate(chromosomes):
span = chromo_num_particles[i]
coords_dict[chromo] = coords[:,j:j+span] # all models, slice
j += span
return coords_dict
def anneal_model(model_data, anneal_schedule, masses, radii, restraints, rep_dists,
ambiguity, temp, time_step, dyn_steps, repulse, dist, bead_size):
import gc
m, model_coords = model_data
# Anneal one model in parallel
time_taken = 0.0
if m == 0:
print_interval = max(1, dyn_steps/2)
else:
print_interval = 0
info(' starting model %d' % m)
nrep_max = np.int32(0)
for temp, repulse in anneal_schedule:
gc.collect() # Try to free some memory
# Update coordinates for this temp
#try:
dt, nrep_max = dyn_util.run_dynamics(model_coords, masses, radii, rep_dists,
restraints['indices'], restraints['dists'],
restraints['weight'], ambiguity,
temp, time_step, dyn_steps, repulse, dist,
bead_size, nrep_max,
print_interval=print_interval)
#except Exception as err:
# return err
nrep_max = np.int32(nrep_max * 1.2)
time_taken += dt
# Center
model_coords -= model_coords.mean(axis=0)
info(' done model %d' % m)
return model_coords
def calc_bins(chromo_limits, particle_size):
from numpy import arange
from math import ceil
bins = {}
for chromo, (start, end) in chromo_limits.items():
first_bin = start // particle_size
last_bin = 1 + (end-1) // particle_size
start = particle_size * first_bin
end = particle_size * last_bin
bins[chromo] = arange(start, end, particle_size, dtype='int32')
return bins
def calc_ambiguity_offsets(groups):
"""
Convert (sorted) ambiguity groups to group-offsets for
annealing calculations.
"""
from numpy import flatnonzero, zeros
group_starts = zeros(len(groups) + 1, dtype='bool')
group_starts[-1] = group_starts[0] = 1
group_starts[:-1] |= groups == 0
group_starts[1:-1] |= groups[1:] != groups[:-1]
return flatnonzero(group_starts).astype('int32')
def backbone_restraints(seq_pos, particle_size, scale=1.0, lower=0.1, upper=1.1, weight=1.0):
from numpy import empty, arange, array
restraints = empty(len(seq_pos) - 1, dtype=Restraint)
offsets = array([0, 1], dtype='int')
restraints['indices'] = arange(len(restraints))[:, None] + offsets
# Normally 1.0 for regular sized particles
bounds = array([lower, upper], dtype='float') * scale
restraints['dists'] = ((seq_pos[1:] - seq_pos[:-1]) / particle_size)[:, None] * bounds
restraints['ambiguity'] = 0 # Use '0' to represent no ambiguity
restraints['weight'] = weight
return restraints
def flatten_dict(d):
r = {}
for key, value in d.items():
if isinstance(value, dict):
r.update({(key,) + k: v for k, v in flatten_dict(value).items()})
else:
r[(key,)] = value
return r
def tree():
from collections import defaultdict
def tree_():
return defaultdict(tree_)
return tree_()
def unflatten_dict(d):
r = tree()
for ks, v in d.items():
d = r