IMP logo
IMP Reference Guide  2.12.0
The Integrative Modeling Platform
SitesPairScore.h
Go to the documentation of this file.
1 /**
2  * \file SitesPairScore.h
3  * \brief A Score on the distance between a pair of particles.
4  *
5  * Copyright 2007-2019 IMP Inventors. All rights reserved.
6  */
7 
8 // TODO: verify if energy units are kcal/mol or KT
9 
10 #ifndef IMPNPCTRANSPORT_SITES_PAIR_SCORE_H
11 #define IMPNPCTRANSPORT_SITES_PAIR_SCORE_H
12 
13 #include "npctransport_config.h"
16 #include "internal/RigidBodyInfo.h"
17 #include "internal/sites.h"
18 #include <IMP/PairScore.h>
19 #include <IMP/UnaryFunction.h>
20 #include <IMP/Pointer.h>
21 #include <IMP/core/XYZR.h>
22 #include <IMP/core/Typed.h>
25 #include <IMP/generic.h>
28 #include <IMP/set_map_macros.h>
30 #include <IMP/atom/estimates.h>
31 #include <boost/unordered_set.hpp>
32 
33 #include <boost/array.hpp>
34 
35 IMPNPCTRANSPORT_BEGIN_NAMESPACE
36 
37 
38 /** \brief Apply a function to the distance between two particles with
39  a set of specific binding sites
40 
41  The sites are expressed in the local reference frame of
42  the two rigid bodies. Care must be taken to pass the bodies
43  in the appropriate order. See construction documentation for more details.
44 */
45 class IMPNPCTRANSPORTEXPORT SitesPairScore
47 {
48  private:
50 
51 
52  /************************* class variables ****************/
53  bool is_orientational_score_; // if true, use orientation-dependent score
56  sites0_, sites1_;
57 
58 
59  //! Maximal square of distance between particles with interacting sites
60  //! based on the interaction range and list of sites
61  double ubound_distance2_;
62 
63  // //! Cache:
64  // typedef IMP_KERNEL_LARGE_UNORDERED_MAP<ParticleIndex,internal::RigidBodyInfo>
65  // t_particles_rb_cache;
66  // mutable t_particles_rb_cache particles_rb_cache_;
67  // mutable bool is_cache_active_;
68  // mutable unsigned int cur_cache_id_; // to keep track of caching rounds
69 
70  public:
71 
72  /**
73  For positive sigmas, this is an orientation dependent score between two
74  spherical particles that
75  contain a fixed set of interaction sites, sites0 and sites1 (for
76  first and second particle, resp.).
77 
78  The interaction is composed of a non-specific interaction term
79  between the bead spheres, and the sum of interactions between
80  specific interactions sites of each bead.
81 
82  If \sigma{0} and \sigma{1} are positive, the attractive force
83  between pairs of sites at an optimal orientation (sites facing
84  each other) depends on their distance x. The attraction force
85  magnitude is k*x when x<=0.5*range and k*(range-x) when x is
86  between 0.5*range and range. When site0 is rotated by \sigma <
87  \sigma{0}, this force decays further by a factor
88  (cos{\sigma}-cos{\sigma}0)/(1.0-cos{\sigma0}). The force decays
89  similarly when site1 is rotated by \sigma < \sigma{1}. The
90  maximal potential energy difference due to such pair of
91  interacting sites is:
92  max\DELTA{U}_{site_site} = 0.25 * k * range^2 [kcal/mol]
93 
94  If sigma0_deg or sigma1_deg are non-positive, the attractive
95  force between pairs of sites within the attraction range is a
96  constant k, and the maximal interaction energy is instead:
97  max\DELTA{U}_{site-site} = k * range [kCal/mol]
98 
99  In addition to site-site interaction, there is a constant
100  attractive force k_nonspec_attraction between the sphere surfaces
101  up to a range range_nonspec_attraction, with maximal energy
102  contribution:
103  max\DELTA{U}_{non-specific} = 0.5 * k_nonspec_attraction * range_nonspec_attraction [kcal/mol]
104 
105  Note that for a specific pair of particles, each particle might have
106  a different reference frame (rigid body translation and rotation),
107  which is applied to the sites list upon score evaluation.
108 
109  @param range Maximal range of site specific attraction in any direction
110  of specific sites placed on particles
111  @param k Maximal site specific attraction coefficient (in
112  kcal/mol/A^2 when sigma0_deg and sigma1_deg are positive,
113  or kcal/mol/A otherwise, i.e., for orientation-independent interactions)
114  @param sigma0_deg, sigma1_deg Maximal rotational range of sites 0 and 1, respectively,
115  on the particle surface, specified in degrees. If either is 0,
116  the pair score between site centers is used with a constant k.
117  @param range_nonspec_attraction Range for non-specific attraction term
118  between particles that contain the sites
119  @param k_nonspec_attraction Non-specific attraction coefficient between particles
120  (constant force in kCal/mol/A within specified range)
121  @param k_nonspec_repulsion Repulsion coefficient between particles (constant force
122  applied when particle spheres overlap in kCal/mol/A)
123  @param sites0 List of sites on the first particle, in its local reference frame
124  @param sites1 List of sites on the second particle, in its local reference frame
125  */
126  SitesPairScore(double range, double k,
127  double sigma0_deg, double sigma1_deg,
128  double range_nonspec_attraction, double k_nonspec_attraction,
129  double k_nonspec_repulsion,
130  const algebra::Sphere3Ds &sites0,
131  const algebra::Sphere3Ds &sites1);
132 
133 
134  public:
135 
136  virtual double evaluate_indexes(Model *m, const ParticleIndexPairs &p,
138  unsigned int lower_bound,
139  unsigned int upper_bound) const IMP_FINAL;
140 
141  //! evaluated indexes for the range from lower_bound to upper_bound
142  //! in p, if score>max then return max value of double
144  ( Model *m, const ParticleIndexPairs &p,
146  double max, unsigned int lower_bound, unsigned int upper_bound) const
147  {
148  // activate_cache();
149  double ret = 0.0;
150  for (unsigned int i = lower_bound; i < upper_bound; ++i) {
151  ret += evaluate_if_good_index(m, p[i], da, max - ret);
152  if (ret > max) return std::numeric_limits<double>::max();
153  }
154  // deactivate_cache();
155  return ret;
156  }
157 
158 
159  /** evaluate the score for the pair of model particle indexes in p,
160  updating score derivatives to da
161  */
162  virtual double evaluate_index(Model *m, const ParticleIndexPair &p,
164 
165 #ifndef SWIG
166  /**
167  EvaluatE all site-site interactions
168  for evaluate_index() for the pair pip in model m. If da is not nullptr,
169  it accumulated appropriate derivatives. If contacts_accumulator is not null,
170  then the number of individual contacts and occupied sites is asscumulated there.
171 
172  @param sphere_table An array storing of sphere coordinates by particle index
173  @param quaternions_tables An array of quaternions by particle index
174  @param sphere_table An array storing of sphere coordinate derivatives by particle index
175  @param torque_tables An array of torques by particle index
176  @param pip the pair of particle indexes in m
177  @param da optional accumulator for force and torque derivatives
178  @param contacts_accumulator A pointer to a tuple of output values
179  [num-contacts, sites0-bound, sites1-bound, is_nonspec].
180  num-contacts is the total number of site-site contacts between pip.
181  sites0-bound and sites1-bound are vectors of contact counts
182  for each site of pip[0] and pip[1], resp, with one entry per site.
183  is_nonspec is true if the molecules have non-zero
184  nonspecific interactions Ignored if Null
185 
186  @return the site-site contributions for the score for the pair
187  pip in model m.
188  */
189  double
190  evaluate_site_contributions_with_internal_tables
191  (algebra::Sphere3D const* spheres_table,
192  double const**quaternions_tables,
193  algebra::Sphere3D *sphere_derivatives_table,
194  double **torques_tables,
195  const ParticleIndexPair &pip,
197  boost::tuple< unsigned int, std::vector<unsigned int>, std::vector<unsigned int>, bool >
198  (*contacts_accumulator) = nullptr
199  ) const;
200 
201  /**
202  EvaluatE all site-site interactions
203  for evaluate_index() for the pair pip in model m. If da is not nullptr,
204  it accumulated appropriate derivatives. If contacts_accumulator is not null,
205  then the number of individual contacts and occupied sites is asscumulated there.
206 
207  @param m the model
208  @param pip the pair of particle indexes in m
209  @param da optional accumulator for force and torque derivatives
210  @param contacts_accumulator A pointer to a tuple of output values
211  [num-contacts, sites0-bound, sites1-bound, is_nonspec].
212  num-contacts is the total number of site-site contacts between pip.
213  sites0-bound and sites1-bound are vectors of contact counts
214  for each site of pip[0] and pip[1], resp. is_nonspec is true if the
215  spheres have non-zero non-specific interactions
216  Ignored if Null
217 
218  @return the site-site contributions for the score for the pair
219  pip in model m.
220  */
221  double
222  evaluate_site_contributions
223  (Model* m,
224  const ParticleIndexPair &pip,
226  boost::tuple< unsigned int, std::vector<unsigned int>, std::vector<unsigned int>, bool >
227  (*contacts_accumulator)
228  ) const;
229 
230 #endif
231 
233  const ParticleIndexes &pis) const;
234 
235  // Restraints do_create_current_decomposition(Model *m,
236  // const ParticleIndexPair &vt)
237  // const IMP_OVERRIDE;
238 
239  //! return the range for site-site attraction
240  double get_sites_range() const { return params_.r; }
241 
242  //! return the k for site-site attraction
243  double get_sites_k() const { return params_.k; }
244 
245  SitesPairScoreParameters get_parameters() const {return params_;}
246 
247  public:
249 
250  private:
251 
252  /** evaluate the score for the pair of model particle indexes in p,
253  updating score derivatives to da, and using internal attribute
254  tables in Model
255  */
256  inline double evaluate_index_with_internal_tables
257  ( Model* m,
258  algebra::Sphere3D const* spheres_table,
259  double const **quaternions_tables,
260  algebra::Sphere3D *sphere_derivatives_table,
261  double **torques_tables,
262  const ParticleIndexPair &p,
263  DerivativeAccumulator *da) const;
264 
265  // gets the rigid body information (e.g., translation, inverse rotation)
266  // associated with particle m.pi, possibly from cache (depending on internal
267  // cache definitions)
268  inline internal::RigidBodyInfo
269  get_rigid_body_info
270  (algebra::Sphere3D const* spheres_table,
271  double const** quaternions_tables,
272  ParticleIndex pi) const;
273 
274  public:
275  // sets the sites associated with each partner to sites0
276  // and sites1, respectively (in local reference frame)
277  void set_sites(const algebra::Sphere3Ds &sites0,
278  const algebra::Sphere3Ds &sites1){
279  sites0_= sites0;
280  sites1_= sites1;
281  }
282 
283  // sets the sites associated with the first partner
284  // (in local reference frame)
285  void set_sites0(const algebra::Sphere3Ds &sites0){
286  sites0_= sites0;
287  }
288 
289  // sets the sites associated with the first partner
290  // (in local reference frame)
291  void set_sites1(const algebra::Sphere3Ds &sites1){
292  sites1_= sites1;
293  }
294 
295  private:
296  /* // maintain a cache for evaluate_index() till call to deactivate_cache() */
297  /* inline void activate_cache() const */
298  /* { */
299  /* // (note: cache is mutable) */
300  /* is_cache_active_ = true; cur_cache_id_++; */
301  /* } */
302 
303  /* inline void deactivate_cache() const */
304  /* { */
305  /* // (note: cache is mutable) */
306  /* is_cache_active_ = false; */
307  /* } */
308 
309 
310 };
311 
312 //!
313 inline double
314 SitesPairScore::evaluate_index
315 (Model *m, const ParticleIndexPair &p,
316  DerivativeAccumulator *da) const{
317  // get internal tables:
318  algebra::Sphere3D const* spheres_table=
319  m->access_spheres_data();
320  double const* quaternions_tables[4];
321  for(unsigned int i = 0; i < 4; i++){
322  quaternions_tables[i]=
323  core::RigidBody::access_quaternion_i_data(m, i);
324  }
325  algebra::Sphere3D* sphere_derivatives_table=
326  m->access_sphere_derivatives_data();
327  double* torques_tables[3];
328  for(unsigned int i = 0; i < 3; i++){
329  torques_tables[i]=
330  core::RigidBody::access_torque_i_data(m, i);
331  }
332  // evaluate:
333  return evaluate_index_with_internal_tables(m,
334  spheres_table,
335  quaternions_tables,
336  sphere_derivatives_table,
337  torques_tables,
338  p,
339  da);
340 }
341 
342 
343 
344 
345 /**
346  the sites of each particle are transformed to a common frame of reference
347  (using the reference frame of each particle), and the site-specific
348  attraction, and the inter-particle non specific attraction and repulsion
349  are evaluated and summed.
350 */
351 inline double
352 SitesPairScore::evaluate_index_with_internal_tables
353 ( Model* m,
354  algebra::Sphere3D const* spheres_table,
355  double const** quaternions_tables,
356  algebra::Sphere3D *sphere_derivatives_table,
357  double **torques_tables,
358  const ParticleIndexPair &pip,
359  DerivativeAccumulator *da) const {
361 
362  // I. evaluate non-specific attraction and repulsion between
363  // parent particles before computing for specific sites :
364  double non_specific_score = P::evaluate_index(m, pip, da);
366  lips_cache= P::get_evaluation_cache();
367 
368  // II. Return if parent particles are out of site-specific interaction range
369  // using cache to avoid some redundant calcs
370  double const& distance2= lips_cache.particles_delta_squared;
371  IMP_LOG(PROGRESS, "distance2 " << distance2
372  << " ; distance upper-bound " << ubound_distance2_ << std::endl);
373  if (distance2 > ubound_distance2_) {
374  IMP_LOG(PROGRESS, "Sites contribution is 0.0 and non-specific score is "
375  << non_specific_score << std::endl);
376  return non_specific_score;
377  }
378 
379  double site_score=evaluate_site_contributions_with_internal_tables
380  (spheres_table,
381  quaternions_tables,
382  sphere_derivatives_table,
383  torques_tables,
384  pip, da);
385  // III. evaluate site-specific contributions :
386  return site_score + non_specific_score;
387 }
388 
389 #ifndef SWIG
390 
391 //!
392 inline double
393 SitesPairScore::evaluate_site_contributions_with_internal_tables
394 ( algebra::Sphere3D const* spheres_table,
395  double const**quaternions_tables,
396  algebra::Sphere3D *sphere_derivatives_table,
397  double **torques_tables,
398  const ParticleIndexPair &pip,
400  boost::tuple<unsigned int,
401  std::vector<unsigned int>,
402  std::vector<unsigned int>,
403  bool>
404  * contacts_accumulator
405  ) const
406 {
408  // interaction statistics variables
409  static unsigned int n_contacts;
410  static std::vector<unsigned int> occupied_sites0; // how many contacts at each site of particle 0
411  static std::vector<unsigned int> occupied_sites1; // how many contacts at each site of particle 1
412  if(contacts_accumulator){
413  n_contacts= 0;
414  occupied_sites0.resize(sites0_.size());
415  occupied_sites1.resize(sites1_.size());
416  std::fill(occupied_sites0.begin(), occupied_sites0.end(),0);
417  std::fill(occupied_sites1.begin(), occupied_sites1.end(),0);
418  }
419 
420  // bring sites_ to the frame of reference of nn_sites_ and nn_
421  ParticleIndex pi0 = pip[0];
422  ParticleIndex pi1 = pip[1];
423  // get rbi0/1 info, update if needed
424  internal::RigidBodyInfo rbi0 = get_rigid_body_info(spheres_table,
425  quaternions_tables,
426  pi0);
427  internal::RigidBodyInfo rbi1 = get_rigid_body_info(spheres_table,
428  quaternions_tables,
429  pi1);
430  IMP_LOG_PROGRESS( "RBI0.pi " << rbi0.pi
431  << " RB0.cache_id " << rbi0.cache_id
432  << "RBI0.tr " << rbi0.tr << std::endl);
433  IMP_LOG_PROGRESS( "RBI1.cache_id " << rbi1.cache_id
434  << "RBI1.pi " << rbi1.pi
435  << "RBI1.tr " << rbi1.tr << std::endl);
436  // sum over specific interactions between all pairs of sites:
437  double sum = 0;
438  if(is_orientational_score_){
439  // Pre-compute a few variables that do not depend on either both sites or on site1
440  algebra::Vector3D const& gRB0= rbi0.tr.get_translation();
441  algebra::Vector3D const& gRB1= rbi1.tr.get_translation();
442  algebra::Vector3D gUnitRB0RB1= gRB1-gRB0;
443  double distRB0RB1= get_magnitude_and_normalize_in_place(gUnitRB0RB1); // distance between centers
444  for (unsigned int i0 = 0; i0 < sites0_.size(); ++i0) {
445  algebra::Vector3D gSite0 = rbi0.tr.get_transformed(sites0_[i0].get_center());
446  algebra::Vector3D gUnitRB0Site0= (gSite0-gRB0)*rbi0.iradius;
447  double cosSigma0 = gUnitRB0Site0*gUnitRB0RB1;
448  if(cosSigma0 < params_.cosSigma1_max) { // not in range... - note the indexing is not an error - sigma0 is equivalent to params_.sigma1
449  continue;
450  }
451  double kFactor0=internal::get_k_factor(cosSigma0, params_.cosSigma1_max); // note the indexing is not an error - sigma0 is equivalent to params_.sigma1
452  algebra::Vector3D gRotSigma0;
453  double dKFactor0;
454  if(da){
455  gRotSigma0 = get_vector_product(gUnitRB0Site0,gUnitRB0RB1);
456  double absSinSigma0 = get_magnitude_and_normalize_in_place(gRotSigma0);
457  dKFactor0=internal::get_derivative_k_factor(absSinSigma0, params_.cosSigma1_max);
458  }
459  for(unsigned int i1 = 0 ; i1 < sites1_.size(); ++i1) {
460  algebra::Vector3D gSite1 = rbi1.tr.get_transformed(sites1_[i1].get_center());
461  IMP_LOG_PROGRESS( "Evaluating sites at global coordinates: " << gSite0
462  << " ; " << gSite1 << std::endl );
463  double cur_score;
464  cur_score =
465  internal::evaluate_pair_of_sites(params_,
466  rbi0, rbi1,
467  gSite1,
468  gUnitRB0RB1, distRB0RB1,
469  gRotSigma0,
470  kFactor0, dKFactor0,
471  da,
472  sphere_derivatives_table,
473  torques_tables);
474  sum += cur_score;
475  if(contacts_accumulator && cur_score!=0.0){
476  n_contacts++;
477  occupied_sites0[i0]++;
478  occupied_sites1[i1]++;
479  }
480  } // j
481  } // i
482  } // is_orientational_score_
483  else
484  {
485  for (unsigned int i = 0; i < sites0_.size(); ++i) {
486  algebra::Vector3D g0 = rbi0.tr.get_transformed(sites0_[i].get_center());
487  for(unsigned int j = 0 ; j < sites1_.size(); ++j) {
488  algebra::Vector3D g1 = rbi1.tr.get_transformed(sites1_[j].get_center());
489  IMP_LOG_PROGRESS( "Evaluating sites at global coordinates: " << g0 << " ; " << g1 << std::endl );
490  double cur_score;
491  // old score
492  cur_score =
493  internal::evaluate_one_site_3(params_.k,
494  params_.r,
495  rbi0, rbi1,
496  sites0_[i], sites1_[j],
497  g0, g1,
498  da,
499  sphere_derivatives_table,
500  torques_tables);
501 
502  sum += cur_score;
503  if(contacts_accumulator && cur_score!=0.0){
504  n_contacts++;
505  occupied_sites0[i]++;
506  occupied_sites1[i]++;
507  }
508  }// j
509  }// i
510  } // else
511  if(contacts_accumulator){
512  double non_specific_range= P::get_range_attraction();
513  double d_spheres= IMP::algebra::get_distance(spheres_table[pip[0].get_index()],
514  spheres_table[pip[1].get_index()]);
515  bool is_nonspecific_interaction= d_spheres < non_specific_range;
516  (*contacts_accumulator)=
517  boost::make_tuple(n_contacts,
518  occupied_sites0,
519  occupied_sites1,
520  is_nonspecific_interaction);
521  }
522 
523  IMP_LOG_PROGRESS( "Sum " << sum << std::endl);
524  return sum;
525 }
526 
527 
528 //!
529 inline double
530 SitesPairScore::evaluate_site_contributions
531 (Model* m,
532  const ParticleIndexPair &pip,
534  boost::tuple< unsigned int, std::vector<unsigned int>, std::vector<unsigned int>, bool >
535  (*contacts_accumulator)
536  ) const
537 {
538  // Get internal tables
539  algebra::Sphere3D const* spheres_table=
540  m->access_spheres_data();
541  double const* quaternions_tables[4];
542  for(unsigned int i = 0; i < 4; i++){
543  quaternions_tables[i]=
544  core::RigidBody::access_quaternion_i_data(m, i);
545  }
546  algebra::Sphere3D* sphere_derivatives_table=
547  m->access_sphere_derivatives_data();
548  double* torques_tables[3];
549  for(unsigned int i = 0; i < 3; i++){
550  torques_tables[i]=
551  core::RigidBody::access_torque_i_data(m, i);
552  }
553  // evaluate:
554  return evaluate_site_contributions_with_internal_tables
555  (spheres_table,
556  quaternions_tables,
557  sphere_derivatives_table,
558  torques_tables,
559  pip,
560  da,
561  contacts_accumulator);
562 }
563 
564 #endif // ifndef SWIG
565 
566 //!
567 inline internal::RigidBodyInfo
568 SitesPairScore::get_rigid_body_info
569 (algebra::Sphere3D const* spheres_table,
570  double const **quaternions_tables,
571  ParticleIndex pi) const
572 {
573  // TODO: add usage check that it has valid quaternions
574  // IMP_USAGE_CHECK(core::RigidBody::get_is_setup(m, pi),
575  // "PI " << pi.get_index() << " not a rigid body");
576  /* if(is_cache_active_){ */
577  /* std::pair<t_particles_rb_cache::iterator, bool> */
578  /* p = particles_rb_cache_.insert */
579  /* (std::make_pair(pi, internal::RigidBodyInfo())); */
580  /* internal::RigidBodyInfo& rbi_cached = p.first->second; */
581  /* bool const rbi_in_cache = !p.second; */
582  /* if(!rbi_in_cache || */
583  /* rbi_cached.cache_id != cur_cache_id_) // = cached version is outdated */
584  /* { */
585  /* rbi_cached.set_particle(spheres_table, */
586  /* quaternions_tables, */
587  /* pi, */
588  /* cur_cache_id_); */
589  /* } */
590  /* return rbi_cached; */
591  /* } */
592  /* else // if is_cache_active_ */
593  /* { */
594  return internal::RigidBodyInfo(spheres_table,
595  quaternions_tables,
596  pi,
597  INVALID_CACHE_ID);
598  /* } */
599 }
600 
601 
602 
603 IMPNPCTRANSPORT_END_NAMESPACE
604 
605 #endif /* IMPNPCTRANSPORT_SITES_PAIR_SCORE_H */
ModelObjectsTemp do_get_inputs(Model *m, const ParticleIndexes &pis) const
Overload this method to specify the inputs.
Apply a PairScore to each Pair in a list.
virtual double evaluate_if_good_indexes(Model *m, const ParticleIndexPairs &o, DerivativeAccumulator *da, double max, unsigned int lower_bound, unsigned int upper_bound) const
#define IMP_FINAL
Have the compiler report an error if anything overrides this method.
#define IMP_OBJECT_METHODS(Name)
Define the basic things needed by any Object.
Definition: object_macros.h:25
#define IMP_OBJECT_LOG
Set the log level to the object's log level.
Definition: log_macros.h:297
Single variable function.
#define IMP_LOG_PROGRESS(expr)
Definition: log_macros.h:105
A particle with a user-defined type.
Macros to choose the best set or map for different purposes.
Class for storing model, its restraints, constraints, and particles.
Definition: Model.h:72
virtual double evaluate_indexes(Model *m, const ParticleIndexPairs &pips, DerivativeAccumulator *da, unsigned int lower_bound, unsigned int upper_bound) const
double get_magnitude_and_normalize_in_place(VT &vt)
Returns the magnitude of vt and turns it to a unit vector in place.
Definition: VectorBaseD.h:260
Apply a function to the distance between two particles with a set of specific binding sites...
Vector3D get_vector_product(const Vector3D &p1, const Vector3D &p2)
Return the vector product (cross product) of two vectors.
Definition: Vector3D.h:31
Ints get_index(const ParticlesTemp &particles, const Subset &subset, const Subsets &excluded)
Define PairScore.
Represent an XYZR particle with a sphere.
Simple 3D transformation class.
A score on the distance between the surfaces of two spheres.
A nullptr-initialized pointer to an IMP Object.
double get_sites_range() const
return the range for site-site attraction
double get_sites_k() const
return the k for site-site attraction
VectorD< 3 > Vector3D
Definition: VectorD.h:421
double evaluate_index(algebra::Sphere3D const &s0, algebra::Sphere3D const &s1, algebra::Sphere3D &ds0, algebra::Sphere3D &ds1, DerivativeAccumulator *da) const
Functions to search over vectors.
A Score on the distance between a pair of particles.
A summary of useful information about rigid bodies and their transformation for eg, caching purposes for SitesPairScore.
Decorator for a sphere-like particle.
double get_distance(const Line3D &s, const Vector3D &p)
Get closest distance between a line and a point.
#define IMP_OVERRIDE
Cause a compile error if this method does not override a parent method.
Class for adding derivatives from restraints to the model.
Estimates of various physical quantities.
Various important functionality for implementing decorators.