Programs

CMDSTAN_HOME

CmdStan.CMDSTAN_HOMEConstant

The directory which contains the cmdstan executables such as bin/stanc and bin/stansummary. Inferred from the environment variable JULIA_CMDSTAN_HOME or ENV["JULIA_CMDSTAN_HOME"] when available.

If these are not available, use set_cmdstan_home! to set the value of CMDSTAN_HOME.

Example: set_cmdstan_home!(homedir() * "/Projects/Stan/cmdstan/")

Executing versioninfo() will display the value of JULIA_CMDSTAN_HOME if defined.

source

set_cmdstan_home!

CmdStan.set_cmdstan_home!Function

Set the path for the CMDSTAN_HOME environment variable.

Example: set_cmdstan_home!(homedir() * "/Projects/Stan/cmdstan/")

source

stanmodel()

CmdStan.StanmodelType

Method Stanmodel

Create a Stanmodel.

Constructors

Stanmodel(
  method=Sample();
  name="noname", 
  nchains=4,
  num_warmup=1000, 
  num_samples=1000,
  thin=1,
  model="",
  monitors=String[],
  random=Random(),
  output=Output(),
  printsummary=false,
  pdir::String=pwd(),
  tmpdir::String=joinpath(pwd(), "tmp"),
  output_format=:array
)

Required arguments

* `method::Method`             : See ?Method

Optional arguments

* `name::String`               : Name for the model
* `nchains::Int`               : Number of chains
* `num_warmup::Int`            : Number of samples used for num_warmup
* `num_samples::Int`           : Sample iterations
* `thin::Int`                  : Stan thinning factor
* `model::String`              : Stan program source
* `monitors::String[] `        : Variables saved for post-processing
* `random::Random`             : Random seed settings
* `output::Output`             : File output options
* `printsummary=true`          : Show computed stan summary
* `pdir::String`               : Working directory
* `tmpdir::String`             : Directory where output files are stored
* `output_format::Symbol `     : Output format

CmdStan.jl supports 3 output_format values:

1. :array           # Returns an array of draws (default)
2. :mcmcchains      # Return an MCMCChains.Chains object
3. :dataframes      # Return an DataFrames.DataFrame object
4. :namedtuple      # Returns a NamedTuple object

The first option (the default) returns an Array{Float64, 3} with ndraws,
nvars, nchains as indices.

The 2nd option returns an MCMCChains.Chains object, the 3rd a DataFrame
object and the final option returns a NamedTuple.

Example

stanmodel = Stanmodel(num_samples=1200, thin=2, name="bernoulli", 
  model=bernoullimodel);

Related help

?stan                                  : Run a Stanmodel
?CmdStan.Sample                        : Sampling settings
?CmdStan.Method                        : List of available methods
?CmdStan.Output                        : Output file settings
?CmdStan.DataDict                      : Input data
?CmdStan.convert_a3d                   : Options for output formats
source
CmdStan.update_model_fileFunction

Method update_model_file

Update Stan language model file if necessary

Method

CmdStan.update_model_file(
  file::String, 
  str::String
)

Required arguments

* `file::AbstractString`                : File holding existing Stan model
* `str::AbstractString`                 : Stan model string

Related help

?CmdStan.Stanmodel                 : Create a StanModel
source

stan()

CmdStan.stanFunction

stan

Execute a Stan model.

Method

rc, sim, cnames = stan(
  model::Stanmodel, 
  data=Nothing, 
  ProjDir=pwd();
  init=Nothing,
  summary=true, 
  diagnostics=false, 
  CmdStanDir=CMDSTAN_HOME,
  file_run_log=true
)

Required arguments

* `model::Stanmodel`              : See ?Method

Optional positional arguments

* `data=Nothing`                  : Observed input data dictionary 
* `ProjDir=pwd()`                 : Project working directory

Keyword arguments

* `init=Nothing`                  : Initial parameter value dictionary
* `summary=true`                  : Use CmdStan's stansummary to display results
* `diagnostics=false`             : Generate diagnostics file
* `CmdStanDir=CMDSTAN_HOME`       : Location of cmdstan directory
* `file_run_log=true`             : Create run log file (if false, write log to stdout)
* `file_make_log=true`            : Create make log file (if false, write log to stdout)

If the type the data or init arguments are an AbstartcString this is interpreted as a path name to an existing .R file. This file is copied to the coresponding .R data and/or init files for for each chain.

Return values

* `rc::Int`                       : Return code from stan(), rc == 0 if all is well
* `samples`                       : Samples
* `cnames`                        : Vector of variable names

Examples

