/
hmc_nuts_diag_e_adapt.hpp
660 lines (639 loc) · 32.7 KB
/
hmc_nuts_diag_e_adapt.hpp
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
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
#ifndef STAN_SERVICES_SAMPLE_HMC_NUTS_DIAG_E_ADAPT_HPP
#define STAN_SERVICES_SAMPLE_HMC_NUTS_DIAG_E_ADAPT_HPP
#include <stan/callbacks/interrupt.hpp>
#include <stan/callbacks/logger.hpp>
#include <stan/callbacks/structured_writer.hpp>
#include <stan/callbacks/writer.hpp>
#include <stan/io/var_context.hpp>
#include <stan/math/prim.hpp>
#include <stan/mcmc/hmc/nuts/adapt_diag_e_nuts.hpp>
#include <stan/services/error_codes.hpp>
#include <stan/services/util/create_rng.hpp>
#include <stan/services/util/inv_metric.hpp>
#include <stan/services/util/initialize.hpp>
#include <stan/services/util/run_adaptive_sampler.hpp>
#include <vector>
namespace stan {
namespace services {
namespace sample {
/**
* Runs HMC with NUTS with adaptation using diagonal Euclidean metric
* with a pre-specified diagonal metric and saves adapted tuning parameters.
*
* @tparam Model Model class
* @param[in] model Input model (with data already instantiated)
* @param[in] init var context for initialization
* @param[in] init_inv_metric var context exposing an initial diagonal
* inverse Euclidean metric (must be positive definite)
* @param[in] random_seed random seed for the random number generator
* @param[in] chain chain id to advance the pseudo random number generator
* @param[in] init_radius radius to initialize
* @param[in] num_warmup Number of warmup samples
* @param[in] num_samples Number of samples
* @param[in] num_thin Number to thin the samples
* @param[in] save_warmup Indicates whether to save the warmup iterations
* @param[in] refresh Controls the output
* @param[in] stepsize initial stepsize for discrete evolution
* @param[in] stepsize_jitter uniform random jitter of stepsize
* @param[in] max_depth Maximum tree depth
* @param[in] delta adaptation target acceptance statistic
* @param[in] gamma adaptation regularization scale
* @param[in] kappa adaptation relaxation exponent
* @param[in] t0 adaptation iteration offset
* @param[in] init_buffer width of initial fast adaptation interval
* @param[in] term_buffer width of final fast adaptation interval
* @param[in] window initial width of slow adaptation interval
* @param[in,out] interrupt Callback for interrupts
* @param[in,out] logger Logger for messages
* @param[in,out] init_writer Writer callback for unconstrained inits
* @param[in,out] sample_writer Writer for draws
* @param[in,out] diagnostic_writer Writer for diagnostic information
* @param[in,out] metric_writer Writer for tuning params
* @return error_codes::OK if successful
*/
template <class Model>
int hmc_nuts_diag_e_adapt(
Model& model, const stan::io::var_context& init,
const stan::io::var_context& init_inv_metric, unsigned int random_seed,
unsigned int chain, double init_radius, int num_warmup, int num_samples,
int num_thin, bool save_warmup, int refresh, double stepsize,
double stepsize_jitter, int max_depth, double delta, double gamma,
double kappa, double t0, unsigned int init_buffer, unsigned int term_buffer,
unsigned int window, callbacks::interrupt& interrupt,
callbacks::logger& logger, callbacks::writer& init_writer,
callbacks::writer& sample_writer, callbacks::writer& diagnostic_writer,
callbacks::structured_writer& metric_writer) {
stan::rng_t rng = util::create_rng(random_seed, chain);
std::vector<double> cont_vector;
Eigen::VectorXd inv_metric;
try {
cont_vector = util::initialize(model, init, rng, init_radius, true, logger,
init_writer);
inv_metric = util::read_diag_inv_metric(init_inv_metric,
model.num_params_r(), logger);
util::validate_diag_inv_metric(inv_metric, logger);
} catch (const std::exception& e) {
logger.error(e.what());
return error_codes::CONFIG;
}
stan::mcmc::adapt_diag_e_nuts<Model, stan::rng_t> sampler(model, rng);
sampler.set_metric(inv_metric);
sampler.set_nominal_stepsize(stepsize);
sampler.set_stepsize_jitter(stepsize_jitter);
sampler.set_max_depth(max_depth);
sampler.get_stepsize_adaptation().set_mu(log(10 * stepsize));
sampler.get_stepsize_adaptation().set_delta(delta);
sampler.get_stepsize_adaptation().set_gamma(gamma);
sampler.get_stepsize_adaptation().set_kappa(kappa);
sampler.get_stepsize_adaptation().set_t0(t0);
sampler.set_window_params(num_warmup, init_buffer, term_buffer, window,
logger);
try {
util::run_adaptive_sampler(sampler, model, cont_vector, num_warmup,
num_samples, num_thin, refresh, save_warmup, rng,
interrupt, logger, sample_writer,
diagnostic_writer, metric_writer);
} catch (const std::exception& e) {
logger.error(e.what());
return error_codes::SOFTWARE;
}
return error_codes::OK;
}
/**
* Runs HMC with NUTS with adaptation using diagonal Euclidean metric
* with a pre-specified diagonal metric.
*
* @tparam Model Model class
* @param[in] model Input model (with data already instantiated)
* @param[in] init var context for initialization
* @param[in] init_inv_metric var context exposing an initial diagonal
* inverse Euclidean metric (must be positive definite)
* @param[in] random_seed random seed for the random number generator
* @param[in] chain chain id to advance the pseudo random number generator
* @param[in] init_radius radius to initialize
* @param[in] num_warmup Number of warmup samples
* @param[in] num_samples Number of samples
* @param[in] num_thin Number to thin the samples
* @param[in] save_warmup Indicates whether to save the warmup iterations
* @param[in] refresh Controls the output
* @param[in] stepsize initial stepsize for discrete evolution
* @param[in] stepsize_jitter uniform random jitter of stepsize
* @param[in] max_depth Maximum tree depth
* @param[in] delta adaptation target acceptance statistic
* @param[in] gamma adaptation regularization scale
* @param[in] kappa adaptation relaxation exponent
* @param[in] t0 adaptation iteration offset
* @param[in] init_buffer width of initial fast adaptation interval
* @param[in] term_buffer width of final fast adaptation interval
* @param[in] window initial width of slow adaptation interval
* @param[in,out] interrupt Callback for interrupts
* @param[in,out] logger Logger for messages
* @param[in,out] init_writer Writer callback for unconstrained inits
* @param[in,out] sample_writer Writer for draws
* @param[in,out] diagnostic_writer Writer for diagnostic information
* @return error_codes::OK if successful
*/
template <class Model>
int hmc_nuts_diag_e_adapt(
Model& model, const stan::io::var_context& init,
const stan::io::var_context& init_inv_metric, unsigned int random_seed,
unsigned int chain, double init_radius, int num_warmup, int num_samples,
int num_thin, bool save_warmup, int refresh, double stepsize,
double stepsize_jitter, int max_depth, double delta, double gamma,
double kappa, double t0, unsigned int init_buffer, unsigned int term_buffer,
unsigned int window, callbacks::interrupt& interrupt,
callbacks::logger& logger, callbacks::writer& init_writer,
callbacks::writer& sample_writer, callbacks::writer& diagnostic_writer) {
callbacks::structured_writer dummy_metric_writer;
return hmc_nuts_diag_e_adapt(
model, init, init_inv_metric, random_seed, chain, init_radius, num_warmup,
num_samples, num_thin, save_warmup, refresh, stepsize, stepsize_jitter,
max_depth, delta, gamma, kappa, t0, init_buffer, term_buffer, window,
interrupt, logger, init_writer, sample_writer, diagnostic_writer,
dummy_metric_writer);
}
/**
* Runs HMC with NUTS with adaptation using diagonal Euclidean metric,
* with identity matrix as initial inv_metric and saves adapted tuning
* parameters stepsize and inverse metric.
*
* @tparam Model Model class
* @param[in] model Input model (with data already instantiated)
* @param[in] init var context for initialization
* @param[in] random_seed random seed for the random number generator
* @param[in] chain chain id to advance the pseudo random number generator
* @param[in] init_radius radius to initialize
* @param[in] num_warmup Number of warmup samples
* @param[in] num_samples Number of samples
* @param[in] num_thin Number to thin the samples
* @param[in] save_warmup Indicates whether to save the warmup iterations
* @param[in] refresh Controls the output
* @param[in] stepsize initial stepsize for discrete evolution
* @param[in] stepsize_jitter uniform random jitter of stepsize
* @param[in] max_depth Maximum tree depth
* @param[in] delta adaptation target acceptance statistic
* @param[in] gamma adaptation regularization scale
* @param[in] kappa adaptation relaxation exponent
* @param[in] t0 adaptation iteration offset
* @param[in] init_buffer width of initial fast adaptation interval
* @param[in] term_buffer width of final fast adaptation interval
* @param[in] window initial width of slow adaptation interval
* @param[in,out] interrupt Callback for interrupts
* @param[in,out] logger Logger for messages
* @param[in,out] init_writer Writer callback for unconstrained inits
* @param[in,out] sample_writer Writer for draws
* @param[in,out] diagnostic_writer Writer for diagnostic information
* @param[in,out] metric_writer Writer for tuning params
* @return error_codes::OK if successful
*/
template <class Model>
int hmc_nuts_diag_e_adapt(
Model& model, const stan::io::var_context& init, unsigned int random_seed,
unsigned int chain, double init_radius, int num_warmup, int num_samples,
int num_thin, bool save_warmup, int refresh, double stepsize,
double stepsize_jitter, int max_depth, double delta, double gamma,
double kappa, double t0, unsigned int init_buffer, unsigned int term_buffer,
unsigned int window, callbacks::interrupt& interrupt,
callbacks::logger& logger, callbacks::writer& init_writer,
callbacks::writer& sample_writer, callbacks::writer& diagnostic_writer,
callbacks::structured_writer& metric_writer) {
auto default_metric
= util::create_unit_e_diag_inv_metric(model.num_params_r());
return hmc_nuts_diag_e_adapt(
model, init, default_metric, random_seed, chain, init_radius, num_warmup,
num_samples, num_thin, save_warmup, refresh, stepsize, stepsize_jitter,
max_depth, delta, gamma, kappa, t0, init_buffer, term_buffer, window,
interrupt, logger, init_writer, sample_writer, diagnostic_writer,
metric_writer);
}
/**
* Runs HMC with NUTS with adaptation using diagonal Euclidean metric,
* with identity matrix as initial inv_metric.
*
* @tparam Model Model class
* @param[in] model Input model (with data already instantiated)
* @param[in] init var context for initialization
* @param[in] random_seed random seed for the random number generator
* @param[in] chain chain id to advance the pseudo random number generator
* @param[in] init_radius radius to initialize
* @param[in] num_warmup Number of warmup samples
* @param[in] num_samples Number of samples
* @param[in] num_thin Number to thin the samples
* @param[in] save_warmup Indicates whether to save the warmup iterations
* @param[in] refresh Controls the output
* @param[in] stepsize initial stepsize for discrete evolution
* @param[in] stepsize_jitter uniform random jitter of stepsize
* @param[in] max_depth Maximum tree depth
* @param[in] delta adaptation target acceptance statistic
* @param[in] gamma adaptation regularization scale
* @param[in] kappa adaptation relaxation exponent
* @param[in] t0 adaptation iteration offset
* @param[in] init_buffer width of initial fast adaptation interval
* @param[in] term_buffer width of final fast adaptation interval
* @param[in] window initial width of slow adaptation interval
* @param[in,out] interrupt Callback for interrupts
* @param[in,out] logger Logger for messages
* @param[in,out] init_writer Writer callback for unconstrained inits
* @param[in,out] sample_writer Writer for draws
* @param[in,out] diagnostic_writer Writer for diagnostic information
* @return error_codes::OK if successful
*/
template <class Model>
int hmc_nuts_diag_e_adapt(
Model& model, const stan::io::var_context& init, unsigned int random_seed,
unsigned int chain, double init_radius, int num_warmup, int num_samples,
int num_thin, bool save_warmup, int refresh, double stepsize,
double stepsize_jitter, int max_depth, double delta, double gamma,
double kappa, double t0, unsigned int init_buffer, unsigned int term_buffer,
unsigned int window, callbacks::interrupt& interrupt,
callbacks::logger& logger, callbacks::writer& init_writer,
callbacks::writer& sample_writer, callbacks::writer& diagnostic_writer) {
auto default_metric
= util::create_unit_e_diag_inv_metric(model.num_params_r());
callbacks::structured_writer dummy_metric_writer;
return hmc_nuts_diag_e_adapt(
model, init, default_metric, random_seed, chain, init_radius, num_warmup,
num_samples, num_thin, save_warmup, refresh, stepsize, stepsize_jitter,
max_depth, delta, gamma, kappa, t0, init_buffer, term_buffer, window,
interrupt, logger, init_writer, sample_writer, diagnostic_writer,
dummy_metric_writer);
}
/**
* Runs multiple chains of HMC with NUTS with adaptation using diagonal
* Euclidean metric with a pre-specified diagonal metric and saves adapted
* tuning parameters stepsize and inverse metric.
*
* @tparam Model Model class
* @tparam InitContextPtr A pointer with underlying type derived from
* `stan::io::var_context`
* @tparam InitInvContextPtr A pointer with underlying type derived from
* `stan::io::var_context`
* @tparam InitWriter A type derived from `stan::callbacks::writer`
* @tparam SamplerWriter A type derived from `stan::callbacks::writer`
* @tparam DiagnosticWriter A type derived from `stan::callbacks::writer`
* @tparam MetricWriter A type derived from `stan::callbacks::structured_writer`
* @param[in] model Input model (with data already instantiated)
* @param[in] num_chains The number of chains to run in parallel. `init`,
* `init_inv_metric`, `init_writer`, `sample_writer`, and `diagnostic_writer`
* must be the same length as this value.
* @param[in] init A std vector of init var contexts for per-chain
* initialization.
* @param[in] init_inv_metric A std vector of var contexts exposing an initial
* diagonal inverse Euclidean metric for each chain (must be positive definite)
* @param[in] random_seed random seed for the random number generator
* @param[in] init_chain_id first chain id. The pseudo random number generator
* will advance for each chain by an integer sequence from `init_chain_id` to
* `init_chain_id + num_chains - 1`
* @param[in] init_radius radius to initialize
* @param[in] num_warmup Number of warmup samples
* @param[in] num_samples Number of samples
* @param[in] num_thin Number to thin the samples
* @param[in] save_warmup Indicates whether to save the warmup iterations
* @param[in] refresh Controls the output
* @param[in] stepsize initial stepsize for discrete evolution
* @param[in] stepsize_jitter uniform random jitter of stepsize
* @param[in] max_depth Maximum tree depth
* @param[in] delta adaptation target acceptance statistic
* @param[in] gamma adaptation regularization scale
* @param[in] kappa adaptation relaxation exponent
* @param[in] t0 adaptation iteration offset
* @param[in] init_buffer width of initial fast adaptation interval
* @param[in] term_buffer width of final fast adaptation interval
* @param[in] window initial width of slow adaptation interval
* @param[in,out] interrupt Callback for interrupts
* @param[in,out] logger Logger for messages
* @param[in,out] init_writer std vector of Writer callbacks for unconstrained
* inits of each chain.
* @param[in,out] sample_writer std vector of Writers for draws of each chain.
* @param[in,out] diagnostic_writer std vector of Writers for diagnostic
* information of each chain.
* @param[in,out] metric_writer std vector of Writers for tuning params
* @return error_codes::OK if successful
*/
template <class Model, typename InitContextPtr, typename InitInvContextPtr,
typename InitWriter, typename SampleWriter, typename DiagnosticWriter,
typename MetricWriter>
int hmc_nuts_diag_e_adapt(
Model& model, size_t num_chains, const std::vector<InitContextPtr>& init,
const std::vector<InitInvContextPtr>& init_inv_metric,
unsigned int random_seed, unsigned int init_chain_id, double init_radius,
int num_warmup, int num_samples, int num_thin, bool save_warmup,
int refresh, double stepsize, double stepsize_jitter, int max_depth,
double delta, double gamma, double kappa, double t0,
unsigned int init_buffer, unsigned int term_buffer, unsigned int window,
callbacks::interrupt& interrupt, callbacks::logger& logger,
std::vector<InitWriter>& init_writer,
std::vector<SampleWriter>& sample_writer,
std::vector<DiagnosticWriter>& diagnostic_writer,
std::vector<MetricWriter>& metric_writer) {
if (num_chains == 1) {
return hmc_nuts_diag_e_adapt(
model, *init[0], *init_inv_metric[0], random_seed, init_chain_id,
init_radius, num_warmup, num_samples, num_thin, save_warmup, refresh,
stepsize, stepsize_jitter, max_depth, delta, gamma, kappa, t0,
init_buffer, term_buffer, window, interrupt, logger, init_writer[0],
sample_writer[0], diagnostic_writer[0], metric_writer[0]);
}
using sample_t = stan::mcmc::adapt_diag_e_nuts<Model, stan::rng_t>;
std::vector<stan::rng_t> rngs;
rngs.reserve(num_chains);
std::vector<std::vector<double>> cont_vectors;
cont_vectors.reserve(num_chains);
std::vector<sample_t> samplers;
samplers.reserve(num_chains);
try {
for (int i = 0; i < num_chains; ++i) {
rngs.emplace_back(util::create_rng(random_seed, init_chain_id + i));
cont_vectors.emplace_back(util::initialize(
model, *init[i], rngs[i], init_radius, true, logger, init_writer[i]));
samplers.emplace_back(model, rngs[i]);
Eigen::VectorXd inv_metric = util::read_diag_inv_metric(
*init_inv_metric[i], model.num_params_r(), logger);
util::validate_diag_inv_metric(inv_metric, logger);
samplers[i].set_metric(inv_metric);
samplers[i].set_nominal_stepsize(stepsize);
samplers[i].set_stepsize_jitter(stepsize_jitter);
samplers[i].set_max_depth(max_depth);
samplers[i].get_stepsize_adaptation().set_mu(log(10 * stepsize));
samplers[i].get_stepsize_adaptation().set_delta(delta);
samplers[i].get_stepsize_adaptation().set_gamma(gamma);
samplers[i].get_stepsize_adaptation().set_kappa(kappa);
samplers[i].get_stepsize_adaptation().set_t0(t0);
samplers[i].set_window_params(num_warmup, init_buffer, term_buffer,
window, logger);
}
} catch (const std::exception& e) {
logger.error(e.what());
return error_codes::CONFIG;
}
try {
tbb::parallel_for(
tbb::blocked_range<size_t>(0, num_chains, 1),
[num_warmup, num_samples, num_thin, refresh, save_warmup, num_chains,
init_chain_id, &samplers, &model, &rngs, &interrupt, &logger,
&sample_writer, &cont_vectors, &diagnostic_writer,
&metric_writer](const tbb::blocked_range<size_t>& r) {
for (size_t i = r.begin(); i != r.end(); ++i) {
util::run_adaptive_sampler(
samplers[i], model, cont_vectors[i], num_warmup, num_samples,
num_thin, refresh, save_warmup, rngs[i], interrupt, logger,
sample_writer[i], diagnostic_writer[i], metric_writer[i],
init_chain_id + i, num_chains);
}
},
tbb::simple_partitioner());
} catch (const std::exception& e) {
logger.error(e.what());
return error_codes::SOFTWARE;
}
return error_codes::OK;
}
/**
* Runs multiple chains of HMC with NUTS with adaptation using diagonal
* Euclidean metric with a pre-specified diagonal metric.
*
* @tparam Model Model class
* @tparam InitContextPtr A pointer with underlying type derived from
* `stan::io::var_context`
* @tparam InitContextPtr A pointer with underlying type derived from
* `stan::io::var_context`
* @tparam InitWriter A type derived from `stan::callbacks::writer`
* @tparam SamplerWriter A type derived from `stan::callbacks::writer`
* @tparam DiagnosticWriter A type derived from `stan::callbacks::writer`
* @param[in] model Input model (with data already instantiated)
* @param[in] num_chains The number of chains to run in parallel. `init`,
* `init_writer`, `sample_writer`, and `diagnostic_writer` must be the same
* length as this value.
* @param[in] init A std vector of init var contexts for initialization
* of each chain.
* @param[in] init_inv_metric A std vector of var contexts exposing an initial
* diagonal inverse Euclidean metric for each chain (must be positive definite)
* @param[in] random_seed random seed for the random number generator
* @param[in] init_chain_id first chain id. The pseudo random number generator
* will advance by for each chain by an integer sequence from `init_chain_id` to
* `init_chain_id+num_chains-1`
* @param[in] init_radius radius to initialize
* @param[in] num_warmup Number of warmup samples
* @param[in] num_samples Number of samples
* @param[in] num_thin Number to thin the samples
* @param[in] save_warmup Indicates whether to save the warmup iterations
* @param[in] refresh Controls the output
* @param[in] stepsize initial stepsize for discrete evolution
* @param[in] stepsize_jitter uniform random jitter of stepsize
* @param[in] max_depth Maximum tree depth
* @param[in] delta adaptation target acceptance statistic
* @param[in] gamma adaptation regularization scale
* @param[in] kappa adaptation relaxation exponent
* @param[in] t0 adaptation iteration offset
* @param[in] init_buffer width of initial fast adaptation interval
* @param[in] term_buffer width of final fast adaptation interval
* @param[in] window initial width of slow adaptation interval
* @param[in,out] interrupt Callback for interrupts
* @param[in,out] logger Logger for messages
* @param[in,out] init_writer std vector of Writer callbacks for unconstrained
* inits of each chain.
* @param[in,out] sample_writer std vector of Writers for draws of each chain.
* @param[in,out] diagnostic_writer std vector of Writers for diagnostic
* information of each chain.
* @return error_codes::OK if successful
*/
template <class Model, typename InitContextPtr, typename InitInvContextPtr,
typename InitWriter, typename SampleWriter, typename DiagnosticWriter>
int hmc_nuts_diag_e_adapt(
Model& model, size_t num_chains, const std::vector<InitContextPtr>& init,
const std::vector<InitInvContextPtr>& init_inv_metric,
unsigned int random_seed, unsigned int init_chain_id, double init_radius,
int num_warmup, int num_samples, int num_thin, bool save_warmup,
int refresh, double stepsize, double stepsize_jitter, int max_depth,
double delta, double gamma, double kappa, double t0,
unsigned int init_buffer, unsigned int term_buffer, unsigned int window,
callbacks::interrupt& interrupt, callbacks::logger& logger,
std::vector<InitWriter>& init_writer,
std::vector<SampleWriter>& sample_writer,
std::vector<DiagnosticWriter>& diagnostic_writer) {
std::vector<stan::callbacks::structured_writer> dummy_metric_writer(
num_chains);
if (num_chains == 1) {
return hmc_nuts_diag_e_adapt(
model, *init[0], *init_inv_metric[0], random_seed, init_chain_id,
init_radius, num_warmup, num_samples, num_thin, save_warmup, refresh,
stepsize, stepsize_jitter, max_depth, delta, gamma, kappa, t0,
init_buffer, term_buffer, window, interrupt, logger, init_writer[0],
sample_writer[0], diagnostic_writer[0], dummy_metric_writer[0]);
}
return hmc_nuts_diag_e_adapt(
model, num_chains, init, init_inv_metric, random_seed, init_chain_id,
init_radius, num_warmup, num_samples, num_thin, save_warmup, refresh,
stepsize, stepsize_jitter, max_depth, delta, gamma, kappa, t0,
init_buffer, term_buffer, window, interrupt, logger, init_writer,
sample_writer, diagnostic_writer, dummy_metric_writer);
}
/**
* Runs multiple chains of HMC with NUTS with adaptation using diagonal
* with identity matrix as initial inv_metric and saves adapted tuning
* parameters stepsize and inverse metric.
*
* @tparam Model Model class
* @tparam InitContextPtr A pointer with underlying type derived from
* `stan::io::var_context`
* @tparam InitWriter A type derived from `stan::callbacks::writer`
* @tparam SamplerWriter A type derived from `stan::callbacks::writer`
* @tparam DiagnosticWriter A type derived from `stan::callbacks::writer`
* @tparam MetricWriter A type derived from `stan::callbacks::structured_writer`
* @param[in] model Input model (with data already instantiated)
* @param[in] num_chains The number of chains to run in parallel. `init`,
* `init_writer`, `sample_writer`, and `diagnostic_writer` must be the same
* length as this value.
* @param[in] init A std vector of init var contexts for initialization of each
* chain.
* @param[in] random_seed random seed for the random number generator
* @param[in] init_chain_id first chain id. The pseudo random number generator
* will advance by for each chain by an integer sequence from `init_chain_id` to
* `init_chain_id+num_chains-1`
* @param[in] init_radius radius to initialize
* @param[in] num_warmup Number of warmup samples
* @param[in] num_samples Number of samples
* @param[in] num_thin Number to thin the samples
* @param[in] save_warmup Indicates whether to save the warmup iterations
* @param[in] refresh Controls the output
* @param[in] stepsize initial stepsize for discrete evolution
* @param[in] stepsize_jitter uniform random jitter of stepsize
* @param[in] max_depth Maximum tree depth
* @param[in] delta adaptation target acceptance statistic
* @param[in] gamma adaptation regularization scale
* @param[in] kappa adaptation relaxation exponent
* @param[in] t0 adaptation iteration offset
* @param[in] init_buffer width of initial fast adaptation interval
* @param[in] term_buffer width of final fast adaptation interval
* @param[in] window initial width of slow adaptation interval
* @param[in,out] interrupt Callback for interrupts
* @param[in,out] logger Logger for messages
* @param[in,out] init_writer std vector of Writer callbacks for unconstrained
* inits of each chain.
* @param[in,out] sample_writer std vector of Writers for draws of each chain.
* @param[in,out] diagnostic_writer std vector of Writers for diagnostic
* information of each chain.
* @param[in,out] metric_writer std vector of Writers for tuning params
* @return error_codes::OK if successful
*/
template <class Model, typename InitContextPtr, typename InitWriter,
typename SampleWriter, typename DiagnosticWriter,
typename MetricWriter>
int hmc_nuts_diag_e_adapt(
Model& model, size_t num_chains, const std::vector<InitContextPtr>& init,
unsigned int random_seed, unsigned int init_chain_id, double init_radius,
int num_warmup, int num_samples, int num_thin, bool save_warmup,
int refresh, double stepsize, double stepsize_jitter, int max_depth,
double delta, double gamma, double kappa, double t0,
unsigned int init_buffer, unsigned int term_buffer, unsigned int window,
callbacks::interrupt& interrupt, callbacks::logger& logger,
std::vector<InitWriter>& init_writer,
std::vector<SampleWriter>& sample_writer,
std::vector<DiagnosticWriter>& diagnostic_writer,
std::vector<MetricWriter>& metric_writer) {
std::vector<std::unique_ptr<stan::io::array_var_context>> unit_e_metric;
unit_e_metric.reserve(num_chains);
for (size_t i = 0; i < num_chains; ++i) {
unit_e_metric.emplace_back(std::make_unique<stan::io::array_var_context>(
util::create_unit_e_diag_inv_metric(model.num_params_r())));
}
if (num_chains == 1) {
return hmc_nuts_diag_e_adapt(
model, *init[0], *unit_e_metric[0], random_seed, init_chain_id,
init_radius, num_warmup, num_samples, num_thin, save_warmup, refresh,
stepsize, stepsize_jitter, max_depth, delta, gamma, kappa, t0,
init_buffer, term_buffer, window, interrupt, logger, init_writer[0],
sample_writer[0], diagnostic_writer[0], metric_writer[0]);
}
return hmc_nuts_diag_e_adapt(
model, num_chains, init, unit_e_metric, random_seed, init_chain_id,
init_radius, num_warmup, num_samples, num_thin, save_warmup, refresh,
stepsize, stepsize_jitter, max_depth, delta, gamma, kappa, t0,
init_buffer, term_buffer, window, interrupt, logger, init_writer,
sample_writer, diagnostic_writer, metric_writer);
}
/**
* Runs multiple chains of HMC with NUTS with adaptation using diagonal
* with identity matrix as initial inv_metric.
*
* @tparam Model Model class
* @tparam InitContextPtr A pointer with underlying type derived from
* `stan::io::var_context`
* @tparam InitWriter A type derived from `stan::callbacks::writer`
* @tparam SamplerWriter A type derived from `stan::callbacks::writer`
* @tparam DiagnosticWriter A type derived from `stan::callbacks::writer`
* @param[in] model Input model (with data already instantiated)
* @param[in] num_chains The number of chains to run in parallel. `init`,
* `init_writer`, `sample_writer`, and `diagnostic_writer` must be the same
* length as this value.
* @param[in] init A std vector of init var contexts for initialization of each
* chain.
* @param[in] random_seed random seed for the random number generator
* @param[in] init_chain_id first chain id. The pseudo random number generator
* will advance by for each chain by an integer sequence from `init_chain_id` to
* `init_chain_id+num_chains-1`
* @param[in] init_radius radius to initialize
* @param[in] num_warmup Number of warmup samples
* @param[in] num_samples Number of samples
* @param[in] num_thin Number to thin the samples
* @param[in] save_warmup Indicates whether to save the warmup iterations
* @param[in] refresh Controls the output
* @param[in] stepsize initial stepsize for discrete evolution
* @param[in] stepsize_jitter uniform random jitter of stepsize
* @param[in] max_depth Maximum tree depth
* @param[in] delta adaptation target acceptance statistic
* @param[in] gamma adaptation regularization scale
* @param[in] kappa adaptation relaxation exponent
* @param[in] t0 adaptation iteration offset
* @param[in] init_buffer width of initial fast adaptation interval
* @param[in] term_buffer width of final fast adaptation interval
* @param[in] window initial width of slow adaptation interval
* @param[in,out] interrupt Callback for interrupts
* @param[in,out] logger Logger for messages
* @param[in,out] init_writer std vector of Writer callbacks for unconstrained
* inits of each chain.
* @param[in,out] sample_writer std vector of Writers for draws of each chain.
* @param[in,out] diagnostic_writer std vector of Writers for diagnostic
* information of each chain.
* @return error_codes::OK if successful
*/
template <class Model, typename InitContextPtr, typename InitWriter,
typename SampleWriter, typename DiagnosticWriter>
int hmc_nuts_diag_e_adapt(
Model& model, size_t num_chains, const std::vector<InitContextPtr>& init,
unsigned int random_seed, unsigned int init_chain_id, double init_radius,
int num_warmup, int num_samples, int num_thin, bool save_warmup,
int refresh, double stepsize, double stepsize_jitter, int max_depth,
double delta, double gamma, double kappa, double t0,
unsigned int init_buffer, unsigned int term_buffer, unsigned int window,
callbacks::interrupt& interrupt, callbacks::logger& logger,
std::vector<InitWriter>& init_writer,
std::vector<SampleWriter>& sample_writer,
std::vector<DiagnosticWriter>& diagnostic_writer) {
std::vector<std::unique_ptr<stan::io::array_var_context>> unit_e_metric;
unit_e_metric.reserve(num_chains);
for (size_t i = 0; i < num_chains; ++i) {
unit_e_metric.emplace_back(std::make_unique<stan::io::array_var_context>(
util::create_unit_e_diag_inv_metric(model.num_params_r())));
}
std::vector<stan::callbacks::structured_writer> dummy_metric_writer(
num_chains);
if (num_chains == 1) {
return hmc_nuts_diag_e_adapt(
model, *init[0], *unit_e_metric[0], random_seed, init_chain_id,
init_radius, num_warmup, num_samples, num_thin, save_warmup, refresh,
stepsize, stepsize_jitter, max_depth, delta, gamma, kappa, t0,
init_buffer, term_buffer, window, interrupt, logger, init_writer[0],
sample_writer[0], diagnostic_writer[0], dummy_metric_writer[0]);
}
return hmc_nuts_diag_e_adapt(
model, num_chains, init, unit_e_metric, random_seed, init_chain_id,
init_radius, num_warmup, num_samples, num_thin, save_warmup, refresh,
stepsize, stepsize_jitter, max_depth, delta, gamma, kappa, t0,
init_buffer, term_buffer, window, interrupt, logger, init_writer,
sample_writer, diagnostic_writer, dummy_metric_writer);
}
} // namespace sample
} // namespace services
} // namespace stan
#endif