Skip to content
Permalink
master
Switch branches/tags
Go to file
 
 
Cannot retrieve contributors at this time
% Generated by roxygen2: do not edit by hand
% Please edit documentation in R/model.R
\name{model-method-sample}
\alias{model-method-sample}
\alias{sample}
\title{Run Stan's MCMC algorithms}
\usage{
sample(
data = NULL,
seed = NULL,
refresh = NULL,
init = NULL,
save_latent_dynamics = FALSE,
output_dir = NULL,
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,
diagnostics = c("divergences", "treedepth", "ebfmi"),
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
)
}
\arguments{
\item{data}{(multiple options) The data to use for the variables specified in
the data block of the Stan program. One of the following:
\itemize{
\item A named list of \R objects with the names corresponding to variables
declared in the data block of the Stan program. Internally this list is then
written to JSON for CmdStan using \code{\link[=write_stan_json]{write_stan_json()}}. See
\code{\link[=write_stan_json]{write_stan_json()}} for details on the conversions performed on \R objects
before they are passed to Stan.
\item A path to a data file compatible with CmdStan (JSON or \R dump). See the
appendices in the CmdStan guide for details on using these formats.
\item \code{NULL} or an empty list if the Stan program has no data block.
}}
\item{seed}{(positive integer(s)) A seed for the (P)RNG to pass to CmdStan.
In the case of multi-chain sampling the single \code{seed} will automatically be
augmented by the the run (chain) ID so that each chain uses a different
seed. The exception is the transformed data block, which defaults to
using same seed for all chains so that the same data is generated for all
chains if RNG functions are used. The only time \code{seed} should be specified
as a vector (one element per chain) is if RNG functions are used in
transformed data and the goal is to generate \emph{different} data for each
chain.}
\item{refresh}{(non-negative integer) The number of iterations between
printed screen updates. If \code{refresh = 0}, only error messages will be
printed.}
\item{init}{(multiple options) The initialization method to use for the
variables declared in the parameters block of the Stan program. One of
the following:
\itemize{
\item A real number \code{x>0}. This initializes \emph{all} parameters randomly between
\verb{[-x,x]} on the \emph{unconstrained} parameter space.;
\item The number \code{0}. This initializes \emph{all} parameters to \code{0};
\item A character vector of paths (one per chain) to JSON or Rdump files
containing initial values for all or some parameters. See
\code{\link[=write_stan_json]{write_stan_json()}} to write \R objects to JSON files compatible with
CmdStan.
\item A list of lists containing initial values for all or some parameters. For
MCMC the list should contain a sublist for each chain. For optimization and
variational inference there should be just one sublist. The sublists should
have named elements corresponding to the parameters for which you are
specifying initial values. See \strong{Examples}.
\item A function that returns a single list with names corresponding to the
parameters for which you are specifying initial values. The function can
take no arguments or a single argument \code{chain_id}. For MCMC, if the function
has argument \code{chain_id} it will be supplied with the chain id (from 1 to
number of chains) when called to generate the initial values. See
\strong{Examples}.
}}
\item{save_latent_dynamics}{(logical) Should auxiliary diagnostic information
about the latent dynamics be written to temporary diagnostic CSV files?
This argument replaces CmdStan's \code{diagnostic_file} argument and the content
written to CSV is controlled by the user's CmdStan installation and not
CmdStanR (for some algorithms no content may be written). The default
is \code{FALSE}, which is appropriate for almost every use case. To save the
temporary files created when \code{save_latent_dynamics=TRUE} see the
\code{\link[=fit-method-save_latent_dynamics_files]{$save_latent_dynamics_files()}}
method.}
\item{output_dir}{(string) A path to a directory where CmdStan should write
its output CSV files. For interactive use this can typically be left at
\code{NULL} (temporary directory) since CmdStanR makes the CmdStan output
(posterior draws and diagnostics) available in \R via methods of the fitted
model objects. The behavior of \code{output_dir} is as follows:
\itemize{
\item If \code{NULL} (the default), then the CSV files are written to a temporary
directory and only saved permanently if the user calls one of the \verb{$save_*}
methods of the fitted model object (e.g.,
\code{\link[=fit-method-save_output_files]{$save_output_files()}}). These temporary
files are removed when the fitted model object is
\link[base:gc]{garbage collected} (manually or automatically).
\item If a path, then the files are created in \code{output_dir} with names
corresponding to the defaults used by \verb{$save_output_files()}.
}}
\item{output_basename}{(string) A string to use as a prefix for the names of
the output CSV files of CmdStan. If \code{NULL} (the default), the basename of
the output CSV files will be comprised from the model name, timestamp, and
5 random characters.}
\item{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 \code{sig_figs} is 18.
Increasing this value will result in larger output CSV files and thus an
increased usage of disk space.}
\item{chains}{(positive integer) The number of Markov chains to run. The
default is 4.}
\item{parallel_chains}{(positive integer) The \emph{maximum} number of MCMC chains
to run in parallel. If \code{parallel_chains} is not specified then the default
is to look for the option \code{"mc.cores"}, which can be set for an entire \R
session by \code{options(mc.cores=value)}. If the \code{"mc.cores"} option has not
been set then the default is \code{1}.}
\item{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 \code{1}).}
\item{threads_per_chain}{(positive integer) If the model was
\link[=model-method-compile]{compiled} with threading support, the number of
threads to use in parallelized sections \emph{within} an MCMC chain (e.g., when
using the Stan functions \code{reduce_sum()} or \code{map_rect()}). This is in
contrast with \code{parallel_chains}, which specifies the number of chains to
run in parallel. The actual number of CPU cores used is
\code{parallel_chains*threads_per_chain}. For an example of using threading see
the Stan case study
\href{https://mc-stan.org/users/documentation/case-studies/reduce_sum_tutorial.html}{Reduce Sum: A Minimal Example}.}
\item{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 \code{cpp_options = list(stan_opencl = TRUE)} for this
argument to have an effect.}
\item{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
\code{num_warmup}.}
\item{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
\code{num_samples}.}
\item{save_warmup}{(logical) Should warmup iterations be saved? The default
is \code{FALSE}.}
\item{thin}{(positive integer) The period between saved samples. This should
typically be left at its default (no thinning) unless memory is a problem.}
\item{max_treedepth}{(positive integer) The maximum allowed tree depth for
the NUTS engine. See the \emph{Tree Depth} section of the CmdStan User's Guide
for more details.}
\item{adapt_engaged}{(logical) Do warmup adaptation? The default is \code{TRUE}.
If a precomputed inverse metric is specified via the \code{inv_metric} argument
(or \code{metric_file}) then, if \code{adapt_engaged=TRUE}, Stan will use the
provided inverse metric just as an initial guess during adaptation. To turn
off adaptation when using a precomputed inverse metric set
\code{adapt_engaged=FALSE}.}
\item{adapt_delta}{(real in \verb{(0,1)}) The adaptation target acceptance
statistic.}
\item{step_size}{(positive real) The \emph{initial} step size for the discrete
approximation to continuous Hamiltonian dynamics. This is further tuned
during warmup.}
\item{metric}{(string) One of \code{"diag_e"}, \code{"dense_e"}, or \code{"unit_e"},
specifying the geometry of the base manifold. See the \emph{Euclidean Metric}
section of the CmdStan User's Guide for more details. To specify a
precomputed (inverse) metric, see the \code{inv_metric} argument below.}
\item{metric_file}{(character vector) The paths to JSON or
Rdump files (one per chain) compatible with CmdStan that contain
precomputed inverse metrics. The \code{metric_file} argument is inherited from
CmdStan but is confusing in that the entry in JSON or Rdump file(s) must be
named \code{inv_metric}, referring to the \emph{inverse} metric. We recommend instead
using CmdStanR's \code{inv_metric} argument (see below) to specify an inverse
metric directly using a vector or matrix from your \R session.}
\item{inv_metric}{(vector, matrix) A vector (if \code{metric='diag_e'}) or a
matrix (if \code{metric='dense_e'}) for initializing the inverse metric. This
can be used as an alternative to the \code{metric_file} argument. A vector is
interpreted as a diagonal metric. The inverse metric is usually set to an
estimate of the posterior covariance. See the \code{adapt_engaged} argument
above for details about (and control over) how specifying a precomputed
inverse metric interacts with adaptation.}
\item{init_buffer}{(nonnegative integer) Width of initial fast timestep
adaptation interval during warmup.}
\item{term_buffer}{(nonnegative integer) Width of final fast timestep
adaptation interval during warmup.}
\item{window}{(nonnegative integer) Initial width of slow timestep/metric
adaptation interval.}
\item{fixed_param}{(logical) When \code{TRUE}, call CmdStan with argument
\code{"algorithm=fixed_param"}. The default is \code{FALSE}. The fixed parameter
sampler generates a new sample without changing the current state of the
Markov chain; only generated quantities may change. This can be useful
when, for example, trying to generate pseudo-data using the generated
quantities block. If the parameters block is empty then using
\code{fixed_param=TRUE} is mandatory. When \code{fixed_param=TRUE} the \code{chains} and
\code{parallel_chains} arguments will be set to \code{1}.}
\item{show_messages}{(logical) When \code{TRUE} (the default), prints all
informational messages, for example rejection of the current proposal.
Disable if you wish to silence these messages, but this is not usually
recommended unless you are very confident that the model is correct up to
numerical error. If the messages are silenced then the
\code{\link[=fit-method-output]{$output()}} method of the resulting fit object can be
used to display the silenced messages.}
\item{diagnostics}{(character vector) The diagnostics to automatically check
and warn about after sampling. Setting this to an empty string \code{""} or
\code{NULL} can be used to prevent CmdStanR from automatically reading in the
sampler diagnostics from CSV if you wish to manually read in the results
and validate them yourself, for example using \code{\link[=read_cmdstan_csv]{read_cmdstan_csv()}}. The
currently available diagnostics are \code{"divergences"}, \code{"treedepth"},
and \code{"ebfmi"} (the default is to check all of them).
These diagnostics are also available after fitting. The
\code{\link[=fit-method-sampler_diagnostics]{$sampler_diagnostics()}} method provides
access the diagnostic values for each iteration and the
\code{\link[=fit-method-diagnostic_summary]{$diagnostic_summary()}} method provides
summaries of the diagnostics and can regenerate the warning messages.
Diagnostics like R-hat and effective sample size are \emph{not} currently
available via the \code{diagnostics} argument but can be checked after fitting
using the \code{\link[=fit-method-summary]{$summary()}} method.}
\item{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.}
}
\value{
A \code{\link{CmdStanMCMC}} object.
}
\description{
The \verb{$sample()} method of a \code{\link{CmdStanModel}} object runs Stan's
main Markov chain Monte Carlo algorithm.
Any argument left as \code{NULL} will default to the default value used by the
installed version of CmdStan. See the
\href{https://mc-stan.org/docs/cmdstan-guide/}{CmdStan User’s Guide}
for more details.
After model fitting any diagnostics specified via the \code{diagnostics}
argument will be checked and warnings will be printed if warranted.
}
\examples{
\dontrun{
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()
# 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()
# 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"))
# Call CmdStan's diagnose and stansummary utilities
fit_mcmc$cmdstan_diagnose()
fit_mcmc$cmdstan_summary()
# For models fit using MCMC, if you like working with RStan's stanfit objects
# then you can create one with rstan::read_stan_csv()
# stanfit <- rstan::read_stan_csv(fit_mcmc$output_files())
# 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 'variational' method to approximate the posterior (default is meanfield ADVI)
fit_vb <- mod$variational(data = stan_data, seed = 123)
fit_vb$summary()
# Plot approximate posterior using bayesplot
mcmc_hist(fit_vb$draws("theta"))
# 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()
}
}
\seealso{
The CmdStanR website
(\href{https://mc-stan.org/cmdstanr/}{mc-stan.org/cmdstanr}) for online
documentation and tutorials.
The Stan and CmdStan documentation:
\itemize{
\item Stan documentation: \href{https://mc-stan.org/users/documentation/}{mc-stan.org/users/documentation}
\item CmdStan Users Guide: \href{https://mc-stan.org/docs/cmdstan-guide/}{mc-stan.org/docs/cmdstan-guide}
}
Other CmdStanModel methods:
\code{\link{model-method-check_syntax}},
\code{\link{model-method-compile}},
\code{\link{model-method-diagnose}},
\code{\link{model-method-format}},
\code{\link{model-method-generate-quantities}},
\code{\link{model-method-optimize}},
\code{\link{model-method-sample_mpi}},
\code{\link{model-method-variables}},
\code{\link{model-method-variational}}
}
\concept{CmdStanModel methods}