6 import modeller.scripts
 
    7 import modeller.optimizers
 
   11 if not hasattr(modeller.terms, 
'EnergyTerm'):
 
   12     modeller.terms.EnergyTerm = modeller.terms.energy_term
 
   13     modeller.Selection = modeller.selection
 
   17     """Make a temporary directory that is deleted when the object is.""" 
   20         self.tmpdir = tempfile.mkdtemp()
 
   23         shutil.rmtree(self.tmpdir, ignore_errors=
True)
 
   27     """A Modeller restraint which evaluates an IMP scoring function. 
   28        This can be used to incorporate IMP Restraints into an existing 
   29        comparative modeling pipeline, or to use Modeller optimizers or 
   33     _physical_type = modeller.physical.absposition
 
   35     def __init__(self, particles, scoring_function=None):
 
   37            @param particles A list of the IMP atoms (as Particle objects), 
   38                             same order as the Modeller atoms. 
   39            @param scoring_function An IMP::ScoringFunction object that will 
   40                             be incorporated into the Modeller score (molpdf). 
   41            @note since Modeller, unlike IMP, is sensitive to the ordering 
   42                  of atoms, it usually makes sense to create the model in 
   43                  Modeller and then use ModelLoader to load it into IMP, 
   44                  since that will preserve the Modeller atom ordering in IMP. 
   46         modeller.terms.EnergyTerm.__init__(self)
 
   47         self._particles = particles
 
   49             self._sf = scoring_function
 
   51             self._sf = particles[0].get_model()
 
   53     def eval(self, mdl, deriv, indats):
 
   54         atoms = self.indices_to_atoms(mdl, indats)
 
   55         _copy_modeller_coords_to_imp(atoms, self._particles)
 
   56         if len(self._particles) == 0:
 
   59             score = self._sf.evaluate(deriv)
 
   61             dvx = [0.] * len(indats)
 
   62             dvy = [0.] * len(indats)
 
   63             dvz = [0.] * len(indats)
 
   64             _get_imp_derivs(self._particles, dvx, dvy, dvz)
 
   65             return (score, dvx, dvy, dvz)
 
   71     """An IMP restraint using all defined Modeller restraints. 
   72        This is useful if you want to use Modeller restraints with an IMP 
   73        optimizer, or in combination with IMP restraints. 
   75        @note Currently only the coordinates of the atoms are translated 
   76              between Modeller and IMP; thus, a Modeller restraint which 
   77              uses any other attribute (e.g. charge) will not react if 
   78              this attribute is changed by IMP. 
   81     def __init__(self, model, modeller_model, particles):
 
   83            @param model The IMP Model object. 
   84            @param modeller_model The Modeller model object. 
   85            @param particles A list of the IMP atoms (as Particle objects), 
   86                             in the same order as the Modeller atoms. 
   87            @note since Modeller, unlike IMP, is sensitive to the ordering 
   88                  of atoms, it usually makes sense to create the model in 
   89                  Modeller and then use ModelLoader to load it into IMP, 
   90                  since that will preserve the Modeller atom ordering in IMP. 
   93             if hasattr(x, 
'get_particle'):
 
   94                 return x.get_particle()
 
   97         IMP.Restraint.__init__(self, model, 
"ModellerRestraints %1%")
 
   98         self._modeller_model = modeller_model
 
   99         self._particles = [get_particle(x) 
for x 
in particles]
 
  102         atoms = self._modeller_model.atoms
 
  103         sel = modeller.Selection(self._modeller_model)
 
  104         _copy_imp_coords_to_modeller(self._particles, atoms)
 
  105         energies = sel.energy()
 
  107             _add_modeller_derivs_to_imp(atoms, self._particles, accum)
 
  114     def do_show(self, fh):
 
  115         fh.write(
"ModellerRestraints")
 
  118         return self._particles
 
  121 def _copy_imp_coords_to_modeller(particles, atoms):
 
  122     """Copy atom coordinates from IMP to Modeller""" 
  126     for (num, at) 
in enumerate(atoms):
 
  127         at.x = particles[num].get_value(xkey)
 
  128         at.y = particles[num].get_value(ykey)
 
  129         at.z = particles[num].get_value(zkey)
 
  132 def _copy_modeller_coords_to_imp(atoms, particles):
 
  133     """Copy atom coordinates from Modeller to IMP""" 
  137     for (num, at) 
