simmr_input object through the main simmr Markov chain Monte
Carlo (MCMC) functionR/simmr_mcmc.R
simmr_mcmc.RdThis is the main function of simmr. It takes a simmr_input object
created via simmr_load, runs an MCMC to determine the dietary
proportions, and then outputs a simmr_output object for further
analysis and plotting via summary.simmr_output and
plot.simmr_output.
An object created via the function simmr_load
A list of values including arguments named: means
and sd which represent the prior means and standard deviations of the
dietary proportions in centralised log-ratio space; shape and
rate which represent the prior distribution on the residual standard
deviation. These can usually be
left at their default values unless you wish to include to include prior
information, in which case you should use the function
simmr_elicit.
A list of values including arguments named iter
(number of iterations), burn (size of burn-in), thin (amount
of thinning), and n.chain (number of MCMC chains).
An object of class simmr_output with two named top-level
components:
The simmr_input object given to the
simmr_mcmc function
A set of MCMC chains of class
mcmc.list from the coda package. These can be analysed using the
summary.simmr_output and plot.simmr_output
functions.
If, after running simmr_mcmc the convergence diagnostics in
summary.simmr_output are not satisfactory, the values of
iter, burn and thin in mcmc_control should be
increased by a factor of 10.
Andrew C. Parnell, Donald L. Phillips, Stuart Bearhop, Brice X. Semmens, Eric J. Ward, Jonathan W. Moore, Andrew L. Jackson, Jonathan Grey, David J. Kelly, and Richard Inger. Bayesian stable isotope mixing models. Environmetrics, 24(6):387–399, 2013.
Andrew C Parnell, Richard Inger, Stuart Bearhop, and Andrew L Jackson. Source partitioning using stable isotopes: coping with too much variation. PLoS ONE, 5(3):5, 2010.
simmr_load for creating objects suitable for this
function, plot.simmr_input for creating isospace plots,
summary.simmr_output for summarising output, and
plot.simmr_output for plotting output.
if (FALSE) {
## See the package vignette for a detailed run through of these 4 examples
# Data set 1: 10 obs on 2 isos, 4 sources, with tefs and concdep
data(geese_data_day1)
simmr_1 <- with(
geese_data_day1,
simmr_load(
mixtures = mixtures,
source_names = source_names,
source_means = source_means,
source_sds = source_sds,
correction_means = correction_means,
correction_sds = correction_sds,
concentration_means = concentration_means
)
)
# Plot
plot(simmr_1)
# Print
simmr_1
# MCMC run
simmr_1_out <- simmr_mcmc(simmr_1)
# Print it
print(simmr_1_out)
# Summary
summary(simmr_1_out, type = "diagnostics")
summary(simmr_1_out, type = "correlations")
summary(simmr_1_out, type = "statistics")
ans <- summary(simmr_1_out, type = c("quantiles", "statistics"))
# Plot
plot(simmr_1_out, type = "boxplot")
plot(simmr_1_out, type = "histogram")
plot(simmr_1_out, type = "density")
plot(simmr_1_out, type = "matrix")
# Compare two sources
compare_sources(simmr_1_out, source_names = c("Zostera", "Enteromorpha"))
# Compare multiple sources
compare_sources(simmr_1_out)
#####################################################################################
# A version with just one observation
data(geese_data_day1)
simmr_2 <- with(
geese_data_day1,
simmr_load(
mixtures = mixtures[1, , drop = FALSE],
source_names = source_names,
source_means = source_means,
source_sds = source_sds,
correction_means = correction_means,
correction_sds = correction_sds,
concentration_means = concentration_means
)
)
# Plot
plot(simmr_2)
# MCMC run - automatically detects the single observation
simmr_2_out <- simmr_mcmc(simmr_2)
# Print it
print(simmr_2_out)
# Summary
summary(simmr_2_out)
summary(simmr_2_out, type = "diagnostics")
ans <- summary(simmr_2_out, type = c("quantiles"))
# Plot
plot(simmr_2_out)
plot(simmr_2_out, type = "boxplot")
plot(simmr_2_out, type = "histogram")
plot(simmr_2_out, type = "density")
plot(simmr_2_out, type = "matrix")
#####################################################################################
# Data set 2: 3 isotopes (d13C, d15N and d34S), 30 observations, 4 sources
data(simmr_data_2)
simmr_3 <- with(
simmr_data_2,
simmr_load(
mixtures = mixtures,
source_names = source_names,
source_means = source_means,
source_sds = source_sds,
correction_means = correction_means,
correction_sds = correction_sds,
concentration_means = concentration_means
)
)
# Get summary
print(simmr_3)
# Plot 3 times
plot(simmr_3)
plot(simmr_3, tracers = c(2, 3))
plot(simmr_3, tracers = c(1, 3))
# See vignette('simmr') for fancier axis labels
# MCMC run
simmr_3_out <- simmr_mcmc(simmr_3)
# Print it
print(simmr_3_out)
# Summary
summary(simmr_3_out)
summary(simmr_3_out, type = "diagnostics")
summary(simmr_3_out, type = "quantiles")
summary(simmr_3_out, type = "correlations")
# Plot
plot(simmr_3_out)
plot(simmr_3_out, type = "boxplot")
plot(simmr_3_out, type = "histogram")
plot(simmr_3_out, type = "density")
plot(simmr_3_out, type = "matrix")
#####################################################################################
# Data set 5 - Multiple groups Geese data from Inger et al 2006
# Do this in raw data format - Note that there's quite a few mixtures!
data(geese_data)
simmr_5 <- with(
geese_data,
simmr_load(
mixtures = mixtures,
source_names = source_names,
source_means = source_means,
source_sds = source_sds,
correction_means = correction_means,
correction_sds = correction_sds,
concentration_means = concentration_means,
group = groups
)
)
# Plot
plot(simmr_5,
xlab = expression(paste(delta^13, "C (per mille)", sep = "")),
ylab = expression(paste(delta^15, "N (per mille)", sep = "")),
title = "Isospace plot of Inger et al Geese data"
)
# Run MCMC for each group
simmr_5_out <- simmr_mcmc(simmr_5)
# Summarise output
summary(simmr_5_out, type = "quantiles", group = 1)
summary(simmr_5_out, type = "quantiles", group = c(1, 3))
summary(simmr_5_out, type = c("quantiles", "statistics"), group = c(1, 3))
# Plot - only a single group allowed
plot(simmr_5_out, type = "boxplot", group = 2, title = "simmr output group 2")
plot(simmr_5_out, type = c("density", "matrix"), grp = 6, title = "simmr output group 6")
# Compare sources within a group
compare_sources(simmr_5_out, source_names = c("Zostera", "U.lactuca"), group = 2)
compare_sources(simmr_5_out, group = 2)
# Compare between groups
compare_groups(simmr_5_out, source = "Zostera", groups = 1:2)
compare_groups(simmr_5_out, source = "Zostera", groups = 1:3)
compare_groups(simmr_5_out, source = "U.lactuca", groups = c(4:5, 7, 2))
}