IMP logo
IMP Reference Guide  2.6.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 
26 class DegreesOfFreedom(object):
27  """A class to simplify create of constraints and movers for an IMP Hierarchy.
28  Call the various create() functions to get started.
29  Can get all enabled movers with get_movers(). Pass this to ReplicaExchange0.
30  """
31  def __init__(self,mdl):
32  self.mdl = mdl
33  self.movers = []
34  self.rigid_bodies = [] #stores rigid body objects
35  self.flexible_beads = [] # stores all beads including nonrigid members of rigid bodies
36  self.nuisances = []
37 
38  #self.particle_map = {} # map from particles/rb objects to relevant movers+constraints
39  # internal mover = [mover obj, list of particles, enabled?] ?
40  # mover map = {particles/rbs : movers}
41 
42  def create_rigid_body(self,
43  rigid_parts,
44  nonrigid_parts=None,
45  max_trans=1.0,
46  max_rot=0.1,
47  nonrigid_max_trans = 0.1,
48  resolution='all',
49  name=None):
50  """Create rigid body constraint and mover
51  @param rigid_parts Can be one of the following inputs:
52  IMP Hierarchy, PMI System/State/Molecule/TempResidue, a list/set (of list/set) of them
53  or a RigidBody object.
54  Must be uniform input, however. No mixing object types.
55  @param nonrigid_parts Same input format as rigid_parts.
56  Must be a subset of rigid_parts particles.
57  @param max_trans Maximum rigid body translation
58  @param max_rot Maximum rigid body rotation
59  @param nonrigid_max_trans Maximum step for the nonrigid (bead) particles
60  @param resolution Only used if you pass PMI objects. Probably you want 'all'.
61  @param name Rigid body name (if None, use IMP default)
62  \note If you want all resolutions, pass PMI objects because this function will get them all.
63  Alternatively you can do your selection elsewhere and just pass hierarchies.
64  Returns tuple (rb_movers,rb_object)
65  """
66 
67  rb_movers = []
68 
69  # ADD CHECK: these particles are not already part of some RB or SRB
70 
71  # First, is this already a rigid body?
72  if type(rigid_parts) is IMP.core.RigidBody:
73  rb = rigid_parts
74  if name is not None:
75  name = rb.get_name()
76  else:
77  ### Otherwise, setup RB
78  hiers = IMP.pmi.tools.input_adaptor(rigid_parts,
79  resolution,
80  flatten=True)
81  if not hiers:
82  print("WARNING: No hierarchies were passed to create_rigid_body()")
83  return []
84  rb = IMP.atom.create_rigid_body(hiers)
85  self.rigid_bodies.append(rb)
86  rb.set_coordinates_are_optimized(True)
87  rb_mover = IMP.core.RigidBodyMover(rb,max_trans,max_rot)
88  if name is not None:
89  rb.set_name(name)
90  rb_mover.set_name(name)
91  rb_movers.append(rb_mover)
92 
93  ### setup nonrigid parts
94  if nonrigid_parts:
95  nr_hiers = IMP.pmi.tools.input_adaptor(nonrigid_parts,
96  resolution,
97  flatten=True)
98  if nr_hiers:
99  floatkeys = [IMP.FloatKey(4), IMP.FloatKey(5), IMP.FloatKey(6)]
100  rb_idxs = set(rb.get_member_indexes())
101  for h in nr_hiers:
102  self.flexible_beads.append(h)
103  p = h.get_particle()
104  if not p.get_index() in rb_idxs:
105  raise Exception("You tried to create nonrigid members from "
106  "particles that aren't in the RigidBody!")
107 
108  rb.set_is_rigid_member(p.get_index(),False)
109  for fk in floatkeys:
110  p.set_is_optimized(fk,True)
111  rb_movers.append(IMP.core.BallMover([p],
112  IMP.FloatKeys(floatkeys),
113  nonrigid_max_trans))
114  self.movers += rb_movers # probably need to store more info
115  return rb_movers,rb
116 
118  srb_parts,
119  max_trans=1.0,
120  max_rot=0.1,
121  chain_min_length=None,
122  chain_max_length=None,
123  resolution='all',
124  name=None):
125  """Create SUPER rigid body mover from one or more hierarchies. Can also create chain of SRBs.
126  @param srb_parts Can be one of the following inputs:
127  IMP Hierarchy, PMI System/State/Molecule/TempResidue, or a list/set (of list/set) of them.
128  Must be uniform input, however. No mixing object types.
129  @param max_trans Maximum super rigid body translation
130  @param max_rot Maximum super rigid body rotation
131  @param chain_min_length Create a CHAIN of super rigid bodies - must provide list
132  This parameter is the minimum chain length.
133  @param chain_max_length Maximum chain length
134  @param resolution Only used if you pass PMI objects. Probably you want 'all'.
135  @param name The name of the SRB (hard to assign a good one automatically)
136  \note If you set the chain parameters, will NOT create an SRB from all of them together,
137  but rather in groups made from the outermost list.
138  """
139 
140  srb_movers = []
141 
142  ## organize hierarchies based on chains
143  if chain_min_length is None and chain_max_length is None:
144  # the "chain" is just everything together
145  h = IMP.pmi.tools.input_adaptor(srb_parts,resolution,flatten=True)
146  if len(h)==0:
147  print('WARNING: No hierarchies were passed to create_super_rigid_body()')
148  return srb_movers
149  srb_groups = [h]
150  else:
151  if not hasattr(srb_parts,'__iter__'):
152  raise Exception("You tried to make a chain without a list!")
153  srb_groups = [IMP.pmi.tools.input_adaptor(h,resolution,flatten=True) for h in srb_parts]
154 
155  ## create SRBs either from all hierarchies or chain
156  if chain_min_length is None and chain_max_length is None:
157  srb_movers.append(self._setup_srb(srb_groups,max_trans,max_rot))
158  elif chain_min_length is not None and chain_max_length is not None:
159  for hs in IMP.pmi.tools.sublist_iterator(srb_groups, chain_min_length, chain_max_length):
160  srb_movers.append(self._setup_srb(hs,max_trans,max_rot))
161  else:
162  raise Exception("DegreesOfFreedom: SetupSuperRigidBody: if you want chain, specify min AND max")
163  self.movers += srb_movers
164  if name is not None:
165  if len(srb_movers)>1:
166  for n,mv in enumerate(srb_movers):
167  mv.set_name(name+'_'+str(n))
168  else:
169  srb_movers[0].set_name(name)
170  return srb_movers
171 
172  def _setup_srb(self,hiers,max_trans,max_rot):
173  srbm = IMP.pmi.TransformMover(hiers[0][0].get_model(), max_trans, max_rot)
174  super_rigid_rbs,super_rigid_xyzs = IMP.pmi.tools.get_rbs_and_beads(hiers)
175  for xyz in super_rigid_xyzs:
176  srbm.add_xyz_particle(xyz)
177  for rb in super_rigid_rbs:
178  srbm.add_rigid_body_particle(rb)
179  return srbm
180 
181 
183  flex_parts,
184  max_trans=0.1,
185  resolution='all'):
186  """Create a chain of flexible beads
187  @param flex_parts Can be one of the following inputs:
188  IMP Hierarchy, PMI System/State/Molecule/TempResidue, or a list/set (of list/set) of them.
189  Must be uniform input, however. No mixing object types.
190  @param max_trans Maximum flexible bead translation
191  @param resolution Only used if you pass PMI objects. Probably you want 'all'.
192  """
193 
194  fb_movers = []
195  hiers = IMP.pmi.tools.input_adaptor(flex_parts,
196  resolution,
197  flatten=True)
198  if not hiers or len(hiers)==0:
199  print('WARNING: No hierarchies were passed to create_flexible_beads()')
200  return fb_movers
201  for h in hiers:
202  p = h.get_particle()
203  IMP.core.XYZ(p).set_coordinates_are_optimized(True)
205  raise Exception("Cannot create flexible beads from members of rigid body")
206  self.flexible_beads.append(h)
207  fb_movers.append(IMP.core.BallMover([p],max_trans))
208  self.movers += fb_movers
209  return fb_movers
210 
212  nuisance_p,
213  step_size,
214  name=None):
215  """Create MC normal mover for nuisance particles.
216  We will add an easier interface to add all of them from a PMI restraint
217  @param nuisance_p The Nuisance particle (an ISD::Scale)
218  @param step_size The maximum step size for Monte Carlo
219  @param name The name of the mover, useful for better output reading.
220  """
221  mv = IMP.core.NormalMover([nuisance_p],
222  IMP.FloatKeys([IMP.FloatKey("nuisance")]),
223  step_size)
224  if name is not None:
225  mv.set_name(name)
226  self.nuisances.append(nuisance_p)
227  self.movers.append(mv)
228  return [mv]
229 
230  def setup_md(self,
231  hspec):
232  """Setup particles for MD simulation. Returns all particles, just
233  pass this to molecular_dynamics_sample_objects in ReplicaExchange0.
234  @param hspec Can be one of the following inputs:
235  IMP Hierarchy, PMI System/State/Molecule/TempResidue, or a list/set (of list/set) of them.
236  Must be uniform input, however. No mixing object types.
237  """
238  vxkey = IMP.FloatKey('vx')
239  vykey = IMP.FloatKey('vy')
240  vzkey = IMP.FloatKey('vz')
241  hiers = IMP.pmi.tools.input_adaptor(hspec,flatten=True)
242  mdl = hiers[0].get_model()
243  all_ps = []
244  for h in hiers:
245  p = h.get_particle()
246  IMP.core.XYZ(mdl,p.get_index()).set_coordinates_are_optimized(True)
247  mdl.add_attribute(vxkey,p.get_index(),0.0)
248  mdl.add_attribute(vykey,p.get_index(),0.0)
249  mdl.add_attribute(vzkey,p.get_index(),0.0)
250  all_ps.append(p)
251  return all_ps
252 
254  references,
255  clones,
256  transform,
257  resolution=1):
258  """Create a symmetry constraint. Checks:
259  same number of particles
260  disable ANY movers involving symmetry copies
261  (later may support moving them WITH references,
262  but requires code to propagate constraint)
263  @param references Can be one of the following inputs:
264  IMP Hierarchy, PMI System/State/Molecule/TempResidue, or a list/set (of list/set) of them
265  @param clones Same format as references
266  @param transform The transform that moves a clone onto a reference
267  IMP.algebra.Transformation3D
268  @param resolution Only used if you pass PMI objects.
269  If you have a multires system, assuming each are rigid
270  bodies you probably only need one resolution.
271  """
272 
273  # get all RBs and particles
274  href = IMP.pmi.tools.input_adaptor(references,resolution,flatten=True)
275  hclones = IMP.pmi.tools.input_adaptor(clones,resolution,flatten=True)
276 
277  ref_rbs,ref_beads = IMP.pmi.tools.get_rbs_and_beads(href)
278  clones_rbs,clones_beads = IMP.pmi.tools.get_rbs_and_beads(hclones)
279 
280  # dumb check for matching numbers of particles
281  if len(ref_rbs)!=len(clones_rbs) or len(ref_beads)!=len(clones_beads):
282  raise Exception("ERROR: Your references don't match your clones")
283 
284  # symmetry RBs
285  for ref,clone in zip(ref_rbs+ref_beads,clones_rbs+clones_beads):
287  sm = IMP.core.TransformationSymmetry(transform.get_inverse())
289  self.mdl,[p.get_particle_index() for p in clones_rbs+clones_beads])
290  c = IMP.container.SingletonsConstraint(sm, None, lsc)
291  self.mdl.add_score_state(c)
292  print('Created symmetry restraint for',len(ref_rbs),'rigid bodies and',
293  len(ref_beads),'flexible beads')
294 
295  def __repr__(self):
296  # would like something fancy like this:
297  #- super-rigid "SRB1"
298  # - rigid "Mol1" (8 rigid, 3 nonrigid)
299  # - rigid "Mol2" (8 rigid, 3 nonrigid)
300  # - rigid "Mol3" (8 rigid, 3 nonrigid)
301  return 'DegreesOfFreedom: ' + \
302  "\n".join(repr(m) for m in self.movers)
303 
304  def optimize_flexible_beads(self, nsteps, temperature=1.0):
305  '''Set up MC run with just flexible beads.
306  Optimization works much better when restraints
307  are already set up.'''
308  import IMP.pmi.samplers
309  pts = IMP.pmi.tools.ParticleToSampleList()
310  for n, fb in enumerate(self.get_flexible_beads()):
311  pts.add_particle(fb, "Floppy_Bodies", 1.0, "Flexible_Bead_" + str(n))
312  if len(pts.get_particles_to_sample()) > 0:
313  mc = IMP.pmi.samplers.MonteCarlo(self.mdl, [pts], temperature)
314  print("optimize_flexible_beads: optimizing %i flexible beads" % len(self.get_flexible_beads()))
315  mc.optimize(nsteps)
316  else:
317  print("optimize_flexible_beads: no particle to optimize")
318 
319  def get_movers(self):
320  """Should only return Enabled movers?"""
321  return self.movers
322 
323  def get_rigid_bodies(self):
324  """Return list of rigid body objects"""
325  return self.rigid_bodies
326 
328  """Return all flexible beads, including nonrigid members of rigid bodies"""
329  return self.flexible_beads
330 
332  """Extract the nuisances from get_particles_to_sample()"""
333  try:
334  pslist = r.get_particles_to_sample()
335  except:
336  raise Exception("dof.get_nuisances_from_restraint(): the passed object does not have a "
337  "get_particles_to_sample() function")
338  for name in pslist:
339  is_sampled = True
340  if len(pslist[name])==3:
341  ps,maxtrans,is_sampled = pslist[name]
342  else:
343  ps,maxtrans = pslist[name]
344  if is_sampled:
345  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.
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
Should only return Enabled movers?
Modify a set of continuous variables by perturbing them within a ball.
static bool get_is_setup(const IMP::ParticleAdaptor &p)
Definition: rigid_bodies.h:469
def setup_md
Setup particles for MD simulation.
def input_adaptor
Adapt things for PMI (degrees of freedom, restraints, ...) Returns list of list of hierarchies...
Definition: tools.py:1455
def create_rigid_body
Create rigid body constraint and mover.
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:1629
def get_flexible_beads
Return all flexible beads, including nonrigid members of rigid bodies.
static bool get_is_setup(const IMP::ParticleAdaptor &p)
Definition: rigid_bodies.h:490
def sublist_iterator
Yield all sublists of length >= lmin and <= lmax.
Definition: tools.py:1096