IMP logo
IMP Reference Guide  develop.d4e9f3251e,2024/04/26
The Integrative Modeling Platform
EMageFit/imp_general/representation.py
1 """@namespace IMP.EMageFit.imp_general.representation
2  Utility functions to handle representation.
3 """
4 
5 import IMP
6 import IMP.atom
7 import IMP.core
8 import IMP.algebra
9 import string
10 import logging
11 
12 log = logging.getLogger("representation")
13 
14 #
15 """
16 
17  Functions to deal with the representation of assemblies and managing
18  rigid bodies.
19 
20 """
21 #
22 
23 
24 def create_assembly_from_pdb(model, fn_pdb, names=False):
25  """
26  Builds the assembly setting the chains in the PDB file as components
27  """
28  temp = read_component(model, fn_pdb)
29  hchains = IMP.atom.get_by_type(temp, IMP.atom.CHAIN_TYPE)
30  ids = [IMP.atom.Chain(h).get_id() for h in hchains]
31  log.debug("Creating assembly from pdb %s,names: %s. Chains %s",
32  fn_pdb, names, ids)
33  IMP.atom.add_radii(temp)
34  if names:
35  for i, h in enumerate(hchains):
36  h.set_name(names[i])
37  assembly = IMP.atom.Molecule.setup_particle(temp)
38  return assembly
39 
40 
41 def create_assembly(model, fn_pdbs, names=False):
42  """ Read all the PDBs given in the list of names fn_pdbs and adds the
43  hierarchies to the model
44  """
46  for i, fn_pdb in enumerate(fn_pdbs):
47  if names:
48  prot = read_component(model, fn_pdb, names[i])
49  else:
50  prot = read_component(model, fn_pdb)
51  assembly.add_child(prot)
52  return assembly
53 
54 
55 def read_component(model, fn_pdb, name=False):
56  """ Read a PDB molecule, add atoms, and set a name
57  """
58  if name:
59  log.debug("reading component %s from %s", name, fn_pdb)
60  else:
61  log.debug("reading component from %s", fn_pdb)
62 
63  hierarchy = IMP.atom.read_pdb(fn_pdb, model,
65  if name:
66  hierarchy.set_name(name)
67  IMP.atom.add_radii(hierarchy)
68  return hierarchy
69 
70 
71 def create_rigid_bodies(assembly):
72  """ set the children of a molecule type hierarchy as rigid bodies
73  In this case, all the children are the components of the complex.
74  I use the function create_rigid_body(), that creates a lot of
75  sub-rigid bodies.
76 
77  I have changed the function and now build the rigid body directly from
78  the leaves of each of the components. With this I guarantee that the
79  number of rigid members is going to be the same if the components have
80  the same number of atoms.
81  """
82  if not IMP.atom.Molecule.get_is_setup(assembly):
83  raise TypeError("create_rigid_bodies(): The argument is not a valid "
84  "hierarchy")
85  molecule = IMP.atom.Molecule(assembly)
86  rbs = []
87  for c in molecule.get_children():
88  p = IMP.Particle(c.get_model())
90  rb = IMP.core.RigidBody(p)
91 # rb = IMP.atom.create_rigid_body(c)
92  rb.set_name(get_rb_name(c.get_name()))
93  rbs.append(rb)
94  return rbs
95 
96 
97 def rename_chains(assembly):
98  """ Rename all the chains of an assembly so there are no conflicts with
99  the ids. The names are added sequentially.
100  """
101  if not IMP.atom.Molecule.get_is_setup(assembly):
102  raise TypeError("The argument is not a valid hierarchy")
103  m = IMP.atom.Molecule(assembly)
104  all_chains_as_hierarchies = get_all_chains(m.get_children())
105  letters = string.ascii_uppercase
106  n_chains = len(all_chains_as_hierarchies)
107  if len(letters) < n_chains:
108  raise ValueError("There are more chains than letter ids")
109  ids = letters[0:n_chains]
110  for h, c_id in zip(all_chains_as_hierarchies, ids):
111  chain = IMP.atom.Chain(h.get_particle())
112  chain.set_id(c_id)
113  chain.set_name("chain %s" % c_id)
114 
115 
116 def create_simplified_dna(dna_hierarchy, n_res):
117  """ Gets a hierarchy containing a molecule of DNA and simplifies it,
118  generating a coarse representation of spheres. The function returns
119  a hierarchy with the spheres.
120  n_res - Number of residues to use per sphere.
121  """
122  if not IMP.atom.Chain.get_is_setup(dna_hierarchy):
123  raise TypeError("create_simplified_dna: the hierarchy provided is "
124  "not a chain.")
125 
126  model = dna_hierarchy.get_model()
127  ph = IMP.Particle(model)
128  simplified_h = IMP.atom.Hierarchy.setup_particle(ph)
130 
131  residues = IMP.atom.get_by_type(dna_hierarchy, IMP.atom.RESIDUE_TYPE)
132  lenr = len(residues)
133  for i in range(0, lenr, n_res):
134  xyzrs = []
135  equivalent_mass = 0.0
136  residues_numbers = []
137  for r in residues[i: i + n_res]:
138  rr = IMP.atom.Residue(r)
139  residues_numbers.append(rr.get_index())
140  # print "residue",rr.get_name(),rr.get_index()
141  residue_xyzrs = [IMP.core.XYZ(a.get_particle())
142  for a in rr.get_children()]
143  xyzrs += residue_xyzrs
144 # print "residue",r,"mass",get_residue_mass(r)
145  equivalent_mass += get_residue_mass(r)
146 
148  p = IMP.Particle(model)
150  xyzr.set_radius(s.get_radius())
151  xyzr.set_coordinates(s.get_center())
153  fragment.set_residue_indexes(residues_numbers)
154  IMP.atom.Mass.setup_particle(p, equivalent_mass)
155  simplified_h.add_child(fragment)
156  simplified_h.set_name("DNA")
157 # print "simplified_h is valid:",simplified_h.get_is_valid(True)
158  return simplified_h
159 
160 
161 def get_residue_mass(residue):
162  if not IMP.atom.Residue.get_is_setup(residue):
163  raise TypeError("The argument is not a residue")
164  r = IMP.atom.Residue(residue)
165  mass = 0.0
166  for leaf in IMP.atom.get_leaves(r):
167  ms = IMP.atom.Mass(leaf)
168  mass += ms.get_residue_mass()
169  return mass
170 
171 
172 def create_simplified_assembly(assembly, components_rbs, n_res):
173  """ Simplifies an assembly, by creating a hierarchy with one ball per
174  n_res residues. Each of the chains in the new hierarchy are added to
175  the rigid bodies for each of the components.
176  There must be correspondence between the children of the assembly
177  (components) and the rigid bodies. I check for the ids.
178  """
179  if not IMP.atom.Molecule.get_is_setup(assembly):
180  raise TypeError("The argument is not a valid hierarchy")
181  molecule = IMP.atom.Molecule(assembly)
182 
183  model = assembly.get_model()
184  n_children = molecule.get_number_of_children()
185 
186  sh = IMP.Particle(model)
187  simplified_hierarchy = IMP.atom.Molecule.setup_particle(sh)
188 
189  for i in range(n_children): # for all members of the assembly
190  component = molecule.get_child(i)
191  name = component.get_name()
192  rb = components_rbs[i]
193  if rb.get_name() != get_rb_name(name):
194  raise ValueError("Rigid body and component do not match")
195 
196  hchains = IMP.atom.get_by_type(component, IMP.atom.CHAIN_TYPE)
197  ch = IMP.Particle(model)
198  coarse_component_h = IMP.atom.Molecule.setup_particle(ch)
199  # simplify all the chains in the member
200  for h in hchains:
201  chain = IMP.atom.Chain(h.get_particle())
202  coarse_h = None
203  if name == "DNA":
204  coarse_h_particle = create_simplified_dna(h, n_res)
205  coarse_h = IMP.atom.Hierarchy(coarse_h_particle)
206  else:
208  n_res)
209 
210  # does not work for DNA
211  chain_rb = IMP.atom.create_rigid_body(coarse_h)
212  chain_rb.set_name("sub_rb" + name)
213  rb.add_member(chain_rb)
214 
215  # are required to have excluded volume
216  coarse_component_h.add_child(IMP.atom.Chain(coarse_h))
217  coarse_component_h.set_name(name)
218  simplified_hierarchy.add_child(coarse_component_h)
219  return simplified_hierarchy
220 
221 
222 def get_component(assembly, name):
223  """ Select a component of the assembly using the name """
224  for c in assembly.get_children():
225  if (c.get_name() == name):
226  return c
227  raise ValueError(
228  "The requested component %s is not in the assembly" %
229  name)
230 
231 
232 def get_rigid_body(rigid_bodies, name):
233  """ Select a rigid body from the rigid_bodies using the name """
234  for rb in rigid_bodies:
235  if (rb.get_name() == name):
236  return rb
237  raise ValueError("This rigid body is not in the set: %s" % name)
238 
239 
240 def get_rb_name(name):
241  """ Name to use for the rigid body of a hierarch"""
242  return "rb_" + name
243 
244 
246  """ Build the rigid body for all the particles in the selection S """
247  ps = S.get_selected_particles()
248  xyzrs = [IMP.core.XYZR(p) for p in ps]
249  p_rbS = IMP.Particle(model)
250  rbS = IMP.core.RigidBody.setup_particle(p_rbS, xyzrs)
251  return rbS
252 
253 
254 def get_selection_as_hierarchy(model, S):
255  ph = IMP.Particle(model)
257  for p in S.get_selected_particles():
259  h.add_child(x)
260  return h
261 
262 
264  """ Gets a selection of particles and decorates them as Atoms.
265  Then all of them are put into a big residue. I have this to use
266  with the multifit.create_coarse_molecule_from_molecule() function
267  """
268  ph = IMP.Particle(model)
270  for p in S.get_selected_particles():
271  h.add_child(IMP.atom.Atom(p))
272  return h
273 
274 
275 def get_coarse_selection(coarse_h, residues_numbers):
276  """ The function returns the particles (fragments) in the coarse hierarchy
277  that were created by summarizing the residues_numbers.
278 
279  Coarse hierarchy - Hierarchy formed by a bunch of
280  fragments. Each fragment must have the residue numbers that it contains
281  residue_numbers - list with the number of the residues that need to
282  be recovered.
283  The function returns the set of particles that are IMP.atom.Fragments
284  """
285  particles = []
286  fragments = IMP.atom.get_by_type(coarse_h, IMP.atom.FRAGMENT_TYPE)
287  for f in fragments:
288  ff = IMP.atom.Fragment(f)
289  residues_in_f = ff.get_residue_indexes()
290  for number in residues_in_f:
291  if number in residues_numbers:
292  particles.append(ff.get_particle())
293  break
294  return particles
295 
296 
298  """
299  Rotates the reference frame of a rigid body around the centroid
300  """
301  c = rb.get_coordinates()
302  ref = rb.get_reference_frame()
303  R = ref.get_transformation_to().get_rotation()
304  R2 = IMP.algebra.compose(rot, R)
307  rb.set_reference_frame(ref)
308 
309 
311  """
312  Applies a transformation around the centroid of a rigid body.
313  First does the rotation around the centroid and
314  then applies the transformation.
315  @param rb A IMP.core.RigidBody object
316  @param T a IMP.algebra.Transformation3D object
317  """
318  apply_rotation_around_centroid(rb, T.get_rotation())
319  rb.set_coordinates(rb.get_coordinates() + T.get_translation())
320 
321 
322 def get_residue_particle(h, chain_id=False, res=1):
323  """
324  Get the particle for a residue in a hierarchy.
325  @param h The hierarchy
326  @param chain_id If chain_id == False, just search for the residue
327  @param res Number of residue in the chain
328  """
329 # log.debug("get_residue_particle: chain_id %s, res %s",chain_id, res)
330  if chain_id:
331  s = IMP.atom.Selection(h, chain=chain_id, residue_index=res)
332  else:
333  s = IMP.atom.Selection(h, residue_index=res)
334  return s.get_selected_particles()[0]
335 
336 
337 def get_residue_coordinates(h, chain_id=False, res=1):
338  """
339  Get the coordinates of a residue (the coordinates of the first
340  particle)
341  @param h Hierarchy
342  @param chain_id See help for get_residue_particle()
343  @param res See help for get_residue_particle()
344  """
345  p = get_residue_particle(h, chain_id, res)
346  return IMP.core.XYZ(p).get_coordinates()
347 
348 
349 def get_residues_distance(hierarchy1, chain_id1, residue1,
350  hierarchy2, chain_id2, residue2):
351  """
352  Distance between two residues. See the help for get_residue_particle()
353  @param hierarchy1
354  @param chain_id1
355  @param residue1
356  @param hierarchy2
357  @param chain_id2
358  @param residue2
359  """
360  coords1 = get_residue_coordinates(hierarchy1, chain_id1, residue1)
361  coords2 = get_residue_coordinates(hierarchy2, chain_id2, residue2)
362  d = IMP.algebra.get_distance(coords1, coords2)
363  return d
364 
365 
366 def get_all_chains(hierarchies):
367  """ Gets all the chains in a set of hierarchies
368  @param hierarchies A set of IMP.atom.Hierarchy objects
369  """
370  chains = []
371  for h in hierarchies:
372  chains_in_h = IMP.atom.get_by_type(h, IMP.atom.CHAIN_TYPE)
373  for ch in chains_in_h:
374  chains.append(IMP.atom.Chain(ch))
375  return chains
376 
377 
378 def set_reference_frames(rigid_bodies, reference_frames):
379  for ref, rb in zip(reference_frames, rigid_bodies):
380  rb.set_reference_frame(ref)
381 
382 
383 def get_nucleic_acid_backbone(hierarchy, backbone='minimal'):
384  """
385  Returns the atoms in the backbone of the nucleic acid contained
386  in the hierarchy.
387  backbone 'minimal' returns the atoms:
388  ["P", "O5'", "C5'", "C4'", "C3'", "O3'"]
389  backbone 'trace' returns the atoms C4'
390  """
391 # log.debug("get_nucleic_acid_backbone")
392  backbone_atoms = []
393  if backbone == 'minimal':
394  backbone_atoms = ["P", "O5'", "C5'", "C4'", "C3'", "O3'"]
395  elif backbone == 'trace':
396  backbone_atoms = ["C4'"]
397  else:
398  raise ValueError("Wrong value for the type of backbone")
399  backbone_atom_types = [IMP.atom.AtomType(t) for t in backbone_atoms]
400  h_chains = IMP.atom.get_by_type(hierarchy, IMP.atom.CHAIN_TYPE)
401  backbone = []
402  if len(h_chains) > 1:
403  raise ValueError("The hierarchy mas more than one chain")
404  h_residues = IMP.atom.get_by_type(hierarchy, IMP.atom.RESIDUE_TYPE)
405  for hr in h_residues:
406  res = IMP.atom.Residue(hr)
407  if not (res.get_is_dna() or res.get_is_rna()):
408  raise ValueError("Residue is not part of a nucleic acid")
409  h_atoms = IMP.atom.get_by_type(hr, IMP.atom.ATOM_TYPE)
410  for at in h_atoms:
411  if IMP.atom.Atom(at).get_atom_type() in backbone_atom_types:
412  backbone.append(at)
413  return backbone
414 
415 
416 def get_calphas(chain_hierarchy):
417  h_residues = IMP.atom.get_by_type(chain_hierarchy, IMP.atom.RESIDUE_TYPE)
419  for r in h_residues]
420  return cas
421 
422 
423 def get_backbone(hierarchy):
424  """
425  Get the backbone atoms for a hierarchy. It can be a protein or a
426  nucleic acid
427  """
428  h_residues = IMP.atom.get_by_type(hierarchy, IMP.atom.RESIDUE_TYPE)
429  if len(h_residues) == 0:
430  raise ValueError("No residues!")
431  atoms = []
432  res = IMP.atom.Residue(h_residues[0])
433  if res.get_is_dna() or res.get_is_rna():
434  atoms = get_nucleic_acid_backbone(hierarchy, 'trace')
435  else:
436  atoms = get_calphas(hierarchy)
437  return atoms
438 
439 
440 def get_all_members(rigid_bodies):
441  """
442  Gets all the members of a set of rigid bodies, removing the subrigid
443  bodies. Returns all the plain atom or bead members
444  """
445 
446  members = []
447  for rb in rigid_bodies:
448  members += get_simple_members(rb)
449  return members
450 
451 
452 def get_simple_members(rb):
453  # Add members if they are not sub-rigid bodies
454  members = [m for m in rb.get_members()
455  if not IMP.core.RigidBody.get_is_setup(m.get_particle())]
456  return members
Select non water and non hydrogen atoms.
Definition: pdb.h:314
static bool get_is_setup(const IMP::ParticleAdaptor &p)
Definition: Molecule.h:35
Add mass to a particle.
Definition: Mass.h:23
Simple 3D transformation class.
def get_selection_as_atom_hierarchy
Gets a selection of particles and decorates them as Atoms.
static bool get_is_setup(const IMP::ParticleAdaptor &p)
Definition: Residue.h:158
A decorator to associate a particle with a part of a protein/DNA/RNA.
Definition: Fragment.h:20
def get_nucleic_acid_backbone
Returns the atoms in the backbone of the nucleic acid contained in the hierarchy. ...
static Fragment setup_particle(Model *m, ParticleIndex pi)
Definition: Fragment.h:67
static XYZR setup_particle(Model *m, ParticleIndex pi)
Definition: XYZR.h:48
def create_simplified_assembly
Simplifies an assembly, by creating a hierarchy with one ball per n_res residues. ...
def get_residue_particle
Get the particle for a residue in a hierarchy.
def create_rigid_bodies
set the children of a molecule type hierarchy as rigid bodies In this case, all the children are the ...
void add_radii(Hierarchy d, const ForceFieldParameters *ffp=get_all_atom_CHARMM_parameters(), FloatKey radius_key=FloatKey("radius"))
Add vdW radius from given force field.
def get_selection_rigid_body
Build the rigid body for all the particles in the selection S.
def get_component
Select a component of the assembly using the name.
def get_backbone
Get the backbone atoms for a hierarchy.
The type of an atom.
def apply_rotation_around_centroid
Rotates the reference frame of a rigid body around the centroid.
static Residue setup_particle(Model *m, ParticleIndex pi, ResidueType t, int index, int insertion_code)
Definition: Residue.h:160
Atom get_atom(Residue rd, AtomType at)
Return a particle atom from the residue.
def get_all_members
Gets all the members of a set of rigid bodies, removing the subrigid bodies.
def get_all_chains
Gets all the chains in a set of hierarchies.
void read_pdb(TextInput input, int model, Hierarchy h)
A reference frame in 3D.
def get_rb_name
Name to use for the rigid body of a hierarch.
static Molecule setup_particle(Model *m, ParticleIndex pi)
Definition: Molecule.h:37
def apply_transformation_around_centroid
Applies a transformation around the centroid of a rigid body.
def create_assembly_from_pdb
Builds the assembly setting the chains in the PDB file as components.
static Hierarchy setup_particle(Model *m, ParticleIndex pi, ParticleIndexesAdaptor children=ParticleIndexesAdaptor())
Create a Hierarchy of level t by adding the needed attributes.
The standard decorator for manipulating molecular structures.
A decorator for a particle representing an atom.
Definition: atom/Atom.h:238
static Mass setup_particle(Model *m, ParticleIndex pi, Float mass)
Definition: Mass.h:48
static Hierarchy setup_particle(Model *m, ParticleIndex pi, DecoratorTraits tr=get_default_decorator_traits())
def create_assembly
Read all the PDBs given in the list of names fn_pdbs and adds the hierarchies to the model...
def get_coarse_selection
The function returns the particles (fragments) in the coarse hierarchy that were created by summarizi...
A decorator for a particle with x,y,z coordinates.
Definition: XYZ.h:30
Transformation3D compose(const Transformation3D &a, const Transformation3D &b)
Compose two transformations.
static bool get_is_setup(const IMP::ParticleAdaptor &p)
Definition: rigid_bodies.h:169
def rename_chains
Rename all the chains of an assembly so there are no conflicts with the ids.
algebra::Sphere3D get_enclosing_sphere(const XYZs &v)
Get a sphere enclosing the set of XYZRs.
A decorator for a residue.
Definition: Residue.h:137
Basic functionality that is expected to be used by a wide variety of IMP users.
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::ParticleAdaptor &p)
Definition: Chain.h:82
def get_rigid_body
Select a rigid body from the rigid_bodies using the name.
Hierarchy create_simplified_along_backbone(Chain input, const IntRanges &residue_segments, bool keep_detailed=false)
IMP::core::RigidBody create_rigid_body(Hierarchy h)
def get_residue_coordinates
Get the coordinates of a residue (the coordinates of the first particle)
Class to handle individual particles of a Model object.
Definition: Particle.h:43
double get_distance(const VectorD< D > &v1, const VectorD< D > &v2)
Compute the distance between two vectors.
Definition: VectorD.h:196
Store info for a chain of a protein.
Definition: Chain.h:61
A decorator for a rigid body.
Definition: rigid_bodies.h:82
def create_simplified_dna
Gets a hierarchy containing a molecule of DNA and simplifies it, generating a coarse representation o...
Functionality for loading, creating, manipulating and scoring atomic structures.
def read_component
Read a PDB molecule, add atoms, and set a name.
static Chain setup_particle(Model *m, ParticleIndex pi, std::string id)
Definition: Chain.h:83
Hierarchies get_leaves(const Selection &h)
A decorator for a molecule.
Definition: Molecule.h:24
Select hierarchy particles identified by the biological name.
Definition: Selection.h:70
static RigidBody setup_particle(Model *m, ParticleIndex pi, ParticleIndexesAdaptor ps)
Definition: rigid_bodies.h:180
A decorator for a particle with x,y,z coordinates and a radius.
Definition: XYZR.h:27