IMP logo
IMP Reference Guide  2.7.0
The Integrative Modeling Platform
pmi/dof/__init__.py
1 """@namespace IMP.pmi.dof
2  Create movers and setup constraints for PMI objects.
3 * Start by creating the DegreesOfFreedom class with `dof = IMP::pmi::dof::DegreesOfFreedom(model)`
4 * The various "create X" functions make movers for system components as well as setup necessary constraints.
5 For each of these functions, you can generally pass PMI objects like [Molecule](@ref IMP::pmi::topology::Molecule) or slices thereof.
6  * DegreesOfFreedom.create_rigid_body() lets you rigidify a molecule (but allows you to also pass "nonrigid" components which move with the body and also independently).
7  * DegreesOfFreedom.create_super_rigid_body() sets up a special "Super Rigid Body" which moves
8 rigidly but is not always constrained to be rigid (so you can later move the parts separately). Good for speeding up sampling.
9  * DegreesOfFreedom.create_flexible_beads() sets up particles to move individually.
10  * DegreesOfFreedom.setup_md() sets up particles to move with molecular dynamics. Note that this is not (yet) compatible with rigid bodies, and only some restraints.
11  * DegreesOfFreedom.constrain_symmetry() makes a symmetry constraint so that clones automatically move with their references. If instead you want a softer restraint, check out the [SymmetryRestraint](@ref IMP::pmi::restraints::stereochemistry::SymmetryRestraint).
12 * When you are done you can access all movers with DegreesOfFreedom.get_movers(). If you have set up rigid, super rigid, or flexible beads, pass the movers to the `monte_carlo_sample_objects` argument of [ReplicaExchange0](@ref IMP::pmi::macros::ReplicaExchange0).
13 * If you are running MD, you have to separately pass the particles (also returned from DegreesOfFreedom.setup_md()) to the `molecular_dynamics_sample_objects` argument of [ReplicaExchange0](@ref IMP::pmi::macros::ReplicaExchange0). Check out an [MD example here](pmi_2atomistic_8py-example.html).
14 """
15 
16 from __future__ import print_function
17 import IMP
18 import IMP.atom
19 import IMP.algebra
20 import IMP.pmi
21 import IMP.pmi.topology
22 import IMP.pmi.samplers
23 import IMP.pmi.tools
24 import itertools
25 import IMP.pmi.samplers
26 
27 class DegreesOfFreedom(object):
28  """A class to simplify create of constraints and movers for an IMP Hierarchy.
29  Call the various create() functions to get started.
30  Can get all enabled movers with get_movers(). Pass this to ReplicaExchange0.
31  """
32  def __init__(self,mdl):
33  self.mdl = mdl
34  self.movers = []
35  self.fb_movers = [] #stores movers corresponding to floppy parts
36  self.rigid_bodies = [] #stores rigid body objects
37  self.flexible_beads = [] # stores all beads including nonrigid members of rigid bodies
38  self.nuisances = []
39  self._rb2mov = {} # Keys are the RigidBody objects, values are list of movers
40  # the following is needed to keep track of disabled movers
41  self.movers_particles_map=IMP.pmi.tools.OrderedDict()
42  self.movers_rb_map={}
43  self.movers_xyz_map={}
44  self.disabled_movers=[]
45 
46  #self.particle_map = {} # map from particles/rb objects to relevant movers+constraints
47  # internal mover = [mover obj, list of particles, enabled?] ?
48  # mover map = {particles/rbs : movers}
49 
51  rigid_parts,
52  nonrigid_parts=None,
53  max_trans=4.0,
54  max_rot=0.5,
55  nonrigid_max_trans = 4.0,
56  resolution='all',
57  name=None):
58  """Create rigid body constraint and mover
59  @param rigid_parts Can be one of the following inputs:
60  IMP Hierarchy, PMI System/State/Molecule/TempResidue, a list/set (of list/set) of them
61  or a RigidBody object.
62  Must be uniform input, however. No mixing object types.
63  @param nonrigid_parts Same input format as rigid_parts.
64  Must be a subset of rigid_parts particles.
65  @param max_trans Maximum rigid body translation
66  @param max_rot Maximum rigid body rotation
67  @param nonrigid_max_trans Maximum step for the nonrigid (bead) particles
68  @param resolution Only used if you pass PMI objects. Probably you want 'all'.
69  @param name Rigid body name (if None, use IMP default)
70  \note If you want all resolutions, pass PMI objects because this function will get them all.
71  Alternatively you can do your selection elsewhere and just pass hierarchies.
72  Returns tuple (rb_movers,rb_object)
73  """
74 
75  rb_movers = []
76 
77  # ADD CHECK: these particles are not already part of some RB or SRB
78 
79  # First, is this already a rigid body?
80  if type(rigid_parts) is IMP.core.RigidBody:
81  rb = rigid_parts
82  model=rb.get_model()
83  if name is not None:
84  name = rb.get_name()
85  hiers= [IMP.atom.get_leaves(IMP.atom.Hierarchy(m.get_particle(i)))[0] for i in rb.get_member_particle_indexes()]
86  else:
87  ### Otherwise, setup RB
88  hiers = IMP.pmi.tools.input_adaptor(rigid_parts,
89  resolution,
90  flatten=True)
91  if not hiers:
92  print("WARNING: No hierarchies were passed to create_rigid_body()")
93  return []
94  rb = IMP.atom.create_rigid_body(hiers)
95  self.rigid_bodies.append(rb)
96  rb.set_coordinates_are_optimized(True)
97  rb_mover = IMP.core.RigidBodyMover(rb.get_model(), rb, max_trans,
98  max_rot)
99  if name is not None:
100  rb.set_name(name)
101  rb_mover.set_name(name)
102  rb_movers.append(rb_mover)
103  self.movers_particles_map[rb_mover]=[]
104  self.movers_rb_map[rb_mover]=[rb]
105  for h in hiers:
106  self.movers_particles_map[rb_mover]+=IMP.atom.get_leaves(h)
107  ### setup nonrigid parts
108  if nonrigid_parts:
109  nr_hiers = IMP.pmi.tools.input_adaptor(nonrigid_parts,
110  resolution,
111  flatten=True)
112  if nr_hiers:
113  floatkeys = [IMP.FloatKey(4), IMP.FloatKey(5), IMP.FloatKey(6)]
114  rb_idxs = set(rb.get_member_indexes())
115  for h in nr_hiers:
116  self.flexible_beads.append(h)
117  p = h.get_particle()
118  if not p.get_index() in rb_idxs:
119  raise Exception("You tried to create nonrigid members from "
120  "particles that aren't in the RigidBody!")
121 
122  rb.set_is_rigid_member(p.get_index(),False)
123  for fk in floatkeys:
124  p.set_is_optimized(fk,True)
125  fbmv=IMP.core.BallMover(p.get_model(), p,
126  IMP.FloatKeys(floatkeys),
127  nonrigid_max_trans)
128  self.fb_movers.append(fbmv)
129  self.movers_particles_map[fbmv]=IMP.atom.get_leaves(h)
130  self.movers_xyz_map[fbmv]=IMP.atom.get_leaves(h)
131  rb_movers.append(fbmv)
132 
133  self.movers += rb_movers # probably need to store more info
134  self._rb2mov[rb] = rb_movers #dictionary relating rb to movers
135 
136  return rb_movers,rb
137 
138  def create_main_chain_mover(self,molecule,resolution=10,lengths=[5,10]):
139  """Create crankshaft moves from a set of SUPER rigid body mover from one molecule.
140  See http://scfbm.biomedcentral.com/articles/10.1186/1751-0473-3-12
141  """
142  hiers=IMP.pmi.tools.input_adaptor(molecule,resolution,flatten=True)
143 
144  for length in lengths:
145  for n in range(len(hiers)-length):
146  hs=hiers[n+1:n+length]
147  #print("MIDDLE",n+1,n+length,hs)
148  self.create_super_rigid_body(hs, max_trans=0.0,max_rot=0.05, axis=(hiers[n].get_particle(),hiers[n+length].get_particle()))
149  #for n in range(1,len(hiers)-1,10):
150  # hs=hiers[n:]
151  # print("END",n,hs)
152  # self.create_super_rigid_body(hs, max_trans=0.01,axis=(hiers[n].get_particle(),hiers[n+1].get_particle()))
153  # hs=hiers[:n+1]
154  # print("BEGIN",n-1,hs)
155  # self.create_super_rigid_body(hs, max_trans=0.01,axis=(hiers[n].get_particle(),hiers[n-1].get_particle()))
156 
157  def create_super_rigid_body(self,
158  srb_parts,
159  max_trans=1.0,
160  max_rot=0.1,
161  chain_min_length=None,
162  chain_max_length=None,
163  resolution='all',
164  name=None,
165  axis=None):
166  """Create SUPER rigid body mover from one or more hierarchies. Can also create chain of SRBs.
167  If you don't pass chain min/max, it'll treat everything you pass as ONE rigid body.
168  If you DO pass chain min/max, it'll expect srb_parts is a list and break it into bits.
169  @param srb_parts Can be one of the following inputs:
170  IMP Hierarchy, PMI System/State/Molecule/TempResidue, or a list/set (of list/set) of them.
171  Must be uniform input, however. No mixing object types.
172  @param max_trans Maximum super rigid body translation
173  @param max_rot Maximum super rigid body rotation
174  @param chain_min_length Create a CHAIN of super rigid bodies - must provide list
175  This parameter is the minimum chain length.
176  @param chain_max_length Maximum chain length
177  @param resolution Only used if you pass PMI objects. Probably you want 'all'.
178  @param name The name of the SRB (hard to assign a good one automatically)
179  @param axis A tuple containing two particles which are used to compute
180  the rotation axis of the SRB. The default is None, meaning that the rotation axis is random.
181  \note If you set the chain parameters, will NOT create an SRB from all of them together,
182  but rather in groups made from the outermost list.
183  """
184 
185  srb_movers = []
186 
187  ## organize hierarchies based on chains
188  if chain_min_length is None and chain_max_length is None:
189  # the "chain" is just everything together
190  h = IMP.pmi.tools.input_adaptor(srb_parts,resolution,flatten=True)
191  if len(h)==0:
192  print('WARNING: No hierarchies were passed to create_super_rigid_body()')
193  return srb_movers
194  srb_groups = [h]
195  else:
196  if not hasattr(srb_parts,'__iter__'):
197  raise Exception("You tried to make a chain without a list!")
198  srb_groups = [IMP.pmi.tools.input_adaptor(h,resolution,flatten=True,
199  warn_about_slices=False) for h in srb_parts]
200 
201  ## create SRBs either from all hierarchies or chain
202  if chain_min_length is None and chain_max_length is None:
203  mv = self._setup_srb(srb_groups,max_trans,max_rot,axis)
204  if mv:
205  srb_movers.append(mv)
206  elif chain_min_length is not None and chain_max_length is not None:
207  for hs in IMP.pmi.tools.sublist_iterator(srb_groups, chain_min_length, chain_max_length):
208  mv = self._setup_srb(hs,max_trans,max_rot,axis)
209  if mv:
210  srb_movers.append(mv)
211  else:
212  raise Exception("DegreesOfFreedom: SetupSuperRigidBody: if you want chain, specify min AND max")
213  self.movers += srb_movers
214  if name is not None:
215  if len(srb_movers)>1:
216  for n,mv in enumerate(srb_movers):
217  mv.set_name(name+'_'+str(n))
218  else:
219  srb_movers[0].set_name(name)
220  return srb_movers
221 
222  def _setup_srb(self,hiers,max_trans,max_rot,axis):
223  if axis is None:
224  srbm = IMP.pmi.TransformMover(hiers[0][0].get_model(), max_trans, max_rot)
225  else:
226  srbm = IMP.pmi.TransformMover(hiers[0][0].get_model(),axis[0],axis[1],max_trans, max_rot)
227  super_rigid_rbs,super_rigid_xyzs = IMP.pmi.tools.get_rbs_and_beads(hiers)
228  ct = 0
229  self.movers_particles_map[srbm]=[]
230  for h in hiers:
231  self.movers_particles_map[srbm]+=IMP.atom.get_leaves(h)
232  for xyz in super_rigid_xyzs:
233  srbm.add_xyz_particle(xyz)
234  ct+=1
235  for rb in super_rigid_rbs:
236  srbm.add_rigid_body_particle(rb)
237  ct+=1
238  if ct>1:
239  return srbm
240  else:
241  return 0
242 
243 
245  flex_parts,
246  max_trans=3.0,
247  resolution='all'):
248  """Create a chain of flexible beads
249  @param flex_parts Can be one of the following inputs:
250  IMP Hierarchy, PMI System/State/Molecule/TempResidue, or a list/set (of list/set) of them.
251  Must be uniform input, however. No mixing object types.
252  @param max_trans Maximum flexible bead translation
253  @param resolution Only used if you pass PMI objects. Probably you want 'all'.
254  """
255 
256  fb_movers = []
257  hiers = IMP.pmi.tools.input_adaptor(flex_parts,
258  resolution,
259  flatten=True)
260  if not hiers or len(hiers)==0:
261  print('WARNING: No hierarchies were passed to create_flexible_beads()')
262  return fb_movers
263  for h in hiers:
264  p = h.get_particle()
265  IMP.core.XYZ(p).set_coordinates_are_optimized(True)
267  raise Exception("Cannot create flexible beads from members of rigid body")
268  self.flexible_beads.append(h)
269  fbmv=IMP.core.BallMover(p.get_model(), p, max_trans)
270  fb_movers.append(fbmv)
271  self.fb_movers.append(fbmv)
272  self.movers_particles_map[fbmv]=IMP.atom.get_leaves(h)
273  self.movers_xyz_map[fbmv]=IMP.atom.get_leaves(h)
274  self.movers += fb_movers
275  return fb_movers
276 
278  nuisance_p,
279  step_size,
280  name=None):
281  """Create MC normal mover for nuisance particles.
282  We will add an easier interface to add all of them from a PMI restraint
283  @param nuisance_p The Nuisance particle (an ISD::Scale)
284  @param step_size The maximum step size for Monte Carlo
285  @param name The name of the mover, useful for better output reading.
286  """
287  mv = IMP.core.NormalMover([nuisance_p],
288  IMP.FloatKeys([IMP.FloatKey("nuisance")]),
289  step_size)
290  if name is not None:
291  mv.set_name(name)
292  self.nuisances.append(nuisance_p)
293  self.movers.append(mv)
294  return [mv]
295 
296  def setup_md(self,
297  hspec):
298  """Setup particles for MD simulation. Returns all particles, just
299  pass this to molecular_dynamics_sample_objects in ReplicaExchange0.
300  @param hspec Can be one of the following inputs:
301  IMP Hierarchy, PMI System/State/Molecule/TempResidue, or a list/set (of list/set) of them.
302  Must be uniform input, however. No mixing object types.
303  """
304  vxkey = IMP.FloatKey('vx')
305  vykey = IMP.FloatKey('vy')
306  vzkey = IMP.FloatKey('vz')
307  hiers = IMP.pmi.tools.input_adaptor(hspec,flatten=True)
308  mdl = hiers[0].get_model()
309  all_ps = []
310  for hl in hiers:
311  for h in IMP.core.get_leaves(hl):
312  p = h.get_particle()
313  IMP.core.XYZ(mdl,p.get_index()).set_coordinates_are_optimized(True)
314  mdl.add_attribute(vxkey,p.get_index(),0.0)
315  mdl.add_attribute(vykey,p.get_index(),0.0)
316  mdl.add_attribute(vzkey,p.get_index(),0.0)
317  all_ps.append(p)
318  return all_ps
319 
321  references,
322  clones,
323  transform,
324  resolution='all'):
325  """Create a symmetry constraint. Checks:
326  same number of particles
327  disable ANY movers involving symmetry copies
328  (later may support moving them WITH references,
329  but requires code to propagate constraint)
330  @param references Can be one of the following inputs:
331  IMP Hierarchy, PMI System/State/Molecule/TempResidue, or a list/set (of list/set) of them
332  @param clones Same format as references
333  @param transform The transform that moves a clone onto a reference
334  IMP.algebra.Transformation3D
335  @param resolution Only used if you pass PMI objects.
336  If you have a multires system, assuming each are rigid
337  bodies you probably only need one resolution.
338  """
339 
340  # get all RBs and particles
341  href = IMP.pmi.tools.input_adaptor(references,resolution,flatten=True)
342  hclones = IMP.pmi.tools.input_adaptor(clones,resolution,flatten=True)
343 
344  ref_rbs,ref_beads = IMP.pmi.tools.get_rbs_and_beads(href)
345  clones_rbs,clones_beads = IMP.pmi.tools.get_rbs_and_beads(hclones)
346 
347  # dumb check for matching numbers of particles
348  #if len(ref_rbs)!=len(clones_rbs) or len(ref_beads)!=len(clones_beads):
349  # raise Exception("ERROR: Your references don't match your clones")
350 
351  # symmetry RBs
352  # this code produces weird results (random flipping of rigid bodies
353  #for ref,clone in zip(ref_rbs+ref_beads,clones_rbs+clones_beads):
354  # print(clone.get_particle().get_index(),ref.get_particle().get_index())
355  # IMP.core.Reference.setup_particle(clone.get_particle(),ref.get_particle())
356  for ref,clone in zip(href,hclones):
357  IMP.core.Reference.setup_particle(clone.get_particle(),ref.get_particle())
358 
359  sm = IMP.core.TransformationSymmetry(transform.get_inverse())
361  self.mdl,[p.get_particle().get_index() for p in hclones])
362  c = IMP.container.SingletonsConstraint(sm, None, lsc)
363  self.mdl.add_score_state(c)
364  print('Created symmetry restraint for',len(ref_rbs),'rigid bodies and',
365  len(ref_beads),'flexible beads')
366 
367  # removing movers involved in clones
368  #sym_movers = []
369 
370  #sym_movers = [m for cl in clones_rbs for m in self._rb2mov[cl]]
371  #self.movers = [m for m in self.movers if m not in sym_movers]
372  self.mdl.update()
373  self.disable_movers(hclones)
374 
375  def __repr__(self):
376  # would like something fancy like this:
377  #- super-rigid "SRB1"
378  # - rigid "Mol1" (8 rigid, 3 nonrigid)
379  # - rigid "Mol2" (8 rigid, 3 nonrigid)
380  # - rigid "Mol3" (8 rigid, 3 nonrigid)
381  return 'DegreesOfFreedom: ' + \
382  "\n".join(repr(m) for m in self.movers)
383 
384  def optimize_flexible_beads(self, nsteps, temperature=1.0):
385  '''Set up MC run with just flexible beads.
386  Optimization works much better when restraints
387  are already set up.'''
388  pts = IMP.pmi.tools.ParticleToSampleList()
389  for n, fb in enumerate(self.get_flexible_beads()):
390  pts.add_particle(fb, "Floppy_Bodies", 1.0, "Flexible_Bead_" + str(n))
391  if len(pts.get_particles_to_sample()) > 0:
392  mc = IMP.pmi.samplers.MonteCarlo(self.mdl, [pts], temperature)
393  print("optimize_flexible_beads: optimizing %i flexible beads" % len(self.get_flexible_beads()))
394  mc.optimize(nsteps)
395  else:
396  print("optimize_flexible_beads: no particle to optimize")
397 
398  def get_movers(self):
399  """Returns Enabled movers"""
400  if self.disabled_movers:
401  filtered_mover_list=[]
402  for mv in self.movers:
403  if not mv in self.disabled_movers:
404  filtered_mover_list.append(mv)
405  return filtered_mover_list
406  else:
407  return self.movers
408 
410  """Return all movers corresponding to individual beads"""
411  return self.fb_movers
412 
413  def get_rigid_bodies(self):
414  """Return list of rigid body objects"""
415  return self.rigid_bodies
416 
418  """Return all flexible beads, including nonrigid members of rigid bodies"""
419  return self.flexible_beads
420 
421  def disable_movers(self,objects,mover_types=None):
422  """Fix the position of the particles by disabling the corresponding movers
423  @param objects Can be one of the following inputs:
424  IMP Hierarchy, PMI System/State/Molecule/TempResidue, or a list/set (of list/set) of them.
425  Must be uniform input, however. No mixing object types.
426  @param mover_types further filter the mover type that will be disabled, it can be a list of IMP.core.RigidBodyMover,
427  IMP.core.BallMover etc etc if one wants to fix the corresponding rigid body, or the floppy bodies.
428  An empty mover_types list is interpreted as all possible movers.
429  It returns the list of fixed xyz particles (ie, floppy bodies/beads) and rigid bodies
430  whose movers were disabled
431  """
432  hierarchies = IMP.pmi.tools.input_adaptor(objects,
433  pmi_resolution='all',
434  flatten=True)
435  tmp_set=set()
436  fixed_rb=set()
437  fixed_xyz=set()
438  if mover_types is None: mover_types=[]
439 
440  for key in self.movers_particles_map:
441  for h in hierarchies:
442  if h in self.movers_particles_map[key] and (type(key) in mover_types or not mover_types):
443  tmp_set.add(key)
444  if key in self.movers_rb_map:
445  fixed_rb|=set(self.movers_rb_map[key])
446  if key in self.movers_xyz_map:
447  fixed_xyz|=set(self.movers_xyz_map[key])
448  print("Fixing %s movers" %(str(len(list(tmp_set)))))
449  self.disabled_movers+=list(tmp_set)
450  return list(fixed_xyz),list(fixed_rb)
451 
452  def enable_all_movers(self):
453  """Reenable all movers: previously fixed particles will be released"""
454  self.disabled_movers=[]
455 
457  """Extract the nuisances from get_particles_to_sample()"""
458  try:
459  pslist = r.get_particles_to_sample()
460  except:
461  raise Exception("dof.get_nuisances_from_restraint(): the passed object does not have a "
462  "get_particles_to_sample() function")
463  for name in pslist:
464  is_sampled = True
465  if len(pslist[name])==3:
466  ps,maxtrans,is_sampled = pslist[name]
467  else:
468  ps,maxtrans = pslist[name]
469  if is_sampled:
470  self.create_nuisance_mover(ps[0],maxtrans,name)
def optimize_flexible_beads
Set up MC run with just flexible beads.
A class to simplify create of constraints and movers for an IMP Hierarchy.
def get_rigid_bodies
Return list of rigid body objects.
Apply a SingletonFunction to a SingletonContainer to maintain an invariant.
Set the coordinates of a particle to be a transformed version of a reference.
Definition: core/symmetry.h:76
Set of python classes to create a multi-state, multi-resolution IMP hierarchy.
def create_flexible_beads
Create a chain of flexible beads.
def enable_all_movers
Reenable all movers: previously fixed particles will be released.
Miscellaneous utilities.
Definition: tools.py:1
def get_nuisances_from_restraint
Extract the nuisances from get_particles_to_sample()
Modify the transformation of a rigid body.
def create_nuisance_mover
Create MC normal mover for nuisance particles.
def get_movers
Returns Enabled movers.
def get_floppy_body_movers
Return all movers corresponding to individual beads.
Move continuous particle variables by perturbing them within a ball.
static bool get_is_setup(const IMP::ParticleAdaptor &p)
Definition: rigid_bodies.h:619
GenericHierarchies get_leaves(Hierarchy mhd)
Get all the leaves of the bit of hierarchy.
def setup_md
Setup particles for MD simulation.
def create_main_chain_mover
Create crankshaft moves from a set of SUPER rigid body mover from one molecule.
def input_adaptor
Adapt things for PMI (degrees of freedom, restraints, ...) Returns list of list of hierarchies...
Definition: tools.py:1608
def create_rigid_body
Create rigid body constraint and mover.
The standard decorator for manipulating molecular structures.
Ints get_index(const ParticlesTemp &particles, const Subset &subset, const Subsets &excluded)
def create_super_rigid_body
Create SUPER rigid body mover from one or more hierarchies.
Store a list of ParticleIndexes.
A decorator for a particle with x,y,z coordinates.
Definition: XYZ.h:30
Modify the transformation of a rigid body.
Modify a set of continuous variables using a normal distribution.
Definition: NormalMover.h:20
Sampling of the system.
Definition: samplers.py:1
General purpose algebraic and geometric methods that are expected to be used by a wide variety of IMP...
Sample using Monte Carlo.
Definition: samplers.py:36
The general base class for IMP exceptions.
Definition: exception.h:49
IMP::core::RigidBody create_rigid_body(Hierarchy h)
static Reference setup_particle(Model *m, ParticleIndex pi, ParticleIndexAdaptor reference)
Definition: core/symmetry.h:32
def constrain_symmetry
Create a symmetry constraint.
Python classes to represent, score, sample and analyze models.
A decorator for a rigid body.
Definition: rigid_bodies.h:75
Functionality for loading, creating, manipulating and scoring atomic structures.
def get_rbs_and_beads
Returns unique objects in original order.
Definition: tools.py:1861
Hierarchies get_leaves(const Selection &h)
def get_flexible_beads
Return all flexible beads, including nonrigid members of rigid bodies.
def disable_movers
Fix the position of the particles by disabling the corresponding movers.
static bool get_is_setup(const IMP::ParticleAdaptor &p)
Definition: rigid_bodies.h:640
def sublist_iterator
Yield all sublists of length >= lmin and <= lmax.
Definition: tools.py:1164