-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathmcmclib_docs_rmhmc.html
331 lines (245 loc) · 11.5 KB
/
mcmclib_docs_rmhmc.html
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
<!DOCTYPE html>
<html lang="en">
<head>
<meta charset="utf-8">
<meta http-equiv="X-UA-Compatible" content="IE=edge">
<meta name="viewport" content="width=device-width, initial-scale=1">
<meta name="description" content="MCMCLib: a C++ MCMC library.">
<meta name="author" content="Keith O'Hara">
<meta name="keywords" content="mcmc, MCMC, Metropolis Hastings, RWMH, Differential Evolution, Metropolis-adjusted Langevin algorithm, MALA, Hamiltonian Monte Carlo, HMC, C++, C++11, Cpp, NYU, New York University, Econometrics, Research" />
<link rel="shortcut icon" type="image/x-icon" href="siteicon.ico">
<title>MCMC: RM-HMC</title>
<!-- Bootstrap Core CSS -->
<link href="css/bootstrap.min.css" rel="stylesheet">
<!-- Custom CSS -->
<link href="css/modern-business.css" rel="stylesheet">
<!-- Custom Fonts -->
<link href="font-awesome/css/font-awesome.min.css" rel="stylesheet" type="text/css">
<!-- Additional Settings -->
<link href="css/yuxuan_settings.css" rel="stylesheet">
<!-- Syntax Highlighter -->
<script type="text/javascript" src="js/syntaxhighlighter.js"></script>
<link type="text/css" rel="stylesheet" href="css/swift_theme.css">
<script>
(function(i,s,o,g,r,a,m){i['GoogleAnalyticsObject']=r;i[r]=i[r]||function(){
(i[r].q=i[r].q||[]).push(arguments)},i[r].l=1*new Date();a=s.createElement(o),
m=s.getElementsByTagName(o)[0];a.async=1;a.src=g;m.parentNode.insertBefore(a,m)
})(window,document,'script','https://www.google-analytics.com/analytics.js','ga');
ga('create', 'UA-93902857-1', 'auto');
ga('send', 'pageview');
</script>
<!-- MathJax -->
<script type="text/x-mathjax-config">
MathJax.Hub.Config({tex2jax: {inlineMath: [['$','$'], ['\\(','\\)']]}});
</script>
<script type="text/javascript" async
src="https://cdn.mathjax.org/mathjax/latest/MathJax.js?config=TeX-AMS_CHTML">
</script>
<script async defer src="https://buttons.github.io/buttons.js"></script>
<script src="js/jquery.js"></script>
<script>
$(function(){
$("#mynavbar").load("navbar.html")
$("#mcmchead").load("mcmclib_header.html")
$("#myfooter").load("footer.html")
});
</script>
</head>
<style>
pre {
display: inline-block;
}
</style>
<body>
<!-- Navigation -->
<div id="mynavbar"></div>
<!-- Page Content -->
<div class="container">
<!-- Page Heading/Breadcrumbs -->
<div id="mcmchead"></div>
<br>
<div class="row">
<div class="col-md-2"></div>
<div class="col-md-7">
<div class="panel panel-default">
<div class="panel-heading">
<a data-toggle="collapse" href="#collapse1"><h4><strong style="font-size: 120%;">MCMC: Riemannian Manifold Hamiltonian Monte Carlo</strong></h4></a>
</div>
<div id="collapse1" class="panel-collapse collapse">
<div class="panel-body">
<a href="#definition">Definition</a> <br>
<a href="#details">Details</a> <br>
<a href="#examples">Examples</a>
</div>
</div>
</div>
</div>
</div>
<!-- -->
<div class="row">
<div class="col-md-2"></div>
<div class="col-md-8">
<p>MCMC via Hamiltonian dynamics.</p>
<hr style="height:2px;border-width:0;background-color:black">
<h3 style="text-align: left;" id="definition"><strong style="font-size: 100%;">Definition and Syntax</strong></h3>
<pre class="brush: cpp;">
bool rmhmc(const arma::vec& initial_vals, arma::mat& draws_out, std::function<double (const arma::vec& vals_inp, arma::vec* grad_out, void* target_data)> target_log_kernel, void* target_data, std::function<arma::mat (const arma::vec& vals_inp, arma::cube* tensor_deriv_out, void* tensor_data)> tensor_fn, void* tensor_data);
bool rmhmc(const arma::vec& initial_vals, arma::mat& draws_out, std::function<double (const arma::vec& vals_inp, arma::vec* grad_out, void* target_data)> target_log_kernel, void* target_data, std::function<arma::mat (const arma::vec& vals_inp, arma::cube* tensor_deriv_out, void* tensor_data)> tensor_fn, void* tensor_data, algo_settings_t& settings);
</pre>
<p><strong>Function arguments:</strong></p>
<ul>
<li><code>initial_vals</code> a column vector of initial values.</li>
<li><code>draws_out</code> a two-dimensional array containing the posterior draws.</li>
<li><code>target_log_kernel</code> the target log-posterior kernel function, taking three arguments:
<ul>
<li><code>vals_inp</code> a vector of input values;</li>
<li><code>grad_out</code> a vector to store the gradient values; and</li>
<li><code>target_data</code> additional parameters passed to the function.</li>
</ul>
<li><code>target_data</code> additional parameters passed to the posterior kernel.</li>
<li><code>tensor_fn</code> metric tensor function, taking three arguments:
<ul>
<li><code>vals_inp</code> a vector of input values;</li>
<li><code>tensor_deriv_out</code> a three-dimensional array to store derivative values of the metric tensor function; and</li>
<li><code>tensor_data</code> additional parameters passed to the function.</li>
</ul>
<li><code>tensor_data</code> additional parameters passed to the tensor function.</li>
<li><code>settings</code> parameters controlling the MCMC routine; see below.</li>
</ul>
<p><strong>MCMC control parameters:</strong></p>
<ul>
<li><code>bool vals_bound</code> whether the search space is bounded. If true, then</li>
<ul>
<li><code>arma::vec lower_bounds</code> this defines the lower bounds.</li>
<li><code>arma::vec upper_bounds</code> this defines the upper bounds.</li>
</ul>
<li><code>int hmc_n_draws</code> number of posterior draws to keep.</li>
<li><code>int hmc_n_burnin</code> number of burnin draws.</li>
<li><code>double hmc_step_size</code> step size $\epsilon$ below.</li>
<li><code>int hmc_leap_steps</code> number of leapfrog steps.</li>
<li><code>int rmhmc_fp_steps</code> number of fixed-point iterations.</li>
</ul>
<hr style="height:2px;border-width:0;background-color:black">
<h3 style="text-align: left;" id="details"><strong style="font-size: 100%;">Details</strong></h3>
<p>Let $\theta^{(i)}$ denote a $d$-dimensional vector of stored values at stage $i$ of the algorithm. The basic HMC algorithm consists of three steps.</p>
<ul>
<li><strong>Initializaton.</strong> Denote the Hamiltonian by
$$
H \left\{ \theta, p \right\} := - \ln K(\theta | X) + \frac{1}{2} \log \left\{ (2 \pi)^d | \mathbf{G}(\theta) | \right\} + \frac{1}{2} p^\top \mathbf{G}^{-1}(\theta) p
$$
where $\mathbf{G}$ is the metric tensor function. Sample $p^{(i)} \sim N(0,\mathbf{G}(\theta^{(i)}))$, and set $\theta^{(*)} = \theta^{(i)}$ and $p^{(*)} = p^{(i)}$.
<li><strong>Leapfrog Steps.</strong> <br> For $k = 1, \ldots,$ <code>hmc_leap_steps</code> do
<ul>
<li><strong>Initializaton.</strong> Set $\theta_o^{(*)} = \theta^{(*)}$ and $p_o^{(*)} = p^{(*)}$.</li>
<li><strong>Momentum Update Step.</strong> <br> For $j = 1, \ldots,$ <code>rmhmc_fp_steps</code> do
$$p_h^{(*)} = p_o^{(*)} - \dfrac{\epsilon}{2} \times \nabla_\theta H \left\{ \theta_o^{(*)},p_h^{(*)} \right\}$$
where the subscript $h$ denotes a half-step. Notice that $p_h^{(*)}$ appears on both sides of this expression.</li>
<li><strong>Position (Proposal) Update Step.</strong> <br> For $j = 1, \ldots,$ <code>rmhmc_fp_steps</code> do
$$\theta^{(*)} = \theta_o^{(*)} + \dfrac{\epsilon}{2} \times \left[ \nabla_p H \left\{ \theta_o^{(*)},p_h^{(*)} \right\} + \nabla_p H \left\{ \theta^{(*)},p_h^{(*)} \right\} \right] $$
Notice that $\theta^{(*)}$ appears on both sides of this expression.</li>
<li><strong>Momentum Update Step.</strong>
$$p^{(*)} = p_h^{(*)} + \epsilon \times \nabla_\theta H \left\{ \theta^{(*)},p_h^{(*)} \right\}$$
</li>
</ul>
<li><strong>Accept/Reject Step.</strong> Define
$$
\alpha = \min \left\{ 1, \exp( H(\theta^{(i)}, p^{(i)}) - H(\theta^{(*)}, p^{(*)}) ) \right\}
$$
Then
$$
\theta^{(i+1)} = \begin{cases} \theta^{(*)} & \text{ if } Z < \alpha \\ \theta^{(i)} & \text{ else } \end{cases}
$$
where $Z \sim \text{Unif}(0,1)$.</li>
</ul>
<hr style="height:2px;border-width:0;background-color:black">
<h3 style="text-align: left;" id="examples"><strong style="font-size: 100%;">Examples</strong></h3>
<p>Normal likelihood.</p>
<br>
<pre class="brush: cpp;">
#include "mcmc.hpp"
struct norm_data {
arma::vec x;
};
double ll_dens(const arma::vec& vals_inp, arma::vec* grad_out, void* ll_data)
{
const double mu = vals_inp(0);
const double sigma = vals_inp(1);
const double pi = arma::datum::pi;
norm_data* dta = reinterpret_cast<norm_data*>(ll_data);
const arma::vec x = dta->x;
const int n_vals = x.n_rows;
//
const double ret = - ((double) n_vals) * (0.5*std::log(2*pi) + std::log(sigma)) - arma::accu( arma::pow(x - mu,2) / (2*sigma*sigma) );
//
if (grad_out)
{
grad_out->set_size(2,1);
//
double m_1 = arma::accu(x - mu);
double m_2 = arma::accu( arma::pow(x - mu,2) );
(*grad_out)(0,0) = m_1 / (sigma*sigma);
(*grad_out)(1,0) = (m_2 / (sigma*sigma*sigma)) - ((double) n_vals) / sigma;
}
//
return ret;
}
arma::mat tensor_fn(const arma::vec& vals_inp, arma::cube* tensor_deriv_out, void* tensor_data)
{
// const double mu = vals_inp(0);
const double sigma = vals_inp(1);
norm_data* dta = reinterpret_cast<norm_data*>(ll_data);
const int n_vals = dta->x.n_rows;
//
const double sigma_sq = sigma*sigma;
arma::mat tensor_out = arma::zeros(2,2);
tensor_out(0,0) = ((double) n_vals) / sigma_sq;
tensor_out(1,1) = 2.0 * ((double) n_vals) / sigma_sq;
//
if (tensor_deriv_out) {
tensor_deriv_out->set_size(2,2,2);
//
tensor_deriv_out->slice(0).zeros();
tensor_deriv_out->slice(1) = -2.0*tensor_out / sigma;
}
//
return tensor_out;
}
double log_target_dens(const arma::vec& vals_inp, arma::vec* grad_out, void* ll_data)
{
return ll_dens(vals_inp,grad_out,ll_data);
}
int main()
{
const int n_data = 200;
const double mu = 2.0;
const double sigma = 2.0;
norm_data dta;
arma::vec x_dta = mu + sigma*arma::randn(n_data,1);
dta.x = x_dta;
arma::vec initial_val(2);
initial_val(0) = mu + 1.0; // mu
initial_val(1) = sigma + 1.0; // sigma
mcmc::algo_settings_t settings;
settings.hmc_step_size = 0.1;
settings.hmc_n_burnin = 1000;
settings.hmc_n_draws = 5000;
printf("begin run\n");
arma::mat draws_out;
mcmc::rmhmc(initial_val,draws_out,log_target_dens,&dta,tensor_fn,&dta,settings);
arma::cout << "draws:\n" << draws_out.rows(0,9) << arma::endl;
std::cout << "acceptance rate = " << settings.hmc_accept_rate << arma::endl;
std::cout << "mcmc mean = " << arma::mean(draws_out) << std::endl;
std::cout << "mcmc var:\n" << arma::cov(draws_out) << std::endl;
return 0;
}
</pre>
</div>
</div>
</div>
<div id="myfooter"></div>
<!-- jQuery -->
<!--<script src="js/jquery.js"></script>-->
<!-- Bootstrap Core JavaScript -->
<script src="js/bootstrap.min.js"></script>
</body>
</html>