# no data, use default ProjDir (pwd)
stan(mymodel)
# default ProjDir (pwd)
stan(mymodel, mydata)
# specify ProjDir
stan(mymodel, mydata, "~/myproject/")
# keyword arguments
stan(mymodel, mydata, "~/myproject/", diagnostics=true, summary=false)
# use default ProjDir (pwd), with keyword arguments
stan(mymodel, mydata, diagnostics=true, summary=false)

Related help

?Stanmodel                         : Create a StanModel
source

Types

CmdStan.MethodType

Available top level Method

Method

*  Sample::Method             : Sampling
*  Optimize::Method           : Optimization
*  Diagnose::Method           : Diagnostics
*  Variational::Method        : Variational Bayes
source
CmdStan.SampleType

Sample type and constructor

Settings for method=Sample() in Stanmodel.

Method

Sample(;
  num_samples=1000,
  num_warmup=1000,
  save_warmup=false,
  thin=1,
  adapt=CmdStan.Adapt(),
  algorithm=SamplingAlgorithm()
)

Optional arguments

* `num_samples::Int64`          : Number of sampling iterations ( >= 0 )
* `num_warmup::Int64`           : Number of warmup iterations ( >= 0 )
* `save_warmup::Bool`           : Include warmup samples in output
* `thin::Int64`                 : Period between saved samples
* `adapt::Adapt`                : Warmup adaptation settings
* `algorithm::SamplingAlgorithm`: Sampling algorithm

Related help

?Stanmodel                      : Create a StanModel
?Adapt
?SamplingAlgorithm
source
CmdStan.AdaptType

Adapt type and constructor

Settings for adapt=CmdStan.Adapt() in Sample().

Method

Adapt(;
  engaged=true,
  gamma=0.05,
  delta=0.8,
  kappa=0.75,
  t0=10.0,
  init_buffer=75,
  term_buffer=50,
  window::25
)

Optional arguments

* `engaged::Bool`              : Adaptation engaged?
* `gamma::Float64`             : Adaptation regularization scale
* `delta::Float64`             : Adaptation target acceptance statistic
* `kappa::Float64`             : Adaptation relaxation exponent
* `t0::Float64`                : Adaptation iteration offset
* `init_buffer::Int64`         : Width of initial adaptation interval
* `term_buffer::Int64`         : Width of final fast adaptation interval
* `window::Int64`              : Initial width of slow adaptation interval

Related help

?Sample                        : Sampling settings
source
CmdStan.HmcType

Hmc type and constructor

Settings for algorithm=CmdStan.Hmc() in Sample().

Method

Hmc(;
  engine=CmdStan.Nuts(),
  metric=CmdStan.diag_e,
  stepsize=1.0,
  stepsize_jitter=1.0
)

Optional arguments

* `engine::Engine`              : Engine for Hamiltonian Monte Carlo
* `metric::Metric`              : Geometry for base manifold
* `stepsize::Float64`           : Stepsize for discrete evolutions
* `stepsize_jitter::Float64`    : Uniform random jitter of the stepsize [%]

Related help

?Sample                         : Sampling settings
?Engine                         : Engine for Hamiltonian Monte Carlo
?Nuts                           : Settings for Nuts
?Static                         : Settings for Static
?Metric                         : Base manifold geometries
source
CmdStan.MetricType

Metric types

Geometry of base manifold

Types defined

* unit_e::Metric      : Euclidean manifold with unit metric
* dense_e::Metric     : Euclidean manifold with dense netric
* diag_e::Metric      : Euclidean manifold with diag netric
source
CmdStan.EngineType

Engine types

Engine for Hamiltonian Monte Carlo

Types defined

* Nuts       : No-U-Turn sampler
* Static     : Static integration time
source
CmdStan.NutsType

Nuts type and constructor

Settings for engine=Nuts() in Hmc().

Method

Nuts(;max_depth=10)

Optional arguments

* `max_depth::Number`           : Maximum tree depth

Related help

?Sample                         : Sampling settings
?Engine                         : Engine for Hamiltonian Monte Carlo
source
CmdStan.StaticType

Static type and constructor

Settings for engine=Static() in Hmc().

Method

Static(;int_time=2 * pi)

Optional arguments

* `;int_time::Number`          : Static integration time

Related help

?Sample                        : Sampling settings
?Engine                        : Engine for Hamiltonian Monte Carlo
source
CmdStan.GradientType

Gradient type and constructor

Settings for diagnostic=Gradient() in Diagnose().

Method

Gradient(;epsilon=1e-6, error=1e-6)

Optional arguments

* `epsilon::Float64`           : Finite difference step size
* `error::Float64`             : Error threshold

Related help

?Diagnose                      : Diagnostic method
source
CmdStan.DiagnoseType

