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