IMP logo
IMP Reference Guide  2.18.0
The Integrative Modeling Platform
restraints/stereochemistry.py
1 """@namespace IMP.pmi.restraints.stereochemistry
2 Restraints for keeping correct stereochemistry.
3 """
4 
5 from __future__ import print_function
6 import IMP
7 import IMP.core
8 import IMP.atom
9 import IMP.container
10 import IMP.isd
11 import IMP.pmi.tools
12 from operator import itemgetter
13 import math
14 
15 
16 class ConnectivityRestraint(IMP.pmi.restraints.RestraintBase):
17  """Create a restraint between consecutive TempResidue objects
18  or an entire PMI Molecule object."""
19 
20  def __init__(self,
21  objects,
22  scale=1.0,
23  disorderedlength=False,
24  upperharmonic=True,
25  resolution=1,
26  label=None):
27  """
28  @param objects - a list of hierarchies, PMI TempResidues OR a
29  single Molecule
30  @param scale Scale the maximal distance between the beads by this
31  factor when disorderedlength is False. The maximal distance
32  is calculated as ((float(residuegap) + 1.0) * 3.6) * scale.
33  @param disorderedlength - This flag uses either disordered length
34  calculated for random coil peptides (True) or zero
35  surface-to-surface distance between beads (False)
36  as optimal distance for the sequence connectivity restraint.
37  @param upperharmonic - This flag uses either harmonic (False)
38  or upperharmonic (True) in the intra-pair
39  connectivity restraint.
40  @param resolution - The resolution to connect things at - only used
41  if you pass PMI objects
42  @param label - A string to identify this restraint in the
43  output/stat file
44  """
45 
46  hiers = IMP.pmi.tools.input_adaptor(objects, resolution)
47  if len(hiers) > 1:
48  raise Exception("ConnectivityRestraint: only pass stuff from "
49  "one Molecule, please")
50  hiers = hiers[0]
51  m = list(hiers)[0].get_model()
52  super(ConnectivityRestraint, self).__init__(m, label=label)
53 
54  self.kappa = 10 # spring constant used for the harmonic restraints
55  SortedSegments = []
56  for h in hiers:
57  try:
58  start = IMP.atom.Hierarchy(h).get_children()[0]
59  except: # noqa: E722
60  start = IMP.atom.Hierarchy(h)
61 
62  try:
63  end = IMP.atom.Hierarchy(h).get_children()[-1]
64  except: # noqa: E722
65  end = IMP.atom.Hierarchy(h)
66 
67  startres = IMP.pmi.tools.get_residue_indexes(start)[0]
68  SortedSegments.append((start, end, startres))
69  SortedSegments = sorted(SortedSegments, key=itemgetter(2))
70 
71  # connect the particles
72  self.particle_pairs = []
73  for x in range(len(SortedSegments) - 1):
74 
75  last = SortedSegments[x][1]
76  first = SortedSegments[x + 1][0]
77 
78  apply_restraint = True
79 
80  # Apply connectivity runless ALL of the following are true:
81  # - first and last both have RigidBodyMember decorators
82  # - first and last are both RigidMembers
83  # - first and last are part of the same RigidBody object
84 
85  # Check for both in a rigid body
90  # Check if the rigid body objects for each particle are
91  # the same object.
92  # if so, skip connectivity restraint
93  if IMP.core.RigidBodyMember(first).get_rigid_body() \
94  == IMP.core.RigidBodyMember(last).get_rigid_body():
95  apply_restraint = False
96 
97  if apply_restraint:
98 
99  nreslast = len(IMP.pmi.tools.get_residue_indexes(last))
100  lastresn = IMP.pmi.tools.get_residue_indexes(last)[-1]
101  nresfirst = len(IMP.pmi.tools.get_residue_indexes(first))
102  firstresn = IMP.pmi.tools.get_residue_indexes(first)[0]
103 
104  residuegap = firstresn - lastresn - 1
105  if disorderedlength and (nreslast / 2 + nresfirst / 2
106  + residuegap) > 20.0:
107  # calculate the distance between the sphere centers
108  # using Kohn PNAS 2004
109  optdist = math.sqrt(5 / 3) * 1.93 * \
110  (nreslast / 2 + nresfirst / 2 + residuegap) ** 0.6
111  if upperharmonic:
112  hu = IMP.core.HarmonicUpperBound(optdist, self.kappa)
113  else:
114  hu = IMP.core.Harmonic(optdist, self.kappa)
116  else: # default
117  optdist = (0.0 + (float(residuegap) + 1.0) * 3.6) * scale
118  if upperharmonic: # default
119  hu = IMP.core.HarmonicUpperBound(optdist, self.kappa)
120  else:
121  hu = IMP.core.Harmonic(optdist, self.kappa)
123 
124  pt0 = last.get_particle()
125  pt1 = first.get_particle()
126  self.particle_pairs.append((pt0, pt1))
128  self.model, dps, (pt0.get_index(), pt1.get_index()))
129 
130  print("Adding sequence connectivity restraint between",
131  pt0.get_name(), " and ", pt1.get_name(), 'of distance',
132  optdist)
133  self.rs.add_restraint(r)
134 
136  """ Returns number of connectivity restraints """
137  return len(self.rs.get_restraints())
138 
140  """ Returns the list of connected particles pairs """
141  return self.particle_pairs
142 
143 
144 class ExcludedVolumeSphere(IMP.pmi.restraints.RestraintBase):
145  """A class to create an excluded volume restraint for a set of particles
146  at a given resolution.
147  Can be initialized as a bipartite restraint between two sets of particles.
148  # Potential additional function: Variable resolution for each PMI object.
149  Perhaps passing selection_tuples with (PMI_object, resolution)
150  """
151  _include_in_rmf = True
152 
153  def __init__(self,
154  included_objects,
155  other_objects=None,
156  resolution=1000,
157  kappa=1.0):
158  """Constructor.
159  @param included_objects Can be one of the following inputs:
160  IMP Hierarchy, PMI System/State/Molecule/TempResidue,
161  or a list/set of them
162  @param other_objects Initializes a bipartite restraint between
163  included_objects and other_objects
164  Same format as included_objects
165  @param resolution The resolution particles at which to impose the
166  restraint. By default, the coarsest particles will be chosen.
167  If a number is chosen, for each particle, the closest
168  resolution will be used (see IMP.atom.Selection).
169  @param kappa Restraint strength
170  """
171 
172  self.kappa = kappa
173  self.cpc = None
174  bipartite = False
175 
176  # gather IMP hierarchies from input objects
177  hierarchies = IMP.pmi.tools.input_adaptor(included_objects,
178  resolution,
179  flatten=True)
180  included_ps = []
181  if other_objects is not None:
182  bipartite = True
183  other_hierarchies = IMP.pmi.tools.input_adaptor(other_objects,
184  resolution,
185  flatten=True)
186  other_ps = []
187 
188  # perform selection
189  if hierarchies is None:
190  raise Exception("Must at least pass included objects")
191  mdl = hierarchies[0].get_model()
192  super(ExcludedVolumeSphere, self).__init__(mdl)
193 
194  included_ps = [h.get_particle() for h in hierarchies]
195  if bipartite:
196  other_ps = [h.get_particle() for h in other_hierarchies]
197 
198  # setup score
199  ssps = IMP.core.SoftSpherePairScore(self.kappa)
200  lsa = IMP.container.ListSingletonContainer(self.model)
201  lsa.add(IMP.get_indexes(included_ps))
202 
203  # setup close pair container
204  if not bipartite:
206  self.cpc = IMP.container.ClosePairContainer(lsa, 0.0, rbcpf, 10.0)
207  evr = IMP.container.PairsRestraint(ssps, self.cpc)
208  else:
209  other_lsa = IMP.container.ListSingletonContainer(self.model)
210  other_lsa.add(IMP.get_indexes(other_ps))
212  lsa,
213  other_lsa,
214  0.0,
215  10.0)
216  evr = IMP.container.PairsRestraint(ssps, self.cpc)
217 
218  self.rs.add_restraint(evr)
219 
220  def add_excluded_particle_pairs(self, excluded_particle_pairs):
221  # add pairs to be filtered when calculating the score
222  inverted = [(p1, p0) for p0, p1 in excluded_particle_pairs]
223  lpc = IMP.container.ListPairContainer(self.model)
224  lpc.add(IMP.get_indexes(excluded_particle_pairs))
225  lpc.add(IMP.get_indexes(inverted))
227  self.cpc.add_pair_filter(icpf)
228 
229 
230 class HelixRestraint(IMP.pmi.restraints.RestraintBase):
231  """Enforce ideal Helix dihedrals and bonds for a selection
232  at resolution 0"""
233  def __init__(self,
234  hierarchy,
235  selection_tuple,
236  weight=1.0,
237  label=None):
238  """Constructor
239  @param hierarchy the root node
240  @param selection_tuple (start, stop, molname, copynum=0)
241  @param weight
242  """
243  m = hierarchy.get_model()
244  super(HelixRestraint, self).__init__(m, weight=weight)
245  start = selection_tuple[0]
246  stop = selection_tuple[1]
247  mol = selection_tuple[2]
248  copy_index = 0
249  if len(selection_tuple) > 3:
250  copy_index = selection_tuple[3]
251 
252  sel = IMP.atom.Selection(
253  hierarchy, molecule=mol, copy_index=copy_index,
254  residue_indexes=range(start, stop+1))
255  ps = sel.get_selected_particles(with_representation=False)
256  res = [IMP.atom.Residue(p) for p in ps]
257  self.r = IMP.atom.HelixRestraint(res)
258  self.rs.add_restraint(self.r)
259  print('Created helix %s.%i.%i-%i with %i dihedrals and %i bonds'
260  % (mol, copy_index, start, stop, self.get_number_of_bonds(),
261  self.get_number_of_dihedrals()))
262 
263  def get_number_of_bonds(self):
264  return self.r.get_number_of_bonds()
265 
266  def get_number_of_dihedrals(self):
267  return self.r.get_number_of_dihedrals()
268 
269 
270 class ResidueBondRestraint(IMP.pmi.restraints.RestraintBase):
271  """ Add bond restraint between pair of consecutive
272  residues/beads to enforce the stereochemistry.
273  """
274  def __init__(self, objects, distance=3.78, strength=10.0, jitter=None):
275  """Constructor
276  @param objects Objects to restrain
277  @param distance Resting distance for restraint
278  @param strength Bond constant
279  @param jitter Defines the +- added to the optimal distance in
280  the harmonic well restraint used to increase the tolerance
281  """
282 
283  particles = IMP.pmi.tools.input_adaptor(objects, 1, flatten=True)
284  m = particles[0].get_model()
285  super(ResidueBondRestraint, self).__init__(m)
286 
287  self.pairslist = []
288 
289  if not jitter:
290  ts = IMP.core.Harmonic(distance, strength)
291  else:
293  (distance - jitter, distance + jitter), strength)
294 
295  for ps in IMP.pmi.tools.sublist_iterator(particles, 2, 2):
296  pair = []
297  if len(ps) != 2:
298  raise ValueError("wrong length of pair")
299  for p in ps:
301  raise TypeError("%s is not a residue" % p)
302  else:
303  pair.append(p)
304  print("ResidueBondRestraint: adding a restraint between %s %s"
305  % (pair[0].get_name(), pair[1].get_name()))
306  self.rs.add_restraint(
307  IMP.core.DistanceRestraint(self.model, ts, pair[0], pair[1]))
308  self.pairslist.append(IMP.ParticlePair(pair[0], pair[1]))
309  self.pairslist.append(IMP.ParticlePair(pair[1], pair[0]))
310 
311  def get_excluded_pairs(self):
312  return self.pairslist
313 
314 
315 class ResidueAngleRestraint(IMP.pmi.restraints.RestraintBase):
316  """Add angular restraint between triplets of consecutive
317  residues/beads to enforce the stereochemistry.
318  """
319  def __init__(self, objects, anglemin=100.0, anglemax=140.0, strength=10.0):
320 
321  particles = IMP.pmi.tools.input_adaptor(objects, 1, flatten=True)
322  m = particles[0].get_model()
323  super(ResidueAngleRestraint, self).__init__(m)
324 
325  self.pairslist = []
326 
328  (math.pi * anglemin / 180.0,
329  math.pi * anglemax / 180.0),
330  strength)
331 
332  for ps in IMP.pmi.tools.sublist_iterator(particles, 3, 3):
333  triplet = []
334  if len(ps) != 3:
335  raise ValueError("wrong length of triplet")
336  for p in ps:
338  raise TypeError("%s is not a residue" % p)
339  else:
340  triplet.append(p)
341  print("ResidueAngleRestraint: adding a restraint between %s %s %s"
342  % (triplet[0].get_name(), triplet[1].get_name(),
343  triplet[2].get_name()))
344  self.rs.add_restraint(
345  IMP.core.AngleRestraint(triplet[0].get_model(), ts,
346  triplet[0],
347  triplet[1],
348  triplet[2]))
349  self.pairslist.append(IMP.ParticlePair(triplet[0], triplet[2]))
350  self.pairslist.append(IMP.ParticlePair(triplet[2], triplet[0]))
351 
352  def get_excluded_pairs(self):
353  return self.pairslist
354 
355 
356 class ResidueDihedralRestraint(IMP.pmi.restraints.RestraintBase):
357  """Add dihedral restraints between quadruplet of consecutive
358  residues/beads to enforce the stereochemistry.
359  Give as input a string of "C" and "T", meaning cys (0+-40)
360  or trans (180+-40) dihedral. The length of the string must be #residue-3.
361  Without the string, the dihedral will be assumed trans.
362  """
363  def __init__(self, objects, stringsequence=None, strength=10.0):
364 
365  particles = IMP.pmi.tools.input_adaptor(objects, 1, flatten=True)
366  m = particles[0].get_model()
367  super(ResidueDihedralRestraint, self).__init__(m)
368 
369  self.pairslist = []
370 
371  if stringsequence is None:
372  stringsequence = "T" * (len(particles) - 3)
373 
374  for n, ps in enumerate(
375  IMP.pmi.tools.sublist_iterator(particles, 4, 4)):
376  quadruplet = []
377  if len(ps) != 4:
378  raise ValueError("wrong length of quadruplet")
379  for p in ps:
381  raise TypeError("%s is not a residue" % p)
382  else:
383  quadruplet.append(p)
384  dihedraltype = stringsequence[n]
385  if dihedraltype == "C":
386  anglemin = -20.0
387  anglemax = 20.0
389  (math.pi * anglemin / 180.0,
390  math.pi * anglemax / 180.0),
391  strength)
392  print("ResidueDihedralRestraint: adding a CYS restraint "
393  "between %s %s %s %s"
394  % (quadruplet[0].get_name(), quadruplet[1].get_name(),
395  quadruplet[2].get_name(), quadruplet[3].get_name()))
396  if dihedraltype == "T":
397  anglemin = 180 - 70.0
398  anglemax = 180 + 70.0
400  (math.pi * anglemin / 180.0,
401  math.pi * anglemax / 180.0),
402  strength)
403  print("ResidueDihedralRestraint: adding a TRANS restraint "
404  "between %s %s %s %s"
405  % (quadruplet[0].get_name(),
406  quadruplet[1].get_name(), quadruplet[2].get_name(),
407  quadruplet[3].get_name()))
408  self.rs.add_restraint(
409  IMP.core.DihedralRestraint(self.model, ts,
410  quadruplet[0],
411  quadruplet[1],
412  quadruplet[2],
413  quadruplet[3]))
414  self.pairslist.append(
415  IMP.ParticlePair(quadruplet[0], quadruplet[3]))
416  self.pairslist.append(
417  IMP.ParticlePair(quadruplet[3], quadruplet[0]))
418 
419 
421  """Add harmonic restraints between all pairs
422  """
423  def __init__(self, hierarchy, selection_tuples=None, resolution=1,
424  strength=0.01, dist_cutoff=10.0, ca_only=True):
425  """Constructor
426  @param hierarchy Root hierarchy to select from
427  @param selection_tuples Selecting regions for the restraint
428  [[start,stop,molname,copy_index=0],...]
429  @param resolution Resolution for applying restraint
430  @param strength Bond strength
431  @param dist_cutoff Cutoff for making restraints
432  @param ca_only Selects only CAlphas. Only matters if resolution=0.
433  """
434 
435  particles = []
436  self.m = hierarchy.get_model()
437  for st in selection_tuples:
438  copy_index = 0
439  if len(st) > 3:
440  copy_index = st[3]
441  if not ca_only:
442  sel = IMP.atom.Selection(
443  hierarchy, molecule=st[2],
444  residue_indexes=range(st[0], st[1]+1),
445  copy_index=copy_index)
446  else:
447  sel = IMP.atom.Selection(
448  hierarchy, molecule=st[2],
449  residue_indexes=range(st[0], st[1]+1),
450  copy_index=copy_index,
451  atom_type=IMP.atom.AtomType("CA"))
452  particles += sel.get_selected_particles()
453 
454  self.weight = 1
455  self.label = "None"
456  self.pairslist = []
457 
458  # create score
460  particles, dist_cutoff, strength)
461  for r in self.rs.get_restraints():
462  a1, a2 = r.get_inputs()
463  self.pairslist.append(IMP.ParticlePair(a1, a2))
464  self.pairslist.append(IMP.ParticlePair(a2, a1))
465  print('ElasticNetwork: created', self.rs.get_number_of_restraints(),
466  'restraints')
467 
468  def set_label(self, label):
469  self.label = label
470  self.rs.set_name(label)
471  for r in self.rs.get_restraints():
472  r.set_name(label)
473 
474  def add_to_model(self):
475  IMP.pmi.tools.add_restraint_to_model(self.m, self.rs)
476 
477  def get_restraint(self):
478  return self.rs
479 
480  def set_weight(self, weight):
481  self.weight = weight
482  self.rs.set_weight(weight)
483 
484  def get_excluded_pairs(self):
485  return self.pairslist
486 
487  def get_output(self):
488  output = {}
489  score = self.weight * self.rs.unprotected_evaluate(None)
490  output["_TotalScore"] = str(score)
491  output["ElasticNetworkRestraint_" + self.label] = str(score)
492  return output
493 
494 
496  """ Enable CHARMM force field """
497  def __init__(self, root, ff_temp=300.0, zone_ps=None, zone_size=10.0,
498  enable_nonbonded=True, enable_bonded=True,
499  zone_nonbonded=False):
500  """Setup the CHARMM restraint on a selection. Expecting atoms.
501  @param root The node at which to apply the restraint
502  @param ff_temp The temperature of the force field
503  @param zone_ps Create a zone around this set of particles
504  Automatically includes the entire residue (incl. backbone)
505  @param zone_size The size for looking for neighbor residues
506  @param enable_nonbonded Allow the repulsive restraint
507  @param enable_bonded Allow the bonded restraint
508  @param zone_nonbonded EXPERIMENTAL: exclude from nonbonded all
509  sidechains that aren't in zone!
510  """
511 
512  kB = (1.381 * 6.02214) / 4184.0
513 
514  self.mdl = root.get_model()
515  self.bonds_rs = IMP.RestraintSet(self.mdl, 1.0 / (kB * ff_temp),
516  'BONDED')
517  self.nonbonded_rs = IMP.RestraintSet(self.mdl, 1.0 / (kB * ff_temp),
518  'NONBONDED')
519  self.weight = 1.0
520  self.label = ""
521 
522  # setup topology and bonds etc
523  if enable_bonded:
525  topology = ff.create_topology(root)
526  topology.apply_default_patches()
527  topology.setup_hierarchy(root)
528  if zone_ps is not None:
529  limit_to_ps = IMP.pmi.topology.get_particles_within_zone(
530  root, zone_ps, zone_size, entire_residues=True,
531  exclude_backbone=False)
533  topology,
534  limit_to_ps)
535  self.ps = limit_to_ps
536  else:
537  r = IMP.atom.CHARMMStereochemistryRestraint(root, topology)
538  self.ps = IMP.core.get_leaves(root)
539  print('init bonds score', r.unprotected_evaluate(None))
540  self.bonds_rs.add_restraint(r)
541  ff.add_radii(root)
542  ff.add_well_depths(root)
543 
544  atoms = IMP.atom.get_by_type(root, IMP.atom.ATOM_TYPE)
545  # non-bonded forces
546  if enable_nonbonded:
547  if (zone_ps is not None) and zone_nonbonded:
548  print('stereochemistry: zone_nonbonded is True')
549  # atoms list should only include backbone plus zone_ps!
550  backbone_types = ['C', 'N', 'CB', 'O']
551  sel = IMP.atom.Selection(
552  root, atom_types=[IMP.atom.AtomType(n)
553  for n in backbone_types])
554  backbone_atoms = sel.get_selected_particles()
555  sel_ps = IMP.pmi.topology.get_particles_within_zone(
556  root, zone_ps, zone_size, entire_residues=True,
557  exclude_backbone=True)
558 
560  IMP.container.ListSingletonContainer(backbone_atoms),
562  4.0)
563  else:
564  cont = IMP.container.ListSingletonContainer(self.mdl, atoms)
565  self.nbl = IMP.container.ClosePairContainer(cont, 4.0)
566  if enable_bonded:
567  self.nbl.add_pair_filter(r.get_full_pair_filter())
568  pairscore = IMP.isd.RepulsiveDistancePairScore(0, 1)
569  pr = IMP.container.PairsRestraint(pairscore, self.nbl)
570  self.nonbonded_rs.add_restraint(pr)
571  print('CHARMM is set up')
572 
573  def set_label(self, label):
574  self.label = label
575  self.rs.set_name(label)
576  for r in self.rs.get_restraints():
577  r.set_name(label)
578 
579  def add_to_model(self):
580  IMP.pmi.tools.add_restraint_to_model(self.mdl, self.bonds_rs)
581  IMP.pmi.tools.add_restraint_to_model(self.mdl, self.nonbonded_rs)
582 
583  def get_restraint(self):
584  return self.rs
585 
586  def get_close_pair_container(self):
587  return self.nbl
588 
589  def set_weight(self, weight):
590  self.weight = weight
591  self.rs.set_weight(weight)
592 
593  def get_output(self):
594  output = {}
595  bonds_score = self.weight * self.bonds_rs.unprotected_evaluate(None)
596  nonbonded_score = \
597  self.weight * self.nonbonded_rs.unprotected_evaluate(None)
598  score = bonds_score+nonbonded_score
599  output["_TotalScore"] = str(score)
600  output["CHARMM_BONDS"] = str(bonds_score)
601  output["CHARMM_NONBONDED"] = str(nonbonded_score)
602  return output
603 
604 
605 class PseudoAtomicRestraint(object):
606  """Add bonds and improper dihedral restraints for the CBs
607  """
608  def __init__(
609  self, rnums, representation, selection_tuple, strength=10.0, kappa=1.0,
610  jitter_angle=0.0, jitter_improper=0.0):
611  '''
612  need to add:
613  ca-ca bond
614  ca-cb is a constraint, no restraint needed
615  ca-ca-ca
616  cb-ca-ca-cb
617  '''
618 
619  self.m = representation.prot.get_model()
620  self.rs = IMP.RestraintSet(self.m, "PseudoAtomic")
621  self.rset_angles = IMP.RestraintSet(self.m, "PseudoAtomic_Angles")
622  self.rset_bonds = IMP.RestraintSet(self.m, "PseudoAtomic_Bonds")
623  self.weight = 1
624  self.label = "None"
625  self.pairslist = []
626 
627  for rnum in rnums:
628  ca, cb = self.get_ca_cb(
629  IMP.pmi.tools.select_by_tuple(representation,
630  (rnum, rnum, 'chainA'),
631  resolution=0))
632  if cb is not None:
633  nter = False
634  cter = False
635  ca_prev, cb_prev = self.get_ca_cb(
636  IMP.pmi.tools.select_by_tuple(representation,
637  (rnum - 1, rnum - 1,
638  'chainA'), resolution=0))
639  ca_next, cb_next = self.get_ca_cb(
640  IMP.pmi.tools.select_by_tuple(representation,
641  (rnum + 1, rnum + 1,
642  'chainA'), resolution=0))
643  if ca_prev is None:
644  nter = True
645  else:
646  if ca_next is None:
647  cter = True
648  else:
649  if (nter and cter):
650  continue
651 
652  # adding a bond restraint between CA and CB
653  # h=IMP.core.Harmonic(6.0,kappa)
654  # dps=IMP.core.DistancePairScore(h)
655  # pr=IMP.core.PairRestraint(dps,IMP.ParticlePair(ca,cb))
656  # self.pairslist.append((ca,cb))
657  # self.rset_bonds.add_restraint(pr)
658 
659  # creating improper dihedral restraint
660  # hus=IMP.core.Harmonic(2.09,kappa)
662  2.09 +
663  jitter_angle /
664  0.5,
665  kappa)
667  2.09 -
668  jitter_angle /
669  0.5,
670  kappa)
671  if not nter:
672  # ar13=IMP.core.AngleRestraint(hus,ca_prev,ca,cb)
673  # self.rset_angles.add_restraint(ar13)
674  ar13u = IMP.core.AngleRestraint(hupp, ca_prev, ca, cb)
675  ar13l = IMP.core.AngleRestraint(hlow, ca_prev, ca, cb)
676  self.rset_angles.add_restraint(ar13u)
677  self.rset_angles.add_restraint(ar13l)
678  if not cter:
679  # ar23=IMP.core.AngleRestraint(hus,ca_next,ca,cb)
680  # self.rset_angles.add_restraint(ar23)
681  ar23u = IMP.core.AngleRestraint(hupp, ca_next, ca, cb)
682  ar23l = IMP.core.AngleRestraint(hlow, ca_next, ca, cb)
683  self.rset_angles.add_restraint(ar23u)
684  self.rset_angles.add_restraint(ar23l)
685  if not nter and not cter:
686  # hus2=IMP.core.Harmonic(0,kappa)
687  # idr=IMP.core.DihedralRestraint(hus2,ca,ca_prev,ca_next,cb)
688  # self.rset_angles.add_restraint(idr)
689 
690  hus2upp = IMP.core.HarmonicUpperBound(
691  jitter_improper,
692  kappa)
693  hus2low = IMP.core.HarmonicLowerBound(
694  -
695  jitter_improper,
696  kappa)
698  hus2upp,
699  ca,
700  ca_prev,
701  ca_next,
702  cb)
704  hus2low,
705  ca,
706  ca_prev,
707  ca_next,
708  cb)
709  self.rset_angles.add_restraint(idru)
710  self.rset_angles.add_restraint(idrl)
711  self.rs.add_restraint(self.rset_bonds)
712  self.rs.add_restraint(self.rset_angles)
713 
714  def get_ca_cb(self, atoms):
715  ca = None
716  cb = None
717  for a in atoms:
718  if IMP.atom.Atom(a).get_atom_type() == IMP.atom.AtomType("CA"):
719  ca = a.get_particle()
720  elif IMP.atom.Atom(a).get_atom_type() == IMP.atom.AtomType("CB"):
721  cb = a.get_particle()
722  return ca, cb
723 
724  def set_label(self, label):
725  self.label = label
726  self.rs.set_name(label)
727  for r in self.rs.get_restraints():
728  r.set_name(label)
729 
730  def add_to_model(self):
731  IMP.pmi.tools.add_restraint_to_model(self.m, self.rs)
732 
733  def get_restraint(self):
734  return self.rs
735 
736  def set_weight(self, weight):
737  self.weight = weight
738  self.rs.set_weight(weight)
739 
740  def get_excluded_pairs(self):
741  return self.pairslist
742 
743  def get_output(self):
744  output = {}
745  score = self.weight * self.rs.unprotected_evaluate(None)
746  output["_TotalScore"] = str(score)
747  output["PseudoAtomicRestraint_" + self.label] = str(score)
748  return output
749 
750 
751 class SymmetryRestraint(object):
752  """Create harmonic restraints between the reference and (transformed)
753  clones.
754 
755  @note Wraps IMP::core::TransformedDistancePairScore with an
756  IMP::core::Harmonic
757  """
758  def __init__(self, references, clones_list, transforms,
759  label='', strength=10.0, ca_only=False):
760  """Constructor
761  @param references Can be one of the following inputs:
762  IMP Hierarchy, PMI System/State/Molecule/TempResidue,
763  or a list/set of them
764  @param clones_list List of lists of the above
765  @param transforms Transforms moving each selection to the first
766  selection
767  @param label Label for output
768  @param strength The elastic bond strength
769  @param ca_only Optionally select so only CAlpha particles are used
770  """
771 
772  refs = IMP.pmi.tools.input_adaptor(references, flatten=True)
773  self.mdl = refs[0].get_model()
774  self.rs = IMP.RestraintSet(self.mdl, "Symmetry")
775  self.weight = 1
776  self.label = label
777  if len(clones_list) != len(transforms):
778  raise Exception(
779  'Error: There should be as many clones as transforms')
780 
781  harmonic = IMP.core.Harmonic(0., strength)
782  for tmp_clones, trans in zip(clones_list, transforms):
783  clones = IMP.pmi.tools.input_adaptor(tmp_clones, flatten=True)
784  if len(clones) != len(refs):
785  raise Exception("Error: len(references)!=len(clones)")
786  pair_score = IMP.core.TransformedDistancePairScore(harmonic, trans)
787  for p0, p1 in zip(refs, clones):
788  if not ca_only or (
789  IMP.atom.Atom(p0).get_atom_type()
790  == IMP.atom.AtomType("CA") and
791  IMP.atom.Atom(p1).get_atom_type()
792  == IMP.atom.AtomType("CA")):
794  self.mdl, pair_score, [p0.get_particle_index(),
795  p1.get_particle_index()])
796  self.rs.add_restraint(r)
797  print('created symmetry network with',
798  self.rs.get_number_of_restraints(), 'restraints')
799 
800  def set_label(self, label):
801  self.label = label
802  self.rs.set_name(label)
803  for r in self.rs.get_restraints():
804  r.set_name(label)
805 
806  def add_to_model(self):
807  IMP.pmi.tools.add_restraint_to_model(self.mdl, self.rs)
808 
809  def get_restraint(self):
810  return self.rs
811 
812  def set_weight(self, weight):
813  self.weight = weight
814  self.rs.set_weight(weight)
815 
816  def get_excluded_pairs(self):
817  return self.pairslist
818 
819  def get_output(self):
820  output = {}
821  score = self.weight * self.rs.unprotected_evaluate(None)
822  output["SymmetryRestraint_" + self.label] = str(score)
823  output["_TotalScore"] = str(score)
824  return output
825 
826 
827 class FusionRestraint(object):
828  """Creates a restraint between the termini two polypeptides, to simulate
829  the sequence connectivity."""
830  def __init__(self, nterminal, cterminal, scale=1.0, disorderedlength=False,
831  upperharmonic=True, resolution=1, label="None"):
832  """
833  @param nterminal - single PMI2 Hierarchy/molecule at the nterminal
834  @param cterminal - single PMI2 Hierarchy/molecule at the cterminal
835  @param scale Scale the maximal distance between the beads by this
836  factor when disorderedlength is False.
837  The maximal distance is calculated as
838  ((float(residuegap) + 1.0) * 3.6) * scale.
839  @param disorderedlength - This flag uses either disordered length
840  calculated for random coil peptides (True) or zero
841  surface-to-surface distance between beads (False)
842  as optimal distance for the sequence connectivity
843  restraint.
844  @param upperharmonic - This flag uses either harmonic (False)
845  or upperharmonic (True) in the intra-pair
846  connectivity restraint.
847  @param resolution - The resolution to connect things at - only used
848  if you pass PMI objects
849  @param label - A string to identify this restraint in the
850  output/stat file
851  """
852  self.label = label
853  self.weight = 1.0
854  ssn = IMP.pmi.tools.get_sorted_segments(nterminal)
855  ssc = IMP.pmi.tools.get_sorted_segments(cterminal)
856  nter_lastres = ssn[-1][1]
857  cter_firstres = ssc[0][0]
858  self.m = nter_lastres.get_model()
859 
860  self.kappa = 10 # spring constant used for the harmonic restraints
861 
862  optdist = (3.6) * scale
863  if upperharmonic: # default
864  hu = IMP.core.HarmonicUpperBound(optdist, self.kappa)
865  else:
866  hu = IMP.core.Harmonic(optdist, self.kappa)
868 
869  pt0 = nter_lastres.get_particle()
870  pt1 = cter_firstres.get_particle()
871  r = IMP.core.PairRestraint(self.m, dps,
872  (pt0.get_index(), pt1.get_index()))
873  self.rs = IMP.RestraintSet(self.m, "fusion_restraint")
874  print("Adding fusion connectivity restraint between", pt0.get_name(),
875  " and ", pt1.get_name(), 'of distance', optdist)
876  self.rs.add_restraint(r)
877 
878  def set_label(self, label):
879  self.label = label
880 
881  def get_weight(self):
882  return self.weight
883 
884  def add_to_model(self):
885  IMP.pmi.tools.add_restraint_to_model(self.m, self.rs)
886 
887  def get_restraint(self):
888  return self.rs
889 
890  def set_weight(self, weight):
891  self.weight = weight
892  self.rs.set_weight(weight)
893 
894  def get_output(self):
895  output = {}
896  score = self.weight * self.rs.unprotected_evaluate(None)
897  output["_TotalScore"] = str(score)
898  output["FusionRestraint_" + self.label] = str(score)
899  return output
900 
901  def evaluate(self):
902  return self.weight * self.rs.unprotected_evaluate(None)
903 
904 
905 class PlaneDihedralRestraint(IMP.pmi.restraints.RestraintBase):
906 
907  """Restrain the dihedral between planes defined by three particles.
908 
909  This restraint is useful for restraining the twist of a string of
910  more or less identical rigid bodies, so long as the curvature is mild.
911  """
912 
913  def __init__(self, particle_triplets, angle=0., k=1., label=None,
914  weight=1.):
915  """Constructor
916  @param particle_triplets List of lists of 3 particles. Each triplet
917  defines a plane. Dihedrals of adjacent planes
918  in list are scored.
919  @param angle Angle of plane dihedral in degrees
920  @param k Strength of restraint
921  @param label Label for output
922  @param weight Weight of restraint
923  @note Particles defining planes should be rigid and more or less
924  parallel for proper behavior
925  """
926  model = particle_triplets[0][0].get_model()
927  super(PlaneDihedralRestraint, self).__init__(model, label=label,
928  weight=weight)
929 
930  angle = math.pi * angle / 180.
931  ds = IMP.core.Cosine(.5 * k, 1, -angle)
932  for i, t1 in enumerate(particle_triplets[:-1]):
933  t2 = particle_triplets[i + 1]
934  q1 = [t1[1], t1[0], t2[0], t2[1]]
935  q2 = [t1[2], t1[0], t2[0], t2[2]]
936  self.rs.add_restraint(
937  IMP.core.DihedralRestraint(self.model, ds, *q1))
938  self.rs.add_restraint(
939  IMP.core.DihedralRestraint(self.model, ds, *q2))
A filter which returns true if a container containers the Pair.
Add dihedral restraints between quadruplet of consecutive residues/beads to enforce the stereochemist...
CHARMMParameters * get_heavy_atom_CHARMM_parameters()
Lower bound harmonic function (non-zero when feature < mean)
static bool get_is_setup(const IMP::ParticleAdaptor &p)
Definition: Residue.h:158
Enforce CHARMM stereochemistry on the given Hierarchy.
A member of a rigid body, it has internal (local) coordinates.
Definition: rigid_bodies.h:616
static bool get_is_setup(const IMP::ParticleAdaptor &p)
Definition: rigid_bodies.h:617
Various classes to hold sets of particles.
Upper bound harmonic function (non-zero when feature > mean)
A class to store an fixed array of same-typed values.
Definition: Array.h:33
Enforce ideal Helix dihedrals and bonds for a selection at resolution 0.
Miscellaneous utilities.
Definition: tools.py:1
Cosine function.
Definition: Cosine.h:22
Creates a restraint between the termini two polypeptides, to simulate the sequence connectivity...
Apply a function to the distance between two particles after transforming the second.
Restrain the dihedral between planes defined by three particles.
Dihedral restraint between four particles.
The type of an atom.
A score on the distance between the surfaces of two spheres.
Return all close unordered pairs of particles taken from the SingletonContainer.
static bool get_is_setup(const IMP::ParticleAdaptor &p)
Definition: rigid_bodies.h:731
Return all spatially-proximals pairs of particles (a,b) from the two SingletonContainers A and B...
GenericHierarchies get_leaves(Hierarchy mhd)
Get all the leaves of the bit of hierarchy.
Distance restraint between two particles.
Object used to hold a set of restraints.
Definition: RestraintSet.h:37
def input_adaptor
Adapt things for PMI (degrees of freedom, restraints, ...) Returns list of list of hierarchies...
Definition: tools.py:946
Store a list of ParticleIndexPairs.
A well with harmonic barriers.
Definition: HarmonicWell.h:25
Angle restraint between three particles.
ParticleIndexPairs get_indexes(const ParticlePairsTemp &ps)
Get the indexes from a list of particle pairs.
Add bond restraint between pair of consecutive residues/beads to enforce the stereochemistry.
The standard decorator for manipulating molecular structures.
RestraintSet * create_elastic_network(const Particles &ps, Float dist_cutoff, Float strength)
Create an elastic network restraint set.
Definition: pmi/utilities.h:26
Add bonds and improper dihedral restraints for the CBs.
Store a list of ParticleIndexes.
A decorator for a particle representing an atom.
Definition: atom/Atom.h:238
def get_particle_pairs
Returns the list of connected particles pairs.
def add_restraint_to_model
Add a PMI restraint to the model.
Definition: tools.py:91
Create harmonic restraints between the reference and (transformed) clones.
Add angular restraint between triplets of consecutive residues/beads to enforce the stereochemistry...
def get_num_restraints
Returns number of connectivity restraints.
A decorator for a residue.
Definition: Residue.h:137
Basic functionality that is expected to be used by a wide variety of IMP users.
Create a restraint between consecutive TempResidue objects or an entire PMI Molecule object...
A class to create an excluded volume restraint for a set of particles at a given resolution.
The general base class for IMP exceptions.
Definition: exception.h:48
def __init__
need to add: ca-ca bond ca-cb is a constraint, no restraint needed ca-ca-ca cb-ca-ca-cb ...
Applies a PairScore to a Pair.
Definition: PairRestraint.h:29
def get_sorted_segments
Returns sequence-sorted segments array, each containing the first particle the last particle and the ...
Definition: tools.py:1081
Functionality for loading, creating, manipulating and scoring atomic structures.
Select hierarchy particles identified by the biological name.
Definition: Selection.h:70
Applies a PairScore to each Pair in a list.
def get_residue_indexes
Retrieve the residue indexes for the given particle.
Definition: tools.py:570
A repulsive potential on the distance between two atoms.
Perform more efficient close pair finding when rigid bodies are involved.
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:645
Harmonic function (symmetric about the mean)
Definition: core/Harmonic.h:24
Restraint a set of residues to use ideal helix dihedrals and bonds.