in enumerate(atoms):
 
  138         particles[num].set_value(xkey, at.x)
 
  139         particles[num].set_value(ykey, at.y)
 
  140         particles[num].set_value(zkey, at.z)
 
  143 def _add_modeller_derivs_to_imp(atoms, particles, accum):
 
  144     """Add atom derivatives from Modeller to IMP""" 
  145     for (num, at) 
in enumerate(atoms):
 
  147         xyz.add_to_derivative(0, at.dvx, accum)
 
  148         xyz.add_to_derivative(1, at.dvy, accum)
 
  149         xyz.add_to_derivative(2, at.dvz, accum)
 
  152 def _get_imp_derivs(particles, dvx, dvy, dvz):
 
  153     """Move atom derivatives from IMP to Modeller""" 
  157     for idx 
in range(0, len(dvx)):
 
  158         dvx[idx] = particles[idx].get_derivative(xkey)
 
  159         dvy[idx] = particles[idx].get_derivative(ykey)
 
  160         dvz[idx] = particles[idx].get_derivative(zkey)
 
  164 def _HarmonicLowerBoundGenerator(parameters, modalities):
 
  165     (mean, stdev) = parameters
 
  170 def _HarmonicUpperBoundGenerator(parameters, modalities):
 
  171     (mean, stdev) = parameters
 
  176 def _HarmonicGenerator(parameters, modalities):
 
  177     (mean, stdev) = parameters
 
  182 def _CosineGenerator(parameters, modalities):
 
  183     (phase, force_constant) = parameters
 
  184     (periodicity,) = modalities
 
  188 def _LinearGenerator(parameters, modalities):
 
  189     (scale,) = parameters
 
  193 def _SplineGenerator(parameters, modalities):
 
  194     (open, low, high, delta, lowderiv, highderiv) = parameters[:6]
 
  196     for v 
in parameters[6:]:
 
  205 _unary_func_generators = {
 
  206     1: _HarmonicLowerBoundGenerator,
 
  207     2: _HarmonicUpperBoundGenerator,
 
  208     3: _HarmonicGenerator,
 
  211     10: _SplineGenerator,
 
  216 def _DistanceRestraintGenerator(form, modalities, atoms, parameters):
 
  217     unary_func_gen = _unary_func_generators[form]
 
  219                                       unary_func_gen(parameters, modalities),
 
  223 def _AngleRestraintGenerator(form, modalities, atoms, parameters):
 
  224     unary_func_gen = _unary_func_generators[form]
 
  226                                    unary_func_gen(parameters, modalities),
 
  227                                    atoms[0], atoms[1], atoms[2])
 
  230 def _MultiBinormalGenerator(form, modalities, atoms, parameters):
 
  231     nterms = modalities[0]
 
  232     if len(parameters) != nterms * 6:
 
  233         raise ValueError(
"Incorrect number of parameters (%d) for multiple " 
  234                          "binormal restraint - expecting %d (%d terms * 6)" 
  235                          % (len(parameters), nterms * 6, nterms))
 
  237                                            atoms[:4], atoms[4:8])
 
  238     for i 
in range(nterms):
 
  240         t.set_weight(parameters[i])
 
  241         t.set_means((parameters[nterms + i * 2],
 
  242                      parameters[nterms + i * 2 + 1]))
 
  243         t.set_standard_deviations((parameters[nterms * 3 + i * 2],
 
  244                                    parameters[nterms * 3 + i * 2 + 1]))
 
  245         t.set_correlation(parameters[nterms * 5 + i])
 
  250 def _DihedralRestraintGenerator(form, modalities, atoms, parameters):
 
  252         return _MultiBinormalGenerator(form, modalities, atoms, parameters)
 
  253     unary_func_gen = _unary_func_generators[form]
 
  255                                       unary_func_gen(parameters, modalities),
 
  256                                       atoms[0], atoms[1], atoms[2], atoms[3])
 
  259 def _get_protein_atom_particles(protein):
 
  260     """Given a protein particle, get the flattened list of all child atoms""" 
  262     for ichain 
