IMP  2.3.1
The Integrative Modeling Platform
vonMises.h
Go to the documentation of this file.
1 /**
2  * \file IMP/isd/vonMises.h \brief Normal distribution of Function
3  *
4  * Copyright 2007-2014 IMP Inventors. All rights reserved.
5  */
6 
7 #ifndef IMPISD_VON_MISES_H
8 #define IMPISD_VON_MISES_H
9 
10 #include <IMP/isd/isd_config.h>
11 #include <IMP/macros.h>
12 #include <IMP/kernel/Model.h>
13 #include <IMP/constants.h>
14 #include <math.h>
15 #include <boost/math/special_functions/bessel.hpp>
16 
17 IMPISD_BEGIN_NAMESPACE
18 
19 //! vonMises
20 /** Probability density function and -log(p) of von Mises distribution
21  \f[ f(x|\mu,\kappa) = \frac{\exp \left(\kappa \cos (x-\mu)\right)}{2\pi
22  I_0(\kappa)} \f]
23  This is the proper treatment for a "normally distributed" angle.
24  When \f$\kappa\f$ becomes infinite, the distribution tends to a gaussian,
25  and \f$\kappa = 1/\sigma^2\f$.
26  */
27 
28 class vonMises : public base::Object {
29  public:
30  vonMises(double x, double mu, double kappa)
31  : base::Object("von Mises %1%"), x_(x), mu_(mu) {
32  force_set_kappa(kappa);
33  }
34 
35  /* energy (score) functions, aka -log(p) */
36  virtual double evaluate() const { return logterm_ - kappa_ * cos(x_ - mu_); }
37 
38  virtual double evaluate_derivative_x() const {
39  return kappa_ * sin(x_ - mu_);
40  }
41 
42  virtual double evaluate_derivative_mu() const {
43  return -kappa_ * sin(x_ - mu_);
44  }
45 
46  virtual double evaluate_derivative_kappa() const {
47  return -cos(x_ - mu_) + I1_ / I0_;
48  }
49 
50  /* probability density function */
51  virtual double density() const {
52  return exp(kappa_ * cos(x_ - mu_)) / (2 * IMP::PI * I0_);
53  }
54 
55  /* change of parameters */
56  void set_x(double x) { x_ = x; }
57  void set_mu(double mu) { mu_ = mu; }
58  void set_kappa(double kappa) {
59  if (kappa_ != kappa) {
60  force_set_kappa(kappa);
61  }
62  }
63 
65  /*IMP_OBJECT_INLINE(vonMises, out << "vonMises: " << x_ << ", " << mu_
66  << ", " << kappa_ <<std::endl, {});*/
67 
68  private:
69  void force_set_kappa(double kappa) {
70  kappa_ = kappa;
71  I0_ = boost::math::cyl_bessel_i(0, kappa);
72  I1_ = boost::math::cyl_bessel_i(1, kappa);
73  logterm_ = log(2 * IMP::PI * I0_);
74  }
75  double x_, mu_, kappa_, I0_, I1_, logterm_;
76 };
77 
78 IMPISD_END_NAMESPACE
79 
80 #endif /* IMPISD_VON_MISES_H */
Declare an efficient stl-compatible map.
static const double PI
the constant pi
Import IMP/kernel/constants.h in the namespace.
#define IMP_OBJECT_METHODS(Name)
Define the basic things needed by any Object.
Definition: object_macros.h:25
Import IMP/kernel/macros.h in the namespace.
Storage of a model, its restraints, constraints and particles.
Common base class for heavy weight IMP objects.
Definition: Object.h:106
vonMises
Definition: vonMises.h:28