[Date Prev][Date Next][Thread Prev][Thread Next][Date Index][Thread Index]

Re: [IMP-dev] SphereRestraint



Daniel Russel wrote:
This diff removes CoordinateRestraint and adds SphereRestraint which restraints a point's distance to a fixed point. It also makes the vectori to Int change. The XYZ decorators patch is also attached again as the previous version didn't quite split things right (I had needed to add methods to Particle.h to changed the optimized status of a float).

The patch breaks 3 of the unit tests, at least on synth, so I'm not committing it. For reference, attached is the "complete" patch against r639.

	Ben
--
                      http://salilab.org/~ben/
"It is a capital mistake to theorize before one has data."
	- Sir Arthur Conan Doyle
Index: test/coordinate/test_coordinate.py
===================================================================
--- test/coordinate/test_coordinate.py	(revision 639)
+++ test/coordinate/test_coordinate.py	(working copy)
@@ -8,167 +8,63 @@
     """Test various absolute position restraints"""
 
     def setUp(self):
-        self.imp_model = IMP.Model()
-        self.particles = []
-        self.restraint_sets = []
-        self.rsrs = []
+        self.model = IMP.Model()
+        p= IMP.Particle()
+        self.pi= self.model.add_particle(p);
+        d= IMP.XYZDecorator.create(p)
+        d.set_coordinates_are_optimized(True)
+        pc= IMP.Particle()
+        self.pic= self.model.add_particle(pc);
+        dc= IMP.XYZDecorator.create(pc)
+        
+        self.opt = IMP.SteepestDescent()
+        self.opt.set_model(self.model)
+        #self.opt.set_threshold(1e-5)
 
-        for p in range(12):
-            self.particles.append(IMP.utils.XYZParticle(self.imp_model,
-                                                        0., 0., 0.))
-        self.opt = IMP.ConjugateGradients()
-        self.opt.set_model(self.imp_model)
-        self.opt.set_threshold(1e-5)
 
-
-    def _do_test_min(self, coord, mask):
+    def _do_test(self, center, sf):
         """All coordinate values should be greater than the set minimum"""
-        self.randomize_particles(self.particles, 50.0)
 
-        rs = IMP.RestraintSet("%s-min" % coord)
-        self.restraint_sets.append(rs)
-        self.imp_model.add_restraint(rs)
+        r= IMP.SphericalRestraint(self.model,
+                                  self.pi,
+                                  center[0], center[1], center[2],
+                                  sf)
+        ri=self.model.add_restraint(r)
 
-        score_func_params = IMP.BasicScoreFuncParams("harmonic_lower_bound",
-                                                     8.0, 0.1)
-        for p in self.imp_model.get_particles():
-            self.rsrs.append(IMP.CoordinateRestraint(self.imp_model,
-                                                     p,
-                                                     "%s_AXIS" % coord.upper(),
-                                                     score_func_params))
-            r = self.rsrs[len(self.rsrs)-1]
-            rs.add_restraint(r)
-
         self.opt.optimize(55)
-        for p in self.particles:
-            self.assert_(self.check_abs_pos(p, '>=', 7.999, *mask),
-                         "%s-min condition" % coord)
-        rs.set_is_active(False)
+        self.model.get_restraint(ri).set_is_active(False)
 
-    def test_min(self):
-        """All coordinate values should be greater than the set minimum"""
-        for (coord, mask) in (('x', (1,0,0)), ('y', (0,1,0)), ('z', (0,0,1))):
-            self._do_test_min(coord, mask)
+    def _get_center(self):
+        v= IMP.Floats()
+        d= IMP.XYZDecorator.cast(self.model.get_particle(self.pic))
+        v.append(d.get_x())
+        v.append(d.get_y())
+        v.append(d.get_z())
+        return v
 