in range(protein.get_number_of_children()):
 
  263         chain = protein.get_child(ichain)
 
  264         for ires 
in range(chain.get_number_of_children()):
 
  265             residue = chain.get_child(ires)
 
  266             for iatom 
in range(residue.get_number_of_children()):
 
  267                 atom = residue.get_child(iatom)
 
  268                 atom_particles.append(atom.get_particle())
 
  269     return atom_particles
 
  272 def _load_restraints_line(line, atom_particles):
 
  273     """Parse a single Modeller restraints file line and return the 
  274        corresponding IMP restraint.""" 
  277     if typ == 
'MODELLER5':
 
  280         raise NotImplementedError(
"Only 'R' lines currently read from " + 
  281                                   "Modeller restraints files")
 
  282     form = int(spl.pop(0))
 
  283     modalities = [int(spl.pop(0))]
 
  284     features = [int(spl.pop(0))]
 
  287     natoms = [int(spl.pop(0))]
 
  288     nparam = int(spl.pop(0))
 
  289     nfeat = int(spl.pop(0))
 
  290     for i 
in range(nfeat - 1):
 
  291         modalities.append(int(spl.pop(0)))
 
  292         features.append(int(spl.pop(0)))
 
  293         natoms.append(int(spl.pop(0)))
 
  294     atoms = [int(spl.pop(0)) 
for x 
in range(natoms[0])]
 
  295     for i 
