IMP logo
IMP Reference Guide  2.7.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 IMP.isd
11 import IMP.pmi
12 import IMP.pmi.topology
13 import collections
14 import itertools
15 from math import log,pi,sqrt,exp
16 import sys,os
17 import random
18 import time
19 import RMF
20 import IMP.rmf
21 from collections import defaultdict
22 try:
23  from collections import OrderedDict
24 except ImportError:
25  from IMP.pmi._compat_collections import OrderedDict
26 
27 def _get_restraint_set_key():
28  if not hasattr(_get_restraint_set_key, 'pmi_rs_key'):
29  _get_restraint_set_key.pmi_rs_key = IMP.ModelKey("PMI restraints")
30  return _get_restraint_set_key.pmi_rs_key
31 
32 def _add_restraint_set(model, mk):
33  rs = IMP.RestraintSet(model, "All PMI restraints")
34  model.add_data(mk, rs)
35  return rs
36 
37 def add_restraint_to_model(model, restraint):
38  """Add a PMI restraint to the model.
39  Since Model.add_restraint() no longer exists (in modern IMP restraints
40  should be added to a ScoringFunction instead) store them instead in
41  a RestraintSet, and keep a reference to it in the Model."""
42  mk = _get_restraint_set_key()
43  if model.get_has_data(mk):
44  rs = IMP.RestraintSet.get_from(model.get_data(mk))
45  else:
46  rs = _add_restraint_set(model, mk)
47  rs.add_restraint(restraint)
48 
49 def get_restraint_set(model):
50  """Get a RestraintSet containing all PMI restraints added to the model"""
51  mk = _get_restraint_set_key()
52  if not model.get_has_data(mk):
53  print("WARNING: no restraints added to model yet")
54  _add_restraint_set(model, mk)
55  return IMP.RestraintSet.get_from(model.get_data(mk))
56 
57 class Stopwatch(object):
58 
59  def __init__(self, isdelta=True):
60 
61  self.starttime = time.clock()
62  self.label = "None"
63  self.isdelta = isdelta
64 
65  def set_label(self, labelstr):
66  self.label = labelstr
67 
68  def get_output(self):
69  output = {}
70  if self.isdelta:
71  newtime = time.clock()
72  output[
73  "Stopwatch_" +
74  self.label +
75  "_delta_seconds"] = str(
76  newtime -
77  self.starttime)
78  self.starttime = newtime
79  else:
80  output["Stopwatch_" + self.label +
81  "_elapsed_seconds"] = str(time.clock() - self.starttime)
82  return output
83 
84 
85 class SetupNuisance(object):
86 
87  def __init__(self, m, initialvalue, minvalue, maxvalue, isoptimized=True):
88 
89  nuisance = IMP.isd.Scale.setup_particle(IMP.Particle(m), initialvalue)
90  if minvalue:
91  nuisance.set_lower(minvalue)
92  if maxvalue:
93  nuisance.set_upper(maxvalue)
94 
95  # m.add_score_state(IMP.core.SingletonConstraint(IMP.isd.NuisanceRangeModifier(),None,nuisance))
96  nuisance.set_is_optimized(nuisance.get_nuisance_key(), isoptimized)
97  self.nuisance = nuisance
98 
99  def get_particle(self):
100  return self.nuisance
101 
102 
103 class SetupWeight(object):
104 
105  def __init__(self, m, isoptimized=True):
106  pw = IMP.Particle(m)
107  self.weight = IMP.isd.Weight.setup_particle(pw)
108  self.weight.set_weights_are_optimized(True)
109 
110  def get_particle(self):
111  return self.weight
112 
113 
114 class SetupSurface(object):
115 
116  def __init__(self, m, center, normal, isoptimized=True):
117  p = IMP.Particle(m)
118  self.surface = IMP.core.Surface.setup_particle(p, center, normal)
119  self.surface.set_coordinates_are_optimized(isoptimized)
120  self.surface.set_normal_is_optimized(isoptimized)
121 
122  def get_particle(self):
123  return self.surface
124 
125 
126 class ParticleToSampleFilter(object):
127  def __init__(self, sampled_objects):
128  self.sampled_objects=sampled_objects
129  self.filters=[]
130 
131  def add_filter(self,filter_string):
132  self.filters.append(filter_string)
133 
134  def get_particles_to_sample(self):
135  particles_to_sample={}
136  for so in self.sampled_objects:
137  ps_dict=so.get_particles_to_sample()
138  for key in ps_dict:
139  for f in self.filters:
140  if f in key:
141  if key not in particles_to_sample:
142  particles_to_sample[key]=ps_dict[key]
143  else:
144  particles_to_sample[key]+=ps_dict[key]
145  return particles_to_sample
146 
147 class ParticleToSampleList(object):
148 
149  def __init__(self, label="None"):
150 
151  self.dictionary_particle_type = {}
152  self.dictionary_particle_transformation = {}
153  self.dictionary_particle_name = {}
154  self.label = label
155 
156  def add_particle(
157  self,
158  particle,
159  particle_type,
160  particle_transformation,
161  name):
162  if not particle_type in ["Rigid_Bodies", "Floppy_Bodies", "Nuisances", "X_coord", "Weights", "Surfaces"]:
163  raise TypeError("not the right particle type")
164  else:
165  self.dictionary_particle_type[particle] = particle_type
166  if particle_type == "Rigid_Bodies":
167  if type(particle_transformation) == tuple and len(particle_transformation) == 2 and type(particle_transformation[0]) == float and type(particle_transformation[1]) == float:
168  self.dictionary_particle_transformation[
169  particle] = particle_transformation
170  self.dictionary_particle_name[particle] = name
171  else:
172  raise TypeError("ParticleToSampleList: not the right transformation format for Rigid_Bodies, should be a tuple of floats")
173  elif particle_type == "Surfaces":
174  if type(particle_transformation) == tuple and len(particle_transformation) == 3 and all(isinstance(x, float) for x in particle_transformation):
175  self.dictionary_particle_transformation[
176  particle] = particle_transformation
177  self.dictionary_particle_name[particle] = name
178  else:
179  raise TypeError("ParticleToSampleList: not the right transformation format for Surfaces, should be a tuple of floats")
180  else:
181  if type(particle_transformation) == float:
182  self.dictionary_particle_transformation[
183  particle] = particle_transformation
184  self.dictionary_particle_name[particle] = name
185  else:
186  raise TypeError("ParticleToSampleList: not the right transformation format, should be a float")
187 
188  def get_particles_to_sample(self):
189  ps = {}
190  for particle in self.dictionary_particle_type:
191  key = self.dictionary_particle_type[
192  particle] + "ParticleToSampleList_" + self.dictionary_particle_name[particle] + "_" + self.label
193  value = (
194  [particle],
195  self.dictionary_particle_transformation[particle])
196  ps[key] = value
197  return ps
198 
199 
200 class Variance(object):
201 
202  def __init__(self, model, tau, niter, prot, th_profile, write_data=False):
203 
204  self.model = model
205  self.write_data = write_data
206  self.tau = tau
207  self.niter = niter
208  #! select particles from the model
209  particles = IMP.atom.get_by_type(prot, IMP.atom.ATOM_TYPE)
210  self.particles = particles
211  # store reference coordinates and theoretical profile
212  self.refpos = [IMP.core.XYZ(p).get_coordinates() for p in particles]
213  self.model_profile = th_profile
214 
215  def perturb_particles(self, perturb=True):
216  for i, p in enumerate(self.particles):
217  newpos = array(self.refpos[i])
218  if perturb:
219  newpos += random.normal(0, self.tau, 3)
220  newpos = IMP.algebra.Vector3D(newpos)
221  IMP.core.XYZ(p).set_coordinates(newpos)
222 
223  def get_profile(self):
224  model_profile = self.model_profile
225  p = model_profile.calculate_profile(self.particles, IMP.saxs.CA_ATOMS)
226  return array([model_profile.get_intensity(i) for i in
227  range(model_profile.size())])
228 
229  def init_variances(self):
230  # create placeholders
231  N = self.model_profile.size()
232  a = self.profiles[0][:]
233  self.m = matrix(a).T # Nx1
234  self.V = self.m * self.m.T
235  self.normm = linalg.norm(self.m)
236  self.normV = linalg.norm(self.V)
237 
238  def update_variances(self):
239  a = matrix(self.profiles[-1]) # 1xN
240  n = float(len(self.profiles))
241  self.m = a.T / n + (n - 1) / n * self.m
242  self.V = a.T * a + self.V
243  self.oldnormm = self.normm
244  self.oldnormV = self.normV
245  self.normm = linalg.norm(self.m)
246  self.normV = linalg.norm(self.V)
247  self.diffm = (self.oldnormm - self.normm) / self.oldnormm
248  self.diffV = (self.oldnormV - self.normV) / self.oldnormV
249 
250  def get_direct_stats(self, a):
251  nq = len(a[0])
252  nprof = len(a)
253  m = [0] * nq
254  for prof in a:
255  for q, I in enumerate(prof):
256  m[q] += I
257  m = array(m) / nprof
258  V = matrix(a)
259  V = V.T * V
260  Sigma = (matrix(a - m))
261  Sigma = Sigma.T * Sigma / (nprof - 1)
262  mi = matrix(diag(1. / m))
263  Sigmarel = mi.T * Sigma * mi
264  return m, V, Sigma, Sigmarel
265 
266  def store_data(self):
267  if not os.path.isdir('data'):
268  os.mkdir('data')
269  profiles = matrix(self.profiles)
270  self.directm, self.directV, self.Sigma, self.Sigmarel = \
271  self.get_direct_stats(array(profiles))
272  directV = self.directV
273  # print "V comparison",(linalg.norm(directV-self.V)/self.normV)
274  save('data/profiles', profiles)
275  # absolute profile differences
276  fl = open('data/profiles.dat', 'w')
277  for i, l in enumerate(array(profiles).T):
278  self.model_profile.get_q(i)
279  fl.write('%s ' % i)
280  for k in l:
281  fl.write('%s ' % (k - self.directm[i]))
282  fl.write('\n')
283  # relative profile differences
284  fl = open('data/profiles_rel.dat', 'w')
285  for i, l in enumerate(array(profiles).T):
286  self.model_profile.get_q(i)
287  fl.write('%s ' % i)
288  for k in l:
289  fl.write('%s ' % ((k - self.directm[i]) / self.directm[i]))
290  fl.write('\n')
291  save('data/m', self.directm)
292  save('data/V', self.directV)
293  Sigma = self.Sigma
294  save('data/Sigma', Sigma)
295  # Sigma matrix
296  fl = open('data/Sigma.dat', 'w')
297  model_profile = self.model_profile
298  for i in range(model_profile.size()):
299  qi = model_profile.get_q(i)
300  for j in range(model_profile.size()):
301  qj = model_profile.get_q(j)
302  vij = self.Sigma[i, j]
303  fl.write('%s %s %s\n' % (qi, qj, vij))
304  fl.write('\n')
305  # Sigma eigenvalues
306  fl = open('data/eigenvals', 'w')
307  for i in linalg.eigvalsh(Sigma):
308  fl.write('%s\n' % i)
309  Sigmarel = self.Sigmarel
310  save('data/Sigmarel', Sigmarel)
311  # Sigmarel matrix
312  fl = open('data/Sigmarel.dat', 'w')
313  model_profile = self.model_profile
314  for i in range(model_profile.size()):
315  qi = model_profile.get_q(i)
316  for j in range(model_profile.size()):
317  qj = model_profile.get_q(j)
318  vij = self.Sigmarel[i, j]
319  fl.write('%s %s %s\n' % (qi, qj, vij))
320  fl.write('\n')
321  # Sigma eigenvalues
322  fl = open('data/eigenvals_rel', 'w')
323  for i in linalg.eigvalsh(Sigmarel):
324  fl.write('%s\n' % i)
325  # mean profile
326  fl = open('data/mean.dat', 'w')
327  for i in range(len(self.directm)):
328  qi = self.model_profile.get_q(i)
329  fl.write('%s ' % qi)
330  fl.write('%s ' % self.directm[i])
331  fl.write('%s ' % sqrt(self.Sigma[i, i]))
332  fl.write('\n')
333 
334  def try_chol(self, jitter):
335  Sigma = self.Sigma
336  try:
337  linalg.cholesky(Sigma + matrix(eye(len(Sigma))) * jitter)
338  except linalg.LinAlgError:
339  print("Decomposition failed with jitter =", jitter)
340  return
341  print("Successful decomposition with jitter =", jitter)
342 
343  def run(self):
344  self.profiles = [self.get_profile()]
345  # self.init_variances()
346  for n in range(self.niter):
347  self.perturb_particles()
348  self.profiles.append(self.get_profile())
349  # self.update_variances()
350  #profiles = matrix(self.profiles)
351  # print n,self.diffm,self.diffV
352  # print
353  #
354  if self.write_data:
355  self.store_data()
356  # self.try_chol(0.)
357  # for i in logspace(-7,0,num=8):
358  # self.try_chol(i)
359 
360  def get_cov(self, relative=True):
361  if not relative:
362  return self.Sigma
363  else:
364  return self.Sigmarel
365 
366  #-------------------------------
367 
368 def get_random_cross_link_dataset(representation,
369  resolution=1.0,
370  number_of_cross_links=10,
371  ambiguity_probability=0.1,
372  confidence_score_range=[0,100],
373  avoid_same_particles=False):
374  '''Return a random cross-link dataset as a string.
375  Every line is a residue pair, together with UniqueIdentifier
376  and XL score.'''
377 
378  residue_pairs=get_random_residue_pairs(representation, resolution, number_of_cross_links, avoid_same_particles=avoid_same_particles)
379 
380  unique_identifier=0
381  cmin=float(min(confidence_score_range))
382  cmax=float(max(confidence_score_range))
383 
384  dataset="#\n"
385 
386  for (name1, r1, name2, r2) in residue_pairs:
387  if random.random() > ambiguity_probability:
388  unique_identifier+=1
389  score=random.random()*(cmax-cmin)+cmin
390  dataset+=str(name1)+" "+str(name2)+" "+str(r1)+" "+str(r2)+" "+str(score)+" "+str(unique_identifier)+"\n"
391 
392  return dataset
393 
394 
395  #-------------------------------
396 
397 def get_cross_link_data(directory, filename, dist, omega, sigma,
398  don=None, doff=None, prior=0, type_of_profile="gofr"):
399 
400  (distmin, distmax, ndist) = dist
401  (omegamin, omegamax, nomega) = omega
402  (sigmamin, sigmamax, nsigma) = sigma
403 
404  filen = IMP.isd.get_data_path("CrossLinkPMFs.dict")
405  xlpot = open(filen)
406 
407  for line in xlpot:
408  dictionary = eval(line)
409  break
410 
411  xpot = dictionary[directory][filename]["distance"]
412  pot = dictionary[directory][filename][type_of_profile]
413 
414  dist_grid = get_grid(distmin, distmax, ndist, False)
415  omega_grid = get_log_grid(omegamin, omegamax, nomega)
416  sigma_grid = get_log_grid(sigmamin, sigmamax, nsigma)
417 
418  if not don is None and not doff is None:
419  xlmsdata = IMP.isd.CrossLinkData(
420  dist_grid,
421  omega_grid,
422  sigma_grid,
423  xpot,
424  pot,
425  don,
426  doff,
427  prior)
428  else:
429  xlmsdata = IMP.isd.CrossLinkData(
430  dist_grid,
431  omega_grid,
432  sigma_grid,
433  xpot,
434  pot)
435  return xlmsdata
436 
437  #-------------------------------
438 
439 
440 def get_cross_link_data_from_length(length, xxx_todo_changeme3, xxx_todo_changeme4, xxx_todo_changeme5):
441  (distmin, distmax, ndist) = xxx_todo_changeme3
442  (omegamin, omegamax, nomega) = xxx_todo_changeme4
443  (sigmamin, sigmamax, nsigma) = xxx_todo_changeme5
444 
445  dist_grid = get_grid(distmin, distmax, ndist, False)
446  omega_grid = get_log_grid(omegamin, omegamax, nomega)
447  sigma_grid = get_log_grid(sigmamin, sigmamax, nsigma)
448 
449  xlmsdata = IMP.isd.CrossLinkData(dist_grid, omega_grid, sigma_grid, length)
450  return xlmsdata
451 
452 
453 def get_grid(gmin, gmax, ngrid, boundaries):
454  grid = []
455  dx = (gmax - gmin) / float(ngrid)
456  for i in range(0, ngrid + 1):
457  if(not boundaries and i == 0):
458  continue
459  if(not boundaries and i == ngrid):
460  continue
461  grid.append(gmin + float(i) * dx)
462  return grid
463 
464  #-------------------------------
465 
466 
467 def get_log_grid(gmin, gmax, ngrid):
468  grid = []
469  for i in range(0, ngrid + 1):
470  grid.append(gmin * exp(float(i) / ngrid * log(gmax / gmin)))
471  return grid
472 
473  #-------------------------------
474 
475 
477  '''
478  example '"{ID_Score}" > 28 AND "{Sample}" ==
479  "%10_1%" OR ":Sample}" == "%10_2%" OR ":Sample}"
480  == "%10_3%" OR ":Sample}" == "%8_1%" OR ":Sample}" == "%8_2%"'
481  '''
482 
483  import pyparsing as pp
484 
485  operator = pp.Regex(">=|<=|!=|>|<|==|in").setName("operator")
486  value = pp.QuotedString(
487  '"') | pp.Regex(
488  r"[+-]?\d+(:?\.\d*)?(:?[eE][+-]?\d+)?")
489  identifier = pp.Word(pp.alphas, pp.alphanums + "_")
490  comparison_term = identifier | value
491  condition = pp.Group(comparison_term + operator + comparison_term)
492 
493  expr = pp.operatorPrecedence(condition, [
494  ("OR", 2, pp.opAssoc.LEFT, ),
495  ("AND", 2, pp.opAssoc.LEFT, ),
496  ])
497 
498  parsedstring = str(expr.parseString(inputstring)) \
499  .replace("[", "(") \
500  .replace("]", ")") \
501  .replace(",", " ") \
502  .replace("'", " ") \
503  .replace("%", "'") \
504  .replace("{", "float(entry['") \
505  .replace("}", "'])") \
506  .replace(":", "str(entry['") \
507  .replace("}", "'])") \
508  .replace("AND", "and") \
509  .replace("OR", "or")
510  return parsedstring
511 
512 
513 def open_file_or_inline_text(filename):
514  try:
515  fl = open(filename, "r")
516  except IOError:
517  fl = filename.split("\n")
518  return fl
519 
520 
521 def get_drmsd(prot0, prot1):
522  drmsd = 0.
523  npairs = 0.
524  for i in range(0, len(prot0) - 1):
525  for j in range(i + 1, len(prot0)):
526  dist0 = IMP.core.get_distance(prot0[i], prot0[j])
527  dist1 = IMP.core.get_distance(prot1[i], prot1[j])
528  drmsd += (dist0 - dist1) ** 2
529  npairs += 1.
530  return math.sqrt(drmsd / npairs)
531 
532  #-------------------------------
533 
534 
535 def get_ids_from_fasta_file(fastafile):
536  ids = []
537  with open(fastafile) as ff:
538  for l in ff:
539  if l[0] == ">":
540  ids.append(l[1:-1])
541  return ids
542 
543 
544 def get_closest_residue_position(hier, resindex, terminus="N"):
545  '''
546  this function works with plain hierarchies, as read from the pdb,
547  no multi-scale hierarchies
548  '''
549  p = []
550  niter = 0
551  while len(p) == 0:
552  niter += 1
553  sel = IMP.atom.Selection(hier, residue_index=resindex,
554  atom_type=IMP.atom.AT_CA)
555 
556  if terminus == "N":
557  resindex += 1
558  if terminus == "C":
559  resindex -= 1
560 
561  if niter >= 10000:
562  print("get_closest_residue_position: exiting while loop without result")
563  break
564  p = sel.get_selected_particles()
565 
566  if len(p) == 1:
567  return IMP.core.XYZ(p[0]).get_coordinates()
568  elif len(p) == 0:
569  print("get_closest_residue_position: got NO residues for hierarchy %s and residue %i" % (hier, resindex))
570  raise Exception("get_closest_residue_position: got NO residues for hierarchy %s and residue %i" % (
571  hier, resindex))
572  else:
573  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])))
574 
575 @IMP.deprecated_function("2.6", "Use get_terminal_residue_position() instead.")
576 def get_position_terminal_residue(hier, terminus="C", resolution=1):
577  '''
578  Get the xyz position of the terminal residue at the given resolution.
579  @param hier hierarchy containing the terminal residue
580  @param terminus either 'N' or 'C'
581  @param resolution resolution to use.
582  '''
583  termresidue = None
584  termparticle = None
585  for p in IMP.atom.get_leaves(hier):
586  if IMP.pmi.Resolution(p).get_resolution() == resolution:
588  if terminus == "C":
589  if max(residues) >= termresidue and not termresidue is None:
590  termresidue = max(residues)
591  termparticle = p
592  elif termresidue is None:
593  termresidue = max(residues)
594  termparticle = p
595  elif terminus == "N":
596  if min(residues) <= termresidue and not termresidue is None:
597  termresidue = min(residues)
598  termparticle = p
599  elif termresidue is None:
600  termresidue = min(residues)
601  termparticle = p
602  else:
603  raise ValueError("terminus argument should be either N or C")
604 
605  return IMP.core.XYZ(termparticle).get_coordinates()
606 
607 def get_terminal_residue(representation, hier, terminus="C", resolution=1):
608  '''
609  Get the particle of the terminal residue at the GIVEN resolution
610  (NOTE: not the closest resolution!).
611  To get the terminal residue at the closest resolution use:
612  particles=IMP.pmi.tools.select_by_tuple(representation,molecule_name)
613  particles[0] and particles[-1] will be the first and last particles
614  corresponding to the two termini.
615  It is needed for instance to determine the last residue of a pdb.
616  @param hier hierarchy containing the terminal residue
617  @param terminus either 'N' or 'C'
618  @param resolution resolution to use.
619  '''
620  termresidue = None
621  termparticle = None
622 
623  ps=select(representation,
624  resolution=resolution,
625  hierarchies=[hier])
626 
627  for p in ps:
628  if IMP.pmi.Resolution(p).get_resolution() == resolution:
630  if terminus == "C":
631  if termresidue is None:
632  termresidue = max(residues)
633  termparticle = p
634  elif max(residues) >= termresidue:
635  termresidue = max(residues)
636  termparticle = p
637  elif terminus == "N":
638  if termresidue is None:
639  termresidue = min(residues)
640  termparticle = p
641  elif min(residues) <= termresidue:
642  termresidue = min(residues)
643  termparticle = p
644  else:
645  raise ValueError("terminus argument should be either N or C")
646  return termparticle
647 
648 def get_terminal_residue_position(representation, hier, terminus="C",
649  resolution=1):
650  """Get XYZ coordinates of the terminal residue at the GIVEN resolution"""
651  p = get_terminal_residue(representation, hier, terminus, resolution)
652  return IMP.core.XYZ(p).get_coordinates()
653 
654 def get_residue_gaps_in_hierarchy(hierarchy, start, end):
655  '''
656  Return the residue index gaps and contiguous segments in the hierarchy.
657 
658  @param hierarchy hierarchy to examine
659  @param start first residue index
660  @param end last residue index
661 
662  @return A list of lists of the form
663  [[1,100,"cont"],[101,120,"gap"],[121,200,"cont"]]
664  '''
665  gaps = []
666  for n, rindex in enumerate(range(start, end + 1)):
667  sel = IMP.atom.Selection(hierarchy, residue_index=rindex,
668  atom_type=IMP.atom.AT_CA)
669 
670  if len(sel.get_selected_particles()) == 0:
671  if n == 0:
672  # set the initial condition
673  rindexgap = start
674  rindexcont = start - 1
675  if rindexgap == rindex - 1:
676  # residue is contiguous with the previously discovered gap
677  gaps[-1][1] += 1
678  else:
679  # residue is not contiguous with the previously discovered gap
680  # hence create a new gap tuple
681  gaps.append([rindex, rindex, "gap"])
682  # update the index of the last residue gap
683  rindexgap = rindex
684  else:
685  if n == 0:
686  # set the initial condition
687  rindexgap = start - 1
688  rindexcont = start
689  if rindexcont == rindex - 1:
690  # residue is contiguous with the previously discovered
691  # continuous part
692  gaps[-1][1] += 1
693  else:
694  # residue is not contiguous with the previously discovered continuous part
695  # hence create a new cont tuple
696  gaps.append([rindex, rindex, "cont"])
697  # update the index of the last residue gap
698  rindexcont = rindex
699  return gaps
700 
701 
702 class map(object):
703 
704  def __init__(self):
705  self.map = {}
706 
707  def set_map_element(self, xvalue, yvalue):
708  self.map[xvalue] = yvalue
709 
710  def get_map_element(self, invalue):
711  if type(invalue) == float:
712  n = 0
713  mindist = 1
714  for x in self.map:
715  dist = (invalue - x) * (invalue - x)
716 
717  if n == 0:
718  mindist = dist
719  minx = x
720  if dist < mindist:
721  mindist = dist
722  minx = x
723  n += 1
724  return self.map[minx]
725  elif type(invalue) == str:
726  return self.map[invalue]
727  else:
728  raise TypeError("wrong type for map")
729 
730 
731 def select(representation,
732  resolution=None,
733  hierarchies=None,
734  selection_arguments=None,
735  name=None,
736  name_is_ambiguous=False,
737  first_residue=None,
738  last_residue=None,
739  residue=None,
740  representation_type=None):
741  '''
742  this function uses representation=SimplifiedModel
743  it returns the corresponding selected particles
744  representation_type="Beads", "Res:X", "Densities", "Representation", "Molecule"
745  '''
746 
747  if resolution is None:
748  allparticles = IMP.atom.get_leaves(representation.prot)
749  resolution_particles = None
750  hierarchies_particles = None
751  names_particles = None
752  residue_range_particles = None
753  residue_particles = None
754  representation_type_particles = None
755 
756  if not resolution is None:
757  resolution_particles = []
758  hs = representation.get_hierarchies_at_given_resolution(resolution)
759  for h in hs:
760  resolution_particles += IMP.atom.get_leaves(h)
761 
762  if not hierarchies is None:
763  hierarchies_particles = []
764  for h in hierarchies:
765  hierarchies_particles += IMP.atom.get_leaves(h)
766 
767  if not name is None:
768  names_particles = []
769  if name_is_ambiguous:
770  for namekey in representation.hier_dict:
771  if name in namekey:
772  names_particles += IMP.atom.get_leaves(
773  representation.hier_dict[namekey])
774  elif name in representation.hier_dict:
775  names_particles += IMP.atom.get_leaves(representation.hier_dict[name])
776  else:
777  print("select: component %s is not there" % name)
778 
779  if not first_residue is None and not last_residue is None:
780  sel = IMP.atom.Selection(representation.prot,
781  residue_indexes=range(first_residue, last_residue + 1))
782  residue_range_particles = [IMP.atom.Hierarchy(p)
783  for p in sel.get_selected_particles()]
784 
785  if not residue is None:
786  sel = IMP.atom.Selection(representation.prot, residue_index=residue)
787  residue_particles = [IMP.atom.Hierarchy(p)
788  for p in sel.get_selected_particles()]
789 
790  if not representation_type is None:
791  representation_type_particles = []
792  if representation_type == "Molecule":
793  for name in representation.hier_representation:
794  for repr_type in representation.hier_representation[name]:
795  if repr_type == "Beads" or "Res:" in repr_type:
796  h = representation.hier_representation[name][repr_type]
797  representation_type_particles += IMP.atom.get_leaves(h)
798 
799  elif representation_type == "PDB":
800  for name in representation.hier_representation:
801  for repr_type in representation.hier_representation[name]:
802  if repr_type == "Res:" in repr_type:
803  h = representation.hier_representation[name][repr_type]
804  representation_type_particles += IMP.atom.get_leaves(h)
805 
806  else:
807  for name in representation.hier_representation:
808  h = representation.hier_representation[
809  name][
810  representation_type]
811  representation_type_particles += IMP.atom.get_leaves(h)
812 
813  selections = [hierarchies_particles, names_particles,
814  residue_range_particles, residue_particles, representation_type_particles]
815 
816  if resolution is None:
817  selected_particles = set(allparticles)
818  else:
819  selected_particles = set(resolution_particles)
820 
821  for s in selections:
822  if not s is None:
823  selected_particles = (set(s) & selected_particles)
824 
825  return list(selected_particles)
826 
827 
828 def select_by_tuple(
829  representation,
830  tupleselection,
831  resolution=None,
832  name_is_ambiguous=False):
833  if isinstance(tupleselection, tuple) and len(tupleselection) == 3:
834  particles = IMP.pmi.tools.select(representation, resolution=resolution,
835  name=tupleselection[2],
836  first_residue=tupleselection[0],
837  last_residue=tupleselection[1],
838  name_is_ambiguous=name_is_ambiguous)
839  elif isinstance(tupleselection, str):
840  particles = IMP.pmi.tools.select(representation, resolution=resolution,
841  name=tupleselection,
842  name_is_ambiguous=name_is_ambiguous)
843  else:
844  raise ValueError('you passed something bad to select_by_tuple()')
845  # now order the result by residue number
846  particles = IMP.pmi.tools.sort_by_residues(particles)
847 
848  return particles
849 
850 def select_by_tuple_2(hier,tuple_selection,resolution):
851  """New tuple format: molname OR (start,stop,molname,copynum,statenum)
852  Copy and state are optional. Can also use 'None' for them which will get all.
853  You can also pass -1 for stop which will go to the end.
854  Returns the particles
855  """
856  kwds = {} # going to accumulate keywords
857  kwds['resolution'] = resolution
858  if type(tuple_selection) is str:
859  kwds['molecule'] = tuple_selection
860  elif type(tuple_selection) is tuple:
861  rbegin = tuple_selection[0]
862  rend = tuple_selection[1]
863  kwds['molecule'] = tuple_selection[2]
864  try:
865  copynum = tuple_selection[3]
866  if copynum is not None:
867  kwds['copy_index'] = copynum
868  except:
869  pass
870  try:
871  statenum = tuple_selection[4]
872  if statenum is not None:
873  kwds['state_index'] = statenum
874  except:
875  pass
876  if rend==-1:
877  if rbegin>1:
878  s = IMP.atom.Selection(hier,**kwds)
879  s -= IMP.atom.Selection(hier,
880  residue_indexes=range(1,rbegin),
881  **kwds)
882  return s.get_selected_particles()
883  else:
884  kwds['residue_indexes'] = range(rbegin,rend+1)
885  s = IMP.atom.Selection(hier,**kwds)
886  return s.get_selected_particles()
887 
888 
889 
890 def get_db_from_csv(csvfilename):
891  import csv
892  outputlist = []
893  with open(csvfilename) as fh:
894  csvr = csv.DictReader(fh)
895  for l in csvr:
896  outputlist.append(l)
897  return outputlist
898 
899 
900 class HierarchyDatabase(object):
901  """Store the representations for a system."""
902 
903  def __init__(self):
904  self.db = {}
905  # this dictionary map a particle to its root hierarchy
906  self.root_hierarchy_dict = {}
907  self.preroot_fragment_hierarchy_dict = {}
908  self.particle_to_name = {}
909  self.model = None
910 
911  def add_name(self, name):
912  if name not in self.db:
913  self.db[name] = {}
914 
915  def add_residue_number(self, name, resn):
916  resn = int(resn)
917  self.add_name(name)
918  if resn not in self.db[name]:
919  self.db[name][resn] = {}
920 
921  def add_resolution(self, name, resn, resolution):
922  resn = int(resn)
923  resolution = float(resolution)
924  self.add_name(name)
925  self.add_residue_number(name, resn)
926  if resolution not in self.db[name][resn]:
927  self.db[name][resn][resolution] = []
928 
929  def add_particles(self, name, resn, resolution, particles):
930  resn = int(resn)
931  resolution = float(resolution)
932  self.add_name(name)
933  self.add_residue_number(name, resn)
934  self.add_resolution(name, resn, resolution)
935  self.db[name][resn][resolution] += particles
936  for p in particles:
937  (rh, prf) = self.get_root_hierarchy(p)
938  self.root_hierarchy_dict[p] = rh
939  self.preroot_fragment_hierarchy_dict[p] = prf
940  self.particle_to_name[p] = name
941  if self.model is None:
942  self.model = particles[0].get_model()
943 
944  def get_model(self):
945  return self.model
946 
947  def get_names(self):
948  names = list(self.db.keys())
949  names.sort()
950  return names
951 
952  def get_particles(self, name, resn, resolution):
953  resn = int(resn)
954  resolution = float(resolution)
955  return self.db[name][resn][resolution]
956 
957  def get_particles_at_closest_resolution(self, name, resn, resolution):
958  resn = int(resn)
959  resolution = float(resolution)
960  closestres = min(self.get_residue_resolutions(name, resn),
961  key=lambda x: abs(float(x) - float(resolution)))
962  return self.get_particles(name, resn, closestres)
963 
964  def get_residue_resolutions(self, name, resn):
965  resn = int(resn)
966  resolutions = list(self.db[name][resn].keys())
967  resolutions.sort()
968  return resolutions
969 
970  def get_molecule_resolutions(self, name):
971  resolutions = set()
972  for resn in self.db[name]:
973  resolutions.update(list(self.db[name][resn].keys()))
974  resolutions.sort()
975  return resolutions
976 
977  def get_residue_numbers(self, name):
978  residue_numbers = list(self.db[name].keys())
979  residue_numbers.sort()
980  return residue_numbers
981 
982  def get_particles_by_resolution(self, name, resolution):
983  resolution = float(resolution)
984  particles = []
985  for resn in self.get_residue_numbers(name):
986  result = self.get_particles_at_closest_resolution(
987  name,
988  resn,
989  resolution)
990  pstemp = [p for p in result if p not in particles]
991  particles += pstemp
992  return particles
993 
994  def get_all_particles_by_resolution(self, resolution):
995  resolution = float(resolution)
996  particles = []
997  for name in self.get_names():
998  particles += self.get_particles_by_resolution(name, resolution)
999  return particles
1000 
1001  def get_root_hierarchy(self, particle):
1002  prerootfragment = particle
1003  while IMP.atom.Atom.get_is_setup(particle) or \
1004  IMP.atom.Residue.get_is_setup(particle) or \
1006  if IMP.atom.Atom.get_is_setup(particle):
1007  p = IMP.atom.Atom(particle).get_parent()
1008  elif IMP.atom.Residue.get_is_setup(particle):
1009  p = IMP.atom.Residue(particle).get_parent()
1010  elif IMP.atom.Fragment.get_is_setup(particle):
1011  p = IMP.atom.Fragment(particle).get_parent()
1012  prerootfragment = particle
1013  particle = p
1014  return (
1015  (IMP.atom.Hierarchy(particle), IMP.atom.Hierarchy(prerootfragment))
1016  )
1017 
1018  def get_all_root_hierarchies_by_resolution(self, resolution):
1019  hierarchies = []
1020  resolution = float(resolution)
1021  particles = self.get_all_particles_by_resolution(resolution)
1022  for p in particles:
1023  rh = self.root_hierarchy_dict[p]
1024  if rh not in hierarchies:
1025  hierarchies.append(IMP.atom.Hierarchy(rh))
1026  return hierarchies
1027 
1028  def get_preroot_fragments_by_resolution(self, name, resolution):
1029  fragments = []
1030  resolution = float(resolution)
1031  particles = self.get_particles_by_resolution(name, resolution)
1032  for p in particles:
1033  fr = self.preroot_fragment_hierarchy_dict[p]
1034  if fr not in fragments:
1035  fragments.append(fr)
1036  return fragments
1037 
1038  def show(self, name):
1039  print(name)
1040  for resn in self.get_residue_numbers(name):
1041  print(resn)
1042  for resolution in self.get_residue_resolutions(name, resn):
1043  print("----", resolution)
1044  for p in self.get_particles(name, resn, resolution):
1045  print("--------", p.get_name())
1046 
1047 
1048 def get_prot_name_from_particle(p, list_of_names):
1049  '''Return the component name provided a particle and a list of names'''
1050  root = p
1051  protname = root.get_name()
1052  is_a_bead = False
1053  while not protname in list_of_names:
1054  root0 = root.get_parent()
1055  if root0 == IMP.atom.Hierarchy():
1056  return (None, None)
1057  protname = root0.get_name()
1058 
1059  # check if that is a bead
1060  # this piece of code might be dangerous if
1061  # the hierarchy was called Bead :)
1062  if "Beads" in protname:
1063  is_a_bead = True
1064  root = root0
1065  return (protname, is_a_bead)
1066 
1067 
1069  '''
1070  Retrieve the residue indexes for the given particle.
1071 
1072  The particle must be an instance of Fragment,Residue or Atom
1073  or else returns an empty list
1074  '''
1075  resind = []
1077  resind = IMP.atom.Fragment(hier).get_residue_indexes()
1078  elif IMP.atom.Residue.get_is_setup(hier):
1079  resind = [IMP.atom.Residue(hier).get_index()]
1080  elif IMP.atom.Atom.get_is_setup(hier):
1081  a = IMP.atom.Atom(hier)
1082  resind = [IMP.atom.Residue(a.get_parent()).get_index()]
1083  else:
1084  resind = []
1085  return resind
1086 
1087 
1088 def sort_by_residues(particles):
1089  particles_residues = [(p, IMP.pmi.tools.get_residue_indexes(p))
1090  for p in particles]
1091  sorted_particles_residues = sorted(
1092  particles_residues,
1093  key=lambda tup: tup[1])
1094  particles = [p[0] for p in sorted_particles_residues]
1095  return particles
1096 
1097 
1098 def get_residue_to_particle_map(particles):
1099  # this function returns a dictionary that map particles to residue indexes
1100  particles = sort_by_residues(particles)
1101  particles_residues = [(p, IMP.pmi.tools.get_residue_indexes(p))
1102  for p in particles]
1103  return dict(zip(particles_residues, particles))
1104 
1105 #
1106 # Parallel Computation
1107 #
1108 
1109 
1111  """Synchronize data over a parallel run"""
1112  from mpi4py import MPI
1113  comm = MPI.COMM_WORLD
1114  rank = comm.Get_rank()
1115  number_of_processes = comm.size
1116  comm.Barrier()
1117  if rank != 0:
1118  comm.send(data, dest=0, tag=11)
1119 
1120  elif rank == 0:
1121  for i in range(1, number_of_processes):
1122  data_tmp = comm.recv(source=i, tag=11)
1123  if type(data) == list:
1124  data += data_tmp
1125  elif type(data) == dict:
1126  data.update(data_tmp)
1127  else:
1128  raise TypeError("data not supported, use list or dictionaries")
1129 
1130  for i in range(1, number_of_processes):
1131  comm.send(data, dest=i, tag=11)
1132 
1133  if rank != 0:
1134  data = comm.recv(source=0, tag=11)
1135  return data
1136 
1138  """Synchronize data over a parallel run"""
1139  from mpi4py import MPI
1140  comm = MPI.COMM_WORLD
1141  rank = comm.Get_rank()
1142  number_of_processes = comm.size
1143  comm.Barrier()
1144  if rank != 0:
1145  comm.send(data, dest=0, tag=11)
1146  elif rank == 0:
1147  for i in range(1, number_of_processes):
1148  data_tmp = comm.recv(source=i, tag=11)
1149  for k in data:
1150  data[k]+=data_tmp[k]
1151 
1152  for i in range(1, number_of_processes):
1153  comm.send(data, dest=i, tag=11)
1154 
1155  if rank != 0:
1156  data = comm.recv(source=0, tag=11)
1157  return data
1158 
1159 
1160 #
1161 ### Lists and iterators
1162 #
1163 
1164 def sublist_iterator(l, lmin=None, lmax=None):
1165  '''
1166  Yield all sublists of length >= lmin and <= lmax
1167  '''
1168  if lmin is None:
1169  lmin = 0
1170  if lmax is None:
1171  lmax = len(l)
1172  n = len(l) + 1
1173  for i in range(n):
1174  for j in range(i + 1, n):
1175  if len(l[i:j]) <= lmax and len(l[i:j]) >= lmin:
1176  yield l[i:j]
1177 
1178 
1179 def flatten_list(l):
1180  return [item for sublist in l for item in sublist]
1181 
1182 
1183 def list_chunks_iterator(list, length):
1184  """ Yield successive length-sized chunks from a list.
1185  """
1186  for i in range(0, len(list), length):
1187  yield list[i:i + length]
1188 
1189 
1190 def chunk_list_into_segments(seq, num):
1191  seq = list(seq)
1192  avg = len(seq) / float(num)
1193  out = []
1194  last = 0.0
1195 
1196  while last < len(seq):
1197  out.append(seq[int(last):int(last + avg)])
1198  last += avg
1199 
1200  return out
1201 
1202 
1203 class Segments(object):
1204 
1205  ''' This class stores integers
1206  in ordered compact lists eg:
1207  [[1,2,3],[6,7,8]]
1208  the methods help splitting and merging the internal lists
1209  Example:
1210  s=Segments([1,2,3]) is [[1,2,3]]
1211  s.add(4) is [[1,2,3,4]] (add right)
1212  s.add(3) is [[1,2,3,4]] (item already existing)
1213  s.add(7) is [[1,2,3,4],[7]] (new list)
1214  s.add([8,9]) is [[1,2,3,4],[7,8,9]] (add item right)
1215  s.add([5,6]) is [[1,2,3,4,5,6,7,8,9]] (merge)
1216  s.remove(3) is [[1,2],[4,5,6,7,8,9]] (split)
1217  etc.
1218  '''
1219 
1220  def __init__(self,index):
1221  '''index can be a integer or a list of integers '''
1222  if type(index) is int:
1223  self.segs=[[index]]
1224  elif type(index) is list:
1225  self.segs=[[index[0]]]
1226  for i in index[1:]:
1227  self.add(i)
1228 
1229  def add(self,index):
1230  '''index can be a integer or a list of integers '''
1231  if type(index) is int:
1232  mergeleft=None
1233  mergeright=None
1234  for n,s in enumerate(self.segs):
1235  if index in s:
1236  return 0
1237  else:
1238  if s[0]-index==1:
1239  mergeleft=n
1240  if index-s[-1]==1:
1241  mergeright=n
1242  if mergeright is None and mergeleft is None:
1243  self.segs.append([index])
1244  if not mergeright is None and mergeleft is None:
1245  self.segs[mergeright].append(index)
1246  if not mergeleft is None and mergeright is None:
1247  self.segs[mergeleft]=[index]+self.segs[mergeleft]
1248  if not mergeleft is None and not mergeright is None:
1249  self.segs[mergeright]=self.segs[mergeright]+[index]+self.segs[mergeleft]
1250  del self.segs[mergeleft]
1251 
1252  for n in range(len(self.segs)):
1253  self.segs[n].sort()
1254 
1255  self.segs.sort(key=lambda tup: tup[0])
1256 
1257  elif type(index) is list:
1258  for i in index:
1259  self.add(i)
1260 
1261  def remove(self,index):
1262  '''index can be a integer'''
1263  for n,s in enumerate(self.segs):
1264  if index in s:
1265  if s[0]==index:
1266  self.segs[n]=s[1:]
1267  elif s[-1]==index:
1268  self.segs[n]=s[:-1]
1269  else:
1270  i=self.segs[n].index(index)
1271  self.segs[n]=s[:i]
1272  self.segs.append(s[i+1:])
1273  for n in range(len(self.segs)):
1274  self.segs[n].sort()
1275  if len(self.segs[n])==0:
1276  del self.segs[n]
1277  self.segs.sort(key=lambda tup: tup[0])
1278 
1279  def get_flatten(self):
1280  ''' Returns a flatten list '''
1281  return [item for sublist in self.segs for item in sublist]
1282 
1283 
1284 
1285 #
1286 # COORDINATE MANIPULATION
1287 #
1288 
1289 
1290 def translate_hierarchy(hierarchy, translation_vector):
1291  '''
1292  Apply a translation to a hierarchy along the input vector.
1293  '''
1294  rbs = set()
1295  xyzs = set()
1296  if type(translation_vector) == list:
1297  transformation = IMP.algebra.Transformation3D(
1298  IMP.algebra.Vector3D(translation_vector))
1299  else:
1300  transformation = IMP.algebra.Transformation3D(translation_vector)
1301  for p in IMP.atom.get_leaves(hierarchy):
1303  rbs.add(IMP.core.RigidBody(p))
1305  rb = IMP.core.RigidMember(p).get_rigid_body()
1306  rbs.add(rb)
1307  else:
1308  xyzs.add(p)
1309  for xyz in xyzs:
1310  IMP.core.transform(IMP.core.XYZ(xyz), transformation)
1311  for rb in rbs:
1312  IMP.core.transform(rb, transformation)
1313 
1314 
1315 def translate_hierarchies(hierarchies, translation_vector):
1316  for h in hierarchies:
1317  IMP.pmi.tools.translate_hierarchy(h, translation_vector)
1318 
1319 
1320 def translate_hierarchies_to_reference_frame(hierarchies):
1321  xc = 0
1322  yc = 0
1323  zc = 0
1324  nc = 0
1325  for h in hierarchies:
1326  for p in IMP.atom.get_leaves(h):
1327  coor = IMP.core.XYZ(p).get_coordinates()
1328  nc += 1
1329  xc += coor[0]
1330  yc += coor[1]
1331  zc += coor[2]
1332  xc = xc / nc
1333  yc = yc / nc
1334  zc = zc / nc
1335  IMP.pmi.tools.translate_hierarchies(hierarchies, (-xc, -yc, -zc))
1336 
1337 
1338 #
1339 # Tools to simulate data
1340 #
1341 
1342 def normal_density_function(expected_value, sigma, x):
1343  return (
1344  1 / math.sqrt(2 * math.pi) / sigma *
1345  math.exp(-(x - expected_value) ** 2 / 2 / sigma / sigma)
1346  )
1347 
1348 
1349 def log_normal_density_function(expected_value, sigma, x):
1350  return (
1351  1 / math.sqrt(2 * math.pi) / sigma / x *
1352  math.exp(-(math.log(x / expected_value) ** 2 / 2 / sigma / sigma))
1353  )
1354 
1355 
1356 def get_random_residue_pairs(representation, resolution,
1357  number,
1358  max_distance=None,
1359  avoid_same_particles=False,
1360  names=None):
1361 
1362  particles = []
1363  if names is None:
1364  names=list(representation.hier_dict.keys())
1365 
1366  for name in names:
1367  prot = representation.hier_dict[name]
1368  particles += select(representation,name=name,resolution=resolution)
1369  random_residue_pairs = []
1370  while len(random_residue_pairs)<=number:
1371  p1 = random.choice(particles)
1372  p2 = random.choice(particles)
1373  if max_distance is not None and \
1374  core.get_distance(core.XYZ(p1), core.XYZ(p2)) > max_distance:
1375  continue
1376  r1 = random.choice(IMP.pmi.tools.get_residue_indexes(p1))
1377  r2 = random.choice(IMP.pmi.tools.get_residue_indexes(p2))
1378  if r1==r2 and avoid_same_particles: continue
1379  name1 = representation.get_prot_name_from_particle(p1)
1380  name2 = representation.get_prot_name_from_particle(p2)
1381  random_residue_pairs.append((name1, r1, name2, r2))
1382 
1383  return random_residue_pairs
1384 
1385 
1386 def get_random_data_point(
1387  expected_value,
1388  ntrials,
1389  sensitivity,
1390  sigma,
1391  outlierprob,
1392  begin_end_nbins_tuple,
1393  log=False,
1394  loggrid=False):
1395  import random
1396 
1397  begin = begin_end_nbins_tuple[0]
1398  end = begin_end_nbins_tuple[1]
1399  nbins = begin_end_nbins_tuple[2]
1400 
1401  if not loggrid:
1402  fmod_grid = get_grid(begin, end, nbins, True)
1403  else:
1404  fmod_grid = get_log_grid(begin, end, nbins)
1405 
1406  norm = 0
1407  cumul = []
1408  cumul.append(0)
1409 
1410  a = []
1411  for i in range(0, ntrials):
1412  a.append([random.random(), True])
1413 
1414  if sigma != 0.0:
1415  for j in range(1, len(fmod_grid)):
1416  fj = fmod_grid[j]
1417  fjm1 = fmod_grid[j - 1]
1418  df = fj - fjm1
1419 
1420  if not log:
1421  pj = normal_density_function(expected_value, sigma, fj)
1422  pjm1 = normal_density_function(expected_value, sigma, fjm1)
1423  else:
1424  pj = log_normal_density_function(expected_value, sigma, fj)
1425  pjm1 = log_normal_density_function(expected_value, sigma, fjm1)
1426 
1427  norm += (pj + pjm1) / 2.0 * df
1428  cumul.append(norm)
1429  # print fj, pj
1430 
1431  random_points = []
1432 
1433  for i in range(len(cumul)):
1434  # print i,a, cumul[i], norm
1435  for aa in a:
1436  if (aa[0] <= cumul[i] / norm and aa[1]):
1437  random_points.append(
1438  int(fmod_grid[i] / sensitivity) * sensitivity)
1439  aa[1] = False
1440 
1441  else:
1442  random_points = [expected_value] * ntrials
1443 
1444  for i in range(len(random_points)):
1445  if random.random() < outlierprob:
1446  a = random.uniform(begin, end)
1447  random_points[i] = int(a / sensitivity) * sensitivity
1448  print(random_points)
1449  '''
1450  for i in range(ntrials):
1451  if random.random() > OUTLIERPROB_:
1452  r=truncnorm.rvs(0.0,1.0,expected_value,BETA_)
1453  if r>1.0: print r,expected_value,BETA_
1454  else:
1455  r=random.random()
1456  random_points.append(int(r/sensitivity)*sensitivity)
1457  '''
1458 
1459  rmean = 0.
1460  rmean2 = 0.
1461  for r in random_points:
1462  rmean += r
1463  rmean2 += r * r
1464 
1465  rmean /= float(ntrials)
1466  rmean2 /= float(ntrials)
1467  stddev = math.sqrt(max(rmean2 - rmean * rmean, 0.))
1468  return rmean, stddev
1469 
1470 def print_multicolumn(list_of_strings, ncolumns=2, truncate=40):
1471 
1472  l = list_of_strings
1473 
1474  cols = ncolumns
1475  # add empty entries after l
1476  for i in range(len(l) % cols):
1477  l.append(" ")
1478 
1479  split = [l[i:i + len(l) / cols] for i in range(0, len(l), len(l) / cols)]
1480  for row in zip(*split):
1481  print("".join(str.ljust(i, truncate) for i in row))
1482 
1483 class ColorChange(object):
1484  '''Change color code to hexadecimal to rgb'''
1485  def __init__(self):
1486  self._NUMERALS = '0123456789abcdefABCDEF'
1487  self._HEXDEC = dict((v, int(v, 16)) for v in (x+y for x in self._NUMERALS for y in self._NUMERALS))
1488  self.LOWERCASE, self.UPPERCASE = 'x', 'X'
1489 
1490  def rgb(self,triplet):
1491  return float(self._HEXDEC[triplet[0:2]]), float(self._HEXDEC[triplet[2:4]]), float(self._HEXDEC[triplet[4:6]])
1492 
1493  def triplet(self,rgb, lettercase=None):
1494  if lettercase is None: lettercase=self.LOWERCASE
1495  return format(rgb[0]<<16 | rgb[1]<<8 | rgb[2], '06'+lettercase)
1496 
1497 # -------------- Collections --------------- #
1498 
1499 class OrderedSet(collections.MutableSet):
1500 
1501  def __init__(self, iterable=None):
1502  self.end = end = []
1503  end += [None, end, end] # sentinel node for doubly linked list
1504  self.map = {} # key --> [key, prev, next]
1505  if iterable is not None:
1506  self |= iterable
1507 
1508  def __len__(self):
1509  return len(self.map)
1510 
1511  def __contains__(self, key):
1512  return key in self.map
1513 
1514  def add(self, key):
1515  if key not in self.map:
1516  end = self.end
1517  curr = end[1]
1518  curr[2] = end[1] = self.map[key] = [key, curr, end]
1519 
1520  def discard(self, key):
1521  if key in self.map:
1522  key, prev, next = self.map.pop(key)
1523  prev[2] = next
1524  next[1] = prev
1525 
1526  def __iter__(self):
1527  end = self.end
1528  curr = end[2]
1529  while curr is not end:
1530  yield curr[0]
1531  curr = curr[2]
1532 
1533  def __reversed__(self):
1534  end = self.end
1535  curr = end[1]
1536  while curr is not end:
1537  yield curr[0]
1538  curr = curr[1]
1539 
1540  def pop(self, last=True):
1541  if not self:
1542  raise KeyError('set is empty')
1543  if last:
1544  key = self.end[1][0]
1545  else:
1546  key = self.end[2][0]
1547  self.discard(key)
1548  return key
1549 
1550  def __repr__(self):
1551  if not self:
1552  return '%s()' % (self.__class__.__name__,)
1553  return '%s(%r)' % (self.__class__.__name__, list(self))
1554 
1555  def __eq__(self, other):
1556  if isinstance(other, OrderedSet):
1557  return len(self) == len(other) and list(self) == list(other)
1558  return set(self) == set(other)
1559 
1560 
1561 class OrderedDefaultDict(OrderedDict):
1562  """Store objects in order they were added, but with default type.
1563  Source: http://stackoverflow.com/a/4127426/2608793
1564  """
1565  def __init__(self, *args, **kwargs):
1566  if not args:
1567  self.default_factory = None
1568  else:
1569  if not (args[0] is None or callable(args[0])):
1570  raise TypeError('first argument must be callable or None')
1571  self.default_factory = args[0]
1572  args = args[1:]
1573  super(OrderedDefaultDict, self).__init__(*args, **kwargs)
1574 
1575  def __missing__ (self, key):
1576  if self.default_factory is None:
1577  raise KeyError(key)
1578  self[key] = default = self.default_factory()
1579  return default
1580 
1581  def __reduce__(self): # optional, for pickle support
1582  args = (self.default_factory,) if self.default_factory else ()
1583  return self.__class__, args, None, None, self.iteritems()
1584 
1585 # -------------- PMI2 Tools --------------- #
1586 
1587 def set_coordinates_from_rmf(hier,rmf_fn,frame_num=0):
1588  """Extract frame from RMF file and fill coordinates. Must be identical topology.
1589  @param hier The (System) hierarchy to fill (e.g. after you've built it)
1590  @param rmf_fn The file to extract from
1591  @param frame_num The frame number to extract
1592  """
1593  rh = RMF.open_rmf_file_read_only(rmf_fn)
1594  IMP.rmf.link_hierarchies(rh,[hier])
1595  IMP.rmf.load_frame(rh, RMF.FrameID(frame_num))
1596  del rh
1597 
1598 def input_adaptor(stuff,
1599  pmi_resolution=0,
1600  flatten=False,
1601  selection_tuple=None,
1602  warn_about_slices=True):
1603  """Adapt things for PMI (degrees of freedom, restraints, ...)
1604  Returns list of list of hierarchies, separated into Molecules if possible.
1605  (iterable of ^2) hierarchy -> returns input as list of list of hierarchies, only one entry.
1606  (iterable of ^2) PMI::System/State/Molecule/TempResidue ->
1607  returns residue hierarchies, grouped in molecules, at requested resolution
1608  @param stuff Can be one of the following inputs:
1609  IMP Hierarchy, PMI System/State/Molecule/TempResidue, or a list/set (of list/set) of them.
1610  Must be uniform input, however. No mixing object types.
1611  @param pmi_resolution For selecting, only does it if you pass PMI objects. Set it to "all"
1612  if you want all resolutions!
1613  @param flatten Set to True if you just want all hierarchies in one list.
1614  @param warn_about_slices Print a warning if you are requesting only part of a bead.
1615  Sometimes you just don't care!
1616  \note since this relies on IMP::atom::Selection, this will not return any objects if they weren't built!
1617  But there should be no problem if you request unbuilt residues, they should be ignored.
1618  """
1619 
1620  # if iterable (of iterable), make sure uniform type
1621  def get_ok_iter(it):
1622  type_set = set(type(a) for a in it)
1623  if len(type_set)!=1:
1624  raise Exception('input_adaptor: can only pass one type of object at a time')
1625  xtp = type_set.pop()
1626  return xtp,list(it)
1627 
1628  if stuff is None:
1629  return stuff
1630  if hasattr(stuff,'__iter__'):
1631  if len(stuff)==0:
1632  return stuff
1633  tp,thelist = get_ok_iter(stuff)
1634 
1635  # iter of iter of should be ok
1636  if hasattr(next(iter(thelist)),'__iter__'):
1637  flatlist = [i for sublist in thelist for i in sublist]
1638  tp,thelist = get_ok_iter(flatlist)
1639 
1640  stuff = thelist
1641  else:
1642  tp = type(stuff)
1643  stuff = [stuff]
1644 
1645  # now that things are ok, do selection if requested
1646  hier_list = []
1647  pmi_input = False
1650  # if PMI, perform selection using gathered indexes
1651  pmi_input = True
1652  indexes_per_mol = OrderedDefaultDict(list) #key is Molecule object, value are residues
1653  if tp==IMP.pmi.topology.System:
1654  for system in stuff:
1655  for state in system.get_states():
1656  mdict = state.get_molecules()
1657  for molname in mdict:
1658  for copy in mdict[molname]:
1659  indexes_per_mol[copy] += [r.get_index() for r in copy.get_residues()]
1660  elif tp==IMP.pmi.topology.State:
1661  for state in stuff:
1662  mdict = state.get_molecules()
1663  for molname in mdict:
1664  for copy in mdict[molname]:
1665  indexes_per_mol[copy] += [r.get_index() for r in copy.get_residues()]
1666  elif tp==IMP.pmi.topology.Molecule:
1667  for molecule in stuff:
1668  indexes_per_mol[molecule] += [r.get_index() for r in molecule.get_residues()]
1670  for tempres in stuff:
1671  indexes_per_mol[tempres.get_molecule()].append(tempres.get_index())
1672  for mol in indexes_per_mol:
1673  if pmi_resolution=='all':
1674  # because you select from the molecule,
1675  # this will start the search from the base resolution
1676  ps = select_at_all_resolutions(mol.get_hierarchy(),
1677  residue_indexes=indexes_per_mol[mol])
1678  else:
1679  sel = IMP.atom.Selection(mol.get_hierarchy(),
1680  resolution=pmi_resolution,
1681  residue_indexes=indexes_per_mol[mol])
1682  ps = sel.get_selected_particles()
1683 
1684  # check that you don't have any incomplete fragments!
1685  if warn_about_slices:
1686  rset = set(indexes_per_mol[mol])
1687  for p in ps:
1689  fset = set(IMP.atom.Fragment(p).get_residue_indexes())
1690  if not fset <= rset:
1691  print('BAD')
1692  minset = min(fset)
1693  maxset = max(fset)
1694  found = fset&rset
1695  minf = min(found)
1696  maxf = max(found)
1697  resbreak = maxf if minf==minset else minset-1
1698  print('WARNING: You are trying to select only part of the bead %s:%i-%i.\n'
1699  'The residues you requested are %i-%i. You can fix this by:\n'
1700  '1) requesting the whole bead/none of it or\n'
1701  '2) break the bead up by passing bead_extra_breaks=[\'%i\'] in '
1702  'molecule.add_representation()'
1703  %(mol.get_name(),minset,maxset,minf,maxf,resbreak))
1704  hier_list.append([IMP.atom.Hierarchy(p) for p in ps])
1705  else:
1706  try:
1707  if IMP.atom.Hierarchy.get_is_setup(stuff[0]):
1708  if flatten:
1709  hier_list = stuff
1710  else:
1711  hier_list = [stuff]
1712  else:
1713  raise Exception('input_adaptor: you passed something of type',tp)
1714  except:
1715  raise Exception('input_adaptor: you passed something of type',tp)
1716 
1717  if flatten and pmi_input:
1718  return [h for sublist in hier_list for h in sublist]
1719  else:
1720  return hier_list
1721 
1722 
1724  """Returns sequence-sorted segments array, each containing the first particle
1725  the last particle and the first residue index."""
1726 
1727  from operator import itemgetter
1728  hiers=IMP.pmi.tools.input_adaptor(mol)
1729  if len(hiers)>1:
1730  raise Exception("IMP.pmi.tools.get_sorted_segments: only pass stuff from one Molecule, please")
1731  hiers = hiers[0]
1732  SortedSegments = []
1733  for h in hiers:
1734  try:
1735  start = IMP.atom.Hierarchy(h).get_children()[0]
1736  except:
1737  start = IMP.atom.Hierarchy(h)
1738 
1739  try:
1740  end = IMP.atom.Hierarchy(h).get_children()[-1]
1741  except:
1742  end = IMP.atom.Hierarchy(h)
1743 
1744  startres = IMP.pmi.tools.get_residue_indexes(start)[0]
1745  endres = IMP.pmi.tools.get_residue_indexes(end)[-1]
1746  SortedSegments.append((start, end, startres))
1747  SortedSegments = sorted(SortedSegments, key=itemgetter(2))
1748  return SortedSegments
1749 
1750 def display_bonds(mol):
1751  """Decorate the sequence-consecutive particles from a PMI2 molecule with a bond,
1752  so that they appear connected in the rmf file"""
1753  SortedSegments=get_sorted_segments(mol)
1754  for x in range(len(SortedSegments) - 1):
1755 
1756  last = SortedSegments[x][1]
1757  first = SortedSegments[x + 1][0]
1758 
1759  p1 = last.get_particle()
1760  p2 = first.get_particle()
1761  if not IMP.atom.Bonded.get_is_setup(p1):
1763  if not IMP.atom.Bonded.get_is_setup(p2):
1765 
1768  IMP.atom.Bonded(p1),
1769  IMP.atom.Bonded(p2),1)
1770 
1771 
1772 def get_residue_type_from_one_letter_code(code,is_nucleic=None):
1773  if not is_nucleic:
1774  threetoone = {'ALA': 'A', 'ARG': 'R', 'ASN': 'N', 'ASP': 'D',
1775  'CYS': 'C', 'GLU': 'E', 'GLN': 'Q', 'GLY': 'G',
1776  'HIS': 'H', 'ILE': 'I', 'LEU': 'L', 'LYS': 'K',
1777  'MET': 'M', 'PHE': 'F', 'PRO': 'P', 'SER': 'S',
1778  'THR': 'T', 'TRP': 'W', 'TYR': 'Y', 'VAL': 'V', 'UNK': 'X'}
1779  else:
1780  threetoone = {'ADE': 'A', 'URA': 'U', 'CYT': 'C', 'GUA': 'G',
1781  'THY': 'T', 'UNK': 'X'}
1782  one_to_three={}
1783  for k in threetoone:
1784  one_to_three[threetoone[k]] = k
1785  return IMP.atom.ResidueType(one_to_three[code])
1786 
1787 
1788 def get_all_leaves(list_of_hs):
1789  """ Just get the leaves from a list of hierarchies """
1790  lvs = list(itertools.chain.from_iterable(IMP.atom.get_leaves(item) for item in list_of_hs))
1791  return lvs
1792 
1793 
1795  hiers=None,
1796  **kwargs):
1797  """Perform selection using the usual keywords but return ALL resolutions (BEADS and GAUSSIANS).
1798  Returns in flat list!
1799  """
1800 
1801  if hiers is None:
1802  hiers = []
1803  if hier is not None:
1804  hiers.append(hier)
1805  if len(hiers)==0:
1806  print("WARNING: You passed nothing to select_at_all_resolutions()")
1807  return []
1808  ret = OrderedSet()
1809  for hsel in hiers:
1810  try:
1811  htest = IMP.atom.Hierarchy.get_is_setup(hsel)
1812  except:
1813  raise Exception('select_at_all_resolutions: you have to pass an IMP Hierarchy')
1814  if not htest:
1815  raise Exception('select_at_all_resolutions: you have to pass an IMP Hierarchy')
1816  if 'resolution' in kwargs or 'representation_type' in kwargs:
1817  raise Exception("don't pass resolution or representation_type to this function")
1818  selB = IMP.atom.Selection(hsel,resolution=IMP.atom.ALL_RESOLUTIONS,
1819  representation_type=IMP.atom.BALLS,**kwargs)
1820  selD = IMP.atom.Selection(hsel,resolution=IMP.atom.ALL_RESOLUTIONS,
1821  representation_type=IMP.atom.DENSITIES,**kwargs)
1822  ret |= OrderedSet(selB.get_selected_particles())
1823  ret |= OrderedSet(selD.get_selected_particles())
1824  return list(ret)
1825 
1826 
1828  target_ps,
1829  sel_zone,
1830  entire_residues,
1831  exclude_backbone):
1832  """Utility to retrieve particles from a hierarchy within a
1833  zone around a set of ps.
1834  @param hier The hierarchy in which to look for neighbors
1835  @param target_ps The particles for zoning
1836  @param sel_zone The maximum distance
1837  @param entire_residues If True, will grab entire residues
1838  @param exclude_backbone If True, will only return sidechain particles
1839  """
1840 
1841  test_sel = IMP.atom.Selection(hier)
1842  backbone_types=['C','N','CB','O']
1843  if exclude_backbone:
1844  test_sel -= IMP.atom.Selection(hier,atom_types=[IMP.atom.AtomType(n)
1845  for n in backbone_types])
1846  test_ps = test_sel.get_selected_particles()
1847  nn = IMP.algebra.NearestNeighbor3D([IMP.core.XYZ(p).get_coordinates()
1848  for p in test_ps])
1849  zone = set()
1850  for target in target_ps:
1851  zone|=set(nn.get_in_ball(IMP.core.XYZ(target).get_coordinates(),sel_zone))
1852  zone_ps = [test_ps[z] for z in zone]
1853  if entire_residues:
1854  final_ps = set()
1855  for z in zone_ps:
1856  final_ps|=set(IMP.atom.Hierarchy(z).get_parent().get_children())
1857  zone_ps = [h.get_particle() for h in final_ps]
1858  return zone_ps
1859 
1860 
1862  """Returns unique objects in original order"""
1863  rbs = set()
1864  beads = []
1865  rbs_ordered = []
1866  if not hasattr(hiers,'__iter__'):
1867  hiers = [hiers]
1868  for p in get_all_leaves(hiers):
1870  rb = IMP.core.RigidMember(p).get_rigid_body()
1871  if rb not in rbs:
1872  rbs.add(rb)
1873  rbs_ordered.append(rb)
1875  rb = IMP.core.NonRigidMember(p).get_rigid_body()
1876  if rb not in rbs:
1877  rbs.add(rb)
1878  rbs_ordered.append(rb)
1879  beads.append(p)
1880  else:
1881  beads.append(p)
1882  return rbs_ordered,beads
1883 
1884 
1885 def get_densities(input_objects):
1886  """Given a list of PMI objects, return all density hierarchies within
1887  these objects. The output of this function can be inputted into
1888  things such as EM restraints.
1889  """
1890  stuff=input_adaptor(input_objects, flatten=True)
1891  densities=[]
1892  try:
1893  for i in iter(stuff):
1894  densities.append(IMP.atom.Selection(i,representation_type=IMP.atom.DENSITIES).get_selected_particles())
1895  except TypeError as te:
1896  densities.append(IMP.atom.Selection(stuff,representation_type=IMP.atom.DENSITIES).get_selected_particles())
1897 
1898  return [h for sublist in densities for h in sublist]
1899 
1900 def shuffle_configuration(objects,
1901  max_translation=300., max_rotation=2.0 * pi,
1902  avoidcollision_rb=True, avoidcollision_fb=False,
1903  cutoff=10.0, niterations=100,
1904  bounding_box=None,
1905  excluded_rigid_bodies=[],
1906  hierarchies_excluded_from_collision=[],
1907  hierarchies_included_in_collision=[],
1908  verbose=False):
1909  """Shuffle particles. Used to restart the optimization.
1910  The configuration of the system is initialized by placing each
1911  rigid body and each bead randomly in a box with a side of
1912  max_translation angstroms, and far enough from each other to
1913  prevent any steric clashes. The rigid bodies are also randomly rotated.
1914  @param objects Can be one of the following inputs:
1915  IMP Hierarchy, PMI System/State/Molecule/TempResidue, or a list/set of them
1916  @param max_translation Max translation (rbs and flexible beads)
1917  @param max_rotation Max rotation (rbs only)
1918  @param avoidcollision_rb check if the particle/rigid body was
1919  placed close to another particle; uses the optional
1920  arguments cutoff and niterations
1921  @param avoidcollision_fb Advanced. Generally you want this False because it's hard to shuffle beads.
1922  @param cutoff Distance less than this is a collision
1923  @param niterations How many times to try avoiding collision
1924  @param bounding_box Only shuffle particles within this box. Defined by ((x1,y1,z1),(x2,y2,z2)).
1925  @param excluded_rigid_bodies Don't shuffle these rigid body objects
1926  @param hierarchies_excluded_from_collision Don't count collision with these bodies
1927  @param hierarchies_included_in_collision Hierarchies that are not shuffled, but should be included in collision calculation (for fixed regions)
1928  @param verbose Give more output
1929  \note Best to only call this function after you've set up degrees of freedom
1930  For debugging purposes, returns: <shuffled indexes>, <collision avoided indexes>
1931  """
1932 
1933  ### checking input
1934  hierarchies = IMP.pmi.tools.input_adaptor(objects,
1935  pmi_resolution='all',
1936  flatten=True)
1937  rigid_bodies,flexible_beads = get_rbs_and_beads(hierarchies)
1938  if len(rigid_bodies)>0:
1939  mdl = rigid_bodies[0].get_model()
1940  elif len(flexible_beads)>0:
1941  mdl = flexible_beads[0].get_model()
1942  else:
1943  raise Exception("Could not find any particles in the hierarchy")
1944  if len(rigid_bodies) == 0:
1945  print("shuffle_configuration: rigid bodies were not intialized")
1946 
1947  ### gather all particles
1949  gcpf.set_distance(cutoff)
1950 
1951  # Add particles from excluded hierarchies to excluded list
1952  collision_excluded_hierarchies = IMP.pmi.tools.input_adaptor(hierarchies_excluded_from_collision,
1953  pmi_resolution='all',
1954  flatten=True)
1955 
1956  collision_included_hierarchies = IMP.pmi.tools.input_adaptor(hierarchies_included_in_collision,
1957  pmi_resolution='all',
1958  flatten=True)
1959 
1960  collision_excluded_idxs = set([l.get_particle().get_index() for h in collision_excluded_hierarchies \
1961  for l in IMP.core.get_leaves(h)])
1962 
1963  collision_included_idxs = set([l.get_particle().get_index() for h in collision_included_hierarchies \
1964  for l in IMP.core.get_leaves(h)])
1965 
1966  # Excluded collision with Gaussians
1967  all_idxs = [] #expand to representations?
1968  for p in IMP.pmi.tools.get_all_leaves(hierarchies):
1970  all_idxs.append(p.get_particle_index())
1972  collision_excluded_idxs.add(p.get_particle_index())
1973 
1974  print(len(all_idxs), len(collision_included_idxs), len(collision_excluded_idxs))
1975 
1976  if bounding_box is not None:
1977  ((x1, y1, z1), (x2, y2, z2)) = bounding_box
1978  ub = IMP.algebra.Vector3D(x1, y1, z1)
1979  lb = IMP.algebra.Vector3D(x2, y2, z2)
1980  bb = IMP.algebra.BoundingBox3D(ub, lb)
1981 
1982  all_idxs = set(all_idxs) | collision_included_idxs
1983  all_idxs = all_idxs - collision_excluded_idxs
1984  print(len(all_idxs), len(collision_included_idxs), len(collision_excluded_idxs))
1985  debug = []
1986  print('shuffling', len(rigid_bodies), 'rigid bodies')
1987  for rb in rigid_bodies:
1988  if rb not in excluded_rigid_bodies:
1989  # gather particles to avoid with this transform
1990  if avoidcollision_rb:
1991  rb_idxs = set(rb.get_member_particle_indexes()) - \
1992  collision_excluded_idxs
1993  other_idxs = all_idxs - rb_idxs
1994  print("----INOI", rb, len(other_idxs), len(rb_idxs))
1995  if not other_idxs:
1996  continue
1997 
1998  # iterate, trying to avoid collisions
1999  niter = 0
2000  while niter < niterations:
2001  rbxyz = (rb.get_x(), rb.get_y(), rb.get_z())
2002 
2003  # local transform
2004  if bounding_box:
2005  translation = IMP.algebra.get_random_vector_in(bb)
2006  # First move to origin
2008  -IMP.core.XYZ(rb).get_coordinates())
2009  IMP.core.transform(rb, transformation_orig)
2010  old_coord = IMP.core.XYZ(rb).get_coordinates()
2012  transformation = IMP.algebra.Transformation3D(rotation,
2013  translation)
2014 
2015  else:
2017  rbxyz,
2018  max_translation,
2019  max_rotation)
2020 
2021  debug.append([rb, other_idxs if avoidcollision_rb else set()])
2022  IMP.core.transform(rb, transformation)
2023  print("TSFM", rb)
2024 
2025  # check collisions
2026  if avoidcollision_rb:
2027  mdl.update()
2028  npairs = len(gcpf.get_close_pairs(mdl,
2029  list(other_idxs),
2030  list(rb_idxs)))
2031  #print("NPAIRS:", npairs)
2032  if npairs==0:
2033  break
2034  else:
2035  niter += 1
2036  if verbose:
2037  print("shuffle_configuration: rigid body placed close to other %d particles, trying again..." % npairs)
2038  print("shuffle_configuration: rigid body name: " + rb.get_name())
2039  if niter == niterations:
2040  raise ValueError("tried the maximum number of iterations to avoid collisions, increase the distance cutoff")
2041  else:
2042  break
2043 
2044  print('shuffling', len(flexible_beads), 'flexible beads')
2045  for fb in flexible_beads:
2046  # gather particles to avoid
2047  if avoidcollision_fb:
2048  fb_idxs = set(IMP.get_indexes([fb]))
2049  other_idxs = all_idxs - fb_idxs
2050  if not other_idxs:
2051  continue
2052 
2053  # iterate, trying to avoid collisions
2054  niter = 0
2055  while niter < niterations:
2056  if bounding_box:
2057  translation = IMP.algebra.get_random_vector_in(bb)
2058  transformation = IMP.algebra.Transformation3D(translation)
2059  else:
2060  fbxyz = IMP.core.XYZ(fb).get_coordinates()
2061  transformation = IMP.algebra.get_random_local_transformation(fbxyz,
2062  max_translation,
2063  max_rotation)
2064 
2065  # For gaussians, treat this fb as an rb
2067  xyz=[fb.get_value(IMP.FloatKey(4)), fb.get_value(IMP.FloatKey(5)), fb.get_value(IMP.FloatKey(6))]
2068  xyz_transformed=transformation.get_transformed(xyz)
2069  if bounding_box:
2070  # Translate to origin
2071  fb.set_value(IMP.FloatKey(4),xyz_transformed[0]-xyz[0])
2072  fb.set_value(IMP.FloatKey(5),xyz_transformed[1]-xyz[1])
2073  fb.set_value(IMP.FloatKey(6),xyz_transformed[2]-xyz[2])
2074  debug.append([xyz,other_idxs if avoidcollision_fb else set()])
2075  xyz2=[fb.get_value(IMP.FloatKey(4)), fb.get_value(IMP.FloatKey(5)), fb.get_value(IMP.FloatKey(6))]
2076  else:
2077  fb.set_value(IMP.FloatKey(4),xyz_transformed[0])
2078  fb.set_value(IMP.FloatKey(5),xyz_transformed[1])
2079  fb.set_value(IMP.FloatKey(6),xyz_transformed[2])
2080  debug.append([xyz,other_idxs if avoidcollision_fb else set()])
2081  else:
2082  d =IMP.core.XYZ(fb)
2083  if bounding_box:
2084  # Translate to origin first
2085  IMP.core.transform(d, -d.get_coordinates())
2086  d =IMP.core.XYZ(fb)
2087  debug.append([d,other_idxs if avoidcollision_fb else set()])
2088  IMP.core.transform(d, transformation)
2089 
2090  if avoidcollision_fb:
2091  mdl.update()
2092  npairs = len(gcpf.get_close_pairs(mdl,
2093  list(other_idxs),
2094  list(fb_idxs)))
2095 
2096  if npairs==0:
2097  break
2098  else:
2099  niter += 1
2100  print("shuffle_configuration: floppy body placed close to other %d particles, trying again..." % npairs)
2101  if niter == niterations:
2102  raise ValueError("tried the maximum number of iterations to avoid collisions, increase the distance cutoff")
2103  else:
2104  break
2105  return debug
2106 
2107 def color2rgb(colorname):
2108  """Given a chimera color name, return RGB"""
2109  d = {'aquamarine': (0.4980392156862745, 1.0, 0.8313725490196079),
2110  'black': (0.0, 0.0, 0.0),
2111  'blue': (0.0, 0.0, 1.0),
2112  'brown': (0.6470588235294118, 0.16470588235294117, 0.16470588235294117),
2113  'chartreuse': (0.4980392156862745, 1.0, 0.0),
2114  'coral': (1.0, 0.4980392156862745, 0.3137254901960784),
2115  'cornflower blue': (0.39215686274509803, 0.5843137254901961, 0.9294117647058824),
2116  'cyan': (0.0, 1.0, 1.0),
2117  'dark cyan': (0.0, 0.5450980392156862, 0.5450980392156862),
2118  'dark gray': (0.6627450980392157, 0.6627450980392157, 0.6627450980392157),
2119  'dark green': (0.0, 0.39215686274509803, 0.0),
2120  'dark khaki': (0.7411764705882353, 0.7176470588235294, 0.4196078431372549),
2121  'dark magenta': (0.5450980392156862, 0.0, 0.5450980392156862),
2122  'dark olive green': (0.3333333333333333, 0.4196078431372549, 0.1843137254901961),
2123  'dark red': (0.5450980392156862, 0.0, 0.0),
2124  'dark slate blue': (0.2823529411764706, 0.23921568627450981, 0.5450980392156862),
2125  'dark slate gray': (0.1843137254901961, 0.30980392156862746, 0.30980392156862746),
2126  'deep pink': (1.0, 0.0784313725490196, 0.5764705882352941),
2127  'deep sky blue': (0.0, 0.7490196078431373, 1.0),
2128  'dim gray': (0.4117647058823529, 0.4117647058823529, 0.4117647058823529),
2129  'dodger blue': (0.11764705882352941, 0.5647058823529412, 1.0),
2130  'firebrick': (0.6980392156862745, 0.13333333333333333, 0.13333333333333333),
2131  'forest green': (0.13333333333333333, 0.5450980392156862, 0.13333333333333333),
2132  'gold': (1.0, 0.8431372549019608, 0.0),
2133  'goldenrod': (0.8549019607843137, 0.6470588235294118, 0.12549019607843137),
2134  'gray': (0.7450980392156863, 0.7450980392156863, 0.7450980392156863),
2135  'green': (0.0, 1.0, 0.0),
2136  'hot pink': (1.0, 0.4117647058823529, 0.7058823529411765),
2137  'khaki': (0.9411764705882353, 0.9019607843137255, 0.5490196078431373),
2138  'light blue': (0.6784313725490196, 0.8470588235294118, 0.9019607843137255),
2139  'light gray': (0.8274509803921568, 0.8274509803921568, 0.8274509803921568),
2140  'light green': (0.5647058823529412, 0.9333333333333333, 0.5647058823529412),
2141  'light sea green': (0.12549019607843137, 0.6980392156862745, 0.6666666666666666),
2142  'lime green': (0.19607843137254902, 0.803921568627451, 0.19607843137254902),
2143  'magenta': (1.0, 0.0, 1.0),
2144  'medium blue': (0.19607843137254902, 0.19607843137254902, 0.803921568627451),
2145  'medium purple': (0.5764705882352941, 0.4392156862745098, 0.8588235294117647),
2146  'navy blue': (0.0, 0.0, 0.5019607843137255),
2147  'olive drab': (0.4196078431372549, 0.5568627450980392, 0.13725490196078433),
2148  'orange red': (1.0, 0.27058823529411763, 0.0),
2149  'orange': (1.0, 0.4980392156862745, 0.0),
2150  'orchid': (0.8549019607843137, 0.4392156862745098, 0.8392156862745098),
2151  'pink': (1.0, 0.7529411764705882, 0.796078431372549),
2152  'plum': (0.8666666666666667, 0.6274509803921569, 0.8666666666666667),
2153  'purple': (0.6274509803921569, 0.12549019607843137, 0.9411764705882353),
2154  'red': (1.0, 0.0, 0.0),
2155  'rosy brown': (0.7372549019607844, 0.5607843137254902, 0.5607843137254902),
2156  'salmon': (0.9803921568627451, 0.5019607843137255, 0.4470588235294118),
2157  'sandy brown': (0.9568627450980393, 0.6431372549019608, 0.3764705882352941),
2158  'sea green': (0.1803921568627451, 0.5450980392156862, 0.3411764705882353),
2159  'sienna': (0.6274509803921569, 0.3215686274509804, 0.17647058823529413),
2160  'sky blue': (0.5294117647058824, 0.807843137254902, 0.9215686274509803),
2161  'slate gray': (0.4392156862745098, 0.5019607843137255, 0.5647058823529412),
2162  'spring green': (0.0, 1.0, 0.4980392156862745),
2163  'steel blue': (0.27450980392156865, 0.5098039215686274, 0.7058823529411765),
2164  'tan': (0.8235294117647058, 0.7058823529411765, 0.5490196078431373),
2165  'turquoise': (0.25098039215686274, 0.8784313725490196, 0.8156862745098039),
2166  'violet red': (0.8156862745098039, 0.12549019607843137, 0.5647058823529412),
2167  'white': (1.0, 1.0, 1.0),
2168  'yellow': (1.0, 1.0, 0.0)}
2169  return d[colorname]
2170 
2171 class Colors(object):
2172  def __init__(self):
2173  self.colors={
2174  "reds":[("maroon","#800000",(128,0,0)),("dark red","#8B0000",(139,0,0)),
2175  ("brown","#A52A2A",(165,42,42)),("firebrick","#B22222",(178,34,34)),
2176  ("crimson","#DC143C",(220,20,60)),("red","#FF0000",(255,0,0)),
2177  ("tomato","#FF6347",(255,99,71)),("coral","#FF7F50",(255,127,80)),
2178  ("indian red","#CD5C5C",(205,92,92)),("light coral","#F08080",(240,128,128)),
2179  ("dark salmon","#E9967A",(233,150,122)),("salmon","#FA8072",(250,128,114)),
2180  ("light salmon","#FFA07A",(255,160,122)),("orange red","#FF4500",(255,69,0)),
2181  ("dark orange","#FF8C00",(255,140,0))],
2182  "yellows":[("orange","#FFA500",(255,165,0)),("gold","#FFD700",(255,215,0)),
2183  ("dark golden rod","#B8860B",(184,134,11)),("golden rod","#DAA520",(218,165,32)),
2184  ("pale golden rod","#EEE8AA",(238,232,170)),("dark khaki","#BDB76B",(189,183,107)),
2185  ("khaki","#F0E68C",(240,230,140)),("olive","#808000",(128,128,0)),
2186  ("yellow","#FFFF00",(255,255,0)),("antique white","#FAEBD7",(250,235,215)),
2187  ("beige","#F5F5DC",(245,245,220)),("bisque","#FFE4C4",(255,228,196)),
2188  ("blanched almond","#FFEBCD",(255,235,205)),("wheat","#F5DEB3",(245,222,179)),
2189  ("corn silk","#FFF8DC",(255,248,220)),("lemon chiffon","#FFFACD",(255,250,205)),
2190  ("light golden rod yellow","#FAFAD2",(250,250,210)),("light yellow","#FFFFE0",(255,255,224))],
2191  "greens":[("yellow green","#9ACD32",(154,205,50)),("dark olive green","#556B2F",(85,107,47)),
2192  ("olive drab","#6B8E23",(107,142,35)),("lawn green","#7CFC00",(124,252,0)),
2193  ("chart reuse","#7FFF00",(127,255,0)),("green yellow","#ADFF2F",(173,255,47)),
2194  ("dark green","#006400",(0,100,0)),("green","#008000",(0,128,0)),
2195  ("forest green","#228B22",(34,139,34)),("lime","#00FF00",(0,255,0)),
2196  ("lime green","#32CD32",(50,205,50)),("light green","#90EE90",(144,238,144)),
2197  ("pale green","#98FB98",(152,251,152)),("dark sea green","#8FBC8F",(143,188,143)),
2198  ("medium spring green","#00FA9A",(0,250,154)),("spring green","#00FF7F",(0,255,127)),
2199  ("sea green","#2E8B57",(46,139,87)),("medium aqua marine","#66CDAA",(102,205,170)),
2200  ("medium sea green","#3CB371",(60,179,113)),("light sea green","#20B2AA",(32,178,170)),
2201  ("dark slate gray","#2F4F4F",(47,79,79)),("teal","#008080",(0,128,128)),
2202  ("dark cyan","#008B8B",(0,139,139))],
2203  "blues":[("dark turquoise","#00CED1",(0,206,209)),
2204  ("turquoise","#40E0D0",(64,224,208)),("medium turquoise","#48D1CC",(72,209,204)),
2205  ("pale turquoise","#AFEEEE",(175,238,238)),("aqua marine","#7FFFD4",(127,255,212)),
2206  ("powder blue","#B0E0E6",(176,224,230)),("cadet blue","#5F9EA0",(95,158,160)),
2207  ("steel blue","#4682B4",(70,130,180)),("corn flower blue","#6495ED",(100,149,237)),
2208  ("deep sky blue","#00BFFF",(0,191,255)),("dodger blue","#1E90FF",(30,144,255)),
2209  ("light blue","#ADD8E6",(173,216,230)),("sky blue","#87CEEB",(135,206,235)),
2210  ("light sky blue","#87CEFA",(135,206,250)),("midnight blue","#191970",(25,25,112)),
2211  ("navy","#000080",(0,0,128)),("dark blue","#00008B",(0,0,139)),
2212  ("medium blue","#0000CD",(0,0,205)),("blue","#0000FF",(0,0,255)),("royal blue","#4169E1",(65,105,225)),
2213  ("aqua","#00FFFF",(0,255,255)),("cyan","#00FFFF",(0,255,255)),("light cyan","#E0FFFF",(224,255,255))],
2214  "violets":[("blue violet","#8A2BE2",(138,43,226)),("indigo","#4B0082",(75,0,130)),
2215  ("dark slate blue","#483D8B",(72,61,139)),("slate blue","#6A5ACD",(106,90,205)),
2216  ("medium slate blue","#7B68EE",(123,104,238)),("medium purple","#9370DB",(147,112,219)),
2217  ("dark magenta","#8B008B",(139,0,139)),("dark violet","#9400D3",(148,0,211)),
2218  ("dark orchid","#9932CC",(153,50,204)),("medium orchid","#BA55D3",(186,85,211)),
2219  ("purple","#800080",(128,0,128)),("thistle","#D8BFD8",(216,191,216)),
2220  ("plum","#DDA0DD",(221,160,221)),("violet","#EE82EE",(238,130,238)),
2221  ("magenta / fuchsia","#FF00FF",(255,0,255)),("orchid","#DA70D6",(218,112,214)),
2222  ("medium violet red","#C71585",(199,21,133)),("pale violet red","#DB7093",(219,112,147)),
2223  ("deep pink","#FF1493",(255,20,147)),("hot pink","#FF69B4",(255,105,180)),
2224  ("light pink","#FFB6C1",(255,182,193)),("pink","#FFC0CB",(255,192,203))],
2225  "browns":[("saddle brown","#8B4513",(139,69,19)),("sienna","#A0522D",(160,82,45)),
2226  ("chocolate","#D2691E",(210,105,30)),("peru","#CD853F",(205,133,63)),
2227  ("sandy brown","#F4A460",(244,164,96)),("burly wood","#DEB887",(222,184,135)),
2228  ("tan","#D2B48C",(210,180,140)),("rosy brown","#BC8F8F",(188,143,143)),
2229  ("moccasin","#FFE4B5",(255,228,181)),("navajo white","#FFDEAD",(255,222,173)),
2230  ("peach puff","#FFDAB9",(255,218,185)),("misty rose","#FFE4E1",(255,228,225)),
2231  ("lavender blush","#FFF0F5",(255,240,245)),("linen","#FAF0E6",(250,240,230)),
2232  ("old lace","#FDF5E6",(253,245,230)),("papaya whip","#FFEFD5",(255,239,213)),
2233  ("sea shell","#FFF5EE",(255,245,238))],
2234  "greys":[("black","#000000",(0,0,0)),("dim gray / dim grey","#696969",(105,105,105)),
2235  ("gray / grey","#808080",(128,128,128)),("dark gray / dark grey","#A9A9A9",(169,169,169)),
2236  ("silver","#C0C0C0",(192,192,192)),("light gray / light grey","#D3D3D3",(211,211,211)),
2237  ("gainsboro","#DCDCDC",(220,220,220)),("white smoke","#F5F5F5",(245,245,245)),
2238  ("white","#FFFFFF",(255,255,255))]}
2239 
2240  def assign_color_group(self,color_group,representation,component_names):
2241  for n,p in enumerate(component_names):
2242  s=IMP.atom.Selection(representation.prot,molecule=p)
2243  psel=s.get_selected_particles()
2244  ctuple=self.colors[color_group][n]
2245  print("Assigning "+p+" to color "+ctuple[0])
2246  c=ctuple[2]
2247  color=IMP.display.Color(float(c[0])/255,float(c[1])/255,float(c[2])/255)
2248  for part in psel:
2250  IMP.display.Colored(part).set_color(color)
2251  else:
2253 
2254  def get_list_distant_colors(self):
2255  cnames = ['#F0F8FF', '#FAEBD7', '#00FFFF', '#7FFFD4', '#F0FFFF', '#F5F5DC',
2256  '#FFE4C4', '#000000', '#FFEBCD', '#0000FF', '#8A2BE2', '#A52A2A', '#DEB887',
2257  '#5F9EA0', '#7FFF00', '#D2691E', '#FF7F50', '#6495ED', '#FFF8DC', '#DC143C',
2258  '#00FFFF', '#00008B', '#008B8B', '#B8860B', '#A9A9A9', '#006400', '#BDB76B',
2259  '#8B008B', '#556B2F', '#FF8C00', '#9932CC', '#8B0000', '#E9967A', '#8FBC8F',
2260  '#483D8B', '#2F4F4F', '#00CED1', '#9400D3', '#FF1493', '#00BFFF', '#696969',
2261  '#1E90FF', '#B22222', '#FFFAF0', '#228B22', '#FF00FF', '#DCDCDC', '#F8F8FF',
2262  '#FFD700', '#DAA520', '#808080', '#008000', '#ADFF2F', '#F0FFF0', '#FF69B4',
2263  '#CD5C5C', '#4B0082', '#FFFFF0', '#F0E68C', '#E6E6FA', '#FFF0F5', '#7CFC00',
2264  '#FFFACD', '#ADD8E6', '#F08080', '#E0FFFF', '#FAFAD2', '#90EE90', '#D3D3D3',
2265  '#FFB6C1', '#FFA07A', '#20B2AA', '#87CEFA', '#778899', '#B0C4DE', '#FFFFE0',
2266  '#00FF00', '#32CD32', '#FAF0E6', '#FF00FF', '#800000', '#66CDAA', '#0000CD',
2267  '#BA55D3', '#9370DB', '#3CB371', '#7B68EE', '#00FA9A', '#48D1CC', '#C71585',
2268  '#191970', '#F5FFFA', '#FFE4E1', '#FFE4B5', '#FFDEAD', '#000080', '#FDF5E6',
2269  '#808000', '#6B8E23', '#FFA500', '#FF4500', '#DA70D6', '#EEE8AA', '#98FB98',
2270  '#AFEEEE', '#DB7093', '#FFEFD5', '#FFDAB9', '#CD853F', '#FFC0CB', '#DDA0DD',
2271  '#B0E0E6', '#800080', '#FF0000', '#BC8F8F', '#4169E1', '#8B4513', '#FA8072',
2272  '#FAA460', '#2E8B57', '#FFF5EE', '#A0522D', '#C0C0C0', '#87CEEB', '#6A5ACD',
2273  '#708090', '#FFFAFA', '#00FF7F', '#4682B4', '#D2B48C', '#008080', '#D8BFD8',
2274  '#FF6347', '#40E0D0', '#EE82EE', '#F5DEB3', '#FFFFFF', '#F5F5F5', '#FFFF00',
2275  '#9ACD32']
2276  return cnames
def list_chunks_iterator
Yield successive length-sized chunks from a list.
Definition: tools.py:1183
Simple 3D transformation class.
static bool get_is_setup(const IMP::ParticleAdaptor &p)
Definition: Residue.h:155
def select_at_all_resolutions
Perform selection using the usual keywords but return ALL resolutions (BEADS and GAUSSIANS).
Definition: tools.py:1794
Represent an RGB color.
Definition: Color.h:24
A decorator to associate a particle with a part of a protein/DNA/RNA.
Definition: Fragment.h:20
def display_bonds
Decorate the sequence-consecutive particles from a PMI2 molecule with a bond, so that they appear con...
Definition: tools.py:1750
def get_restraint_set
Get a RestraintSet containing all PMI restraints added to the model.
Definition: tools.py:49
Store the representations for a system.
Definition: tools.py:900
def remove
index can be a integer
Definition: tools.py:1261
def get_random_cross_link_dataset
Return a random cross-link dataset as a string.
Definition: tools.py:368
def shuffle_configuration
Shuffle particles.
Definition: tools.py:1913
static bool get_is_setup(const IMP::ParticleAdaptor &p)
Definition: atom/Atom.h:241
void add_particles(RMF::FileHandle fh, const ParticlesTemp &hs)
def get_particles_within_zone
Utility to retrieve particles from a hierarchy within a zone around a set of ps.
Definition: tools.py:1827
Set of python classes to create a multi-state, multi-resolution IMP hierarchy.
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:575
A decorator for a particle which has bonds.
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
Rotation3D get_random_rotation_3d(const Rotation3D &center, double distance)
Pick a rotation at random near the provided one.
def scatter_and_gather_dict_append
Synchronize data over a parallel run.
Definition: tools.py:1137
def deprecated_function
Python decorator to mark a function as deprecated.
Definition: __init__.py:9647
Change color code to hexadecimal to rgb.
Definition: tools.py:1483
IMP::Vector< Color > Colors
Definition: Color.h:188
ParticlesTemp get_particles(Model *m, const ParticleIndexes &ps)
This class stores integers in ordered compact lists eg: [[1,2,3],[6,7,8]] the methods help splitting ...
Definition: tools.py:1203
The type of an atom.
static Surface setup_particle(Model *m, ParticleIndex pi)
Definition: Surface.h:43
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:619
static bool get_is_setup(const IMP::ParticleAdaptor &p)
Definition: XYZ.h:49
GenericHierarchies get_leaves(Hierarchy mhd)
Get all the leaves of the bit of hierarchy.
def get_prot_name_from_particle
Return the component name provided a particle and a list of names.
Definition: tools.py:1048
This class initializes the root node of the global IMP.atom.Hierarchy.
double get_distance(XYZR a, XYZR b)
Compute the sphere distance between a and b.
Definition: XYZR.h:89
def get_residue_gaps_in_hierarchy
Return the residue index gaps and contiguous segments in the hierarchy.
Definition: tools.py:660
def add
index can be a integer or a list of integers
Definition: tools.py:1229
Vector3D get_random_vector_in(const Cylinder3D &c)
Generate a random vector in a cylinder with uniform density.
def get_all_leaves
Just get the leaves from a list of hierarchies.
Definition: tools.py:1788
Add resolution to a particle.
Definition: Resolution.h:24
Bond create_bond(Bonded a, Bonded b, Bond o)
Connect the two wrapped particles by a custom bond.
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)
Return a dense grid containing the voxels of the passed density map.
Stores a named protein chain.
def color2rgb
Given a chimera color name, return RGB.
Definition: tools.py:2107
def input_adaptor
Adapt things for PMI (degrees of freedom, restraints, ...) Returns list of list of hierarchies...
Definition: tools.py:1608
def get_terminal_residue
Get the particle of the terminal residue at the GIVEN resolution (NOTE: not the closest resolution!)...
Definition: tools.py:607
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:1290
ParticleIndexPairs get_indexes(const ParticlePairsTemp &ps)
static bool get_is_setup(Model *m, ParticleIndex pi)
Definition: Colored.h:48
The standard decorator for manipulating molecular structures.
Ints get_index(const ParticlesTemp &particles, const Subset &subset, const Subsets &excluded)
def select_by_tuple_2
New tuple format: molname OR (start,stop,molname,copynum,statenum) Copy and state are optional...
Definition: tools.py:850
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:1110
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:37
static Bonded setup_particle(Model *m, ParticleIndex pi)
The type for a residue.
def get_terminal_residue_position
Get XYZ coordinates of the terminal residue at the GIVEN resolution.
Definition: tools.py:648
def __init__
index can be a integer or a list of integers
Definition: tools.py:1220
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.
Definition: XYZ.h:30
def get_flatten
Returns a flatten list.
Definition: tools.py:1279
static Scale setup_particle(Model *m, ParticleIndex pi)
Definition: Scale.h:30
A base class for Keys.
Definition: Key.h:44
static Colored setup_particle(Model *m, ParticleIndex pi, Color color)
Definition: Colored.h:62
A decorator for a particle that is part of a rigid body but not rigid.
Definition: rigid_bodies.h:638
def get_closest_residue_position
this function works with plain hierarchies, as read from the pdb, no multi-scale hierarchies ...
Definition: tools.py:544
std::ostream & show(Hierarchy h, std::ostream &out=std::cout)
Print the hierarchy using a given decorator to display each node.
Find all nearby pairs by testing all pairs.
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...
static bool get_is_setup(const IMP::ParticleAdaptor &p)
static bool get_is_setup(Model *m, ParticleIndex p)
Check if the particle has the needed attributes for a cast to succeed.
The general base class for IMP exceptions.
Definition: exception.h:49
VectorD< 3 > Vector3D
Definition: VectorD.h:395
Rotation3D get_identity_rotation_3d()
Return a rotation that does not do anything.
Definition: Rotation3D.h:244
def get_densities
Given a list of PMI objects, return all density hierarchies within these objects. ...
Definition: tools.py:1885
void link_hierarchies(RMF::FileConstHandle fh, const atom::Hierarchies &hs)
Class to handle individual particles of a Model object.
Definition: Particle.h:41
Bond get_bond(Bonded a, Bonded b)
Get the bond between two particles.
Stores a list of Molecules all with the same State index.
std::string get_data_path(std::string file_name)
Return the full path to one of this module's data files.
def set_coordinates_from_rmf
Extract frame from RMF file and fill coordinates.
Definition: tools.py:1587
def select
this function uses representation=SimplifiedModel it returns the corresponding selected particles rep...
Definition: tools.py:731
Python classes to represent, score, sample and analyze models.
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:476
static bool get_is_setup(const IMP::ParticleAdaptor &p)
Definition: Gaussian.h:47
def get_sorted_segments
Returns sequence-sorted segments array, each containing the first particle the last particle and the ...
Definition: tools.py:1723
def get_rbs_and_beads
Returns unique objects in original order.
Definition: tools.py:1861
Hierarchies get_leaves(const Selection &h)
Select hierarchy particles identified by the biological name.
Definition: Selection.h:66
Support for the RMF file format for storing hierarchical molecular data and markup.
def get_residue_indexes
Retrieve the residue indexes for the given particle.
Definition: tools.py:1068
static bool get_is_setup(const IMP::ParticleAdaptor &p)
Definition: rigid_bodies.h:640
Transformation3D get_random_local_transformation(Vector3D origin, double max_translation=5., double max_angle_in_rad=0.26)
Get a local transformation.
Temporarily stores residue information, even without structure available.
Store objects in order they were added, but with default type.
Definition: tools.py:1561
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:1164
A particle with a color.
Definition: Colored.h:23