-    def _do_test_max(self, coord, mask):
-        """All coordinate values should be less than the set maximum"""
-        self.randomize_particles(self.particles, 50.0)
+    def test_in_ball(self):
+        """Testing that restraint keeps point in ball"""
+        self.randomize_particles(self.model.get_particles(), 50.0)
+        f= IMP.HarmonicUpperBound(10,.1)
+        c= self._get_center()
+        self._do_test(c, f)
+        pd= IMP.XYZDecorator.cast(self.model.get_particle(self.pi))
+        cd= IMP.XYZDecorator.cast(self.model.get_particle(self.pic))
+        d= IMP.distance(pd, cd)
+        self.assert_(d < 11)
 
-        rs = IMP.RestraintSet("%s-max" % coord)
-        self.restraint_sets.append(rs)
-        self.imp_model.add_restraint(rs)
+    def test_on_ball(self):
+        """Testing that restraint keeps point on sphere"""
+        self.randomize_particles(self.model.get_particles(), 50.0)
+        f= IMP.Harmonic(10,.1)
+        c= self._get_center()
+        self._do_test( c, f)
+        pd= IMP.XYZDecorator.cast(self.model.get_particle(self.pi))
+        cd= IMP.XYZDecorator.cast(self.model.get_particle(self.pic))
+        d= IMP.distance(pd, cd)
+        self.assert_(d < 11)
+        self.assert_(d > 9)
+  
 
-        score_func_params = IMP.BasicScoreFuncParams("harmonic_upper_bound",
-                                                     -8.0, 0.1)
-        for p in self.imp_model.get_particles():
-            self.rsrs.append(IMP.CoordinateRestraint(self.imp_model,
-                                                     p,
-                                                     "%s_AXIS" % coord.upper(),
-                                                     score_func_params))
-            r = self.rsrs[len(self.rsrs)-1]
-            rs.add_restraint(r)
-
-        self.opt.optimize(55)
-        for p in self.particles:
-            self.assert_(self.check_abs_pos(p, '<=', -7.9999, *mask),
-                         "%s-max condition" % coord)
-        rs.set_is_active(False)
-
-    def test_max(self):
-        """All coordinate values should be less than the set maximum"""
-        for (coord, mask) in (('x', (1,0,0)), ('y', (0,1,0)), ('z', (0,0,1))):
-            self._do_test_max(coord, mask)
-
-    def test_xy_radial(self):
-        """All xy distances should be less than the set maximum"""
-        self.randomize_particles(self.particles, 50.0)
-
-        rs = IMP.RestraintSet("xy_radial")
-        self.restraint_sets.append(rs)
-        self.imp_model.add_restraint(rs)
-
-        score_func_params = IMP.BasicScoreFuncParams("harmonic_upper_bound",
-                                                     8.0, 0.1)
-        for p in self.imp_model.get_particles():
-            r = IMP.CoordinateRestraint(self.imp_model, p,
-                                        "XY_RADIAL", score_func_params)
-            self.rsrs.append(r)
-            rs.add_restraint(r)
-
-        self.opt.optimize(55)
-
-        for p in self.particles:
-            self.assert_(self.check_abs_pos(p, '<=', 8.001, 1, 1, 0),
-                         "xy_radial condition")
-        rs.set_is_active(False)
-
-    def test_xz_radial(self):
-        """All xz distances should be less than the set maximum"""
-        self.randomize_particles(self.particles, 50.0)
-
-        rs = IMP.RestraintSet("xz_radial")
-        self.restraint_sets.append(rs)
-        self.imp_model.add_restraint(rs)
-
-        score_func_params = IMP.BasicScoreFuncParams("harmonic_upper_bound",
-                                                     8.0, 0.1)
-        for p in self.imp_model.get_particles():
-            r = IMP.CoordinateRestraint(self.imp_model, p,
-                                        "XZ_RADIAL", score_func_params)
-            self.rsrs.append(r)
-            rs.add_restraint(r)
-
-        self.opt.optimize(55)
-
-        for p in self.particles:
-            self.assert_(self.check_abs_pos(p, '<=', 8.001, 1, 0, 1),
-                         "xz_radial condition")
-        rs.set_is_active(False)
-
-    def test_yz_radial(self):
-        """All yz distances should be less than the set maximum"""
-        self.randomize_particles(self.particles, 50.0)
-
-        self.restraint_sets.append(IMP.RestraintSet("yz_radial"))
-        rs = self.restraint_sets[len(self.restraint_sets)-1]
-        self.imp_model.add_restraint(rs)
-
-        score_func_params = IMP.BasicScoreFuncParams("harmonic_upper_bound",
-                                                     8.0, 0.1)
-        for p in self.imp_model.get_particles():
-            r = IMP.CoordinateRestraint(self.imp_model, p,
-                                        "YZ_RADIAL", score_func_params)
-            self.rsrs.append(r)
-            rs.add_restraint(r)
-
-        self.opt.optimize(55)
-        for p in self.particles:
-            self.assert_(self.check_abs_pos(p, '<=', 8.001, 0, 1, 1),
-                         "yz_radial condition")
-        rs.set_is_active(False)
-
-    def test_xyz_sphere(self):
-        """All xyz distances should be less than the set maximum"""
-        self.randomize_particles(self.particles, 50.0)
-
-        self.restraint_sets.append(IMP.RestraintSet("xyz_sphere"))
-        rs = self.restraint_sets[len(self.restraint_sets)-1]
-        self.imp_model.add_restraint(rs)
-
-        score_func_params = IMP.BasicScoreFuncParams("harmonic_upper_bound",
-                                                     8.0, 0.1)
-        for p in self.imp_model.get_particles():
-            r = IMP.CoordinateRestraint(self.imp_model, p,
-                                        "XYZ_SPHERE", score_func_params)
-            self.rsrs.append(r)
-            rs.add_restraint(r)
-
-        self.opt.optimize(55)
-
-        for p in self.particles:
-            self.assert_(self.check_abs_pos(p, '<=', 8.001, 1, 1, 1),
-                         "xyz_sphere condition")
-        rs.set_is_active(False)
-
 if __name__ == '__main__':
     unittest.main()
