1 """@namespace IMP.pmi.io.xltable
2 Tools to plot a contact map overlaid with cross-links.
5 from __future__
import print_function
7 from scipy.spatial.distance
import cdist
9 import matplotlib
as mpl
11 import matplotlib.cm
as cm
12 import matplotlib.pyplot
as plt
14 from collections
import defaultdict
30 """ class to read, analyze, and plot xlink data on contact maps
31 Canonical way to read the data:
32 1) load sequences and name them
33 2) load coordinates for those sequences from PDB file
39 def __init__(self,contact_threshold):
41 self.cross_link_db =
None
42 self.residue_pair_list = []
43 self.distance_maps = []
44 self.contact_freqs =
None
47 self.index_dict = defaultdict(list)
48 self.contact_threshold = contact_threshold
60 def _colormap_distance(self, dist, threshold=35, tolerance=0):
61 if dist < threshold - tolerance:
63 elif dist >= threshold + tolerance:
68 def _colormap_satisfaction(self, sat, threshold=0.5, tolerance=0.1):
69 if sat >= threshold + tolerance:
72 elif sat < threshold + tolerance
and sat >= threshold - tolerance :
79 def _get_percentage_satisfaction(self,r1,c1,r2,c2,threshold=35):
81 idx1=self.index_dict[c1][r1]
85 idx2=self.index_dict[c2][r2]
89 for dists
in self.dist_maps:
91 if dist<threshold: nsatisfied+=1
92 return float(nsatisfied)/len(self.dist_maps)
94 def _get_distance(self,r1,c1,r2,c2):
95 if self.index_dict
is not None:
97 idx1=self.index_dict[c1][r1]
101 idx2=self.index_dict[c2][r2]
104 return self.av_dist_map[idx1,idx2]
106 if (r1,c1,r2,c2)
not in self.stored_dists.keys():
108 selpart=sel.get_selected_particles()
109 selpart_res_one=list(set(self.particles_resolution_one) & set(selpart))
110 if len(selpart_res_one)>1:
return None
111 if len(selpart_res_one)==0:
return None
112 selpart_res_one_1=selpart_res_one[0]
114 selpart=sel.get_selected_particles()
115 selpart_res_one=list(set(self.particles_resolution_one) & set(selpart))
116 if len(selpart_res_one)>1:
return None
117 if len(selpart_res_one)==0:
return None
118 selpart_res_one_2=selpart_res_one[0]
122 self.stored_dists[(r1,c1,r2,c2)]=dist
124 dist=self.stored_dists[(r1,c1,r2,c2)]
130 def _internal_load_maps(self,maps_fn):
131 npzfile = np.load(maps_fn)
132 cname_array=npzfile[
'cname_array']
133 idx_array=npzfile[
'idx_array']
135 for cname,idxs
in zip(cname_array,idx_array):
138 index_dict[cname]=tmp[0:tmp.index(-1)]
140 index_dict[cname]=tmp
141 av_dist_map = npzfile[
'av_dist_map']
142 contact_map = npzfile[
'contact_map']
143 return index_dict,av_dist_map,contact_map
146 """ read sequence. structures are displayed in the same order as sequences are read.
147 fasta_file: file to read
148 id_in_fasta_file: id of desired sequence
149 protein_name: identifier for this sequence (use same name when handling coordinates)
150 can provide the fasta name (for retrieval) and the protein name (for storage) """
153 if id_in_fasta_file
is None:
154 id_in_fasta_file = name
155 if id_in_fasta_file
not in record_dict:
156 raise KeyError(
"id %s not found in fasta file" % id_in_fasta_file)
157 length = len(record_dict[id_in_fasta_file])
159 self.sequence_dict[protein_name] = str(record_dict[id_in_fasta_file])
162 """ read coordinates from a pdb file. also appends to distance maps
163 @param pdbfile file for reading coords
164 @param chain_to_name_map correspond chain ID with protein name (will ONLY read these chains)
165 Key: PDB chain ID. Value: Protein name (set in sequence reading)
166 @note This function returns an error if the sequence for each chain has NOT been read
169 total_len = sum(len(self.sequence_dict[s])
for s
in self.sequence_dict)
170 coords = np.ones((total_len,3)) * 1e5
172 for cid
in chain_to_name_map:
173 cname = chain_to_name_map[cid]
174 if cname
not in self.sequence_dict:
175 print(
"ERROR: chain",cname,
'has not been read or has a naming mismatch')
178 self.index_dict[cname]=range(prev_stop,prev_stop+len(self.sequence_dict[cname]))
180 for p
in sel.get_selected_particles():
182 coords[rnum+prev_stop-1,:] = list(
IMP.core.XYZ(p).get_coordinates())
183 prev_stop+=len(self.sequence_dict[cname])
184 dists = cdist(coords, coords)
185 binary_dists = np.where((dists <= self.contact_threshold) & (dists >= 1.0), 1.0, 0.0)
187 self.dist_maps= [dists]
188 self.av_dist_map = dists
189 self.contact_freqs = binary_dists
192 self.dist_maps.append(dists)
193 self.av_dist_map += dists
194 self.contact_freqs += binary_dists
200 """ read coordinates from a rmf file. It needs IMP to run.
201 rmf has been created using IMP.pmi conventions. It gets the
202 highest resolution automatically. Also appends to distance maps
203 @param rmf_name file for reading coords
204 @param rmf_frame_index frame index from the rmf
205 @param nomap Default False, if True it will not calculate the contact map
207 (particles_resolution_one, prots)=self._get_rmf_structure(rmf_name,rmf_frame_index)
209 total_len = sum(len(self.sequence_dict[s])
for s
in self.sequence_dict)
210 print(self.sequence_dict,total_len)
212 coords = np.ones((total_len,3)) * 1e6
214 sorted_particles=IMP.pmi.tools.sort_by_residues(particles_resolution_one)
217 self.particles_resolution_one=particles_resolution_one
222 for cname
in chain_names:
225 self.index_dict[cname]=range(prev_stop,prev_stop+len(self.sequence_dict[cname]))
226 rindexes=range(1,len(self.sequence_dict[cname])+1)
227 for rnum
in rindexes:
229 selpart=sel.get_selected_particles()
230 selpart_res_one=list(set(particles_resolution_one) & set(selpart))
231 if len(selpart_res_one)>1:
continue
232 if len(selpart_res_one)==0:
continue
233 selpart_res_one=selpart_res_one[0]
235 coords[rnum+prev_stop-1,:]=
IMP.core.XYZ(selpart_res_one).get_coordinates()
237 print(
"Error: exceed max size",prev_stop,total_len,cname,rnum)
239 prev_stop+=len(self.sequence_dict[cname])
240 dists = cdist(coords, coords)
241 binary_dists = np.where((dists <= self.contact_threshold) & (dists >= 1.0), 1.0, 0.0)
243 self.dist_maps= [dists]
244 self.av_dist_map = dists
245 self.contact_freqs = binary_dists
248 self.dist_maps.append(dists)
249 self.av_dist_map += dists
250 self.contact_freqs += binary_dists
254 def _get_rmf_structure(self,rmf_name,rmf_frame_index):
255 rh= RMF.open_rmf_file_read_only(rmf_name)
258 print(
"getting coordinates for frame %i rmf file %s" % (rmf_frame_index, rmf_name))
263 protein_names=particle_dict.keys()
264 particles_resolution_one=[]
265 for k
in particle_dict:
266 particles_resolution_one+=(particle_dict[k])
268 return particles_resolution_one, prots
271 def save_maps(self,maps_fn):
272 maxlen=max(len(self.index_dict[key])
for key
in self.index_dict)
275 for cname,idx
in self.index_dict.items():
277 idxs.append(idx+[-1]*(maxlen-len(idx)))
278 idx_array=np.array(idxs)
279 cname_array=np.array(cnames)
281 cname_array=cname_array,
283 av_dist_map=self.av_dist_map,
284 contact_map=self.contact_freqs)
286 def load_maps(self,maps_fn):
287 self.index_dict,self.av_dist_map,self.contact_freqs=self._internal_load_maps(maps_fn)
290 """ read crosslinks from a CSV file.
291 provide a CrossLinkDataBaseKeywordsConverter to explain the columns
292 @distance_field is the optional keyword for the distance to be read form the file.
293 This can skip the rmf reading to calculate the distance of cross-links if
294 already provided in the csv file."""
295 if type(CrossLinkDataBase)
is not IMP.pmi.io.crosslink.CrossLinkDataBase:
296 raise TypeError(
"Crosslink database must be a IMP.pmi.io.CrossLinkDataBase type")
297 self.cross_link_db=CrossLinkDataBase
298 if distance_field
is not None:
299 total_len = sum(len(self.sequence_dict[s])
for s
in self.sequence_dict)
300 zeros = np.zeros((total_len,3))
301 self.av_dist_map = cdist(zeros,zeros)
302 for xl
in self.cross_link_db:
303 c1=xl[self.cross_link_db.protein1_key]
304 c2=xl[self.cross_link_db.protein2_key]
305 r1=xl[self.cross_link_db.residue1_key]
306 r2=xl[self.cross_link_db.residue2_key]
308 self.stored_dists[(r1,c1,r2,c2)]=float(xl[distance_field])
309 self.stored_dists[(r2,c2,r1,c1)]=float(xl[distance_field])
311 self.stored_dists[(r1,c1,r2,c2)]=10e6
312 self.stored_dists[(r2,c2,r1,c1)]=10e6
314 for xl
in self.cross_link_db:
315 c1=xl[self.cross_link_db.protein1_key]
316 c2=xl[self.cross_link_db.protein2_key]
317 r1=xl[self.cross_link_db.residue1_key]
318 r2=xl[self.cross_link_db.residue2_key]
319 self.stored_dists[(r1,c1,r2,c2)]=self._get_distance(r1,c1,r2,c2)
320 self.stored_dists[(r2,c2,r1,c1)]=self._get_distance(r2,c2,r1,c1)
323 """ select the atom names of residue pairs to plot on the contact map
324 list of residues types must be single letter code
325 e.g. residue_type_pair=("K","K")
327 rtp=sorted(residue_type_pair)
328 for prot1
in self.sequence_dict:
329 seq1=self.sequence_dict[prot1]
330 for nres1,res1
in enumerate(seq1):
331 for prot2
in self.sequence_dict:
332 seq2=self.sequence_dict[prot2]
333 for nres2,res2
in enumerate(seq2):
334 if sorted((res1,res2))==rtp:
335 self.residue_pair_list.append((nres1+1,prot1,nres2+1,prot2))
338 """ loop through each distance map and get frequency of contacts
340 if self.num_pdbs!=0
and self.num_rmfs==0:
341 self.av_dist_map = 1.0/self.num_pdbs * self.av_dist_map
342 self.contact_freqs = 1.0/self.num_pdbs * self.contact_freqs
343 if self.num_pdbs==0
and self.num_rmfs!=0:
344 self.av_dist_map = 1.0/self.num_rmfs * self.av_dist_map
345 self.contact_freqs = 1.0/self.num_rmfs * self.contact_freqs
347 def setup_difference_map(self,maps_fn1,maps_fn2,thresh):
348 idx1,av1,contact1=self._internal_load_maps(maps_fn1)
349 idx2,av2,contact2=self._internal_load_maps(maps_fn2)
351 print(
"UH OH: index dictionaries do not match!")
359 elif c1>thresh
and c2<thresh:
361 elif c1<thresh
and c2>thresh:
365 f = np.vectorize(logic,otypes=[np.int])
366 print(
'computing contact map')
367 self.contact_freqs = f(contact1,contact2)
372 def spring_layout(self,ax,plt,data, annotations, iterations = 100, k=1):
374 import networkx
as nx
382 Author: G. Bouvier, Pasteur Institute, Paris
383 Website: http://bloggb.fr/2015/10/19/spring_layout_for_annotating_plot_in_matplotlib.html
384 - data: coordinates of your points [(x1,y1), (x2,y2), ..., (xn,yn)]
385 - annotations: text for your annotation ['str_1', 'str_2', ..., 'str_n']
386 - iterations: number of iterations for spring layout
387 - k: optimal distance between nodes
391 x, y = [e[0]
for e
in data], [e[1]
for e
in data]
392 xmin, xmax = min(x), max(x)
393 ymin, ymax = min(y), max(y)
402 for i,xy
in enumerate(data):
429 ps = particles.values()+particles_fixed.values()
436 rs.add_restraint(evr)
440 #G.add_edge(xy, text)
443 delTri = scipy.spatial.Delaunay(data)
444 #plt.scatter([e[0] for e in data], [e[1] for e in data])
445 # create a set for edges that are indexes of the points
447 # for each Delaunay triangle
448 for n in xrange(delTri.nsimplex):
449 # for each edge of the triangle
451 # (sorting avoids duplicated edges being added to the set)
452 # and add to the edges set
453 edge = sorted([delTri.vertices[n,0], delTri.vertices[n,1]])
454 edges.add((edge[0], edge[1]))
455 edge = sorted([delTri.vertices[n,0], delTri.vertices[n,2]])
456 edges.add((edge[0], edge[1]))
457 edge = sorted([delTri.vertices[n,1], delTri.vertices[n,2]])
458 edges.add((edge[0], edge[1]))
462 G.add_edge(e[0],e[1])
464 d1=IMP.core.XYZ(particles[e[0]])
465 d2=IMP.core.XYZ(particles[e[1]])
466 dist=IMP.core.get_distance(d1,d2)
467 ts1 = IMP.core.Harmonic(dist+150, 0.001)
469 IMP.core.DistanceRestraint(m, ts1,
477 mc.set_scoring_function(rs)
478 mc.set_return_best(
False)
484 pos = nx.spring_layout(G ,pos=init_pos)
487 for i,xy in enumerate(data):
488 xynew = pos[i] * [xmax-xmin, ymax-ymin] + [xmin, ymin]
489 plt.plot((xy[0],xynew[0]), (xy[1],xynew[1]), 'k-')
490 x_list.append(xynew[0])
491 y_list.append(xynew[1])
492 #plt.annotate(name, xy, size=5, xycoords='data', xytext=xytext, textcoords='data', bbox=dict(boxstyle='round,pad=0.2', fc='yellow', alpha=0.3),\
493 # arrowprops=dict(arrowstyle="-", connectionstyle="arc3", color='gray'))
495 points=ax.scatter(x_list,y_list,alpha=1)
496 #pos*=[numpy.ptp([e[0] for e in data]), numpy.ptp([e[1] for e in data])]
502 mc.optimize(1000*len(particles.keys())*2)
503 print(rs.evaluate(
False))
506 for i,xy
in enumerate(data):
509 xynew = (coord[0],coord[1])
510 plt.plot((xy[0],xynew[0]), (xy[1],xynew[1]),
'k-')
511 x_list.append(xynew[0])
512 y_list.append(xynew[1])
513 points=ax.scatter(x_list,y_list,alpha=1,facecolors=
'none', edgecolors=
'k')
518 def show_mpld3(self,fig,ax,points,xl_list,xl_labels):
520 from mpld3
import plugins
527 border-collapse: collapse;
532 background-color: #ffffff;
536 background-color: #cccccc;
540 font-family:Arial, Helvetica, sans-serif;
541 border: 1px solid black;
546 df = pd.DataFrame(index=xl_labels)
548 sorted_keys=sorted(xl_list[0].keys())
550 for k
in sorted_keys:
551 df[k] = np.array([xl[k]
for xl
in xl_list])
554 for i
in range(len(xl_labels)):
555 label = df.ix[[i], :].T
557 labels.append(str(label.to_html()))
559 tooltip = plugins.PointHTMLTooltip(points, labels,
560 voffset=10, hoffset=10, css=css)
561 plugins.connect(fig, tooltip)
562 mpld3.save_html(fig,
"output.html")
565 def get_residue_contacts(self,prot_listx,prot_listy):
566 for px
in prot_listx:
567 for py
in prot_listy:
568 indexes_x = self.index_dict[px]
569 minx = min(indexes_x)
570 maxx = max(indexes_x)
571 indexes_y = self.index_dict[py]
572 miny = min(indexes_y)
573 maxy = max(indexes_y)
574 array = self.contact_freqs[minx:maxx,miny:maxy]
575 (xresidues,yresidues)=np.where(array>0)
576 for n,xr
in enumerate(xresidues):
577 print(xr,yresidues[n],px,py,array[xr,yresidues[n]])
583 confidence_info=
False,
585 display_residue_pairs=
False,
588 confidence_classes=
None,
590 scale_symbol_size=1.0,
591 gap_between_components=0,
592 dictionary_of_gaps={},
594 crosslink_threshold=
None,
599 color_crosslinks_by_distance=
True):
600 """ plot the xlink table with optional contact map.
601 prot_listx: list of protein names on the x-axis
602 prot_listy: list of protein names on the y-axis
603 no_dist_info: plot only the cross-links as grey spots
605 filter: list of tuples to filter on. each one contains:
606 keyword in the database to be filtered on
607 relationship ">","==","<"
609 example ("ID_Score",">",40)
610 display_residue_pairs: display all pairs defined in self.residue_pair_list
611 contactmap: display the contact map
612 filename: save to file (adds .pdf extension)
615 scale_symbol_size: rescale the symbol for the crosslink
616 gap_between_components: the numbeber of residues to leave blannk between each component
617 dictionary_of_gaps: add extra space after the given protein. dictionary_of_gaps={prot_name:extra_gap}
622 fig = plt.figure(figsize=(100, 100))
623 ax = fig.add_subplot(111)
627 if cbar_labels
is not None:
628 if len(cbar_labels)!=4:
629 print(
"to provide cbar labels, give 3 fields (first=first input file, last=last input) in oppose order of input contact maps")
632 if prot_listx
is None:
633 prot_listx = self.sequence_dict.keys()
635 nresx = gap_between_components + \
636 sum([len(self.sequence_dict[name])
637 + gap_between_components
for name
in prot_listx]) + \
638 sum([dictionary_of_gaps[prot]
for prot
in dictionary_of_gaps.keys()])
641 if prot_listy
is None:
642 prot_listy = self.sequence_dict.keys()
644 nresy = gap_between_components + \
645 sum([len(self.sequence_dict[name])
646 + gap_between_components
for name
in prot_listy]) + \
647 sum([dictionary_of_gaps[prot]
for prot
in dictionary_of_gaps.keys()])
652 res = gap_between_components
653 for prot
in prot_listx:
654 resoffsetx[prot] = res
655 res += len(self.sequence_dict[prot])
657 res += gap_between_components
658 if prot
in dictionary_of_gaps.keys(): res+= dictionary_of_gaps[prot]
662 res = gap_between_components
663 for prot
in prot_listy:
664 resoffsety[prot] = res
665 res += len(self.sequence_dict[prot])
667 res += gap_between_components
668 if prot
in dictionary_of_gaps.keys(): res+= dictionary_of_gaps[prot]
670 resoffsetdiagonal = {}
671 res = gap_between_components
672 for prot
in IMP.pmi.io.utilities.OrderedSet(prot_listx + prot_listy):
673 resoffsetdiagonal[prot] = res
674 res += len(self.sequence_dict[prot])
675 res += gap_between_components
676 if prot
in dictionary_of_gaps.keys(): res+= dictionary_of_gaps[prot]
681 for n, prot
in enumerate(prot_listx):
682 res = resoffsetx[prot]
684 for proty
in prot_listy:
685 resy = resoffsety[proty]
686 endy = resendy[proty]
687 ax.plot([res, res], [resy, endy], linestyle=
'-',color=
'gray', lw=0.4)
688 ax.plot([end, end], [resy, endy], linestyle=
'-',color=
'gray', lw=0.4)
689 xticks.append((float(res) + float(end)) / 2)
694 for n, prot
in enumerate(prot_listy):
695 res = resoffsety[prot]
697 for protx
in prot_listx:
698 resx = resoffsetx[protx]
699 endx = resendx[protx]
700 ax.plot([resx, endx], [res, res], linestyle=
'-',color=
'gray', lw=0.4)
701 ax.plot([resx, endx], [end, end], linestyle=
'-',color=
'gray', lw=0.4)
702 yticks.append((float(res) + float(end)) / 2)
707 tmp_array = np.zeros((nresx, nresy))
708 for px
in prot_listx:
709 for py
in prot_listy:
710 resx = resoffsetx[px]
711 lengx = resendx[px] - 1
712 resy = resoffsety[py]
713 lengy = resendy[py] - 1
714 indexes_x = self.index_dict[px]
715 minx = min(indexes_x)
716 maxx = max(indexes_x)
717 indexes_y = self.index_dict[py]
718 miny = min(indexes_y)
719 maxy = max(indexes_y)
720 tmp_array[resx:lengx,resy:lengy] = self.contact_freqs[minx:maxx,miny:maxy]
722 cax = ax.imshow(tmp_array,
727 interpolation=
'nearest')
729 ax.set_xticks(xticks)
730 ax.set_xticklabels(xlabels, rotation=90)
731 ax.set_yticks(yticks)
732 ax.set_yticklabels(ylabels)
733 plt.setp(ax.get_xticklabels(), fontsize=6)
734 plt.setp(ax.get_yticklabels(), fontsize=6)
737 already_added_xls = []
738 xl_coordinates_tuple_list = []
745 markersize = 5 * scale_symbol_size
746 if self.cross_link_db:
747 for xl
in self.cross_link_db:
749 (c1,c2,r1,r2)=IMP.pmi.io.crosslink._ProteinsResiduesArray(xl)
750 label=xl[self.cross_link_db.unique_sub_id_key]
751 if color_crosslinks_by_distance:
754 mdist=self._get_distance(r1,c1,r2,c2)
755 if mdist
is None:
continue
756 color = self._colormap_distance(mdist,threshold=crosslink_threshold)
763 ps=self._get_percentage_satisfaction(r1,c1,r2,c2)
764 if ps
is None:
continue
765 color = self._colormap_satisfaction(ps,threshold=0.2,tolerance=0.1)
770 pos1 = r1 + resoffsetx[c1]
774 pos2 = r2 + resoffsety[c2]
781 pos_for_diagonal1 = r1 + resoffsetdiagonal[c1]
782 pos_for_diagonal2 = r2 + resoffsetdiagonal[c2]
784 if confidence ==
'0.01':
785 markersize = 14 * scale_symbol_size
786 elif confidence ==
'0.05':
787 markersize = 9 * scale_symbol_size
788 elif confidence ==
'0.1':
789 markersize = 6 * scale_symbol_size
791 markersize = 15 * scale_symbol_size
793 markersize = 5 * scale_symbol_size
800 markersize=markersize)
807 markersize=markersize)
813 color_list.append(color)
814 color_list.append(color)
819 xl_labels.append(label)
820 xl_coordinates_tuple_list.append((float(pos1),float(pos2)))
821 xl_labels.append(label+
"*")
822 xl_coordinates_tuple_list.append((float(pos2),float(pos1)))
824 points=ax.scatter(x_list,y_list,s=markersize,c=color_list,alpha=alphablend)
827 if display_residue_pairs:
828 for rp
in self.residue_pair_list:
835 dist=self._get_distance(r1,c1,r2,c2)
842 pos1 = r1 + resoffsetx[c1]
846 pos2 = r2 + resoffsety[c2]
855 markersize=markersize)
858 fig.set_size_inches(0.002 * nresx, 0.002 * nresy)
859 [i.set_linewidth(2.0)
for i
in ax.spines.values()]
860 if cbar_labels
is not None:
861 cbar = fig.colorbar(cax, ticks=[0.5,1.5,2.5,3.5])
862 cbar.ax.set_yticklabels(cbar_labels)
866 points=self.spring_layout(ax,plt,xl_coordinates_tuple_list, xl_labels)
870 plt.savefig(filename, dpi=300,transparent=
"False")
874 self.show_mpld3(fig,ax,points,xl_list,xl_labels)
Utility classes and functions for IO.
atom::Hierarchies create_hierarchies(RMF::FileConstHandle fh, Model *m)
Set of Python classes to create a multi-state, multi-resolution IMP hierarchy.
Upper bound harmonic function (non-zero when feature > mean)
class to read, analyze, and plot xlink data on contact maps Canonical way to read the data: 1) load s...
static XYZR setup_particle(Model *m, ParticleIndex pi)
Utility classes and functions for reading and storing PMI files.
def plot_table
plot the xlink table with optional contact map.
def load_rmf_coordinates
read coordinates from a rmf file.
Handles cross-link data sets.
Return all close unordered pairs of particles taken from the SingletonContainer.
def load_pdb_coordinates
read coordinates from a pdb file.
Distance restraint between two particles.
double get_distance(XYZR a, XYZR b)
Compute the sphere distance between a and b.
void read_pdb(TextInput input, int model, Hierarchy h)
Object used to hold a set of restraints.
def setup_contact_map
loop through each distance map and get frequency of contacts
Class for storing model, its restraints, constraints, and particles.
def load_sequence_from_fasta_file
read sequence.
ParticleIndexPairs get_indexes(const ParticlePairsTemp &ps)
def set_residue_pairs_to_display
select the atom names of residue pairs to plot on the contact map list of residues types must be sing...
Ints get_index(const ParticlesTemp &particles, const Subset &subset, const Subsets &excluded)
Store a list of ParticleIndexes.
def deprecated_method
Python decorator to mark a method as deprecated.
A decorator for a particle representing an atom.
def load_crosslinks
read crosslinks from a CSV file.
void load_frame(RMF::FileConstHandle file, RMF::FrameID frame)
Load the given RMF frame into the state of the linked objects.
A decorator for a particle with x,y,z coordinates.
def get_particles_at_resolution_one
Get particles at res 1, or any beads, based on the name.
Tools for clustering and cluster analysis.
Modify a set of continuous variables using a normal distribution.
def deprecated_object
Python decorator to mark a class as deprecated.
Residue get_residue(Atom d, bool nothrow=false)
Return the Residue containing this atom.
Class to handle individual particles of a Model object.
Select all CA ATOM records.
static const FloatKeys & get_xyz_keys()
Get a vector containing the keys for x,y,z.
Python classes to represent, score, sample and analyze models.
A dictionary-like wrapper for reading and storing sequence data.
Functionality for loading, creating, manipulating and scoring atomic structures.
Select hierarchy particles identified by the biological name.
Support for the RMF file format for storing hierarchical molecular data and markup.
Applies a list of movers one at a time.
Applies a PairScore to each Pair in a list.
Perform more efficient close pair finding when rigid bodies are involved.
A decorator for a particle with x,y,z coordinates and a radius.