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