Index: test/pair_connectivity/test_pair_connectivity.py
===================================================================
--- test/pair_connectivity/test_pair_connectivity.py	(revision 639)
+++ test/pair_connectivity/test_pair_connectivity.py	(working copy)
@@ -96,8 +96,8 @@
 
         # add connectivity restraints
 
-        particle_indexes1 = IMP.vectori()
-        particle_indexes2 = IMP.vectori()
+        particle_indexes1 = IMP.Ints()
+        particle_indexes2 = IMP.Ints()
         rsrs = []
         # connect 2 proteins together by two beads
         particle_indexes1.clear()
Index: test/modeller/test_exclusion_volumes.py
===================================================================
--- test/modeller/test_exclusion_volumes.py	(revision 639)
+++ test/modeller/test_exclusion_volumes.py	(working copy)
@@ -108,8 +108,8 @@
 
         # add connectivity restraints
 
-        particle_indexes1 = IMP.vectori()
-        particle_indexes2 = IMP.vectori()
+        particle_indexes1 = IMP.Ints()
+        particle_indexes2 = IMP.Ints()
         rsrs = []
 
         # connect 2 proteins together by two beads
@@ -168,8 +168,8 @@
 
         # add connectivity restraints
 
-        particle_indexes1 = IMP.vectori()
-        particle_indexes2 = IMP.vectori()
+        particle_indexes1 = IMP.Ints()
+        particle_indexes2 = IMP.Ints()
         rsrs = []
 
         # connect 2 proteins together by two beads
Index: test/modeller/test_proximity.py
===================================================================
--- test/modeller/test_proximity.py	(revision 639)
+++ test/modeller/test_proximity.py	(working copy)
@@ -79,7 +79,7 @@
 
         # add proximity restraints
 
-        particle_indexes = IMP.vectori()
+        particle_indexes = IMP.Ints()
         rsrs = []
 
         # set up exclusion volumes
Index: test/modeller/test_easy_proximity.py
===================================================================
--- test/modeller/test_easy_proximity.py	(revision 639)
+++ test/modeller/test_easy_proximity.py	(working copy)
@@ -69,7 +69,7 @@
 
         # add proximity restraints
 
-        particle_indexes = IMP.vectori()
+        particle_indexes = IMP.Ints()
         rsrs = []
 
         # set up exclusion volumes