in range(len(atoms)):
 
  296         atoms[i] = atom_particles[atoms[i] - 1]
 
  297     parameters = [float(spl.pop(0)) 
for x 
in range(nparam)]
 
  298     restraint_generators = {
 
  299         1: _DistanceRestraintGenerator,
 
  300         2: _AngleRestraintGenerator,
 
  301         3: _DihedralRestraintGenerator,
 
  302         4: _DihedralRestraintGenerator,
 
  304     restraint_gen = restraint_generators[features[0]]
 
  305     return restraint_gen(form, modalities, atoms, parameters)
 
  308 def _load_entire_restraints_file(filename, protein):
 
  309     """Yield a set of IMP restraints from a Modeller restraints file.""" 
  310     atoms = _get_protein_atom_particles(protein)
 
  311     with open(filename, 
'r') as fh: 
  314                 rsr = _load_restraints_line(line, atoms)
 
  318                 print(
"Cannot read restraints file line:\n" + line)
 
  322 def _copy_residue(r, model):
 
  323     """Copy residue information from modeller to imp""" 
  327     p.set_name(str(
"residue "+r.num))
 
  331 def _copy_atom(a, model):
 
  332     """Copy atom information from modeller""" 
  337     if hasattr(a, 
'charge'):
 
  339     if hasattr(a, 
'type'):
 
  341     ap.set_input_index(a.index)
 
  345 def _copy_chain(c, model):
 
  346     """Copy chain information from modeller""" 
  353 def _get_forcefield(submodel):
 
  366     """Add radii to the hierarchy using the Modeller radius library, radii.lib. 
  367        Each radius is scaled by the given scale (Modeller usually scales radii 
  368        by a factor of 0.82). submodel specifies the topology submodel, which is 
  369        the column in radii.lib to use.""" 
  373     with open(filename) 
as fh:
 
  375             if line.startswith(
'#'):
 
  379                 radii[spl[0]] = float(spl[submodel])
 
  380     atoms = IMP.atom.get_by_type(hierarchy, IMP.atom.ATOM_TYPE)
 
  385             radius = radii[ct] * scale
 
  393     """Read a Modeller model into IMP. After creating this object, the atoms 
  394        in the Modeller model can be loaded into IMP using the load_atoms() 
  395        method, then optionally any Modeller static restraints can be read in 
  396        with load_static_restraints() or load_static_restraints_file(). 
  398        This class can also be used to read Modeller alignment structures; 
  399        however, only load_atoms() will be useful in such a case (since 
  400        alignment structures don't have restraints or other information). 
  406            @param modeller_model The Modeller model or alignment structure 
  409         self._modeller_model = modeller_model
 
  412         """Construct an IMP::atom::Hierarchy that contains the same atoms as 
  413            the Modeller model or alignment structure. 
  415            IMP atoms created from a Modeller model will be given charges and 
  416            CHARMM types, extracted from the model. Alignment structures don't 
  417            contain this information, so the IMP atoms won't either. 
  419            @param model The IMP::Model object in which the hierarchy will be 
  420                         created. The highest level hierarchy node is a PROTEIN. 
  421            @return the newly-created root IMP::atom::Hierarchy. 
  426         for chain 
in self._modeller_model.chains:
 
  431             for residue 
in chain.residues:
 
  432                 rp = _copy_residue(residue, model)
 
  435                 for atom 
in residue.atoms:
 
  436                     ap = _copy_atom(atom, model)
 
  439                     self._atoms[atom.index] = ap
 
  440         self._modeller_hierarchy = hpp
 
  443     def _get_nonbonded_list(self, atoms, pair_filter, edat, distance):
 
  448         if pair_filter 
is None:
 
  450             if edat.excl_local[0]:
 
  451                 pair_filter.set_bonds(list(self.
load_bonds()))
 
  452             if edat.excl_local[1]:
 
  454             if edat.excl_local[2]:
 
  456         nbl.add_pair_filter(pair_filter)
 
  460         """Load the Modeller bond topology into the IMP model. Each bond is 
  461            represented in IMP as an IMP::atom::Bond, with no defined length 
  462            or stiffness. These bonds are primarily useful as input to 
  463            IMP::atom::StereochemistryPairFilter, to exclude bond interactions 
  464            from the nonbonded list. Typically the contribution to the scoring 
  465            function from the bonds is included in the Modeller static 
  466            restraints (use load_static_restraints() or 
  467            load_static_restraints_file() to load these). If you want to 
  468            regenerate the stereochemistry in IMP, do not use these functions 
  469            (as then stereochemistry scoring terms and exclusions would be 
  470            double-counted) and instead use the IMP::atom::CHARMMTopology class. 
  472            You must call load_atoms() prior to using this function. 
  473            @see load_angles(), load_dihedrals(), load_impropers() 
  474            @return A generator listing all of the bonds. 
  476         if not hasattr(self, 
'_modeller_hierarchy'):
 
  477             raise ValueError(
"Call load_atoms() first.")
 
  478         for (maa, mab) 
in self._modeller_model.bonds:
 
  479             pa = self._atoms[maa.index]
 
  480             pb = self._atoms[mab.index]
 
  490                                        IMP.atom.Bond.SINGLE).get_particle()
 
  493         """Load the Modeller angle topology into the IMP model. 
  494            See load_bonds() for more details.""" 
  495         return self._internal_load_angles(self._modeller_model.angles,
 
  499         """Load the Modeller dihedral topology into the IMP model. 
  500            See load_bonds() for more details.""" 
  501         return self._internal_load_angles(self._modeller_model.dihedrals,
 
  505         """Load the Modeller improper topology into the IMP model. 
  506            See load_bonds() for more details.""" 
  507         return self._internal_load_angles(self._modeller_model.impropers,
 
  510     def _internal_load_angles(self, angles, angle_class):
 
  511         if not hasattr(self, 
'_modeller_hierarchy'):
 
  512             raise ValueError(
"Call load_atoms() first.")
 
  513         for modeller_atoms 
in angles:
 
  514             imp_particles = [self._atoms[x.index] 
for x 
in modeller_atoms]
 
  516             a = angle_class.setup_particle(
 
  518             yield a.get_particle()
 
  521         """Convert a Modeller static restraints file into equivalent 
  522            IMP::Restraints. load_atoms() must have been called first to read 
  523            in the atoms that the restraints will act upon. 
  524            @param filename Name of the Modeller restraints file. The restraints 
  525                   in this file are assumed to act upon the model read in by 
  526                   load_atoms(); no checking is done to enforce this. 
  527            @return A Python generator of the newly-created IMP::Restraint 
  530         if not hasattr(self, 
'_modeller_hierarchy'):
 
  531             raise ValueError(
"Call load_atoms() first.")
 
  532         return _load_entire_restraints_file(filename, self._modeller_hierarchy)
 
  535         """Convert the current set of Modeller static restraints into 
  536            equivalent IMP::Restraints. load_atoms() must have been called 
  537            first to read in the atoms that the restraints will act upon. 
  538            @return A Python generator of the newly-created IMP::Restraint 
  541         class _RestraintGenerator:
 
  542             """Simple generator wrapper""" 
  546             def __iter__(self, *args, **keys):
 
  549             def close(self, *args, **keys):
 
  550                 return self._gen.close(*args, **keys)
 
  553                 return next(self._gen)
 
  557             def send(self, *args, **keys):
 
  558                 return self._gen.send(*args, **keys)
 
  560             def throw(self, *args, **keys):
 
  561                 return self._gen.throw(*args, **keys)
 
  565         rsrfile = os.path.join(t.tmpdir, 
'restraints.rsr')
 
  566         self._modeller_model.restraints.write(file=rsrfile)
 
  569         wrap = _RestraintGenerator(gen)
 
  575         """Convert Modeller dynamic restraints into IMP::Restraint objects. 
  577            For each currently active Modeller dynamic restraint 
  578            (e.g. soft-sphere, electrostatics) an equivalent IMP::Restraint 
  580            load_atoms() must have been called first to read 
  581            in the atoms that the restraints will act upon. 
  583            If pair_filter is given, it is an IMP::PairFilter object to exclude 
  584            pairs from the nonbonded lists used by the dynamic restraints. 
  585            Otherwise, an IMP::atom::StereochemistryPairFilter object is created 
  586            to exclude Modeller bonds, angles and dihedrals, as specified by 
  587            edat.excl_local. (Note that this calls load_bonds(), load_angles() 
  588            and load_dihedrals(), so will create duplicate lists of bonds if 
  589            those methods are called manually as well.) 
  591            @note Currently only soft-sphere, electrostatic and Lennard-Jones 
  592                  restraints are loaded. 
  593            @return A Python generator of the newly-created IMP::Restraint 
  596         if not hasattr(self, 
'_modeller_hierarchy'):
 
  597             raise ValueError(
"Call load_atoms() first.")
 
  598         edat = self._modeller_model.env.edat
 
  599         libs = self._modeller_model.env.libs
 
  601         m = atoms[0].get_model()
 
  604         if edat.dynamic_sphere:
 
  607             nbl = self._get_nonbonded_list(atoms, pair_filter, edat, 0.)
 
  610                                   libs.topology.submodel, edat.radii_factor)
 
  617         if edat.dynamic_lennard 
