model-method-sample | R Documentation |
The $sample()
method of a CmdStanModel
object runs Stan's
main Markov chain Monte Carlo algorithm.
Any argument left as NULL
will default to the default value used by the
installed version of CmdStan. See the
CmdStan User’s Guide
for more details.
After model fitting any diagnostics specified via the diagnostics
argument will be checked and warnings will be printed if warranted.
sample(
data = NULL,
seed = NULL,
refresh = NULL,
init = NULL,
save_latent_dynamics = FALSE,
output_dir = getOption("cmdstanr_output_dir"),
output_basename = NULL,
sig_figs = NULL,
chains = 4,
parallel_chains = getOption("mc.cores", 1),
chain_ids = seq_len(chains),
threads_per_chain = NULL,
opencl_ids = NULL,
iter_warmup = NULL,
iter_sampling = NULL,
save_warmup = FALSE,
thin = NULL,
max_treedepth = NULL,
adapt_engaged = TRUE,
adapt_delta = NULL,
step_size = NULL,
metric = NULL,
metric_file = NULL,
inv_metric = NULL,
init_buffer = NULL,
term_buffer = NULL,
window = NULL,
fixed_param = FALSE,
show_messages = TRUE,
show_exceptions = TRUE,
diagnostics = c("divergences", "treedepth", "ebfmi"),
save_metric = NULL,
save_cmdstan_config = NULL,
cores = NULL,
num_cores = NULL,
num_chains = NULL,
num_warmup = NULL,
num_samples = NULL,
validate_csv = NULL,
save_extra_diagnostics = NULL,
max_depth = NULL,
stepsize = NULL
)
data |
(multiple options) The data to use for the variables specified in the data block of the Stan program. One of the following:
|
seed |
(positive integer(s)) A seed for the (P)RNG to pass to CmdStan.
In the case of multi-chain sampling the single |
refresh |
(non-negative integer) The number of iterations between
printed screen updates. If |
init |
(multiple options) The initialization method to use for the variables declared in the parameters block of the Stan program. One of the following:
|
save_latent_dynamics |
(logical) Should auxiliary diagnostic information
about the latent dynamics be written to temporary diagnostic CSV files?
This argument replaces CmdStan's |
output_dir |
(string) A path to a directory where CmdStan should write
its output CSV files. For MCMC there will be one file per chain; for other
methods there will be a single file. For interactive use this can typically
be left at
|
output_basename |
(string) A string to use as a prefix for the names of
the output CSV files of CmdStan. If |
sig_figs |
(positive integer) The number of significant figures used
when storing the output values. By default, CmdStan represent the output
values with 6 significant figures. The upper limit for |
chains |
(positive integer) The number of Markov chains to run. The default is 4. |
parallel_chains |
(positive integer) The maximum number of MCMC chains
to run in parallel. If |
chain_ids |
(integer vector) A vector of chain IDs. Must contain as many
unique positive integers as the number of chains. If not set, the default
chain IDs are used (integers starting from |
threads_per_chain |
(positive integer) If the model was
compiled with threading support, the number of
threads to use in parallelized sections within an MCMC chain (e.g., when
using the Stan functions |
opencl_ids |
(integer vector of length 2) The platform and device IDs of
the OpenCL device to use for fitting. The model must be compiled with
|
iter_warmup |
(positive integer) The number of warmup iterations to run
per chain. Note: in the CmdStan User's Guide this is referred to as
|
iter_sampling |
(positive integer) The number of post-warmup iterations
to run per chain. Note: in the CmdStan User's Guide this is referred to as
|
save_warmup |
(logical) Should warmup iterations be saved? The default
is |
thin |
(positive integer) The period between saved samples. This should typically be left at its default (no thinning) unless memory is a problem. |
max_treedepth |
(positive integer) The maximum allowed tree depth for the NUTS engine. See the Tree Depth section of the CmdStan User's Guide for more details. |
adapt_engaged |
(logical) Do warmup adaptation? The default is |
adapt_delta |
(real in |
step_size |
(positive real) The initial step size for the discrete approximation to continuous Hamiltonian dynamics. This is further tuned during warmup. |
metric |
(string) One of |
metric_file |
(character vector) The paths to JSON or Rdump files (one
per chain) compatible with CmdStan that contain precomputed inverse
metrics. The |
inv_metric |
(vector, matrix) A vector (if |
init_buffer |
(nonnegative integer) Width of initial fast timestep adaptation interval during warmup. |
term_buffer |
(nonnegative integer) Width of final fast timestep adaptation interval during warmup. |
window |
(nonnegative integer) Initial width of slow timestep/metric adaptation interval. |
fixed_param |
(logical) When |
show_messages |
(logical) When |
show_exceptions |
(logical) When |
diagnostics |
(character vector) The diagnostics to automatically check
and warn about after sampling. Setting this to an empty string These diagnostics are also available after fitting. The
Diagnostics like R-hat and effective sample size are not currently
available via the |
save_metric |
(logical) When |
save_cmdstan_config |
(logical) When |
cores , num_cores , num_chains , num_warmup , num_samples , save_extra_diagnostics , max_depth , stepsize , validate_csv |
Deprecated and will be removed in a future release. |
A CmdStanMCMC
object.
The CmdStanR website (mc-stan.org/cmdstanr) for online documentation and tutorials.
The Stan and CmdStan documentation:
Stan documentation: mc-stan.org/users/documentation
CmdStan User’s Guide: mc-stan.org/docs/cmdstan-guide
Other CmdStanModel methods:
model-method-check_syntax
,
model-method-compile
,
model-method-diagnose
,
model-method-expose_functions
,
model-method-format
,
model-method-generate-quantities
,
model-method-laplace
,
model-method-optimize
,
model-method-pathfinder
,
model-method-sample_mpi
,
model-method-variables
,
model-method-variational
## Not run:
library(cmdstanr)
library(posterior)
library(bayesplot)
color_scheme_set("brightblue")
# Set path to CmdStan
# (Note: if you installed CmdStan via install_cmdstan() with default settings
# then setting the path is unnecessary but the default below should still work.
# Otherwise use the `path` argument to specify the location of your
# CmdStan installation.)
set_cmdstan_path(path = NULL)
# Create a CmdStanModel object from a Stan program,
# here using the example model that comes with CmdStan
file <- file.path(cmdstan_path(), "examples/bernoulli/bernoulli.stan")
mod <- cmdstan_model(file)
mod$print()
# Print with line numbers. This can be set globally using the
# `cmdstanr_print_line_numbers` option.
mod$print(line_numbers = TRUE)
# Data as a named list (like RStan)
stan_data <- list(N = 10, y = c(0,1,0,0,0,0,0,0,0,1))
# Run MCMC using the 'sample' method
fit_mcmc <- mod$sample(
data = stan_data,
seed = 123,
chains = 2,
parallel_chains = 2
)
# Use 'posterior' package for summaries
fit_mcmc$summary()
# Check sampling diagnostics
fit_mcmc$diagnostic_summary()
# Get posterior draws
draws <- fit_mcmc$draws()
print(draws)
# Convert to data frame using posterior::as_draws_df
as_draws_df(draws)
# Plot posterior using bayesplot (ggplot2)
mcmc_hist(fit_mcmc$draws("theta"))
# Run 'optimize' method to get a point estimate (default is Stan's LBFGS algorithm)
# and also demonstrate specifying data as a path to a file instead of a list
my_data_file <- file.path(cmdstan_path(), "examples/bernoulli/bernoulli.data.json")
fit_optim <- mod$optimize(data = my_data_file, seed = 123)
fit_optim$summary()
# Run 'optimize' again with 'jacobian=TRUE' and then draw from Laplace approximation
# to the posterior
fit_optim <- mod$optimize(data = my_data_file, jacobian = TRUE)
fit_laplace <- mod$laplace(data = my_data_file, mode = fit_optim, draws = 2000)
fit_laplace$summary()
# Run 'variational' method to use ADVI to approximate posterior
fit_vb <- mod$variational(data = stan_data, seed = 123)
fit_vb$summary()
mcmc_hist(fit_vb$draws("theta"))
# Run 'pathfinder' method, a new alternative to the variational method
fit_pf <- mod$pathfinder(data = stan_data, seed = 123)
fit_pf$summary()
mcmc_hist(fit_pf$draws("theta"))
# Run 'pathfinder' again with more paths, fewer draws per path,
# better covariance approximation, and fewer LBFGSs iterations
fit_pf <- mod$pathfinder(data = stan_data, num_paths=10, single_path_draws=40,
history_size=50, max_lbfgs_iters=100)
# Specifying initial values as a function
fit_mcmc_w_init_fun <- mod$sample(
data = stan_data,
seed = 123,
chains = 2,
refresh = 0,
init = function() list(theta = runif(1))
)
fit_mcmc_w_init_fun_2 <- mod$sample(
data = stan_data,
seed = 123,
chains = 2,
refresh = 0,
init = function(chain_id) {
# silly but demonstrates optional use of chain_id
list(theta = 1 / (chain_id + 1))
}
)
fit_mcmc_w_init_fun_2$init()
# Specifying initial values as a list of lists
fit_mcmc_w_init_list <- mod$sample(
data = stan_data,
seed = 123,
chains = 2,
refresh = 0,
init = list(
list(theta = 0.75), # chain 1
list(theta = 0.25) # chain 2
)
)
fit_optim_w_init_list <- mod$optimize(
data = stan_data,
seed = 123,
init = list(
list(theta = 0.75)
)
)
fit_optim_w_init_list$init()
## End(Not run)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.