Index: test/connectivity/test_connectivity.py
===================================================================
--- test/connectivity/test_connectivity.py	(revision 639)
+++ test/connectivity/test_connectivity.py	(working copy)
@@ -95,7 +95,7 @@
 
         # add connectivity restraints
 
-        particle_indexes = IMP.vectori()
+        particle_indexes = IMP.Ints()
         rsrs = []
 
         score_func_params_ub = IMP.BasicScoreFuncParams("harmonic_upper_bound",
Index: include/IMP/ModelData.h
===================================================================
--- include/IMP/ModelData.h	(revision 639)
+++ include/IMP/ModelData.h	(working copy)
@@ -147,6 +147,8 @@
     return string_data_[idx.get_index()];
   }
 
+
+  void show(std::ostream &out=std::cout) const;
 protected:
   ModelData();
   ~ModelData();
@@ -187,6 +189,8 @@
   bool check_particles_active_;
 };
 
+IMP_OUTPUT_OPERATOR(ModelData);
+
 //! Optimizable variable iterator
 /** Returns all optimizable Floats in the ModelData.
  */
Index: include/IMP/base_types.h
===================================================================
--- include/IMP/base_types.h	(revision 639)
+++ include/IMP/base_types.h	(working copy)
@@ -28,6 +28,12 @@
 
 
 
+typedef std::vector<Float> Floats;
+typedef std::vector<Int> Ints;
+typedef std::vector<String> Strings;
+
+
+
 struct ParticleTag {};
 struct RestraintTag {};
 struct StateTag {};
@@ -45,6 +51,11 @@
    We need this to have a uniform return type for python.
  */
 typedef std::vector<Particle*> Particles; 
+
+  typedef std::vector<ParticleIndex> ParticleIndexes;
+  typedef std::vector<FloatIndex> FloatIndexes;
+  typedef std::vector<IntIndex> IntIndexes;
+  typedef std::vector<StringIndex> StringIndexes;
 // typedefs for the particle variable and attribute indexes
 // typedef DataIndex<Float> FloatIndex;
 // typedef DataIndex<Int> IntIndex;
@@ -80,7 +91,8 @@
     inline unsigned int attribute_table_index(Tag) {                    \
       return Name##_attribute_table_index_;                             \
     }                                                                   \
-  }
+  }                                                                     \
+  typedef std::vector<Name> Name##s;
 
 /**
    Declare static data necessary for a new key type.
Index: include/IMP/restraints/SConscript
===================================================================
--- include/IMP/restraints/SConscript	(revision 639)
+++ include/IMP/restraints/SConscript	(working copy)
@@ -2,7 +2,7 @@
 Import('env')
 
 files = ['ConnectivityRestraint.h', 'ProximityRestraint.h',
-         'CoordinateRestraint.h', 'DistanceRestraint.h',
+         'SphericalRestraint.h', 'DistanceRestraint.h',
          'RestraintSet.h', 'ExclusionVolumeRestraint.h', 'TorusRestraint.h',
          'PairConnectivityRestraint.h']
 
Index: include/IMP/restraints/SphericalRestraint.h
===================================================================
--- include/IMP/restraints/SphericalRestraint.h	(revision 0)
+++ include/IMP/restraints/SphericalRestraint.h	(revision 0)
@@ -0,0 +1,47 @@
+/**
+ *  \file SphericalRestraint.h   \brief Absolute position restraint.
+ *
+ *  Optimize based on distance from an absolute position.
+ *
+ *  Copyright 2007 Sali Lab. All rights reserved.
+ *
+ */
+
+#ifndef __IMP_SPHERICAL_RESTRAINT_H
+#define __IMP_SPHERICAL_RESTRAINT_H
+
+#include <list>
+
+#include "../IMP_config.h"
+#include "../Restraint.h"
+
+namespace IMP
+{
+
+//! Restrict particle position based on its distance to a point.
+class IMPDLLEXPORT SphericalRestraint : public Restraint
+{
+public:
+  /**
+     \param[in] model
+     \param[in] p The particle to restrict.
+     \param[in] x The x coordinate to take distance to.
+     \param[in] y The y coordinate to take distance to.
+     \param[in] z The z coordinate to take distance to.
+     \param[in] score_func The scoring function. It is deleted in the destructor.
+   */
+  SphericalRestraint(Model* model, ParticleIndex p,
+                     Float x, Float y, Float z,
+                     ScoreFunc* score_func);
+  virtual ~SphericalRestraint();
+
+  IMP_RESTRAINT("0.5", "Daniel Russel");
+
+protected:
+  Float center_[3];
+  ScoreFunc* score_func_;
+};
+
+} // namespace IMP
+
+#endif /* __IMP_COORDINATE_RESTRAINT_H */
Index: include/IMP.h
===================================================================
--- include/IMP.h	(revision 639)
+++ include/IMP.h	(working copy)
@@ -31,7 +31,7 @@
 #include "IMP/restraints/RestraintSet.h"
 #include "IMP/restraints/DistanceRestraint.h"
 #include "IMP/restraints/TorusRestraint.h"
