View source: R/stanfit-class.R
| sflist2stanfit | R Documentation |
This function takes a list of stanfit objects and returns a
consolidated stanfit object. The stanfit objects to be merged
need to have the same configuration of iteration, warmup, and thin, besides
being from the same model. This could facilitate some parallel usage of RStan.
For example, if we call stan by parallel and it returns a list of
stanfit objects, this function can be used to create one stanfit
object from the list.
sflist2stanfit(sflist)
sflist |
A list of |
An S4 object of stanfit consolidated from all the input stanfit
objects.
This function should be called in rare circumstances because
sampling has a cores argument that allows multiple
chains to be executed in parallel. However, if you need to depart from that,
the best practice is to use sflist2stanfit on a list of stanfit
objects created with the same seed but different chain_id (see
example below). Using the same seed but different chain_id can make sure
the random number generations for all chains are not correlated.
This function would do some check to see if the stanfit objects in the input list
can be merged. But the check is not sufficient. So generally, it is the
user's responsibility to make sure the input is correct so that the merging
makes sense.
The date in the new stanfit object is when it is merged.
get_seed function for the new consolidated stanfit object only returns
the seed used in the first chain of the new object.
The sampler such as NUTS2 that is displayed in the printout by print
is the sampler used for the first chain. The print method assumes the samplers
are the same for all chains.
The included stanmodel object, which includes the compiled model,
in the new stanfit object is from the first element of the input list.
The Stan Development Team Stan Modeling Language User's Guide and Reference Manual. https://mc-stan.org/.
stan
## Not run:
library(rstan)
scode <- "
data {
int<lower=1> N;
}
parameters {
real y1[N];
real y2[N];
}
model {
y1 ~ normal(0, 1);
y2 ~ double_exponential(0, 2);
}
"
seed <- 123 # or any other integer
foo_data <- list(N = 2)
foo <- stan(model_code = scode, data = foo_data, chains = 1, iter = 1)
f1 <- stan(fit = foo, data = foo_data, chains = 1, seed = seed, chain_id = 1)
f2 <- stan(fit = foo, data = foo_data, chains = 2, seed = seed, chain_id = 2:3)
f12 <- sflist2stanfit(list(f1, f2))
## parallel stan call for unix-like OS
library(parallel)
if (.Platform$OS.type == "unix") {
sflist1 <-
mclapply(1:4, mc.cores = 2,
function(i) stan(fit = foo, data = foo_data, seed = seed,
chains = 1, chain_id = i, refresh = -1))
f3 <- sflist2stanfit(sflist1)
}
if (.Platform$OS.type == "windows") { # also works on non-Windows
CL <- makeCluster(2)
clusterExport(cl = CL, c("foo_data", "foo", "seed"))
sflist1 <- parLapply(CL, 1:4, fun = function(cid) {
require(rstan)
stan(fit = foo, data = foo_data, chains = 1,
iter = 2000, seed = seed, chain_id = cid)
})
fit <- sflist2stanfit(sflist1)
print(fit)
stopCluster(CL)
} # end example for Windows
## End(Not run)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.