Diagnose type and constructor

Method

Diagnose(;d=Gradient())

Optional arguments

* `d::Diagnostics`            : Finite difference step size

Related help

?Diagnostics                  : Diagnostic methods
?Gradient                     : Currently sole diagnostic method
source
CmdStan.OptimizeAlgorithmType

OptimizeAlgorithm types

Types defined

* Lbfgs::OptimizeAlgorithm   : Euclidean manifold with unit metric
* Bfgs::OptimizeAlgorithm    : Euclidean manifold with unit metric
* Newton::OptimizeAlgorithm  : Euclidean manifold with unit metric
source
CmdStan.OptimizeType

Optimize type and constructor

Settings for Optimize top level method.

Method

Optimize(;
  method=Lbfgs(),
  iter=2000,
  save_iterations=false
)

Optional arguments

* `method::OptimizeMethod`      : Optimization algorithm
* `iter::Int`                   : Total number of iterations
* `save_iterations::Bool`       : Stream optimization progress to output

Related help

?Stanmodel                      : Create a StanModel
?OptimizeAlgorithm              : Available algorithms
source
CmdStan.LbfgsType

Lbfgs type and constructor

Settings for method=Lbfgs() in Optimize().

Method

Lbfgs(;init_alpha=0.001, tol_obj=1e-8, tol_grad=1e-8, tol_param=1e-8, history_size=5)

Optional arguments

* `init_alpha::Float64`        : Linear search step size for first iteration
* `tol_obj::Float64`           : Convergence tolerance for objective function
* `tol_rel_obj::Float64`       : Relative change tolerance in objective function
* `tol_grad::Float64`          : Convergence tolerance on norm of gradient
* `tol_rel_grad::Float64`      : Relative change tolerance on norm of gradient
* `tol_param::Float64`         : Convergence tolerance on parameter values
* `history_size::Int`          : No of update vectors to use in Hessian approx

Related help

?OptimizeMethod               : List of available optimize methods
?Optimize                      : Optimize arguments
source
CmdStan.BfgsType

Bfgs type and constructor

Settings for method=Bfgs() in Optimize().

Method

Bfgs(;init_alpha=0.001, tol_obj=1e-8, tol_rel_obj=1e4, 
  tol_grad=1e-8, tol_rel_grad=1e7, tol_param=1e-8)

Optional arguments

* `init_alpha::Float64`        : Linear search step size for first iteration
* `tol_obj::Float64`           : Convergence tolerance for objective function
* `tol_rel_obj::Float64`       : Relative change tolerance in objective function
* `tol_grad::Float64`          : Convergence tolerance on norm of gradient
* `tol_rel_grad::Float64`      : Relative change tolerance on norm of gradient
* `tol_param::Float64`         : Convergence tolerance on parameter values

Related help

?OptimizeMethod                : List of available optimize methods
?Optimize                      : Optimize arguments
source
CmdStan.NewtonType

Newton type and constructor

Settings for method=Newton() in Optimize().

Method

Newton()

Related help

?OptimizeMethod                : List of available optimize methods
?Optimize                      : Optimize arguments
source
CmdStan.VariationalType

Variational type and constructor

Settings for method=Variational() in Stanmodel.

Method

Variational(;
  grad_samples=1,
  elbo_samples=100,
  eta_adagrad=0.1,
  iter=10000,
  tol_rel_obj=0.01,
  eval_elbo=100,
  algorithm=:meanfield,          
  output_samples=10000
)

Optional arguments

* `algorithm::Symbol`             : Variational inference algorithm
                                  : meanfield
                                  : fullrank
* `iter::Int64`                   : Maximum number of iterations
* `grad_samples::Int`             : No of samples for MC estimate of gradients
* `elbo_samples::Int`             : No of samples for MC estimate of ELBO
* `eta::Float64`                  : Step size weighting parameter for adaptive sequence
* `adapt::Adapt`                  : Warmup adaptations
* `tol_rel_obj::Float64`          : Tolerance on the relative norm of objective
* `eval_elbo::Int`                : Tolerance on the relative norm of objective
* `output_samples::Int`           : Number of posterior samples to draw and save

Related help

?Stanmodel                        : Create a StanModel
?Stan.Method                      : Available top level methods
?Stan.Adapt                       : Adaptation settings
source

Utilities

CmdStan.cmdlineFunction

cmdline

Recursively parse the model to construct command line.

Method

cmdline(m)

Required arguments

* `m::Stanmodel`                : Stanmodel

Related help

?Stanmodel                      : Create a StanModel
source
CmdStan.check_dct_typeFunction