-#include "IMP/restraints/CoordinateRestraint.h"
+#include "IMP/restraints/SphericalRestraint.h"
 #include "IMP/restraints/ProximityRestraint.h"
 #include "IMP/restraints/ConnectivityRestraint.h"
 #include "IMP/restraints/PairConnectivityRestraint.h"
Index: doc/examples/misc/energy_demo.py
===================================================================
--- doc/examples/misc/energy_demo.py	(revision 639)
+++ doc/examples/misc/energy_demo.py	(working copy)
@@ -93,7 +93,7 @@
 
 # add connectivity restraints
 
-particle_indexes = IMP.vectori()
+particle_indexes = IMP.Ints()
 rsrs = []
 
 # connect 3 proteins together
Index: impEM/test/test_derivs.py
===================================================================
--- impEM/test/test_derivs.py	(revision 639)
+++ impEM/test/test_derivs.py	(working copy)
@@ -59,7 +59,7 @@
         init_particle(self.particles,2,3.0,12.0,12.0,rad,wei,1)
         IMP.modeller_intf.copy_imp_coords_to_modeller(self.particles,self.modeller_model.atoms)
         #modeller_model.write(file='xxx.pdb')
-        self.particle_indexes = IMP.vectori()
+        self.particle_indexes = IMP.Ints()
         self.particle_indexes.clear()
         for i in xrange(npart):
             self.particle_indexes.push_back(i)
Index: impEM/test/test_em_fit.py
===================================================================
--- impEM/test/test_em_fit.py	(revision 639)
+++ impEM/test/test_em_fit.py	(working copy)
@@ -28,7 +28,7 @@
 
 
         self.particles = []