or edat.dynamic_coulomb:
 
  619             d = max(edat.contact_shell - 3.0, 0.0)
 
  620             nbl = self._get_nonbonded_list(atoms, pair_filter, edat, d)
 
  621             ff = _get_forcefield(libs.topology.submodel)
 
  622             ff.add_radii(self._modeller_hierarchy)
 
  624             if edat.dynamic_lennard:
 
  625                 ff.add_lennard_jones_types(self._modeller_hierarchy)
 
  627                                           edat.lennard_jones_switch[1])
 
  631             if edat.dynamic_coulomb:
 
  633                                           edat.coulomb_switch[1])
 
  635                 ps.set_relative_dielectric(edat.relative_dielectric)
 
  639 __version__ = 
"2.23.0" 
  642     '''Return the version of this module, as a string''' 
  646     '''Return the fully-qualified name of this module''' 
  647     return "IMP::modeller" 
  650     '''Return the full path to one of this module's data files''' 
  652     return IMP._get_module_data_path(
"modeller", fname)
 
  655     '''Return the full path to one of this module's example files''' 
  657     return IMP._get_module_example_path(
"modeller", fname)
 
static CHARMMAtom setup_particle(Model *m, ParticleIndex pi, String charmm_type)
Create a decorator with the passed CHARMM type. 
 
static Charged setup_particle(Model *m, ParticleIndex pi, Float charge)
 
Lower bound harmonic function (non-zero when feature < mean) 
 
def load_static_restraints
Convert the current set of Modeller static restraints into equivalent IMP::Restraints. 
 
def get_module_version
Return the version of this module, as a string. 
 
def load_bonds
Load the Modeller bond topology into the IMP model. 
 
A Modeller restraint which evaluates an IMP scoring function. 
 