check_dct_type(

Check if dct == Dict[] and has length > 0.

Method

check_dct_type((dct)

Required arguments

* `dct::Dict{String, Any}`      : Observed data or parameter init data
source
CmdStan.convert_a3dFunction

convert_a3d

Convert the output file created by cmdstan to the shape of choice.

Method

convert_a3d(a3d_array, cnames, ::Val{Symbol}; start=1)

Required arguments

* `a3d_array::Array{Float64, 3},`      : Read in from output files created by cmdstan                                   
* `cnames::Vector{AbstractString}`     : Monitored variable names

Optional arguments

* `::Val{Symbol}`                      : Output format, default is :mcmcchains
* `::start=1`                          : First draw for MCMCChains.Chains

Method called is based on the output_format defined in the stanmodel, e.g.:

stanmodel = Stanmodel(num_samples=1200, thin=2, name="bernoulli", model=bernoullimodel, output_format=:mcmcchains);

Current formats supported are:

  1. :array (a3d_array format, the default for CmdStan)
  2. :mcmcchains (MCMCChains.Chains object)
  3. :dataframes (Array of DataFrames.DataFrame objects)
  4. :namedtuple (NamedTuple object)

Option 2 is available if MCMCChains is loaded.


### Return values

julia

  • res : Draws converted to the specified format.

```

source

convert_a3d

Convert the output file created by cmdstan to an Array of DataFrames.

convert_a3d(a3d_array, cnames, ; start)
source

convert_a3d

Convert the output file created by cmdstan to NamedTuples.

convert_a3d(a3d_array, cnames, ; start)
source
CmdStan.Fixed_paramType

Fixed_param type and constructor

Settings for algorithm=CmdStan.Fixed_param() in Sample().

Method

Fixed_param()

Related help

?Sample                        : Sampling settings
?Engine                        : Engine for Hamiltonian Monte Carlo
?Nuts                          : Settings for Nuts
?Static                        : Settings for Static
?Metric                        : Base manifold geometries
source
CmdStan.parFunction

par

Rewrite dct to R format in file.

Method

par(cmds)

Required arguments

* `cmds::Array{Base.AbstractCmd,1}`    : Multiple commands to concatenate

or

* `cmd::Base.AbstractCmd`              : Single command to be
* `n::Number`                            inserted n times into cmd


or
* `cmd::Array{String, 1}`              : Array of cmds as Strings
source
CmdStan.read_optimizeFunction

read_optimize

Read optimize output file created by cmdstan.

Method

read_optimize(m::Stanmodel)

Required arguments

* `m::Stanmodel`    : Stanmodel object
source
CmdStan.read_samplesFunction

read_samples

Read sample output files created by cmdstan.

Method

read_samples(m::Stanmodel)

Required arguments

* `m::Stanmodel`    : Stanmodel object
source
CmdStan.read_variationalFunction

read_variational

Read variational sample output files created by cmdstan.

Method

read_variational(m::Stanmodel)

Required arguments

* `m::Stanmodel`    : Stanmodel object
source
CmdStan.read_diagnoseFunction

read_diagnose

Read diagnose output file created by cmdstan.

Method

read_diagnose(m::Stanmodel)

Required arguments

* `m::Stanmodel`    : Stanmodel object
source
CmdStan.stan_summaryFunction

Method stan_summary

Display cmdstan summary

Method

stan_summary(
  model::StanModel,
  file::String; 
  CmdStanDir=CMDSTAN_HOME
)

Required arguments

* `model::Stanmodel             : Stanmodel
* `file::String`                : Name of file with samples

Optional arguments

* CmdStanDir=CMDSTAN_HOME       : cmdstan directory for stansummary program

Related help

?Stan.stan                      : Execute a StanModel
source

Method stan_summary

Display cmdstan summary

Method

stan_summary(
  model::Stanmodel
  filecmd::Cmd; 
  CmdStanDir=CMDSTAN_HOME
)

Required arguments

* `model::Stanmodel`            : Stanmodel
* `filecmd::Cmd`                : Run command containing names of sample files

Optional arguments

* CmdStanDir=CMDSTAN_HOME       : cmdstan directory for stansummary program

Related help

?Stan.stan                      : Create a StanModel
source
CmdStan.read_summaryFunction

read_summary

Read summary output file created by stansummary.

Method

read_summary(m::Stanmodel)

Required arguments

* `m::Stanmodel`    : Stanmodel object
source
CmdStan.update_R_fileFunction

updateRfile

Rewrite a dictionary of observed data or initial parameter values in R dump file format to a file.

Method

update_R_file{T<:Any}(file, dct)

Required arguments

* `file::String`                : R file name
* `dct::Dict{String, T}`        : Dictionary to format in R
source

Index