-        self.particle_indexes = IMP.vectori()
+        self.particle_indexes = IMP.Ints()
 
         origin =  3.0
         self.particles.append(IMP.utils.XYZParticle(self.imp_model, 9.+origin,
Index: impEM/test/test_sample_particles.py
===================================================================
--- impEM/test/test_sample_particles.py	(revision 639)
+++ impEM/test/test_sample_particles.py	(working copy)
@@ -45,7 +45,7 @@
         init_particle(self.particles,0,9.0,9.0,9.0,rad,wei,1)
         init_particle(self.particles,1,12.0,3.0,3.0,rad,wei,1)
         init_particle(self.particles,2,3.0,12.0,12.0,rad,wei,1)
-        self.particle_indexes = IMP.vectori()
+        self.particle_indexes = IMP.Ints()
         self.particle_indexes.clear()
         for i in xrange(npart):
             self.particle_indexes.push_back(i)
Index: src/restraints/SConscript
===================================================================
--- src/restraints/SConscript	(revision 639)
+++ src/restraints/SConscript	(working copy)
@@ -2,7 +2,7 @@
 
 files = ['PairConnectivityRestraint.cpp',
          'ConnectivityRestraint.cpp', 'ProximityRestraint.cpp',
-         'CoordinateRestraint.cpp', 'RestraintSet.cpp',
+         'SphericalRestraint.cpp', 'RestraintSet.cpp',
          'DistanceRestraint.cpp', 'TorusRestraint.cpp',
          'ExclusionVolumeRestraint.cpp']
 
Index: src/restraints/SphericalRestraint.cpp
===================================================================
--- src/restraints/SphericalRestraint.cpp	(revision 0)
+++ src/restraints/SphericalRestraint.cpp	(revision 0)
@@ -0,0 +1,88 @@
+/**
+ *  \file SphericalRestraint.cpp \brief Absolute position restraint.
+ *
+ *  Optimize based on distance from an absolute position.
+ *
+ *  Copyright 2007 Sali Lab. All rights reserved.
+ *
+ */
+
+#include <cmath>
+
+#include "IMP/Model.h"
+#include "IMP/Particle.h"
+#include "IMP/log.h"
+#include "IMP/restraints/SphericalRestraint.h"
+#include "IMP/decorators/XYZDecorator.h"
+
+namespace IMP
+{
+
+static const float MIN_DISTANCE_SQUARED=.001;
+
+
+
+SphericalRestraint::SphericalRestraint(Model* m, ParticleIndex p,
+                                       Float x, Float y, Float z,
+                                       ScoreFunc* score_func)
+{
+  
+  add_particle(m->get_particle(p));
+  center_[0]=x;
+  center_[1]=y;
+  center_[2]=z;
+  score_func_ =score_func;
+}
+
+SphericalRestraint::~SphericalRestraint()
+{
+  delete score_func_;
+}
+
+
+
+Float SphericalRestraint::evaluate(DerivativeAccumulator *accum)
+{
+
+  Float d2=0;
+  Float diff[3];
+  XYZDecorator xyzd= XYZDecorator::cast(get_particle(0));
+  for (unsigned int i=0; i< 3; ++i) {
+    diff[i] =xyzd.get_coordinate(i) - center_[i];
+    d2+= diff[i]*diff[i];
+  }
+  Float d= std::sqrt(d2);
+  if (d2 < MIN_DISTANCE_SQUARED) {
+    return 0;
+  } 
+  
+  Float ret=0;
+  if (accum) {
+    Float deriv;
+    ret= (*score_func_)(d, deriv);
+    for (unsigned int i=0; i< 3; ++i) {
+      Float dd= deriv*diff[i]/d;
+      xyzd.add_to_coordinate_derivative(i, dd, *accum);
+    }
+  } else {
+    ret= (*score_func_)(d);
+  }
+  return ret;
+}
+
+
+void SphericalRestraint::show(std::ostream& out) const
+{
+  if (get_is_active()) {
+    out << "Spherical restraint (active):" << std::endl;
+  } else {
+    out << "Spherical restraint (inactive):" << std::endl;
+  }
+
+  out << "version: " << version() << "  ";
+  out << "last_modified_by: " << last_modified_by() << std::endl;
+  out << "  particle: " << get_particle(0)->get_index();
+  out << std::endl;
+}
+
+}  // namespace IMP
Index: src/ModelData.cpp
===================================================================
--- src/ModelData.cpp	(revision 639)
+++ src/ModelData.cpp	(working copy)
@@ -160,6 +160,15 @@
 }
 
 
+void ModelData::show(std::ostream &out) const {
+  out << "Float variables: " << std::endl;
+  for (unsigned int i=0; i< float_data_.size(); ++i) {
+    out << "  " << float_data_[i].value_ << " " << float_data_[i].deriv_ 
+        << (float_data_[i].is_optimized_?" optimized" : " not optimized")
+        << std::endl;
+  }
+}
+
 //! Reset the iterator.
 /** After the next call to next(), get() will return the first optimizable
     Float variable.
Index: pyext/IMP/test.py
===================================================================
--- pyext/IMP/test.py	(revision 639)
+++ pyext/IMP/test.py	(working copy)
@@ -1,5 +1,6 @@
 import re, math, unittest
 import random
+import IMP
 
 
 class IMPTestCase(unittest.TestCase):
@@ -8,9 +9,10 @@
     def randomize_particles(self, particles, deviation):
         """Randomize the xyz coordinates of a list of particles"""
         for p in particles:
-            p.set_x(random.uniform(-deviation, deviation))
-            p.set_y(random.uniform(-deviation, deviation))
-            p.set_z(random.uniform(-deviation, deviation))
+            d= IMP.XYZDecorator.cast(p)
+            d.set_x(random.uniform(-deviation, deviation))
+            d.set_y(random.uniform(-deviation, deviation))
+            d.set_z(random.uniform(-deviation, deviation))
 
     def load_coordinates(self, pdb_file):
         """Load coordinates from a PDB file"""
Index: pyext/IMP/xml_loader.py
===================================================================
--- pyext/IMP/xml_loader.py	(revision 639)
+++ pyext/IMP/xml_loader.py	(working copy)
@@ -329,7 +329,7 @@
         return
 
     # create a vector of particle indexes
-    particle_indexes = IMP.vectori()
+    particle_indexes = IMP.Ints()
     for p_idx in particle_list:
         particle_indexes.push_back(int(p_idx))
 
@@ -393,12 +393,12 @@
         return
 
     # create a vector of particle indexes
-    particle_indexes1 = IMP.vectori()
+    particle_indexes1 = IMP.Ints()
     particle_indexes1.clear()
     for p_idx in particle_list1:
         particle_indexes1.push_back(int(p_idx))
 
-    particle_indexes2 = IMP.vectori()
+    particle_indexes2 = IMP.Ints()
     particle_indexes2.clear()
     for p_idx in particle_list2:
         particle_indexes2.push_back(int(p_idx))
@@ -447,7 +447,7 @@
         return
 
     # create a vector of particle indexes
-    particle_indexes = IMP.vectori()
+    particle_indexes = IMP.Ints()
     for p_idx in particle_list:
         particle_indexes.push_back(int(p_idx))
 
@@ -489,14 +489,14 @@
         return
 
     # create a vector of particle indexes
-    particle_indexes1 = IMP.vectori()
+    particle_indexes1 = IMP.Ints()
     particle_indexes1.clear()
     for p_idx in particle_list1:
         particle_indexes1.push_back(int(p_idx))
 
     score_func_params = get_basic_score_func_params("harmonic_lower_bound", 0.0, sd)
     if (num_lists == 2):
-        particle_indexes2 = IMP.vectori()
+        particle_indexes2 = IMP.Ints()
         particle_indexes2.clear()
         for p_idx in particle_list2:
             particle_indexes2.push_back(int(p_idx))
Index: pyext/IMP.i
===================================================================
--- pyext/IMP.i	(revision 639)
+++ pyext/IMP.i	(working copy)
@@ -10,10 +10,6 @@
 %include "IMP_macros.i"
 %include "IMP_exceptions.i"
 
-namespace std {
-  %template(vectori) vector<int>;
-  %template(vectorf) vector<float>;
-}
 
 namespace IMP {
   %pythonprepend Model::add_particle %{
@@ -28,6 +24,9 @@
   %pythonprepend RestraintSet::add_restraint %{
         args[1].thisown=0
   %}
+  %pythonprepend SphericalRestraint::SphericalRestraint %{
+        args[5].thisown=0
+  %}
 }
 
 %feature("director");
@@ -61,7 +60,7 @@
 %include "IMP/optimizers/MolecularDynamics.h"
 %include "IMP/restraints/DistanceRestraint.h"
 %include "IMP/restraints/TorusRestraint.h"
-%include "IMP/restraints/CoordinateRestraint.h"
+%include "IMP/restraints/SphericalRestraint.h"
 %include "IMP/restraints/ProximityRestraint.h"
 %include "IMP/restraints/ConnectivityRestraint.h"
 %include "IMP/restraints/PairConnectivityRestraint.h"
@@ -86,4 +85,7 @@
   %template(FloatKeys) ::std::vector<FloatKey>;
   %template(StringKeys) ::std::vector<StringKey>;
   %template(IntKeys) ::std::vector<IntKey>;
+  %template(Floats) ::std::vector<Float>;
+  %template(Strings) ::std::vector<String>;
+  %template(Ints) ::std::vector<Int>;
 }