Coulomb (electrostatic) score between a pair of particles. 
 
static Atom setup_particle(Model *m, ParticleIndex pi, Atom other)
 
Various classes to hold sets of particles. 
 
Upper bound harmonic function (non-zero when feature > mean) 
 
static XYZR setup_particle(Model *m, ParticleIndex pi)
 
A decorator for a particle which has bonds. 
 
def load_impropers
Load the Modeller improper topology into the IMP model. 
 
virtual double unprotected_evaluate(DerivativeAccumulator *da) const 
Return the unweighted score for the restraint. 
 
def get_example_path
Return the full path to one of this module's example files. 
 
Dihedral restraint between four particles. 
 
def load_static_restraints_file
Convert a Modeller static restraints file into equivalent IMP::Restraints. 
 
A score on the distance between the surfaces of two spheres. 
 
Return all close unordered pairs of particles taken from the SingletonContainer. 
 
static Residue setup_particle(Model *m, ParticleIndex pi, ResidueType t, int index, int insertion_code)
 
A single binormal term in a MultipleBinormalRestraint. 
 
def load_dynamic_restraints
Convert Modeller dynamic restraints into IMP::Restraint objects. 
 
Distance restraint between two particles. 
 
static XYZ setup_particle(Model *m, ParticleIndex pi)
 
An IMP restraint using all defined Modeller restraints. 
 
static bool get_is_setup(const IMP::ParticleAdaptor &p)
 
CHARMM force field parameters. 
 
Bond create_bond(Bonded a, Bonded b, Bond o)
Connect the two wrapped particles by a custom bond. 
 
def load_dihedrals
Load the Modeller dihedral topology into the IMP model. 
 
static Float get_k_from_standard_deviation(Float sd, Float t=297.15)
Return the k to use for a given Gaussian standard deviation. 
 
def get_data_path
Return the full path to one of this module's data files. 
 
Angle restraint between three particles. 
 
static Hierarchy setup_particle(Model *m, ParticleIndex pi, ParticleIndexesAdaptor children=ParticleIndexesAdaptor())
Create a Hierarchy of level t by adding the needed attributes. 
 
ParticleIndexPairs get_indexes(const ParticlePairsTemp &ps)
Get the indexes from a list of particle pairs. 
 
A particle that describes an angle between three particles. 
 
The standard decorator for manipulating molecular structures. 
 
Store a list of ParticleIndexes. 
 
Lennard-Jones score between a pair of particles. 
 
static Bonded setup_particle(Model *m, ParticleIndex pi)
 
def load_atoms
Construct an IMP::atom::Hierarchy that contains the same atoms as the Modeller model or alignment str...
 
Version and module information for Objects. 
 
A decorator for a particle with x,y,z coordinates. 
 
Modeller-style multiple binormal (phi/psi) restraint. 
 
def add_soft_sphere_radii
Add radii to the hierarchy using the Modeller radius library, radii.lib. 
 
A particle that describes a dihedral angle between four particles. 
 
def load_angles
Load the Modeller angle topology into the IMP model. 
 
def get_module_name
Return the fully-qualified name of this module. 
 
static bool get_is_setup(const IMP::ParticleAdaptor &p)
 
virtual VersionInfo get_version_info() const 
Get information about the module and version of the object. 
 
A filter that excludes bonds, angles and dihedrals. 
 
Read a Modeller model into IMP. 
 
Class to handle individual particles of a Model object. 
 
std::string get_data_path(std::string file_name)
Return the full path to one of this module's data files. 
 
Smooth interaction scores by switching the derivatives (force switch). 
 
Functionality for loading, creating, manipulating and scoring atomic structures. 
 
static Chain setup_particle(Model *m, ParticleIndex pi, std::string id)
 
Hierarchies get_leaves(const Selection &h)
 
Applies a PairScore to each Pair in a list. 
 
virtual ModelObjectsTemp do_get_inputs() const =0
 
Closed cubic spline function. 
 
A decorator for an atom that has a defined CHARMM type. 
 
A restraint is a term in an IMP ScoringFunction. 
 
Harmonic function (symmetric about the mean) 
 
A decorator for a particle with x,y,z coordinates and a radius.