IMP logo
IMP Reference Guide  develop.63b38c487d,2024/12/22
The Integrative Modeling Platform
nude_modeling.py
1 ##########################################################################################
2 ######################### IMP Modeling Script for NuDe Complex ###########################
3 ##########################################################################################
4 
5 
6 # Imports
7 from __future__ import print_function
8 import IMP
9 import RMF
10 import IMP.rmf
11 import IMP.pmi
12 import IMP.pmi.io
14 import IMP.pmi.tools
15 import IMP.pmi.topology
16 import IMP.pmi.macros
17 import IMP.pmi.restraints
20 import ihm.cross_linkers
23 import IMP.pmi.dof
24 from IMP.nestor import NestedSampling
25 import IMP.atom
26 import contextlib
27 
28 # import IMP.saxs
29 import os
30 import sys
31 
32 IMP.setup_from_argv(sys.argv, "Application of NestOR to the NuDe subcomplex")
33 
34 dat_dir = IMP.nestor.get_example_path("input")
35 
36 if not os.path.exists(os.path.join(dat_dir, "gmm/emd_22904.mrc")):
37  print("To run this example, first download the EM map from EMD22094,")
38  print("extract it, rename it as `emd_22904.mrc` and place it in")
39  print("the `input/gmm/` directory.")
40  sys.exit(0)
41 
42 if len(sys.argv) < 2:
43  print(
44  "Command-line arguments needed to run this script. See the readme for more details"
45  )
46  sys.exit(0)
47 
48 run_output_dir = "run_" + sys.argv[1]
49 topology_file = os.path.join(dat_dir, sys.argv[2])
50 h_param_file = os.path.join(dat_dir, sys.argv[3])
51 
52 max_shuffle_core = 5
53 max_shuffle_set2 = 50
54 rex_max_temp = 2.4
55 
56 # Identify data files
57 sampling_adh_xl_data = os.path.join(dat_dir, "xlms/sampling_filtered_adh.dat")
58 sampling_bs3dss_xl_data = os.path.join(dat_dir, "xlms/sampling_filtered_bs3dss.dat")
59 sampling_dmtmm_xl_data = os.path.join(dat_dir, "xlms/sampling_filtered_dmtmm.dat")
60 
61 evi_adh_xl_data = os.path.join(dat_dir, "xlms/evicalc_filtered_adh.dat")
62 evi_bs3dss_xl_data = os.path.join(dat_dir, "xlms/evicalc_filtered_adh.dat")
63 evi_dmtmm_xl_data = os.path.join(dat_dir, "xlms/evicalc_filtered_adh.dat")
64 
65 gmm_data = os.path.join(dat_dir, "gmm/emd_22904.gmm.40.txt")
66 # Restraint weights
67 intra_xl_weight = 1.0
68 inter_xl_weight = 10.0
69 xl_weight = 10
70 em_weight = 1000.0
71 
72 
73 def modeling(output_dir, topology_file, h_param_file):
74  mdl = IMP.Model()
76  topology_file,
77  pdb_dir=os.path.join(dat_dir, "pdb"),
78  fasta_dir=os.path.join(dat_dir, "fasta"),
79  gmm_dir=os.path.join(dat_dir, "gmm"),
80  )
82  bs.add_state(t)
83 
84  root_hier, dof = bs.execute_macro(
85  max_rb_trans=1,
86  max_rb_rot=0.1,
87  max_bead_trans=3.2,
88  max_srb_trans=0.01,
89  max_srb_rot=0.04,
90  )
91 
92  # Fixing particles
93  set1_core = []
94  for prot in ["MTA1"]:
95  for cp in [0, 1]:
96  set1_core += IMP.atom.Selection(
97  root_hier, molecule=prot, copy_index=cp, residue_indexes=range(165, 334)
98  ).get_selected_particles()
99  for prot in ["HDAC1"]:
100  for cp in [0, 1]:
101  set1_core += IMP.atom.Selection(
102  root_hier, molecule=prot, copy_index=cp, residue_indexes=range(8, 377)
103  ).get_selected_particles()
104 
105  set2 = []
106  for prot in ["MTA1"]:
107  for cp in [0, 1]:
108  set2 += IMP.atom.Selection(
109  root_hier, molecule=prot, copy_index=cp, residue_indexes=range(1, 165)
110  ).get_selected_particles()
111  for prot in ["MTA1"]:
112  for cp in [0, 1]:
113  set2 += IMP.atom.Selection(
114  root_hier, molecule=prot, copy_index=cp, residue_indexes=range(334, 432)
115  ).get_selected_particles()
116  for prot in ["HDAC1"]:
117  for cp in [0, 1]:
118  set2 += IMP.atom.Selection(
119  root_hier, molecule=prot, copy_index=cp, residue_indexes=range(1, 8)
120  ).get_selected_particles()
121  for prot in ["HDAC1"]:
122  for cp in [0, 1]:
123  set2 += IMP.atom.Selection(
124  root_hier, molecule=prot, copy_index=cp, residue_indexes=range(377, 483)
125  ).get_selected_particles()
126  for prot in ["MBD3"]:
127  set2 += IMP.atom.Selection(
128  root_hier, molecule=prot, copy_index=0, residue_indexes=range(1, 296)
129  ).get_selected_particles()
130  for prot in ["P66A"]:
131  set2 += IMP.atom.Selection(
132  root_hier, molecule=prot, copy_index=0, residue_indexes=range(136, 179)
133  ).get_selected_particles()
134  for prot in ["RBBP4"]:
135  for cp in [0, 1, 2, 3]:
136  set2 += IMP.atom.Selection(
137  root_hier, molecule=prot, copy_index=cp, residue_indexes=range(1, 425)
138  ).get_selected_particles()
139 
140  fixed_set1_core_beads, fixed_set1_core = dof.disable_movers(
142  )
143 
144  fixed_set2_beads, fixed_set2 = dof.disable_movers(
146  )
147  molecules = t.get_components()
148 
149  #####################################################
150  ##################### RESTRAINTS ####################
151  #####################################################
152 
153  output_objects = []
154  nestor_restraints = []
155 
156  # -----------------------------
157  # CONNECTIVITY RESTRAINT
158 
159  for m in root_hier.get_children()[0].get_children():
161  cr.add_to_model()
162  output_objects.append(cr)
163 
164  print("Connectivity restraint applied")
165 
166  # -----------------------------
167  # EXCLUDED VOLUME RESTRAINT
168 
170  included_objects=[root_hier], resolution=1000
171  )
172  output_objects.append(evr)
173 
174  print("Excluded volume restraint applied")
175 
176  # -------------------------
177  # CROSSLINKING RESTRAINT
178  # Sampling XL restraints
179 
181  xldbkc.set_standard_keys()
182 
183  xldb_adh_sampling = IMP.pmi.io.crosslink.CrossLinkDataBase()
184  xldb_adh_sampling.create_set_from_file(
185  file_name=sampling_adh_xl_data, converter=xldbkc
186  )
187  xlr_adh_sampling = (
189  root_hier=root_hier,
190  database=xldb_adh_sampling,
191  length=25,
192  resolution=1,
193  slope=0.0001,
194  label="adh",
195  weight=xl_weight,
196  linker=ihm.cross_linkers.edc,
197  )
198  )
199 
200  xldb_bs3dss_sampling = IMP.pmi.io.crosslink.CrossLinkDataBase()
201  xldb_bs3dss_sampling.create_set_from_file(
202  file_name=sampling_bs3dss_xl_data, converter=xldbkc
203  )
204  xlr_bs3dss_sampling = (
206  root_hier=root_hier,
207  database=xldb_bs3dss_sampling,
208  length=25,
209  resolution=1,
210  slope=0.0001,
211  label="bs3dss",
212  weight=xl_weight,
213  linker=ihm.cross_linkers.bs3,
214  )
215  )
216 
217  xldb_dmtmm_sampling = IMP.pmi.io.crosslink.CrossLinkDataBase()
218  xldb_dmtmm_sampling.create_set_from_file(
219  file_name=sampling_dmtmm_xl_data, converter=xldbkc
220  )
221  xlr_dmtmm_sampling = (
223  root_hier=root_hier,
224  database=xldb_dmtmm_sampling,
225  length=16,
226  resolution=1,
227  slope=0.0001,
228  label="dmtmm",
229  weight=xl_weight,
230  linker=ihm.cross_linkers.dsso,
231  )
232  )
233 
234  output_objects.append(xlr_adh_sampling)
235  output_objects.append(xlr_bs3dss_sampling)
236  output_objects.append(xlr_dmtmm_sampling)
237 
238  # ----------------------------------------------------
239  # Nesting XL restraints
240 
242  xldbkc.set_standard_keys()
243 
244  xldb_adh_evicalc = IMP.pmi.io.crosslink.CrossLinkDataBase()
245  xldb_adh_evicalc.create_set_from_file(file_name=evi_adh_xl_data, converter=xldbkc)
246  xlr_adh_evicalc = (
248  root_hier=root_hier,
249  database=xldb_adh_evicalc,
250  length=25,
251  resolution=1,
252  slope=0.0001,
253  label="adh",
254  weight=0,
255  linker=ihm.cross_linkers.edc,
256  )
257  )
258 
259  xldb_bs3dss_evicalc = IMP.pmi.io.crosslink.CrossLinkDataBase()
260  xldb_bs3dss_evicalc.create_set_from_file(
261  file_name=evi_bs3dss_xl_data, converter=xldbkc
262  )
263  xlr_bs3dss_evicalc = (
265  root_hier=root_hier,
266  database=xldb_bs3dss_evicalc,
267  length=25,
268  resolution=1,
269  slope=0.0001,
270  label="bs3dss",
271  weight=0,
272  linker=ihm.cross_linkers.bs3,
273  )
274  )
275 
276  xldb_dmtmm_evicalc = IMP.pmi.io.crosslink.CrossLinkDataBase()
277  xldb_dmtmm_evicalc.create_set_from_file(
278  file_name=evi_dmtmm_xl_data, converter=xldbkc
279  )
280  xlr_dmtmm_evicalc = (
282  root_hier=root_hier,
283  database=xldb_dmtmm_evicalc,
284  length=16,
285  resolution=1,
286  slope=0.0001,
287  label="dmtmm",
288  weight=0,
289  linker=ihm.cross_linkers.dsso,
290  )
291  )
292 
293  output_objects.append(xlr_adh_evicalc)
294  output_objects.append(xlr_bs3dss_evicalc)
295  output_objects.append(xlr_dmtmm_evicalc)
296 
297  nestor_restraints.append(xlr_adh_evicalc)
298  nestor_restraints.append(xlr_bs3dss_evicalc)
299  nestor_restraints.append(xlr_dmtmm_evicalc)
300 
301  print("Cross-linking restraint applied")
302 
303  # EM RESTRAINT
304  densities = IMP.atom.Selection(
305  root_hier, representation_type=IMP.atom.DENSITIES
306  ).get_selected_particles()
307 
309  densities,
310  target_fn=gmm_data,
311  slope=0.00000001,
312  scale_target_to_mass=True,
313  weight=0,
314  )
315 
316  output_objects.append(emr)
317  nestor_restraints.append(emr)
318  print("EM Restraint Applied")
319 
320  #####################################################
321  ###################### SAMPLING #####################
322  #####################################################
323 
324  exit_code = None
325  try:
327  root_hier,
328  max_translation=max_shuffle_set2,
329  excluded_rigid_bodies=fixed_set1_core,
330  )
331 
333  root_hier,
334  max_translation=max_shuffle_core,
335  excluded_rigid_bodies=fixed_set2,
336  )
337 
338  except ValueError:
339  exit_code = 11
340  with open("shuffle_config.err", "w") as shufferr:
341  shufferr.write(str(11))
342  shufferr.flush()
343 
344  print(f"{'-'*50}\nExit code: {exit_code}\n{'-'*50}")
345  dof.optimize_flexible_beads(50)
346 
348 
349  # Now, add all of the other restraints to the scoring function to start sampling
350  evr.add_to_model()
351  emr.add_to_model()
352  xlr_adh_sampling.add_to_model()
353  xlr_bs3dss_sampling.add_to_model()
354  xlr_dmtmm_sampling.add_to_model()
355  xlr_adh_evicalc.add_to_model()
356  xlr_bs3dss_evicalc.add_to_model()
357  xlr_dmtmm_evicalc.add_to_model()
358 
359  print("Replica Exchange Maximum Temperature : " + str(rex_max_temp))
360 
362  mdl,
363  root_hier=root_hier,
364  monte_carlo_temperature=1.0,
365  replica_exchange_minimum_temperature=1.0,
366  replica_exchange_maximum_temperature=rex_max_temp,
367  monte_carlo_sample_objects=dof.get_movers(),
368  global_output_directory=run_output_dir,
369  output_objects=output_objects,
370  monte_carlo_steps=10,
371  number_of_best_scoring_models=0,
372  number_of_frames=0,
373  use_nestor=True,
374  )
375 
376  ns = NestedSampling(
377  rex_macro=rex,
378  nestor_restraints=nestor_restraints,
379  h_param_file=h_param_file,
380  exit_code=exit_code,
381  )
382 
383  ns_return, ns_exit_code = ns.execute_nested_sampling2()
384  return ns_return, ns_exit_code
385 
386 
387 with open("errors.log", "w") as errf:
388  with contextlib.redirect_stdout(None):
389  with contextlib.redirect_stderr(errf):
390  nested_sampling_output, nested_sampling_exitcode = modeling(
391  run_output_dir, topology_file, h_param_file
392  )
393 
394 
395 if nested_sampling_output is not None:
396  print(f"////{nested_sampling_output}")
397 
398 sys.exit(nested_sampling_exitcode)
def get_example_path
Return the full path to one of this module's example files.
Setup cross-link distance restraints from mass spectrometry data.
Definition: crosslinking.py:22
Restraints for keeping correct stereochemistry.
Fit Gaussian-decorated particles to an EM map (also represented with a set of Gaussians) ...
def shuffle_configuration
Shuffle particles.
Definition: tools.py:1242
A macro to help setup and run replica exchange.
Definition: macros.py:65
Set of Python classes to create a multi-state, multi-resolution IMP hierarchy.
Strings setup_from_argv(const Strings &argv, std::string description, std::string positional_description, int num_positional)
Some miscellaneous simple restraints.
Definition: basic.py:1
Utility classes and functions for reading and storing PMI files.
def enable_all_movers
Re-enable all movers: previously fixed particles will be released.
Miscellaneous utilities.
Definition: tools.py:1
Modify the transformation of a rigid body.
A macro to build a IMP::pmi::topology::System based on a TopologyReader object.
Definition: macros.py:630
Protocols for sampling structures and analyzing them.
Definition: macros.py:1
Class for storing model, its restraints, constraints, and particles.
Definition: Model.h:86
Classes to handle different kinds of restraints.
Restraints for handling electron microscopy maps.
Modify the transformation of a rigid body.
Nested sampling-based optimization of representation.
Create a restraint between consecutive TempResidue objects or an entire PMI Molecule object...
Create movers and set up constraints for PMI objects.
A class to create an excluded volume restraint for a set of particles at a given resolution.
Automatically setup System and Degrees of Freedom with a formatted text file.
Restraints for handling cross-linking data.
Definition: crosslinking.py:1
Python classes to represent, score, sample and analyze models.
Functionality for loading, creating, manipulating and scoring atomic structures.
Select hierarchy particles identified by the biological name.
Definition: Selection.h:70
Support for the RMF file format for storing hierarchical molecular data and markup.