Stan  2.5.0
probability, sampling & optimization
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
dense_e_metric.hpp
Go to the documentation of this file.
1 #ifndef STAN__MCMC__DENSE__E__METRIC__BETA
2 #define STAN__MCMC__DENSE__E__METRIC__BETA
3 
4 #include <boost/random/variate_generator.hpp>
5 #include <boost/random/normal_distribution.hpp>
6 
8 #include <Eigen/Cholesky>
9 
13 
14 namespace stan {
15 
16  namespace mcmc {
17 
18  // Euclidean manifold with dense metric
19  template <typename M, typename BaseRNG>
20  class dense_e_metric: public base_hamiltonian<M, dense_e_point, BaseRNG> {
21 
22  public:
23 
24  dense_e_metric(M& m, std::ostream* e):
25  base_hamiltonian<M, dense_e_point, BaseRNG>(m, e) {};
27 
28  double T(dense_e_point& z) {
29  return 0.5 * z.p.transpose() * z.mInv * z.p;
30  }
31 
32  double tau(dense_e_point& z) { return T(z); }
33  double phi(dense_e_point& z) { return this->V(z); }
34 
35  const Eigen::VectorXd dtau_dq(dense_e_point& z) {
36  return Eigen::VectorXd::Zero(this->model_.num_params_r());
37  }
38 
39  const Eigen::VectorXd dtau_dp(dense_e_point& z) {
40  return z.mInv * z.p;
41  }
42 
43  const Eigen::VectorXd dphi_dq(dense_e_point& z) {
44  return z.g;
45  }
46 
47  void sample_p(dense_e_point& z, BaseRNG& rng) {
48  typedef typename stan::math::index_type<Eigen::VectorXd>::type idx_t;
49  boost::variate_generator<BaseRNG&, boost::normal_distribution<> >
50  rand_dense_gaus(rng, boost::normal_distribution<>());
51 
52  Eigen::VectorXd u(z.p.size());
53 
54  for (idx_t i = 0; i < u.size(); ++i)
55  u(i) = rand_dense_gaus();
56 
57  z.p = z.mInv.llt().matrixL().solve(u);
58 
59  }
60 
61  };
62 
63  } // mcmc
64 
65 } // stan
66 
67 
68 #endif
const Eigen::VectorXd dtau_dp(dense_e_point &z)
void sample_p(dense_e_point &z, BaseRNG &rng)
dense_e_metric(M &m, std::ostream *e)
Primary template class for the metaprogram to compute the index type of a container.
Definition: index_type.hpp:21
Eigen::VectorXd p
Definition: ps_point.hpp:48
const Eigen::VectorXd dphi_dq(dense_e_point &z)
Eigen::VectorXd g
Definition: ps_point.hpp:51
double tau(dense_e_point &z)
double e()
Return the base of the natural logarithm.
Definition: constants.hpp:86
double phi(dense_e_point &z)
const Eigen::VectorXd dtau_dq(dense_e_point &z)
double T(dense_e_point &z)

     [ Stan Home Page ] © 2011–2014, Stan Development Team.