IMP logo
IMP Reference Guide  2.5.0
The Integrative Modeling Platform
tools.py
1 #!/usr/bin/env python
2 
3 """@namespace IMP.pmi.tools
4  Miscellaneous utilities.
5 """
6 
7 from __future__ import print_function
8 import IMP
9 import IMP.algebra
10 import collections
11 
12 def _get_restraint_set_key():
13  if not hasattr(_get_restraint_set_key, 'pmi_rs_key'):
14  _get_restraint_set_key.pmi_rs_key = IMP.ModelKey("PMI restraints")
15  return _get_restraint_set_key.pmi_rs_key
16 
17 def _add_restraint_set(model, mk):
18  rs = IMP.RestraintSet(model, "All PMI restraints")
19  model.add_data(mk, rs)
20  return rs
21 
22 def add_restraint_to_model(model, restraint):
23  """Add a PMI restraint to the model.
24  Since Model.add_restraint() no longer exists (in modern IMP restraints
25  should be added to a ScoringFunction instead) store them instead in
26  a RestraintSet, and keep a reference to it in the Model."""
27  mk = _get_restraint_set_key()
28  if model.get_has_data(mk):
29  rs = IMP.RestraintSet.get_from(model.get_data(mk))
30  else:
31  rs = _add_restraint_set(model, mk)
32  rs.add_restraint(restraint)
33 
34 def get_restraint_set(model):
35  """Get a RestraintSet containing all PMI restraints added to the model"""
36  mk = _get_restraint_set_key()
37  if not model.get_has_data(mk):
38  print("WARNING: no restraints added to model yet")
39  _add_restraint_set(model, mk)
40  return IMP.RestraintSet.get_from(model.get_data(mk))
41 
42 class Stopwatch(object):
43 
44  def __init__(self, isdelta=True):
45  global time
46  import time
47 
48  self.starttime = time.clock()
49  self.label = "None"
50  self.isdelta = isdelta
51 
52  def set_label(self, labelstr):
53  self.label = labelstr
54 
55  def get_output(self):
56  output = {}
57  if self.isdelta:
58  newtime = time.clock()
59  output[
60  "Stopwatch_" +
61  self.label +
62  "_delta_seconds"] = str(
63  newtime -
64  self.starttime)
65  self.starttime = newtime
66  else:
67  output["Stopwatch_" + self.label +
68  "_elapsed_seconds"] = str(time.clock() - self.starttime)
69  return output
70 
71 
72 class SetupNuisance(object):
73 
74  def __init__(self, m, initialvalue, minvalue, maxvalue, isoptimized=True):
75  import IMP.isd
76 
77  nuisance = IMP.isd.Scale.setup_particle(IMP.Particle(m), initialvalue)
78  if minvalue:
79  nuisance.set_lower(minvalue)
80  if maxvalue:
81  nuisance.set_upper(maxvalue)
82 
83  # m.add_score_state(IMP.core.SingletonConstraint(IMP.isd.NuisanceRangeModifier(),None,nuisance))
84  nuisance.set_is_optimized(nuisance.get_nuisance_key(), isoptimized)
85  self.nuisance = nuisance
86 
87  def get_particle(self):
88  return self.nuisance
89 
90 
91 class SetupWeight(object):
92 
93  def __init__(self, m, isoptimized=True):
94  import IMP.isd
95  pw = IMP.Particle(m)
96  self.weight = IMP.isd.Weight.setup_particle(pw)
97  self.weight.set_weights_are_optimized(True)
98 
99  def get_particle(self):
100  return self.weight
101 
102 
103 class ParticleToSampleFilter(object):
104  def __init__(self, sampled_objects):
105  self.sampled_objects=sampled_objects
106  self.filters=[]
107 
108  def add_filter(self,filter_string):
109  self.filters.append(filter_string)
110 
111  def get_particles_to_sample(self):
112  particles_to_sample={}
113  for so in self.sampled_objects:
114  ps_dict=so.get_particles_to_sample()
115  for key in ps_dict:
116  for f in self.filters:
117  if f in key:
118  if key not in particles_to_sample:
119  particles_to_sample[key]=ps_dict[key]
120  else:
121  particles_to_sample[key]+=ps_dict[key]
122  return particles_to_sample
123 
124 class ParticleToSampleList(object):
125 
126  def __init__(self, label="None"):
127 
128  self.dictionary_particle_type = {}
129  self.dictionary_particle_transformation = {}
130  self.dictionary_particle_name = {}
131  self.label = label
132 
133  def add_particle(
134  self,
135  particle,
136  particle_type,
137  particle_transformation,
138  name):
139  if not particle_type in ["Rigid_Bodies", "Floppy_Bodies", "Nuisances", "X_coord", "Weights"]:
140  raise TypeError("not the right particle type")
141  else:
142  self.dictionary_particle_type[particle] = particle_type
143  if particle_type == "Rigid_Bodies":
144  if type(particle_transformation) == tuple and len(particle_transformation) == 2 and type(particle_transformation[0]) == float and type(particle_transformation[1]) == float:
145  self.dictionary_particle_transformation[
146  particle] = particle_transformation
147  self.dictionary_particle_name[particle] = name
148  else:
149  raise TypeError("ParticleToSampleList: not the right transformation format for Rigid_Bodies, should be a tuple a floats")
150  else:
151  if type(particle_transformation) == float:
152  self.dictionary_particle_transformation[
153  particle] = particle_transformation
154  self.dictionary_particle_name[particle] = name
155  else:
156  raise TypeError("ParticleToSampleList: not the right transformation format sould be a float")
157 
158  def get_particles_to_sample(self):
159  ps = {}
160  for particle in self.dictionary_particle_type:
161  key = self.dictionary_particle_type[
162  particle] + "ParticleToSampleList_" + self.dictionary_particle_name[particle] + "_" + self.label
163  value = (
164  [particle],
165  self.dictionary_particle_transformation[particle])
166  ps[key] = value
167  return ps
168 
169 
170 class Variance(object):
171 
172  def __init__(self, model, tau, niter, prot, th_profile, write_data=False):
173  global sqrt, os, random
174  from math import sqrt
175  import os
176  import random
177 
178  self.model = model
179  self.write_data = write_data
180  self.tau = tau
181  self.niter = niter
182  #! select particles from the model
183  particles = IMP.atom.get_by_type(prot, IMP.atom.ATOM_TYPE)
184  self.particles = particles
185  # store reference coordinates and theoretical profile
186  self.refpos = [IMP.core.XYZ(p).get_coordinates() for p in particles]
187  self.model_profile = th_profile
188 
189  def perturb_particles(self, perturb=True):
190  for i, p in enumerate(self.particles):
191  newpos = array(self.refpos[i])
192  if perturb:
193  newpos += random.normal(0, self.tau, 3)
194  newpos = IMP.algebra.Vector3D(newpos)
195  IMP.core.XYZ(p).set_coordinates(newpos)
196 
197  def get_profile(self):
198  model_profile = self.model_profile
199  p = model_profile.calculate_profile(self.particles, IMP.saxs.CA_ATOMS)
200  return array([model_profile.get_intensity(i) for i in
201  range(model_profile.size())])
202 
203  def init_variances(self):
204  # create placeholders
205  N = self.model_profile.size()
206  a = self.profiles[0][:]
207  self.m = matrix(a).T # Nx1
208  self.V = self.m * self.m.T
209  self.normm = linalg.norm(self.m)
210  self.normV = linalg.norm(self.V)
211 
212  def update_variances(self):
213  a = matrix(self.profiles[-1]) # 1xN
214  n = float(len(self.profiles))
215  self.m = a.T / n + (n - 1) / n * self.m
216  self.V = a.T * a + self.V
217  self.oldnormm = self.normm
218  self.oldnormV = self.normV
219  self.normm = linalg.norm(self.m)
220  self.normV = linalg.norm(self.V)
221  self.diffm = (self.oldnormm - self.normm) / self.oldnormm
222  self.diffV = (self.oldnormV - self.normV) / self.oldnormV
223 
224  def get_direct_stats(self, a):
225  nq = len(a[0])
226  nprof = len(a)
227  m = [0] * nq
228  for prof in a:
229  for q, I in enumerate(prof):
230  m[q] += I
231  m = array(m) / nprof
232  V = matrix(a)
233  V = V.T * V
234  Sigma = (matrix(a - m))
235  Sigma = Sigma.T * Sigma / (nprof - 1)
236  mi = matrix(diag(1. / m))
237  Sigmarel = mi.T * Sigma * mi
238  return m, V, Sigma, Sigmarel
239 
240  def store_data(self):
241  if not os.path.isdir('data'):
242  os.mkdir('data')
243  profiles = matrix(self.profiles)
244  self.directm, self.directV, self.Sigma, self.Sigmarel = \
245  self.get_direct_stats(array(profiles))
246  directV = self.directV
247  # print "V comparison",(linalg.norm(directV-self.V)/self.normV)
248  save('data/profiles', profiles)
249  # absolute profile differences
250  fl = open('data/profiles.dat', 'w')
251  for i, l in enumerate(array(profiles).T):
252  self.model_profile.get_q(i)
253  fl.write('%s ' % i)
254  for k in l:
255  fl.write('%s ' % (k - self.directm[i]))
256  fl.write('\n')
257  # relative profile differences
258  fl = open('data/profiles_rel.dat', 'w')
259  for i, l in enumerate(array(profiles).T):
260  self.model_profile.get_q(i)
261  fl.write('%s ' % i)
262  for k in l:
263  fl.write('%s ' % ((k - self.directm[i]) / self.directm[i]))
264  fl.write('\n')
265  save('data/m', self.directm)
266  save('data/V', self.directV)
267  Sigma = self.Sigma
268  save('data/Sigma', Sigma)
269  # Sigma matrix
270  fl = open('data/Sigma.dat', 'w')
271  model_profile = self.model_profile
272  for i in range(model_profile.size()):
273  qi = model_profile.get_q(i)
274  for j in range(model_profile.size()):
275  qj = model_profile.get_q(j)
276  vij = self.Sigma[i, j]
277  fl.write('%s %s %s\n' % (qi, qj, vij))
278  fl.write('\n')
279  # Sigma eigenvalues
280  fl = open('data/eigenvals', 'w')
281  for i in linalg.eigvalsh(Sigma):
282  fl.write('%s\n' % i)
283  Sigmarel = self.Sigmarel
284  save('data/Sigmarel', Sigmarel)
285  # Sigmarel matrix
286  fl = open('data/Sigmarel.dat', 'w')
287  model_profile = self.model_profile
288  for i in range(model_profile.size()):
289  qi = model_profile.get_q(i)
290  for j in range(model_profile.size()):
291  qj = model_profile.get_q(j)
292  vij = self.Sigmarel[i, j]
293  fl.write('%s %s %s\n' % (qi, qj, vij))
294  fl.write('\n')
295  # Sigma eigenvalues
296  fl = open('data/eigenvals_rel', 'w')
297  for i in linalg.eigvalsh(Sigmarel):
298  fl.write('%s\n' % i)
299  # mean profile
300  fl = open('data/mean.dat', 'w')
301  for i in range(len(self.directm)):
302  qi = self.model_profile.get_q(i)
303  fl.write('%s ' % qi)
304  fl.write('%s ' % self.directm[i])
305  fl.write('%s ' % sqrt(self.Sigma[i, i]))
306  fl.write('\n')
307 
308  def try_chol(self, jitter):
309  Sigma = self.Sigma
310  try:
311  linalg.cholesky(Sigma + matrix(eye(len(Sigma))) * jitter)
312  except linalg.LinAlgError:
313  print("Decomposition failed with jitter =", jitter)
314  return
315  print("Successful decomposition with jitter =", jitter)
316 
317  def run(self):
318  self.profiles = [self.get_profile()]
319  # self.init_variances()
320  for n in range(self.niter):
321  self.perturb_particles()
322  self.profiles.append(self.get_profile())
323  # self.update_variances()
324  #profiles = matrix(self.profiles)
325  # print n,self.diffm,self.diffV
326  # print
327  #
328  if self.write_data:
329  self.store_data()
330  # self.try_chol(0.)
331  # for i in logspace(-7,0,num=8):
332  # self.try_chol(i)
333 
334  def get_cov(self, relative=True):
335  if not relative:
336  return self.Sigma
337  else:
338  return self.Sigmarel
339 
340  #-------------------------------
341 
342 def get_random_cross_link_dataset(representation,
343  resolution=1.0,
344  number_of_cross_links=10,
345  ambiguity_probability=0.1,
346  confidence_score_range=[0,100],
347  avoid_same_particles=False):
348  '''Return a random cross-link dataset as a string.
349  Every line is a residue pair, together with UniqueIdentifier
350  and XL score.'''
351 
352  residue_pairs=get_random_residue_pairs(representation, resolution, number_of_cross_links, avoid_same_particles)
353 
354  from random import random
355  unique_identifier=0
356  cmin=float(min(confidence_score_range))
357  cmax=float(max(confidence_score_range))
358 
359  dataset="#\n"
360 
361  for (name1, r1, name2, r2) in residue_pairs:
362  if random() > ambiguity_probability:
363  unique_identifier+=1
364  score=random()*(cmax-cmin)+cmin
365  dataset+=str(name1)+" "+str(name2)+" "+str(r1)+" "+str(r2)+" "+str(score)+" "+str(unique_identifier)+"\n"
366 
367  return dataset
368 
369 
370  #-------------------------------
371 
372 def get_cross_link_data(directory, filename, dist, omega, sigma,
373  don=None, doff=None, prior=0, type_of_profile="gofr"):
374 
375  (distmin, distmax, ndist) = dist
376  (omegamin, omegamax, nomega) = omega
377  (sigmamin, sigmamax, nsigma) = sigma
378  import IMP.isd
379 
380  filen = IMP.isd.get_data_path("CrossLinkPMFs.dict")
381  xlpot = open(filen)
382 
383  for line in xlpot:
384  dictionary = eval(line)
385  break
386 
387  xpot = dictionary[directory][filename]["distance"]
388  pot = dictionary[directory][filename][type_of_profile]
389 
390  dist_grid = get_grid(distmin, distmax, ndist, False)
391  omega_grid = get_log_grid(omegamin, omegamax, nomega)
392  sigma_grid = get_log_grid(sigmamin, sigmamax, nsigma)
393 
394  if not don is None and not doff is None:
395  xlmsdata = IMP.isd.CrossLinkData(
396  dist_grid,
397  omega_grid,
398  sigma_grid,
399  xpot,
400  pot,
401  don,
402  doff,
403  prior)
404  else:
405  xlmsdata = IMP.isd.CrossLinkData(
406  dist_grid,
407  omega_grid,
408  sigma_grid,
409  xpot,
410  pot)
411  return xlmsdata
412 
413  #-------------------------------
414 
415 
416 def get_cross_link_data_from_length(length, xxx_todo_changeme3, xxx_todo_changeme4, xxx_todo_changeme5):
417  (distmin, distmax, ndist) = xxx_todo_changeme3
418  (omegamin, omegamax, nomega) = xxx_todo_changeme4
419  (sigmamin, sigmamax, nsigma) = xxx_todo_changeme5
420  import IMP.isd
421 
422  dist_grid = get_grid(distmin, distmax, ndist, False)
423  omega_grid = get_log_grid(omegamin, omegamax, nomega)
424  sigma_grid = get_log_grid(sigmamin, sigmamax, nsigma)
425 
426  xlmsdata = IMP.isd.CrossLinkData(dist_grid, omega_grid, sigma_grid, length)
427  return xlmsdata
428 
429 
430 def get_grid(gmin, gmax, ngrid, boundaries):
431  grid = []
432  dx = (gmax - gmin) / float(ngrid)
433  for i in range(0, ngrid + 1):
434  if(not boundaries and i == 0):
435  continue
436  if(not boundaries and i == ngrid):
437  continue
438  grid.append(gmin + float(i) * dx)
439  return grid
440 
441  #-------------------------------
442 
443 
444 def get_log_grid(gmin, gmax, ngrid):
445  from math import exp, log
446  grid = []
447  for i in range(0, ngrid + 1):
448  grid.append(gmin * exp(float(i) / ngrid * log(gmax / gmin)))
449  return grid
450 
451  #-------------------------------
452 
453 
455  '''
456  example '"{ID_Score}" > 28 AND "{Sample}" ==
457  "%10_1%" OR ":Sample}" == "%10_2%" OR ":Sample}"
458  == "%10_3%" OR ":Sample}" == "%8_1%" OR ":Sample}" == "%8_2%"'
459  '''
460 
461  import pyparsing as pp
462 
463  operator = pp.Regex(">=|<=|!=|>|<|==|in").setName("operator")
464  value = pp.QuotedString(
465  '"') | pp.Regex(
466  r"[+-]?\d+(:?\.\d*)?(:?[eE][+-]?\d+)?")
467  identifier = pp.Word(pp.alphas, pp.alphanums + "_")
468  comparison_term = identifier | value
469  condition = pp.Group(comparison_term + operator + comparison_term)
470 
471  expr = pp.operatorPrecedence(condition, [
472  ("OR", 2, pp.opAssoc.LEFT, ),
473  ("AND", 2, pp.opAssoc.LEFT, ),
474  ])
475 
476  parsedstring = str(expr.parseString(inputstring)) \
477  .replace("[", "(") \
478  .replace("]", ")") \
479  .replace(",", " ") \
480  .replace("'", " ") \
481  .replace("%", "'") \
482  .replace("{", "float(entry['") \
483  .replace("}", "'])") \
484  .replace(":", "str(entry['") \
485  .replace("}", "'])") \
486  .replace("AND", "and") \
487  .replace("OR", "or")
488  return parsedstring
489 
490 
491 def open_file_or_inline_text(filename):
492  try:
493  fl = open(filename, "r")
494  except IOError:
495  fl = filename.split("\n")
496  return fl
497 
498 
499 def get_drmsd(prot0, prot1):
500  drmsd = 0.
501  npairs = 0.
502  for i in range(0, len(prot0) - 1):
503  for j in range(i + 1, len(prot0)):
504  dist0 = IMP.core.get_distance(prot0[i], prot0[j])
505  dist1 = IMP.core.get_distance(prot1[i], prot1[j])
506  drmsd += (dist0 - dist1) ** 2
507  npairs += 1.
508  return math.sqrt(drmsd / npairs)
509 
510  #-------------------------------
511 
512 
513 def get_ids_from_fasta_file(fastafile):
514  ids = []
515  with open(fastafile) as ff:
516  for l in ff:
517  if l[0] == ">":
518  ids.append(l[1:-1])
519  return ids
520 
521 
522 def get_closest_residue_position(hier, resindex, terminus="N"):
523  '''
524  this function works with plain hierarchies, as read from the pdb,
525  no multi-scale hierarchies
526  '''
527  p = []
528  niter = 0
529  while len(p) == 0:
530  niter += 1
531  sel = IMP.atom.Selection(hier, residue_index=resindex,
532  atom_type=IMP.atom.AT_CA)
533 
534  if terminus == "N":
535  resindex += 1
536  if terminus == "C":
537  resindex -= 1
538 
539  if niter >= 10000:
540  print("get_closest_residue_position: exiting while loop without result")
541  break
542  p = sel.get_selected_particles()
543 
544  if len(p) == 1:
545  return IMP.core.XYZ(p[0]).get_coordinates()
546  elif len(p) == 0:
547  print("get_closest_residue_position: got NO residues for hierarchy %s and residue %i" % (hier, resindex))
548  raise Exception("get_closest_residue_position: got NO residues for hierarchy %s and residue %i" % (
549  hier, resindex))
550  else:
551  raise ValueError("got multiple residues for hierarchy %s and residue %i; the list of particles is %s" % (hier, resindex, str([pp.get_name() for pp in p])))
552 
553 def get_position_terminal_residue(hier, terminus="C", resolution=1):
554  '''
555  Get the xyz position of the terminal residue at the given resolution.
556  @param hier hierarchy containing the terminal residue
557  @param terminus either 'N' or 'C'
558  @param resolution resolution to use.
559  '''
560  termresidue = None
561  termparticle = None
562  for p in IMP.atom.get_leaves(hier):
563  if IMP.pmi.Resolution(p).get_resolution() == resolution:
565  if terminus == "C":
566  if max(residues) >= termresidue and not termresidue is None:
567  termresidue = max(residues)
568  termparticle = p
569  elif termresidue is None:
570  termresidue = max(residues)
571  termparticle = p
572  elif terminus == "N":
573  if min(residues) <= termresidue and not termresidue is None:
574  termresidue = min(residues)
575  termparticle = p
576  elif termresidue is None:
577  termresidue = min(residues)
578  termparticle = p
579  else:
580  raise ValueError("terminus argument should be either N or C")
581 
582  return IMP.core.XYZ(termparticle).get_coordinates()
583 
584 
585 def get_residue_gaps_in_hierarchy(hierarchy, start, end):
586  '''
587  Return the residue index gaps and contiguous segments in the hierarchy.
588 
589  @param hierarchy hierarchy to examine
590  @param start first residue index
591  @param end last residue index
592 
593  @return A list of lists of the form
594  [[1,100,"cont"],[101,120,"gap"],[121,200,"cont"]]
595  '''
596  gaps = []
597  for n, rindex in enumerate(range(start, end + 1)):
598  sel = IMP.atom.Selection(hierarchy, residue_index=rindex,
599  atom_type=IMP.atom.AT_CA)
600 
601  if len(sel.get_selected_particles()) == 0:
602  if n == 0:
603  # set the initial condition
604  rindexgap = start
605  rindexcont = start - 1
606  if rindexgap == rindex - 1:
607  # residue is contiguous with the previously discovered gap
608  gaps[-1][1] += 1
609  else:
610  # residue is not contiguous with the previously discovered gap
611  # hence create a new gap tuple
612  gaps.append([rindex, rindex, "gap"])
613  # update the index of the last residue gap
614  rindexgap = rindex
615  else:
616  if n == 0:
617  # set the initial condition
618  rindexgap = start - 1
619  rindexcont = start
620  if rindexcont == rindex - 1:
621  # residue is contiguous with the previously discovered
622  # continuous part
623  gaps[-1][1] += 1
624  else:
625  # residue is not contiguous with the previously discovered continuous part
626  # hence create a new cont tuple
627  gaps.append([rindex, rindex, "cont"])
628  # update the index of the last residue gap
629  rindexcont = rindex
630  return gaps
631 
632 
633 class map(object):
634 
635  def __init__(self):
636  self.map = {}
637 
638  def set_map_element(self, xvalue, yvalue):
639  self.map[xvalue] = yvalue
640 
641  def get_map_element(self, invalue):
642  if type(invalue) == float:
643  n = 0
644  mindist = 1
645  for x in self.map:
646  dist = (invalue - x) * (invalue - x)
647 
648  if n == 0:
649  mindist = dist
650  minx = x
651  if dist < mindist:
652  mindist = dist
653  minx = x
654  n += 1
655  return self.map[minx]
656  elif type(invalue) == str:
657  return self.map[invalue]
658  else:
659  raise TypeError("wrong type for map")
660 
661 
662 def select(representation,
663  resolution=None,
664  hierarchies=None,
665  selection_arguments=None,
666  name=None,
667  name_is_ambiguous=False,
668  first_residue=None,
669  last_residue=None,
670  residue=None,
671  representation_type=None):
672  '''
673  this function uses representation=SimplifiedModel
674  it returns the corresponding selected particles
675  representation_type="Beads", "Res:X", "Densities", "Representation", "Molecule"
676  '''
677 
678  if resolution is None:
679  allparticles = IMP.atom.get_leaves(representation.prot)
680  resolution_particles = None
681  hierarchies_particles = None
682  names_particles = None
683  residue_range_particles = None
684  residue_particles = None
685  representation_type_particles = None
686 
687  if not resolution is None:
688  resolution_particles = []
689  hs = representation.get_hierarchies_at_given_resolution(resolution)
690  for h in hs:
691  resolution_particles += IMP.atom.get_leaves(h)
692 
693  if not hierarchies is None:
694  hierarchies_particles = []
695  for h in hierarchies:
696  hierarchies_particles += IMP.atom.get_leaves(h)
697 
698  if not name is None:
699  names_particles = []
700  if name_is_ambiguous:
701  for namekey in representation.hier_dict:
702  if name in namekey:
703  names_particles += IMP.atom.get_leaves(
704  representation.hier_dict[namekey])
705  elif name in representation.hier_dict:
706  names_particles += IMP.atom.get_leaves(representation.hier_dict[name])
707  else:
708  print("select: component %s is not there" % name)
709 
710  if not first_residue is None and not last_residue is None:
711  sel = IMP.atom.Selection(representation.prot,
712  residue_indexes=range(first_residue, last_residue + 1))
713  residue_range_particles = [IMP.atom.Hierarchy(p)
714  for p in sel.get_selected_particles()]
715 
716  if not residue is None:
717  sel = IMP.atom.Selection(representation.prot, residue_index=residue)
718  residue_particles = [IMP.atom.Hierarchy(p)
719  for p in sel.get_selected_particles()]
720 
721  if not representation_type is None:
722  representation_type_particles = []
723  if representation_type == "Molecule":
724  for name in representation.hier_representation:
725  for repr_type in representation.hier_representation[name]:
726  if repr_type == "Beads" or "Res:" in repr_type:
727  h = representation.hier_representation[name][repr_type]
728  representation_type_particles += IMP.atom.get_leaves(h)
729 
730  elif representation_type == "PDB":
731  for name in representation.hier_representation:
732  for repr_type in representation.hier_representation[name]:
733  if repr_type == "Res:" in repr_type:
734  h = representation.hier_representation[name][repr_type]
735  representation_type_particles += IMP.atom.get_leaves(h)
736 
737  else:
738  for name in representation.hier_representation:
739  h = representation.hier_representation[
740  name][
741  representation_type]
742  representation_type_particles += IMP.atom.get_leaves(h)
743 
744  selections = [hierarchies_particles, names_particles,
745  residue_range_particles, residue_particles, representation_type_particles]
746 
747  if resolution is None:
748  selected_particles = set(allparticles)
749  else:
750  selected_particles = set(resolution_particles)
751 
752  for s in selections:
753  if not s is None:
754  selected_particles = (set(s) & selected_particles)
755 
756  return list(selected_particles)
757 
758 
759 def select_by_tuple(
760  representation,
761  tupleselection,
762  resolution=None,
763  name_is_ambiguous=False):
764  if isinstance(tupleselection, tuple) and len(tupleselection) == 3:
765  particles = IMP.pmi.tools.select(representation, resolution=resolution,
766  name=tupleselection[2],
767  first_residue=tupleselection[0],
768  last_residue=tupleselection[1],
769  name_is_ambiguous=name_is_ambiguous)
770  elif isinstance(tupleselection, str):
771  particles = IMP.pmi.tools.select(representation, resolution=resolution,
772  name=tupleselection,
773  name_is_ambiguous=name_is_ambiguous)
774  else:
775  raise ValueError('you passed something bad to select_by_tuple()')
776  # now order the result by residue number
777  particles = IMP.pmi.tools.sort_by_residues(particles)
778 
779  return particles
780 
781 
782 def get_db_from_csv(csvfilename):
783  import csv
784  outputlist = []
785  csvr = csv.DictReader(open(csvfilename, "rU"))
786  for l in csvr:
787  outputlist.append(l)
788  return outputlist
789 
790 
791 class HierarchyDatabase(object):
792  """Store the representations for a system."""
793 
794  def __init__(self):
795  self.db = {}
796  # this dictionary map a particle to its root hierarchy
797  self.root_hierarchy_dict = {}
798  self.preroot_fragment_hierarchy_dict = {}
799  self.particle_to_name = {}
800  self.model = None
801 
802  def add_name(self, name):
803  if name not in self.db:
804  self.db[name] = {}
805 
806  def add_residue_number(self, name, resn):
807  resn = int(resn)
808  self.add_name(name)
809  if resn not in self.db[name]:
810  self.db[name][resn] = {}
811 
812  def add_resolution(self, name, resn, resolution):
813  resn = int(resn)
814  resolution = float(resolution)
815  self.add_name(name)
816  self.add_residue_number(name, resn)
817  if resolution not in self.db[name][resn]:
818  self.db[name][resn][resolution] = []
819 
820  def add_particles(self, name, resn, resolution, particles):
821  resn = int(resn)
822  resolution = float(resolution)
823  self.add_name(name)
824  self.add_residue_number(name, resn)
825  self.add_resolution(name, resn, resolution)
826  self.db[name][resn][resolution] += particles
827  for p in particles:
828  (rh, prf) = self.get_root_hierarchy(p)
829  self.root_hierarchy_dict[p] = rh
830  self.preroot_fragment_hierarchy_dict[p] = prf
831  self.particle_to_name[p] = name
832  if self.model is None:
833  self.model = particles[0].get_model()
834 
835  def get_model(self):
836  return self.model
837 
838  def get_names(self):
839  names = list(self.db.keys())
840  names.sort()
841  return names
842 
843  def get_particles(self, name, resn, resolution):
844  resn = int(resn)
845  resolution = float(resolution)
846  return self.db[name][resn][resolution]
847 
848  def get_particles_at_closest_resolution(self, name, resn, resolution):
849  resn = int(resn)
850  resolution = float(resolution)
851  closestres = min(self.get_residue_resolutions(name, resn),
852  key=lambda x: abs(float(x) - float(resolution)))
853  return self.get_particles(name, resn, closestres)
854 
855  def get_residue_resolutions(self, name, resn):
856  resn = int(resn)
857  resolutions = list(self.db[name][resn].keys())
858  resolutions.sort()
859  return resolutions
860 
861  def get_molecule_resolutions(self, name):
862  resolutions = set()
863  for resn in self.db[name]:
864  resolutions.update(list(self.db[name][resn].keys()))
865  resolutions.sort()
866  return resolutions
867 
868  def get_residue_numbers(self, name):
869  residue_numbers = list(self.db[name].keys())
870  residue_numbers.sort()
871  return residue_numbers
872 
873  def get_particles_by_resolution(self, name, resolution):
874  resolution = float(resolution)
875  particles = []
876  for resn in self.get_residue_numbers(name):
877  result = self.get_particles_at_closest_resolution(
878  name,
879  resn,
880  resolution)
881  pstemp = [p for p in result if p not in particles]
882  particles += pstemp
883  return particles
884 
885  def get_all_particles_by_resolution(self, resolution):
886  resolution = float(resolution)
887  particles = []
888  for name in self.get_names():
889  particles += self.get_particles_by_resolution(name, resolution)
890  return particles
891 
892  def get_root_hierarchy(self, particle):
893  prerootfragment = particle
894  while IMP.atom.Atom.get_is_setup(particle) or \
895  IMP.atom.Residue.get_is_setup(particle) or \
897  if IMP.atom.Atom.get_is_setup(particle):
898  p = IMP.atom.Atom(particle).get_parent()
899  elif IMP.atom.Residue.get_is_setup(particle):
900  p = IMP.atom.Residue(particle).get_parent()
901  elif IMP.atom.Fragment.get_is_setup(particle):
902  p = IMP.atom.Fragment(particle).get_parent()
903  prerootfragment = particle
904  particle = p
905  return (
906  (IMP.atom.Hierarchy(particle), IMP.atom.Hierarchy(prerootfragment))
907  )
908 
909  def get_all_root_hierarchies_by_resolution(self, resolution):
910  hierarchies = []
911  resolution = float(resolution)
912  particles = self.get_all_particles_by_resolution(resolution)
913  for p in particles:
914  rh = self.root_hierarchy_dict[p]
915  if rh not in hierarchies:
916  hierarchies.append(IMP.atom.Hierarchy(rh))
917  return hierarchies
918 
919  def get_preroot_fragments_by_resolution(self, name, resolution):
920  fragments = []
921  resolution = float(resolution)
922  particles = self.get_particles_by_resolution(name, resolution)
923  for p in particles:
924  fr = self.preroot_fragment_hierarchy_dict[p]
925  if fr not in fragments:
926  fragments.append(fr)
927  return fragments
928 
929  def show(self, name):
930  print(name)
931  for resn in self.get_residue_numbers(name):
932  print(resn)
933  for resolution in self.get_residue_resolutions(name, resn):
934  print("----", resolution)
935  for p in self.get_particles(name, resn, resolution):
936  print("--------", p.get_name())
937 
938 
939 def get_prot_name_from_particle(p, list_of_names):
940  '''Return the component name provided a particle and a list of names'''
941  root = p
942  protname = root.get_name()
943  is_a_bead = False
944  while not protname in list_of_names:
945  root0 = root.get_parent()
946  if root0 == IMP.atom.Hierarchy():
947  return (None, None)
948  protname = root0.get_name()
949 
950  # check if that is a bead
951  # this piece of code might be dangerous if
952  # the hierarchy was called Bead :)
953  if "Beads" in protname:
954  is_a_bead = True
955  root = root0
956  return (protname, is_a_bead)
957 
958 
960  '''
961  Retrieve the residue indexes for the given particle.
962 
963  The particle must be an instance of Fragment,Residue or Atom
964  '''
965  resind = []
967  resind = IMP.atom.Fragment(hier).get_residue_indexes()
969  resind = [IMP.atom.Residue(hier).get_index()]
970  elif IMP.atom.Atom.get_is_setup(hier):
971  a = IMP.atom.Atom(hier)
972  resind = [IMP.atom.Residue(a.get_parent()).get_index()]
973  else:
974  print("get_residue_indexes> input is not Fragment, Residue or Atom")
975  resind = []
976  return resind
977 
978 
979 def sort_by_residues(particles):
980  particles_residues = [(p, IMP.pmi.tools.get_residue_indexes(p))
981  for p in particles]
982  sorted_particles_residues = sorted(
983  particles_residues,
984  key=lambda tup: tup[1])
985  particles = [p[0] for p in sorted_particles_residues]
986  return particles
987 
988 
989 def get_residue_to_particle_map(particles):
990  # this function returns a dictionary that map particles to residue indexes
991  particles = sort_by_residues(particles)
992  particles_residues = [(p, IMP.pmi.tools.get_residue_indexes(p))
993  for p in particles]
994  return dict(zip(particles_residues, particles))
995 
996 #
997 # Parallel Computation
998 #
999 
1000 
1002  """Synchronize data over a parallel run"""
1003  from mpi4py import MPI
1004  comm = MPI.COMM_WORLD
1005  rank = comm.Get_rank()
1006  number_of_processes = comm.size
1007  comm.Barrier()
1008  if rank != 0:
1009  comm.send(data, dest=0, tag=11)
1010 
1011  elif rank == 0:
1012  for i in range(1, number_of_processes):
1013  data_tmp = comm.recv(source=i, tag=11)
1014  if type(data) == list:
1015  data += data_tmp
1016  elif type(data) == dict:
1017  data.update(data_tmp)
1018  else:
1019  raise TypeError("data not supported, use list or dictionaries")
1020 
1021  for i in range(1, number_of_processes):
1022  comm.send(data, dest=i, tag=11)
1023 
1024  if rank != 0:
1025  data = comm.recv(source=0, tag=11)
1026  return data
1027 
1029  """Synchronize data over a parallel run"""
1030  from mpi4py import MPI
1031  comm = MPI.COMM_WORLD
1032  rank = comm.Get_rank()
1033  number_of_processes = comm.size
1034  comm.Barrier()
1035  if rank != 0:
1036  comm.send(data, dest=0, tag=11)
1037  elif rank == 0:
1038  for i in range(1, number_of_processes):
1039  data_tmp = comm.recv(source=i, tag=11)
1040  for k in data:
1041  data[k]+=data_tmp[k]
1042 
1043  for i in range(1, number_of_processes):
1044  comm.send(data, dest=i, tag=11)
1045 
1046  if rank != 0:
1047  data = comm.recv(source=0, tag=11)
1048  return data
1049 
1050 
1051 #
1052 ### Lists and iterators
1053 #
1054 
1055 def sublist_iterator(l, lmin=None, lmax=None):
1056  '''
1057  Yield all sublists of length >= lmin and <= lmax
1058  '''
1059  if lmin is None:
1060  lmin = 0
1061  if lmax is None:
1062  lmax = len(l)
1063  n = len(l) + 1
1064  for i in range(n):
1065  for j in range(i + 1, n):
1066  if len(l[i:j]) <= lmax and len(l[i:j]) >= lmin:
1067  yield l[i:j]
1068 
1069 
1070 def flatten_list(l):
1071  return [item for sublist in l for item in sublist]
1072 
1073 
1074 def list_chunks_iterator(list, length):
1075  """ Yield successive length-sized chunks from a list.
1076  """
1077  for i in range(0, len(list), length):
1078  yield list[i:i + length]
1079 
1080 
1081 def chunk_list_into_segments(seq, num):
1082  seq = list(seq)
1083  avg = len(seq) / float(num)
1084  out = []
1085  last = 0.0
1086 
1087  while last < len(seq):
1088  out.append(seq[int(last):int(last + avg)])
1089  last += avg
1090 
1091  return out
1092 
1093 #
1094 # COORDINATE MANIPULATION
1095 #
1096 
1097 
1098 def translate_hierarchy(hierarchy, translation_vector):
1099  '''
1100  Apply a translation to a hierarchy along the input vector.
1101  '''
1102  rbs = set()
1103  xyzs = set()
1104  if type(translation_vector) == list:
1105  transformation = IMP.algebra.Transformation3D(
1106  IMP.algebra.Vector3D(translation_vector))
1107  else:
1108  transformation = IMP.algebra.Transformation3D(translation_vector)
1109  for p in IMP.atom.get_leaves(hierarchy):
1111  rbs.add(IMP.core.RigidBody(p))
1113  rb = IMP.core.RigidMember(p).get_rigid_body()
1114  rbs.add(rb)
1115  else:
1116  xyzs.add(p)
1117  for xyz in xyzs:
1118  IMP.core.transform(IMP.core.XYZ(xyz), transformation)
1119  for rb in rbs:
1120  IMP.core.transform(rb, transformation)
1121 
1122 
1123 def translate_hierarchies(hierarchies, translation_vector):
1124  for h in hierarchies:
1125  IMP.pmi.tools.translate_hierarchy(h, translation_vector)
1126 
1127 
1128 def translate_hierarchies_to_reference_frame(hierarchies):
1129  xc = 0
1130  yc = 0
1131  zc = 0
1132  nc = 0
1133  for h in hierarchies:
1134  for p in IMP.atom.get_leaves(h):
1135  coor = IMP.core.XYZ(p).get_coordinates()
1136  nc += 1
1137  xc += coor[0]
1138  yc += coor[1]
1139  zc += coor[2]
1140  xc = xc / nc
1141  yc = yc / nc
1142  zc = zc / nc
1143  IMP.pmi.tools.translate_hierarchies(hierarchies, (-xc, -yc, -zc))
1144 
1145 
1146 #
1147 # Tools to simulate data
1148 #
1149 
1150 def normal_density_function(expected_value, sigma, x):
1151  return (
1152  1 / math.sqrt(2 * math.pi) / sigma *
1153  math.exp(-(x - expected_value) ** 2 / 2 / sigma / sigma)
1154  )
1155 
1156 
1157 def log_normal_density_function(expected_value, sigma, x):
1158  return (
1159  1 / math.sqrt(2 * math.pi) / sigma / x *
1160  math.exp(-(math.log(x / expected_value) ** 2 / 2 / sigma / sigma))
1161  )
1162 
1163 
1164 def get_random_residue_pairs(representation, resolution,
1165  number,
1166  avoid_same_particles=False,
1167  names=None):
1168 
1169  from random import choice
1170  particles = []
1171  if names is None:
1172  names=list(representation.hier_dict.keys())
1173 
1174  for name in names:
1175  prot = representation.hier_dict[name]
1176  particles += select(representation,name=name,resolution=resolution)
1177  random_residue_pairs = []
1178  while len(random_residue_pairs)<=number:
1179  p1 = choice(particles)
1180  p2 = choice(particles)
1181  if p1==p2 and avoid_same_particles: continue
1182  r1 = choice(IMP.pmi.tools.get_residue_indexes(p1))
1183  r2 = choice(IMP.pmi.tools.get_residue_indexes(p2))
1184  name1 = representation.get_prot_name_from_particle(p1)
1185  name2 = representation.get_prot_name_from_particle(p2)
1186  random_residue_pairs.append((name1, r1, name2, r2))
1187 
1188  return random_residue_pairs
1189 
1190 
1191 def get_random_data_point(
1192  expected_value,
1193  ntrials,
1194  sensitivity,
1195  sigma,
1196  outlierprob,
1197  begin_end_nbins_tuple,
1198  log=False,
1199  loggrid=False):
1200  import random
1201 
1202  begin = begin_end_nbins_tuple[0]
1203  end = begin_end_nbins_tuple[1]
1204  nbins = begin_end_nbins_tuple[2]
1205 
1206  if not loggrid:
1207  fmod_grid = get_grid(begin, end, nbins, True)
1208  else:
1209  fmod_grid = get_log_grid(begin, end, nbins)
1210 
1211  norm = 0
1212  cumul = []
1213  cumul.append(0)
1214 
1215  a = []
1216  for i in range(0, ntrials):
1217  a.append([random.random(), True])
1218 
1219  if sigma != 0.0:
1220  for j in range(1, len(fmod_grid)):
1221  fj = fmod_grid[j]
1222  fjm1 = fmod_grid[j - 1]
1223  df = fj - fjm1
1224 
1225  if not log:
1226  pj = normal_density_function(expected_value, sigma, fj)
1227  pjm1 = normal_density_function(expected_value, sigma, fjm1)
1228  else:
1229  pj = log_normal_density_function(expected_value, sigma, fj)
1230  pjm1 = log_normal_density_function(expected_value, sigma, fjm1)
1231 
1232  norm += (pj + pjm1) / 2.0 * df
1233  cumul.append(norm)
1234  # print fj, pj
1235 
1236  random_points = []
1237 
1238  for i in range(len(cumul)):
1239  # print i,a, cumul[i], norm
1240  for aa in a:
1241  if (aa[0] <= cumul[i] / norm and aa[1]):
1242  random_points.append(
1243  int(fmod_grid[i] / sensitivity) * sensitivity)
1244  aa[1] = False
1245 
1246  else:
1247  random_points = [expected_value] * ntrials
1248 
1249  for i in range(len(random_points)):
1250  if random.random() < outlierprob:
1251  a = random.uniform(begin, end)
1252  random_points[i] = int(a / sensitivity) * sensitivity
1253  print(random_points)
1254  '''
1255  for i in range(ntrials):
1256  if random.random() > OUTLIERPROB_:
1257  r=truncnorm.rvs(0.0,1.0,expected_value,BETA_)
1258  if r>1.0: print r,expected_value,BETA_
1259  else:
1260  r=random.random()
1261  random_points.append(int(r/sensitivity)*sensitivity)
1262  '''
1263 
1264  rmean = 0.
1265  rmean2 = 0.
1266  for r in random_points:
1267  rmean += r
1268  rmean2 += r * r
1269 
1270  rmean /= float(ntrials)
1271  rmean2 /= float(ntrials)
1272  stddev = math.sqrt(max(rmean2 - rmean * rmean, 0.))
1273  return rmean, stddev
1274 
1275 is_already_printed = {}
1276 
1277 
1278 def print_deprecation_warning(old_name, new_name):
1279  if old_name not in is_already_printed:
1280  print("WARNING: " + old_name + " is deprecated, use " + new_name + " instead")
1281  is_already_printed[old_name] = True
1282 
1283 
1284 def print_multicolumn(list_of_strings, ncolumns=2, truncate=40):
1285 
1286  l = list_of_strings
1287 
1288  cols = ncolumns
1289  # add empty entries after l
1290  for i in range(len(l) % cols):
1291  l.append(" ")
1292 
1293  split = [l[i:i + len(l) / cols] for i in range(0, len(l), len(l) / cols)]
1294  for row in zip(*split):
1295  print("".join(str.ljust(i, truncate) for i in row))
1296 
1297 
1298 def parse_dssp(dssp_fn, limit_to_chains=''):
1299  '''Read dssp file and get secondary structure information.
1300  Values are all PDB residue numbering.
1301  @return dict of sel tuples
1302 helix : [ [ ['A',5,7] ] , [['B',15,17]] , ...] two helices A:5-7,B:15-17
1303 beta : [ [ ['A',1,3] , ['A',100,102] ] , ...] one sheet: A:1-3 & A:100-102
1304 loop : same format as helix, it's the contiguous loops
1305 '''
1306 
1307  from collections import defaultdict
1308 
1309  # setup
1310  sses = {'helix': [],
1311  'beta': [],
1312  'loop': []}
1313  helix_classes = 'GHI'
1314  strand_classes = 'EB'
1315  loop_classes = [' ', '', 'T', 'S']
1316  sse_dict = {}
1317  for h in helix_classes:
1318  sse_dict[h] = 'helix'
1319  for s in strand_classes:
1320  sse_dict[s] = 'beta'
1321  for l in loop_classes:
1322  sse_dict[l] = 'loop'
1323 
1324  # read file and parse
1325  start = False
1326  # temporary beta dictionary indexed by DSSP's ID
1327  beta_dict = defaultdict(list)
1328  prev_sstype = None
1329  cur_sse = []
1330  prev_beta_id = None
1331  for line in open(dssp_fn, 'r'):
1332  fields = line.split()
1333  chain_break = False
1334  if len(fields) < 2:
1335  continue
1336  if fields[1] == "RESIDUE":
1337  # Start parsing from here
1338  start = True
1339  continue
1340  if not start:
1341  continue
1342  if line[9] == " ":
1343  chain_break = True
1344  elif limit_to_chains != '' and line[11] not in limit_to_chains:
1345  break
1346 
1347  # gather line info
1348  if not chain_break:
1349  pdb_res_num = int(line[5:10])
1350  chain = 'chain' + line[11]
1351  sstype = line[16]
1352  beta_id = line[33]
1353 
1354  # decide whether to extend or store the SSE
1355  if prev_sstype is None:
1356  cur_sse = [pdb_res_num, pdb_res_num, chain]
1357  elif sstype != prev_sstype or chain_break:
1358  # add cur_sse to the right place
1359  if sse_dict[prev_sstype] in ['helix', 'loop']:
1360  sses[sse_dict[prev_sstype]].append([cur_sse])
1361  if sse_dict[prev_sstype] == 'beta':
1362  beta_dict[prev_beta_id].append(cur_sse)
1363  cur_sse = [pdb_res_num, pdb_res_num, chain]
1364  else:
1365  cur_sse[1] = pdb_res_num
1366  if chain_break:
1367  prev_sstype = None
1368  prev_beta_id = None
1369  else:
1370  prev_sstype = sstype
1371  prev_beta_id = beta_id
1372 
1373  # final SSE processing
1374  if not prev_sstype is None:
1375  if sse_dict[prev_sstype] in ['helix', 'loop']:
1376  sses[sse_dict[prev_sstype]].append([cur_sse])
1377  if sse_dict[prev_sstype] == 'beta':
1378  beta_dict[prev_beta_id].append(cur_sse)
1379 
1380  # gather betas
1381  for beta_sheet in beta_dict:
1382  sses['beta'].append(beta_dict[beta_sheet])
1383 
1384  return sses
1385 
1386 
1387 def sse_selections_to_chimera_colors(dssp_dict, chimera_model_num=0):
1388  ''' get chimera command to check if you've correctly made the dssp dictionary
1389  colors each helix and beta sheet'''
1390  cmds = {
1391  'helix': 'color green ',
1392  'beta': 'color blue ',
1393  'loop': 'color red '}
1394  for skey in dssp_dict.keys():
1395  for sgroup in dssp_dict[skey]:
1396  for sse in sgroup:
1397  start, stop, chain = sse
1398  chain = chain.strip('chain')
1399  cmds[
1400  skey] += '#%i:%s-%s.%s ' % (chimera_model_num, start, stop, chain)
1401  print('; '.join([cmds[k] for k in cmds]))
1402 
1403 
1404 
1405 
1406 class ColorChange(object):
1407  '''Change color code to hexadecimal to rgb'''
1408  def __init__(self):
1409  self._NUMERALS = '0123456789abcdefABCDEF'
1410  self._HEXDEC = dict((v, int(v, 16)) for v in (x+y for x in self._NUMERALS for y in self._NUMERALS))
1411  self.LOWERCASE, self.UPPERCASE = 'x', 'X'
1412 
1413  def rgb(self,triplet):
1414  return float(self._HEXDEC[triplet[0:2]]), float(self._HEXDEC[triplet[2:4]]), float(self._HEXDEC[triplet[4:6]])
1415 
1416  def triplet(self,rgb, lettercase=None):
1417  if lettercase is None: lettercase=self.LOWERCASE
1418  return format(rgb[0]<<16 | rgb[1]<<8 | rgb[2], '06'+lettercase)
1419 
1420 
1421 class OrderedSet(collections.MutableSet):
1422 
1423  def __init__(self, iterable=None):
1424  self.end = end = []
1425  end += [None, end, end] # sentinel node for doubly linked list
1426  self.map = {} # key --> [key, prev, next]
1427  if iterable is not None:
1428  self |= iterable
1429 
1430  def __len__(self):
1431  return len(self.map)
1432 
1433  def __contains__(self, key):
1434  return key in self.map
1435 
1436  def add(self, key):
1437  if key not in self.map:
1438  end = self.end
1439  curr = end[1]
1440  curr[2] = end[1] = self.map[key] = [key, curr, end]
1441 
1442  def discard(self, key):
1443  if key in self.map:
1444  key, prev, next = self.map.pop(key)
1445  prev[2] = next
1446  next[1] = prev
1447 
1448  def __iter__(self):
1449  end = self.end
1450  curr = end[2]
1451  while curr is not end:
1452  yield curr[0]
1453  curr = curr[2]
1454 
1455  def __reversed__(self):
1456  end = self.end
1457  curr = end[1]
1458  while curr is not end:
1459  yield curr[0]
1460  curr = curr[1]
1461 
1462  def pop(self, last=True):
1463  if not self:
1464  raise KeyError('set is empty')
1465  if last:
1466  key = self.end[1][0]
1467  else:
1468  key = self.end[2][0]
1469  self.discard(key)
1470  return key
1471 
1472  def __repr__(self):
1473  if not self:
1474  return '%s()' % (self.__class__.__name__,)
1475  return '%s(%r)' % (self.__class__.__name__, list(self))
1476 
1477  def __eq__(self, other):
1478  if isinstance(other, OrderedSet):
1479  return len(self) == len(other) and list(self) == list(other)
1480  return set(self) == set(other)
def list_chunks_iterator
Yield successive length-sized chunks from a list.
Definition: tools.py:1074
Simple 3D transformation class.
static bool get_is_setup(const IMP::ParticleAdaptor &p)
Definition: Residue.h:155
A decorator to associate a particle with a part of a protein/DNA/RNA.
Definition: Fragment.h:20
def get_restraint_set
Get a RestraintSet containing all PMI restraints added to the model.
Definition: tools.py:34
Store the representations for a system.
Definition: tools.py:791
def get_random_cross_link_dataset
Return a random cross-link dataset as a string.
Definition: tools.py:342
static bool get_is_setup(const IMP::ParticleAdaptor &p)
Definition: atom/Atom.h:241
void add_particles(RMF::FileHandle fh, const ParticlesTemp &hs)
static Weight setup_particle(Model *m, ParticleIndex pi)
Definition: Weight.h:31
def get_position_terminal_residue
Get the xyz position of the terminal residue at the given resolution.
Definition: tools.py:553
double get_drmsd(const Vector3DsOrXYZs0 &m0, const Vector3DsOrXYZs1 &m1)
Calculate distance the root mean square deviation between two sets of 3D points.
Definition: atom/distance.h:49
def scatter_and_gather_dict_append
Synchronize data over a parallel run.
Definition: tools.py:1028
Change color code to hexadecimal to rgb.
Definition: tools.py:1406
ParticlesTemp get_particles(Model *m, const ParticleIndexes &ps)
void add_particle(RMF::FileHandle fh, Particle *hs)
log
Definition: log.py:1
static bool get_is_setup(const IMP::ParticleAdaptor &p)
Definition: rigid_bodies.h:469
def get_prot_name_from_particle
Return the component name provided a particle and a list of names.
Definition: tools.py:939
double get_distance(XYZR a, XYZR b)
Compute the sphere distance between a and b.
Definition: XYZR.h:87
def get_residue_gaps_in_hierarchy
Return the residue index gaps and contiguous segments in the hierarchy.
Definition: tools.py:585
Add resolution to a particle.
Definition: Resolution.h:24
def sse_selections_to_chimera_colors
get chimera command to check if you've correctly made the dssp dictionary colors each helix and beta ...
Definition: tools.py:1387
Object used to hold a set of restraints.
Definition: RestraintSet.h:36
algebra::GridD< 3, algebra::DenseGridStorageD< 3, float >, float > get_grid(DensityMap *in_map)
static bool get_is_setup(Model *m, ParticleIndex pi)
Definition: Fragment.h:46
def translate_hierarchy
Apply a translation to a hierarchy along the input vector.
Definition: tools.py:1098
The standard decorator for manipulating molecular structures.
Ints get_index(const ParticlesTemp &particles, const Subset &subset, const Subsets &excluded)
A decorator for a particle representing an atom.
Definition: atom/Atom.h:234
def scatter_and_gather
Synchronize data over a parallel run.
Definition: tools.py:1001
void transform(XYZ a, const algebra::Transformation3D &tr)
Apply a transformation to the particle.
def add_restraint_to_model
Add a PMI restraint to the model.
Definition: tools.py:22
Hierarchies get_by_type(Hierarchy mhd, GetByType t)
Gather all the molecular particles of a certain level in the hierarchy.
A decorator for a particle with x,y,z coordinates.
Definition: XYZ.h:30
static Scale setup_particle(Model *m, ParticleIndex pi)
Definition: Scale.h:30
A base class for Keys.
Definition: Key.h:46
def get_closest_residue_position
this function works with plain hierarchies, as read from the pdb, no multi-scale hierarchies ...
Definition: tools.py:522
static bool get_is_setup(const IMP::ParticleAdaptor &p)
Definition: rigid_bodies.h:147
A decorator for a residue.
Definition: Residue.h:134
General purpose algebraic and geometric methods that are expected to be used by a wide variety of IMP...
void show(Hierarchy h, std::ostream &out=std::cout)
Print out a molecular hierarchy.
The general base class for IMP exceptions.
Definition: exception.h:49
def parse_dssp
Read dssp file and get secondary structure information.
Definition: tools.py:1298
VectorD< 3 > Vector3D
Definition: VectorD.h:395
Class to handle individual model particles.
Definition: Particle.h:37
std::string get_data_path(std::string file_name)
Return the full path to one of this module's data files.
def select
this function uses representation=SimplifiedModel it returns the corresponding selected particles rep...
Definition: tools.py:662
double get_resolution(Model *m, ParticleIndex pi)
Estimate the resolution of the hierarchy as used by Representation.
A decorator for a rigid body.
Definition: rigid_bodies.h:75
def cross_link_db_filter_parser
example '"{ID_Score}" > 28 AND "{Sample}" == "%10_1%" OR ":Sample}" == "%10_2%" OR ":Sample...
Definition: tools.py:454
Hierarchies get_leaves(const Selection &h)
Select hierarchy particles identified by the biological name.
Definition: Selection.h:65
def get_residue_indexes
Retrieve the residue indexes for the given particle.
Definition: tools.py:959
Inferential scoring building on methods developed as part of the Inferential Structure Determination ...
def sublist_iterator
Yield all sublists of length >= lmin and <= lmax.
Definition: tools.py:1055