1 #ifndef STAN__COMMON__COMMAND_HPP
2 #define STAN__COMMON__COMMAND_HPP
7 #include <boost/date_time/posix_time/posix_time_types.hpp>
8 #include <boost/math/special_functions/fpclassify.hpp>
9 #include <boost/random/additive_combine.hpp>
10 #include <boost/random/uniform_real_distribution.hpp>
66 template <
class Model>
67 int command(
int argc,
const char* argv[]) {
69 std::vector<stan::gm::argument*> valid_arguments;
77 int err_code = parser.
parse_args(argc, argv, &std::cout, &std::cout);
80 std::cout <<
"Failed to parse arguments, terminating Stan" << std::endl;
94 unsigned int random_seed = 0;
99 random_seed = (boost::posix_time::microsec_clock::universal_time() -
100 boost::posix_time::ptime(boost::posix_time::min_date_time))
101 .total_milliseconds();
106 random_seed = random_arg->
value();
109 typedef boost::ecuyer1988 rng_t;
110 rng_t base_rng(random_seed);
113 static boost::uintmax_t DISCARD_STRIDE =
static_cast<boost::uintmax_t
>(1) << 50;
114 base_rng.discard(DISCARD_STRIDE * (
id - 1));
121 std::string data_file
124 std::fstream data_stream(data_file.c_str(),
131 parser.
arg(
"output")->
arg(
"file"))->value();
132 std::fstream* output_stream = 0;
133 if (output_file !=
"") {
134 output_stream =
new std::fstream(output_file.c_str(),
140 parser.
arg(
"output")->
arg(
"diagnostic_file"))->value();
142 std::fstream* diagnostic_stream = 0;
143 if (diagnostic_file !=
"") {
144 diagnostic_stream =
new std::fstream(diagnostic_file.c_str(),
150 parser.
arg(
"output")->
arg(
"refresh"))->value();
156 Model model(data_var_context, &std::cout);
158 Eigen::VectorXd cont_params = Eigen::VectorXd::Zero(model.num_params_r());
160 parser.
print(&std::cout);
161 std::cout << std::endl;
165 write_model(output_stream, model.model_name(),
"#");
166 parser.
print(output_stream,
"#");
169 if (diagnostic_stream) {
171 write_model(diagnostic_stream, model.model_name(),
"#");
172 parser.
print(diagnostic_stream,
"#");
176 parser.
arg(
"init"))->value();
179 if (!initialize_state<dump_factory>
180 (init, cont_params, model, base_rng, &std::cout,
181 var_context_factory))
188 if (parser.
arg(
"method")->
arg(
"diagnose")) {
190 std::vector<double> cont_vector(cont_params.size());
191 for (
int i = 0; i < cont_params.size(); ++i)
192 cont_vector.at(i) = cont_params(i);
193 std::vector<int> disc_vector;
196 (parser.
arg(
"method")->
arg(
"diagnose")->
arg(
"test"));
198 if (test->
value() ==
"gradient") {
199 std::cout << std::endl <<
"TEST GRADIENT MODE" << std::endl;
202 (test->
arg(
"gradient")->
arg(
"epsilon"))->value();
205 (test->
arg(
"gradient")->
arg(
"error"))->value();
208 = stan::model::test_gradients<true,true>(model,cont_vector, disc_vector,
209 epsilon, error, std::cout);
213 = stan::model::test_gradients<true,true>(model,cont_vector, disc_vector,
214 epsilon, error, *output_stream);
217 if (diagnostic_stream) {
219 = stan::model::test_gradients<true,true>(model,cont_vector, disc_vector,
220 epsilon, error, *diagnostic_stream);
234 if (parser.
arg(
"method")->
arg(
"optimize")) {
235 std::vector<double> cont_vector(cont_params.size());
236 for (
int i = 0; i < cont_params.size(); ++i)
237 cont_vector.at(i) = cont_params(i);
238 std::vector<int> disc_vector;
241 (parser.
arg(
"method")->
arg(
"optimize")->
arg(
"algorithm"));
244 parser.
arg(
"method")->
arg(
"optimize")->
arg(
"iter"))->value();
249 ->
arg(
"save_iterations"))->value();
251 std::vector<std::string> names;
252 names.push_back(
"lp__");
253 model.constrained_param_names(names,
true,
true);
255 (*output_stream) << names.at(0);
256 for (
size_t i = 1; i < names.size(); ++i) {
257 (*output_stream) <<
"," << names.at(i);
259 (*output_stream) << std::endl;
264 if (algo->
value() ==
"newton") {
267 lp = model.template log_prob<false, false>(cont_vector, disc_vector, &std::cout);
268 }
catch (
const std::exception&
e) {
270 lp = -std::numeric_limits<double>::infinity();
273 std::cout <<
"initial log joint probability = " << lp << std::endl;
274 if (output_stream && save_iterations) {
276 lp, cont_vector, disc_vector);
279 double lastlp = lp * 1.1;
281 std::cout <<
"(lp - lastlp) / lp > 1e-8: " << ((lp - lastlp) /
fabs(lp)) << std::endl;
282 while ((lp - lastlp) /
fabs(lp) > 1e-8) {
286 std::cout <<
"Iteration ";
287 std::cout << std::setw(2) << (m + 1) <<
". ";
288 std::cout <<
"Log joint probability = " << std::setw(10) << lp;
289 std::cout <<
". Improved by " << (lp - lastlp) <<
".";
290 std::cout << std::endl;
294 if (output_stream && save_iterations) {
296 lp, cont_vector, disc_vector);
301 }
else if (algo->
value() ==
"bfgs") {
304 Optimizer bfgs(model, cont_vector, disc_vector, &std::cout);
307 algo->
arg(
"bfgs")->
arg(
"init_alpha"))->value();
309 algo->
arg(
"bfgs")->
arg(
"tol_obj"))->value();
311 algo->
arg(
"bfgs")->
arg(
"tol_rel_obj"))->value();
313 algo->
arg(
"bfgs")->
arg(
"tol_grad"))->value();
315 algo->
arg(
"bfgs")->
arg(
"tol_rel_grad"))->value();
317 algo->
arg(
"bfgs")->
arg(
"tol_param"))->value();
318 bfgs._conv_opts.maxIts = num_iterations;
321 lp, cont_vector, disc_vector,
322 output_stream, &std::cout,
323 save_iterations, refresh,
325 }
else if (algo->
value() ==
"lbfgs") {
328 Optimizer bfgs(model, cont_vector, disc_vector, &std::cout);
330 bfgs.get_qnupdate().set_history_size(dynamic_cast<gm::int_argument*>(
331 algo->
arg(
"lbfgs")->
arg(
"history_size"))->value());
333 algo->
arg(
"lbfgs")->
arg(
"init_alpha"))->value();
335 algo->
arg(
"lbfgs")->
arg(
"tol_obj"))->value();
337 algo->
arg(
"lbfgs")->
arg(
"tol_rel_obj"))->value();
339 algo->
arg(
"lbfgs")->
arg(
"tol_grad"))->value();
341 algo->
arg(
"lbfgs")->
arg(
"tol_rel_grad"))->value();
343 algo->
arg(
"lbfgs")->
arg(
"tol_param"))->value();
344 bfgs._conv_opts.maxIts = num_iterations;
347 lp, cont_vector, disc_vector,
348 output_stream, &std::cout,
349 save_iterations, refresh,
357 lp, cont_vector, disc_vector);
358 output_stream->close();
359 delete output_stream;
368 if (parser.
arg(
"method")->
arg(
"sample")) {
371 clock_t start_check = clock();
373 double init_log_prob;
374 Eigen::VectorXd init_grad = Eigen::VectorXd::Zero(model.num_params_r());
378 clock_t end_check = clock();
379 double deltaT = (double)(end_check - start_check) / CLOCKS_PER_SEC;
381 std::cout << std::endl;
382 std::cout <<
"Gradient evaluation took " << deltaT <<
" seconds" << std::endl;
383 std::cout <<
"1000 transitions using 10 leapfrog steps per transition would take "
384 << 1e4 * deltaT <<
" seconds." << std::endl;
385 std::cout <<
"Adjust your expectations accordingly!" << std::endl << std::endl;
386 std::cout << std::endl;
395 writer(sample_recorder, diagnostic_recorder, message_recorder, &std::cout);
399 parser.
arg(
"method")->
arg(
"sample")->
arg(
"num_warmup"))->value();
402 parser.
arg(
"method")->
arg(
"sample")->
arg(
"num_samples"))->value();
405 parser.
arg(
"method")->
arg(
"sample")->
arg(
"thin"))->value();
408 parser.
arg(
"method")->
arg(
"sample")->
arg(
"save_warmup"))->value();
419 (parser.
arg(
"method")->
arg(
"sample")->
arg(
"algorithm"));
422 parser.
arg(
"method")->
arg(
"sample")->
arg(
"adapt"));
425 if (model.num_params_r() == 0 && algo->
value() !=
"fixed_param") {
427 <<
"Must use algorithm=fixed_param for model that has no parameters."
432 if (algo->
value() ==
"fixed_param") {
436 adapt_engaged =
false;
438 if (num_warmup != 0) {
439 std::cout <<
"Warning: warmup will be skipped for the fixed parameter sampler!" << std::endl;
443 }
else if (algo->
value() ==
"rwm") {
448 }
else if (algo->
value() ==
"hmc") {
450 int engine_index = 0;
453 if (engine->
value() ==
"static") {
455 }
else if (engine->
value() ==
"nuts") {
459 int metric_index = 0;
462 if (metric->
value() ==
"unit_e") {
464 }
else if (metric->
value() ==
"diag_e") {
466 }
else if (metric->
value() ==
"dense_e") {
470 int sampler_select = engine_index
472 + 100 *
static_cast<int>(adapt_engaged);
474 switch (sampler_select) {
478 sampler_ptr =
new sampler(model, base_rng, &std::cout, &std::cout);
479 if (!init_static_hmc<sampler>(sampler_ptr, algo))
return 0;
485 sampler_ptr =
new sampler(model, base_rng, &std::cout, &std::cout);
486 if (!init_nuts<sampler>(sampler_ptr, algo))
return 0;
492 sampler_ptr =
new sampler(model, base_rng, &std::cout, &std::cout);
493 if (!init_static_hmc<sampler>(sampler_ptr, algo))
return 0;
499 sampler_ptr =
new sampler(model, base_rng, &std::cout, &std::cout);
500 if (!init_nuts<sampler>(sampler_ptr, algo))
return 0;
506 sampler_ptr =
new sampler(model, base_rng, &std::cout, &std::cout);
507 if (!init_static_hmc<sampler>(sampler_ptr, algo))
return 0;
513 sampler_ptr =
new sampler(model, base_rng, &std::cout, &std::cout);
514 if (!init_nuts<sampler>(sampler_ptr, algo))
return 0;
520 sampler_ptr =
new sampler(model, base_rng, &std::cout, &std::cout);
521 if (!init_static_hmc<sampler>(sampler_ptr, algo))
return 0;
522 if (!init_adapt<sampler>(sampler_ptr, adapt, cont_params))
return 0;
528 sampler_ptr =
new sampler(model, base_rng, &std::cout, &std::cout);
529 if (!init_nuts<sampler>(sampler_ptr, algo))
return 0;
530 if (!init_adapt<sampler>(sampler_ptr, adapt, cont_params))
return 0;
536 sampler_ptr =
new sampler(model, base_rng, &std::cout, &std::cout);
537 if (!init_static_hmc<sampler>(sampler_ptr, algo))
return 0;
538 if (!init_windowed_adapt<sampler>(sampler_ptr, adapt, num_warmup, cont_params))
545 sampler_ptr =
new sampler(model, base_rng, &std::cout, &std::cout);
546 if (!init_nuts<sampler>(sampler_ptr, algo))
return 0;
547 if (!init_windowed_adapt<sampler>(sampler_ptr, adapt, num_warmup, cont_params))
554 sampler_ptr =
new sampler(model, base_rng, &std::cout, &std::cout);
555 if (!init_static_hmc<sampler>(sampler_ptr, algo))
return 0;
556 if (!init_windowed_adapt<sampler>(sampler_ptr, adapt, num_warmup, cont_params))
563 sampler_ptr =
new sampler(model, base_rng, &std::cout, &std::cout);
564 if (!init_nuts<sampler>(sampler_ptr, algo))
return 0;
565 if (!init_windowed_adapt<sampler>(sampler_ptr, adapt, num_warmup, cont_params))
571 std::cout <<
"No sampler matching HMC specification!" << std::endl;
578 writer.write_sample_names(s, sampler_ptr, model);
579 writer.write_diagnostic_names(s, sampler_ptr, model);
581 std::string prefix =
"";
582 std::string suffix =
"\n";
586 clock_t start = clock();
588 warmup<Model, rng_t>(sampler_ptr, num_warmup, num_samples, num_thin,
589 refresh, save_warmup,
592 prefix, suffix, std::cout,
593 startTransitionCallback);
595 clock_t end = clock();
596 warmDeltaT = (double)(end - start) / CLOCKS_PER_SEC;
600 writer.write_adapt_finish(sampler_ptr);
606 stan::common::sample<Model, rng_t>(sampler_ptr, num_warmup, num_samples, num_thin,
610 prefix, suffix, std::cout,
611 startTransitionCallback);
614 sampleDeltaT = (double)(end - start) / CLOCKS_PER_SEC;
616 writer.write_timing(warmDeltaT, sampleDeltaT);
618 if (sampler_ptr)
delete sampler_ptr;
623 output_stream->close();
624 delete output_stream;
627 if (diagnostic_stream) {
628 diagnostic_stream->close();
629 delete diagnostic_stream;
632 for (
size_t i = 0; i < valid_arguments.size(); ++i)
633 delete valid_arguments.at(i);
bool set_value(const T &value)
int parse_args(int argc, const char *argv[], std::ostream *out=0, std::ostream *err=0)
void gradient(const F &f, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, double &fx, Eigen::Matrix< double, Eigen::Dynamic, 1 > &grad_fx)
Calculate the value and the gradient of the specified function at the specified argument.
fvar< T > fabs(const fvar< T > &x)
void write_stan(std::ostream *s, const std::string prefix="")
argument * arg(std::string name)
std::string description() const
double newton_step(M &model, std::vector< double > ¶ms_r, std::vector< int > ¶ms_i, std::ostream *output_stream=0)
argument * arg(const std::string name)
Writes out a vector as string.
void gradient(const M &model, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, double &f, Eigen::Matrix< double, Eigen::Dynamic, 1 > &grad_f, std::ostream *msgs=0)
argument * arg(std::string name)
void write_error_msg(std::ostream *error_stream, const std::exception &e)
void write_iteration(std::ostream &output_stream, Model &model, RNG &base_rng, double lp, std::vector< double > &cont_vector, std::vector< int > &disc_vector)
int do_bfgs_optimize(ModelT &model, BFGSOptimizerT &bfgs, RNGT &base_rng, double &lp, std::vector< double > &cont_vector, std::vector< int > &disc_vector, std::fstream *output_stream, std::ostream *notice_stream, bool save_iterations, int refresh, StartIterationCallback &callback)
void write_model(std::ostream *s, const std::string model_name, const std::string prefix="")
double e()
Return the base of the natural logarithm.
mcmc_writer writes out headers and samples
int command(int argc, const char *argv[])
virtual argument * arg(const std::string name)
Represents named arrays with dimensions.
void print(std::ostream *s, const std::string prefix="")