Title: | Capture-Mark-Recapture Analysis using Multiple Non-Invasive Marks |
---|---|
Description: | Traditional and spatial capture-mark-recapture analysis with multiple non-invasive marks. The models implemented in 'multimark' combine encounter history data arising from two different non-invasive "marks", such as images of left-sided and right-sided pelage patterns of bilaterally asymmetrical species, to estimate abundance and related demographic parameters while accounting for imperfect detection. Bayesian models are specified using simple formulae and fitted using Markov chain Monte Carlo. Addressing deficiencies in currently available software, 'multimark' also provides a user-friendly interface for performing Bayesian multimodel inference using non-spatial or spatial capture-recapture data consisting of a single conventional mark or multiple non-invasive marks. See McClintock (2015) <doi:10.1002/ece3.1676> and Maronde et al. (2020) <doi:10.1002/ece3.6990>. |
Authors: | Brett T. McClintock [aut, cre], Acho Arnold [ctb, cph] (C original matrix library, https://github.com/najela/matrix.h), Barry Brown [ctb] (Fortran original ranlib library), James Lovato [ctb] (Fortran original ranlib library), John Burkardt [ctb] (C original ranlib library, http://people.sc.fsu.edu/~jburkardt/c_src/ranlib), Cleve Moler [ctb] (C original linpack library, http://www.kkant.net/geist/ranlib/), Arjun Gopalaswamy [ctb] (modified snippets of R package SPACECAP code) |
Maintainer: | Brett T. McClintock <[email protected]> |
License: | GPL-2 |
Version: | 2.1.6 |
Built: | 2024-11-21 04:20:29 UTC |
Source: | https://github.com/bmcclintock/multimark |
Example bobcat data for multimark
package.
The data are summarized in a 46x8 matrix containing observed encounter histories for 46 bobcats across 8 sampling occasions. Bobcats are bilaterially asymmetrical, and sampling was conducted using camera stations consisting of a single camera.
Because the left-side cannot be reconciled with the right-side, the two types of “marks” in this case are the pelage patterns on the left- and right-side of each individual. Encounter type 0 corresponds to non-detection, encounter type 1 corresponds to left-sided detection, encounter type 2 corresponds to right-sided detection.
Both-sided encounters were never observed in this dataset, hence the most appropriate multimark
data type is data.type="never".
McClintock, B. T., Conn, P. B., Alonso, R. S., and Crooks, K. R. 2013. Integrated modeling of bilateral photo-identification data in mark-recapture analyses. Ecology 94: 1464-1471.
data(bobcat)
data(bobcat)
Example spatial bobcat data for multimark
package.
These spatial capture-recapture data with multiple mark types are summarized in a list of length 3 containing the following objects:
Enc.Mat
is a 42 x (noccas*ntraps) matrix containing observed encounter histories for 42 bobcats across noccas=187
sampling occasions and ntraps=30
traps. The first 187 columns correspond to trap 1, the second 187 columns corresopond to trap 2, etc.
trapCoords
is a matrix of dimension ntraps
x (2 + noccas
) indicating the Cartesian coordinates and operating occasions for the traps, where rows correspond to trap, the first column the x-coordinate, and the second column the y-coordinate. The last noccas
columns indicate whether or not the trap was operating on each of the occasions, where ‘1’ indciates the trap was operating and ‘0’ indicates the trap was not operating.
studyArea
is a 3-column matrix containing the coordinates for the centroids of the contiguous grid of 1023 cells that define the study area and available habitat. Each row corresponds to a grid cell. The first 2 columns indicate the Cartesian x- and y-coordinate for the centroid of each grid cell, and the third column indicates whether the cell is available habitat (=1) or not (=0). The grid cells are 0.65x0.65km resolution.
Bobcats are bilaterially asymmetrical, and sampling was conducted using camera stations consisting of a single camera. Because the left-side cannot be reconciled with the right-side, the two types of “marks” in this case are the pelage patterns on the left- and right-side of each individual. Encounter type 0 corresponds to non-detection, encounter type 1 corresponds to left-sided detection, encounter type 2 corresponds to right-sided detection.
Both-sided encounters were never observed in this dataset, hence the most appropriate multimark
data type is data.type="never".
The first 15 rows of bobcatSCR$Enc.Mat
correspond to individuals for which both the left and right sides were known because they were physically captured for telemetry deployments prior to sampling surveys. The encounter histories for these 15 individuals are therefore known with certainty and should be specified as such using the known
argument in processdataSCR
and/or multimarkClosedSCR
(see example below).
These data were obtained from the R package SPIM
(Augustine et al. 2017) and modified by projecting onto a regular rectangular grid consisting of square grid cells (as is required by the spatial capture-recapture models in multimark
).
We thank B. Augustine and co-authors for making these data publicly available in the SPIM
package (Augustine et al. 2017).
Augustine, B., Royle, J.A., Kelly, M., Satter, C., Alonso, R., Boydston, E. and Crooks, K. 2017. Spatial capture-recapture with partial identity: an application to camera traps. bioRxiv doi: https://doi.org/10.1101/056804
multimarkClosedSCR
, processdataSCR
data(bobcatSCR) #plot the traps and available habitat within the study area plotSpatialData(trapCoords=bobcatSCR$trapCoords,studyArea=bobcatSCR$studyArea) # This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin # Fit spatial model to tiger data Enc.Mat <- bobcatSCR$Enc.Mat trapCoords <- bobcatSCR$trapCoords studyArea <- bobcatSCR$studyArea # specify known encounter histories known <- c(rep(1,15),rep(0,nrow(Enc.Mat)-15)) # specify prior bounds for sigma2_scr sig_bounds <- c(0.1,max(diff(range(studyArea[,"x"])),diff(range(studyArea[,"y"])))) mmsSCR <- processdataSCR(Enc.Mat,trapCoords,studyArea,known=known) bobcatSCR.dot.type <- multimarkClosedSCR(mms=mmsSCR,iter=200,adapt=100,burnin=100, sigma_bounds=sig_bounds) summary(bobcatSCR.dot.type$mcmc)
data(bobcatSCR) #plot the traps and available habitat within the study area plotSpatialData(trapCoords=bobcatSCR$trapCoords,studyArea=bobcatSCR$studyArea) # This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin # Fit spatial model to tiger data Enc.Mat <- bobcatSCR$Enc.Mat trapCoords <- bobcatSCR$trapCoords studyArea <- bobcatSCR$studyArea # specify known encounter histories known <- c(rep(1,15),rep(0,nrow(Enc.Mat)-15)) # specify prior bounds for sigma2_scr sig_bounds <- c(0.1,max(diff(range(studyArea[,"x"])),diff(range(studyArea[,"y"])))) mmsSCR <- processdataSCR(Enc.Mat,trapCoords,studyArea,known=known) bobcatSCR.dot.type <- multimarkClosedSCR(mms=mmsSCR,iter=200,adapt=100,burnin=100, sigma_bounds=sig_bounds) summary(bobcatSCR.dot.type$mcmc)
This function calculates posterior population density estimates from multimarkClosedSCR
output as D = N/A, where D is density, N is abundance, and A is the area of available habitat within the study area.
getdensityClosedSCR(out)
getdensityClosedSCR(out)
out |
List of output returned by |
An object of class mcmc.list
containing the following:
D |
Posterior samples for density. |
Brett T. McClintock
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Run behavior model for simulated data with constant detection probability (i.e., mod.p=~c) sim.data<-simdataClosedSCR() Enc.Mat<-sim.data$Enc.Mat trapCoords<-sim.data$spatialInputs$trapCoords studyArea<-sim.data$spatialInputs$studyArea example.dot <- multimarkClosedSCR(Enc.Mat,trapCoords,studyArea,mod.p=~1) #Calculate capture and recapture probabilities D <- getdensityClosedSCR(example.dot) summary(D)
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Run behavior model for simulated data with constant detection probability (i.e., mod.p=~c) sim.data<-simdataClosedSCR() Enc.Mat<-sim.data$Enc.Mat trapCoords<-sim.data$spatialInputs$trapCoords studyArea<-sim.data$spatialInputs$studyArea example.dot <- multimarkClosedSCR(Enc.Mat,trapCoords,studyArea,mod.p=~1) #Calculate capture and recapture probabilities D <- getdensityClosedSCR(example.dot) summary(D)
This function calculates posterior capture () and survival (
) probabilities for each sampling occasion from
multimarkCJS
output.
getprobsCJS(out, link = "probit")
getprobsCJS(out, link = "probit")
out |
List of output returned by |
link |
Link function for |
An object of class mcmc.list
containing the following:
p |
Posterior samples for capture probability ( |
phi |
Posterior samples for survival probability ( |
Brett T. McClintock
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Simulate open population data with temporal variation in survival noccas <- 5 data <- simdataCJS(noccas=noccas, phibeta=rnorm(noccas-1,1.6,0.1)) #Fit open population model with temporal variation in survival sim.time <- multimarkCJS(data$Enc.Mat,mod.phi=~time) #Calculate capture and survival probabilities for each cohort and time pphi <- getprobsCJS(sim.time) summary(pphi)
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Simulate open population data with temporal variation in survival noccas <- 5 data <- simdataCJS(noccas=noccas, phibeta=rnorm(noccas-1,1.6,0.1)) #Fit open population model with temporal variation in survival sim.time <- multimarkCJS(data$Enc.Mat,mod.phi=~time) #Calculate capture and survival probabilities for each cohort and time pphi <- getprobsCJS(sim.time) summary(pphi)
This function calculates posterior capture () and recapture (
) probabilities for each sampling occasion from
multimarkClosed
output.
getprobsClosed(out, link = "logit")
getprobsClosed(out, link = "logit")
out |
List of output returned by |
link |
Link function for detection probability. Must be " |
An object of class mcmc.list
containing the following:
p |
Posterior samples for capture probability ( |
c |
Posterior samples for recapture probability ( |
Brett T. McClintock
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Run behavior model for bobcat data with constant detection probability (i.e., mod.p=~c) bobcat.c <- multimarkClosed(bobcat,mod.p=~c) #Calculate capture and recapture probabilities pc <- getprobsClosed(bobcat.c) summary(pc)
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Run behavior model for bobcat data with constant detection probability (i.e., mod.p=~c) bobcat.c <- multimarkClosed(bobcat,mod.p=~c) #Calculate capture and recapture probabilities pc <- getprobsClosed(bobcat.c) summary(pc)
This function calculates posterior spatial capture () and recapture (
) probabilities (at zero distance from an activity center) for each sampling occasion from
multimarkClosedSCR
output.
getprobsClosedSCR(out, link = "cloglog")
getprobsClosedSCR(out, link = "cloglog")
out |
List of output returned by |
link |
Link function for detection probability. Must be " |
An object of class mcmc.list
containing the following:
p |
Posterior samples for capture probability ( |
c |
Posterior samples for recapture probability ( |
Brett T. McClintock
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Run behavior model for simulated data with constant detection probability (i.e., mod.p=~c) sim.data<-simdataClosedSCR() Enc.Mat<-sim.data$Enc.Mat trapCoords<-sim.data$spatialInputs$trapCoords studyArea<-sim.data$spatialInputs$studyArea example.c <- multimarkClosedSCR(Enc.Mat,trapCoords,studyArea,mod.p=~c, iter=1000,adapt=500,burnin=500) #Calculate capture and recapture probabilities pc <- getprobsClosedSCR(example.c) summary(pc)
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Run behavior model for simulated data with constant detection probability (i.e., mod.p=~c) sim.data<-simdataClosedSCR() Enc.Mat<-sim.data$Enc.Mat trapCoords<-sim.data$spatialInputs$trapCoords studyArea<-sim.data$spatialInputs$studyArea example.c <- multimarkClosedSCR(Enc.Mat,trapCoords,studyArea,mod.p=~c, iter=1000,adapt=500,burnin=500) #Calculate capture and recapture probabilities pc <- getprobsClosedSCR(example.c) summary(pc)
This function fits Cormack-Jolly-Seber (CJS) open population models for survival probability () and capture probability (
) for “traditional” capture-mark-recapture data consisting of a single mark type. Using Bayesian analysis methods, Markov chain Monte Carlo (MCMC) is used to draw samples from the joint posterior distribution.
markCJS( Enc.Mat, covs = data.frame(), mod.p = ~1, mod.phi = ~1, parms = c("pbeta", "phibeta"), nchains = 1, iter = 12000, adapt = 1000, bin = 50, thin = 1, burnin = 2000, taccept = 0.44, tuneadjust = 0.95, proppbeta = 0.1, propzp = 1, propsigmap = 1, propphibeta = 0.1, propzphi = 1, propsigmaphi = 1, pbeta0 = 0, pSigma0 = 1, phibeta0 = 0, phiSigma0 = 1, l0p = 1, d0p = 0.01, l0phi = 1, d0phi = 0.01, initial.values = NULL, link = "probit", printlog = FALSE, ... )
markCJS( Enc.Mat, covs = data.frame(), mod.p = ~1, mod.phi = ~1, parms = c("pbeta", "phibeta"), nchains = 1, iter = 12000, adapt = 1000, bin = 50, thin = 1, burnin = 2000, taccept = 0.44, tuneadjust = 0.95, proppbeta = 0.1, propzp = 1, propsigmap = 1, propphibeta = 0.1, propzphi = 1, propsigmaphi = 1, pbeta0 = 0, pSigma0 = 1, phibeta0 = 0, phiSigma0 = 1, l0p = 1, d0p = 0.01, l0phi = 1, d0phi = 0.01, initial.values = NULL, link = "probit", printlog = FALSE, ... )
Enc.Mat |
A matrix of observed encounter histories with rows corresponding to individuals and columns corresponding to sampling occasions. With a single mark type, encounter histories consist of only non-detections (0) and type 1 encounters (1). |
covs |
A data frame of temporal covariates for detection probabilities (ignored unless |
mod.p |
Model formula for detection probability ( |
mod.phi |
Model formula for survival probability ( |
parms |
A character vector giving the names of the parameters and latent variables to monitor. Possible parameters are probit-scale detection probability parameters (" |
nchains |
The number of parallel MCMC chains for the model. |
iter |
The number of MCMC iterations. |
adapt |
Ignored; no adaptive phase is needed for "probit" link. |
bin |
Ignored; no adaptive phase is needed for "probit" link. |
thin |
Thinning interval for monitored parameters. |
burnin |
Number of burn-in iterations ( |
taccept |
Ignored; no adaptive phase is needed for "probit" link. |
tuneadjust |
Ignored; no adaptive phase is needed for "probit" link. |
proppbeta |
Ignored; no adaptive phase is needed for "probit" link. |
propzp |
Ignored; no adaptive phase is needed for "probit" link. |
propsigmap |
Ignored; no adaptive phase is needed for "probit" link. |
propphibeta |
Ignored; no adaptive phase is needed for "probit" link. |
propzphi |
Ignored; no adaptive phase is needed for "probit" link. |
propsigmaphi |
Ignored; no adaptive phase is needed for "probit" link. |
pbeta0 |
Scaler or vector (of length k) specifying mean of pbeta ~ multivariateNormal(pbeta0, pSigma0) prior. If |
pSigma0 |
Scaler or k x k matrix specifying covariance matrix of pbeta ~ multivariateNormal(pbeta0, pSigma0) prior. If |
phibeta0 |
Scaler or vector (of length k) specifying mean of phibeta ~ multivariateNormal(phibeta0, phiSigma0) prior. If |
phiSigma0 |
Scaler or k x k matrix specifying covariance matrix of phibeta ~ multivariateNormal(phibeta0, phiSigma0) prior. If |
l0p |
Specifies "shape" parameter for [sigma2_zp] ~ invGamma(l0p,d0p) prior. Default is |
d0p |
Specifies "scale" parameter for [sigma2_zp] ~ invGamma(l0p,d0p) prior. Default is |
l0phi |
Specifies "shape" parameter for [sigma2_zphi] ~ invGamma(l0phi,d0phi) prior. Default is |
d0phi |
Specifies "scale" parameter for [sigma2_zphi] ~ invGamma(l0phi,d0phi) prior. Default is |
initial.values |
OOptional list of |
link |
Link function for survival and capture probabilities. Only probit link is currently implemented. |
printlog |
Logical indicating whether to print the progress of chains and any errors to a log file in the working directory. Ignored when |
... |
Additional " |
The first time markCJS
(or markClosed
) is called, it will likely produce a firewall warning alerting users that R has requested the ability to accept incoming network connections. Incoming network connections are required to use parallel processing as implemented in multimarkCJS
. Note that setting parms="all"
is required for any markCJS
model output to be used in multimodelCJS
.
A list containing the following:
mcmc |
Markov chain Monte Carlo object of class |
mod.p |
Model formula for detection probability (as specified by |
mod.phi |
Model formula for survival probability (as specified by |
mod.delta |
Formula always |
DM |
A list of design matrices for detection and survival probability respectively generated by |
initial.values |
A list containing the parameter and latent variable values at iteration |
mms |
An object of class |
Brett T. McClintock
# These examples are excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Simulate open population data using defaults data <- simdataCJS(delta_1=1,delta_2=0)$Enc.Mat #Fit default open population model sim.dot <- markCJS(data) #Posterior summary for monitored parameters summary(sim.dot$mcmc) plot(sim.dot$mcmc) #Fit ``age'' model with 2 age classes (e.g., juvenile and adult) for survival #using 'parameters' and 'right' arguments from RMark::make.design.data sim.age <- markCJS(data,mod.phi=~age, parameters=list(Phi=list(age.bins=c(0,1,4))),right=FALSE) summary(getprobsCJS(sim.age))
# These examples are excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Simulate open population data using defaults data <- simdataCJS(delta_1=1,delta_2=0)$Enc.Mat #Fit default open population model sim.dot <- markCJS(data) #Posterior summary for monitored parameters summary(sim.dot$mcmc) plot(sim.dot$mcmc) #Fit ``age'' model with 2 age classes (e.g., juvenile and adult) for survival #using 'parameters' and 'right' arguments from RMark::make.design.data sim.age <- markCJS(data,mod.phi=~age, parameters=list(Phi=list(age.bins=c(0,1,4))),right=FALSE) summary(getprobsCJS(sim.age))
This function fits closed population abundance models for “traditional” capture-mark-recapture data consisting of a single mark type using Bayesian analysis methods. Markov chain Monte Carlo (MCMC) is used to draw samples from the joint posterior distribution.
markClosed( Enc.Mat, covs = data.frame(), mod.p = ~1, parms = c("pbeta", "N"), nchains = 1, iter = 12000, adapt = 1000, bin = 50, thin = 1, burnin = 2000, taccept = 0.44, tuneadjust = 0.95, proppbeta = 0.1, propzp = 1, propsigmap = 1, npoints = 500, a = 25, mu0 = 0, sigma2_mu0 = 1.75, initial.values = NULL, printlog = FALSE, ... )
markClosed( Enc.Mat, covs = data.frame(), mod.p = ~1, parms = c("pbeta", "N"), nchains = 1, iter = 12000, adapt = 1000, bin = 50, thin = 1, burnin = 2000, taccept = 0.44, tuneadjust = 0.95, proppbeta = 0.1, propzp = 1, propsigmap = 1, npoints = 500, a = 25, mu0 = 0, sigma2_mu0 = 1.75, initial.values = NULL, printlog = FALSE, ... )
Enc.Mat |
A matrix of observed encounter histories with rows corresponding to individuals and columns corresponding to sampling occasions. With a single mark type, encounter histories consist of only non-detections (0) and type 1 encounters (1). |
covs |
A data frame of temporal covariates for detection probabilities (ignored unless |
mod.p |
Model formula for detection probability. For example, |
parms |
A character vector giving the names of the parameters and latent variables to monitor. Possible parameters are logit-scale detection probability parameters (" |
nchains |
The number of parallel MCMC chains for the model. |
iter |
The number of MCMC iterations. |
adapt |
The number of iterations for proposal distribution adaptation. If |
bin |
Bin length for calculating acceptance rates during adaptive phase ( |
thin |
Thinning interval for monitored parameters. |
burnin |
Number of burn-in iterations ( |
taccept |
Target acceptance rate during adaptive phase ( |
tuneadjust |
Adjustment term during adaptive phase ( |
proppbeta |
Scaler or vector (of length k) specifying the initial standard deviation of the Normal(pbeta[j], proppbeta[j]) proposal distribution. If |
propzp |
Scaler or vector (of length M) specifying the initial standard deviation of the Normal(zp[i], propzp[i]) proposal distribution. If |
propsigmap |
Scaler specifying the initial Gamma(shape = 1/ |
npoints |
Number of Gauss-Hermite quadrature points to use for numerical integration. Accuracy increases with number of points, but so does computation time. |
a |
Scale parameter for [sigma_z] ~ half-Cauchy(a) prior for the individual hetegeneity term sigma_zp = sqrt(sigma2_zp). Default is “uninformative” |
mu0 |
Scaler or vector (of length k) specifying mean of pbeta[j] ~ Normal(mu0[j], sigma2_mu0[j]) prior. If |
sigma2_mu0 |
Scaler or vector (of length k) specifying variance of pbeta[j] ~ Normal(mu0[j], sigma2_mu0[j]) prior. If |
initial.values |
Optional list of |
printlog |
Logical indicating whether to print the progress of chains and any errors to a log file in the working directory. Ignored when |
... |
Additional " |
The first time markClosed
(or markCJS
) is called, it will likely produce a firewall warning alerting users that R has requested the ability to accept incoming network connections. Incoming network connections are required to use parallel processing as implemented in markClosed
. Note that setting parms="all"
is required for any markClosed
model output to be used in multimodelClosed
.
A list containing the following:
mcmc |
Markov chain Monte Carlo object of class |
mod.p |
Model formula for detection probability (as specified by |
mod.delta |
Formula always |
DM |
A list of design matrices for detection probability generated for model |
initial.values |
A list containing the parameter and latent variable values at iteration |
mms |
An object of class |
Brett T. McClintock
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Run single chain using the default model for simulated ``traditional'' data data<-simdataClosed(delta_1=1,delta_2=0)$Enc.Mat sim.dot<-markClosed(data) #Posterior summary for monitored parameters summary(sim.dot$mcmc) plot(sim.dot$mcmc)
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Run single chain using the default model for simulated ``traditional'' data data<-simdataClosed(delta_1=1,delta_2=0)$Enc.Mat sim.dot<-markClosed(data) #Posterior summary for monitored parameters summary(sim.dot$mcmc) plot(sim.dot$mcmc)
This function fits spatial population abundance models for “traditional” capture-mark-recapture data consisting of a single mark type using Bayesian analysis methods. Markov chain Monte Carlo (MCMC) is used to draw samples from the joint posterior distribution.
markClosedSCR( Enc.Mat, trapCoords, studyArea = NULL, buffer = NULL, ncells = 1024, covs = data.frame(), mod.p = ~1, detection = "half-normal", parms = c("pbeta", "N"), nchains = 1, iter = 12000, adapt = 1000, bin = 50, thin = 1, burnin = 2000, taccept = 0.44, tuneadjust = 0.95, proppbeta = 0.1, propsigma = 1, propcenter = NULL, sigma_bounds = NULL, mu0 = 0, sigma2_mu0 = 1.75, initial.values = NULL, scalemax = 10, printlog = FALSE, ... )
markClosedSCR( Enc.Mat, trapCoords, studyArea = NULL, buffer = NULL, ncells = 1024, covs = data.frame(), mod.p = ~1, detection = "half-normal", parms = c("pbeta", "N"), nchains = 1, iter = 12000, adapt = 1000, bin = 50, thin = 1, burnin = 2000, taccept = 0.44, tuneadjust = 0.95, proppbeta = 0.1, propsigma = 1, propcenter = NULL, sigma_bounds = NULL, mu0 = 0, sigma2_mu0 = 1.75, initial.values = NULL, scalemax = 10, printlog = FALSE, ... )
Enc.Mat |
A matrix containing the observed encounter histories with rows corresponding to individuals and ( |
trapCoords |
A matrix of dimension |
studyArea |
is a 3-column matrix containing the coordinates for the centroids a contiguous grid of cells that define the study area and available habitat. Each row corresponds to a grid cell. The first 2 columns (“x” and “y”) indicate the Cartesian x- and y-coordinate for the centroid of each grid cell, and the third column (“avail”) indicates whether the cell is available habitat (=1) or not (=0). All cells must have the same resolution. If |
buffer |
A scaler in same units as |
ncells |
The number of grid cells in the study area when |
covs |
A data frame of time- and/or trap-dependent covariates for detection probabilities (ignored unless |
mod.p |
Model formula for detection probability. For example, |
detection |
Model for detection probability as a function of distance from activity centers . Must be " |
parms |
A character vector giving the names of the parameters and latent variables to monitor. Possible parameters are cloglog-scale detection probability parameters (" |
nchains |
The number of parallel MCMC chains for the model. |
iter |
The number of MCMC iterations. |
adapt |
The number of iterations for proposal distribution adaptation. If |
bin |
Bin length for calculating acceptance rates during adaptive phase ( |
thin |
Thinning interval for monitored parameters. |
burnin |
Number of burn-in iterations ( |
taccept |
Target acceptance rate during adaptive phase ( |
tuneadjust |
Adjustment term during adaptive phase ( |
proppbeta |
Scaler or vector (of length k) specifying the initial standard deviation of the Normal(pbeta[j], proppbeta[j]) proposal distribution. If |
propsigma |
Scaler specifying the initial Gamma(shape = 1/ |
propcenter |
Scaler specifying the neighborhood distance when proposing updates to activity centers. When |
sigma_bounds |
Positive vector of length 2 for the lower and upper bounds for the [sigma_scr] ~ Uniform(sigma_bounds[1], sigma_bounds[2]) (or [sqrt(lambda)] when |
mu0 |
Scaler or vector (of length k) specifying mean of pbeta[j] ~ Normal(mu0[j], sigma2_mu0[j]) prior. If |
sigma2_mu0 |
Scaler or vector (of length k) specifying variance of pbeta[j] ~ Normal(mu0[j], sigma2_mu0[j]) prior. If |
initial.values |
Optional list of |
scalemax |
Upper bound for internal re-scaling of grid cell centroid coordinates. Default is |
printlog |
Logical indicating whether to print the progress of chains and any errors to a log file in the working directory. Ignored when |
... |
Additional " |
The first time markClosedSCR
is called, it will likely produce a firewall warning alerting users that R has requested the ability to accept incoming network connections. Incoming network connections are required to use parallel processing as implemented in markClosed
. Note that setting parms="all"
is required for any markClosed
model output to be used in multimodelClosed
.
A list containing the following:
mcmc |
Markov chain Monte Carlo object of class |
mod.p |
Model formula for detection probability (as specified by |
mod.delta |
Formula always |
mod.det |
Model formula for detection function (as specified by |
DM |
A list of design matrices for detection probability generated for model |
initial.values |
A list containing the parameter and latent variable values at iteration |
mms |
An object of class |
Brett T. McClintock
Gopalaswamy, A.M., Royle, J.A., Hines, J.E., Singh, P., Jathanna, D., Kumar, N. and Karanth, K.U. 2012. Program SPACECAP: software for estimating animal density using spatially explicit capture-recapture models. Methods in Ecology and Evolution 3:1067-1072.
King, R., McClintock, B. T., Kidney, D., and Borchers, D. L. 2016. Capture-recapture abundance estimation using a semi-complete data likelihood approach. The Annals of Applied Statistics 10: 264-285
Royle, J.A., Karanth, K.U., Gopalaswamy, A.M. and Kumar, N.S. 2009. Bayesian inference in camera trapping studies for a class of spatial capture-recapture models. Ecology 90: 3233-3244.
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Run single chain using the default model for ``traditional'' tiger data of Royle et al (2009) Enc.Mat<-tiger$Enc.Mat trapCoords<-tiger$trapCoords studyArea<-tiger$studyArea tiger.dot<-markClosedSCR(Enc.Mat,trapCoords,studyArea,iter=100,adapt=50,burnin=50) #Posterior summary for monitored parameters summary(tiger.dot$mcmc) plot(tiger.dot$mcmc)
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Run single chain using the default model for ``traditional'' tiger data of Royle et al (2009) Enc.Mat<-tiger$Enc.Mat trapCoords<-tiger$trapCoords studyArea<-tiger$studyArea tiger.dot<-markClosedSCR(Enc.Mat,trapCoords,studyArea,iter=100,adapt=50,burnin=50) #Posterior summary for monitored parameters summary(tiger.dot$mcmc) plot(tiger.dot$mcmc)
This function fits Cormack-Jolly-Seber (CJS) open population models for survival probability () and capture probability (
) from capture-mark-recapture data consisting of multiple non-invasive marks. Using Bayesian analysis methods, Markov chain Monte Carlo (MCMC) is used to draw samples from the joint posterior distribution.
multimarkCJS( Enc.Mat, data.type = "never", covs = data.frame(), mms = NULL, mod.p = ~1, mod.phi = ~1, mod.delta = ~type, parms = c("pbeta", "phibeta", "delta"), nchains = 1, iter = 12000, adapt = 1000, bin = 50, thin = 1, burnin = 2000, taccept = 0.44, tuneadjust = 0.95, proppbeta = 0.1, propzp = 1, propsigmap = 1, propphibeta = 0.1, propzphi = 1, propsigmaphi = 1, maxnumbasis = 1, pbeta0 = 0, pSigma0 = 1, phibeta0 = 0, phiSigma0 = 1, l0p = 1, d0p = 0.01, l0phi = 1, d0phi = 0.01, a0delta = 1, a0alpha = 1, b0alpha = 1, a0psi = 1, b0psi = 1, initial.values = NULL, known = integer(), link = "probit", printlog = FALSE, ... )
multimarkCJS( Enc.Mat, data.type = "never", covs = data.frame(), mms = NULL, mod.p = ~1, mod.phi = ~1, mod.delta = ~type, parms = c("pbeta", "phibeta", "delta"), nchains = 1, iter = 12000, adapt = 1000, bin = 50, thin = 1, burnin = 2000, taccept = 0.44, tuneadjust = 0.95, proppbeta = 0.1, propzp = 1, propsigmap = 1, propphibeta = 0.1, propzphi = 1, propsigmaphi = 1, maxnumbasis = 1, pbeta0 = 0, pSigma0 = 1, phibeta0 = 0, phiSigma0 = 1, l0p = 1, d0p = 0.01, l0phi = 1, d0phi = 0.01, a0delta = 1, a0alpha = 1, b0alpha = 1, a0psi = 1, b0psi = 1, initial.values = NULL, known = integer(), link = "probit", printlog = FALSE, ... )
Enc.Mat |
A matrix of observed encounter histories with rows corresponding to individuals and columns corresponding to sampling occasions (ignored unless |
data.type |
Specifies the encounter history data type. All data types include non-detections (type 0 encounter), type 1 encounter (e.g., left-side), and type 2 encounters (e.g., right-side). When both type 1 and type 2 encounters occur for the same individual within a sampling occasion, these can either be "non-simultaneous" (type 3 encounter) or "simultaneous" (type 4 encounter). Three data types are currently permitted:
|
covs |
A data frame of temporal covariates for detection probabilities (ignored unless |
mms |
An optional object of class |
mod.p |
Model formula for detection probability ( |
mod.phi |
Model formula for survival probability ( |
mod.delta |
Model formula for conditional probabilities of type 1 (delta_1) and type 2 (delta_2) encounters, given detection. Currently only |
parms |
A character vector giving the names of the parameters and latent variables to monitor. Possible parameters are probit-scale detection probability parameters (" |
nchains |
The number of parallel MCMC chains for the model. |
iter |
The number of MCMC iterations. |
adapt |
Ignored; no adaptive phase is needed for "probit" link. |
bin |
Ignored; no adaptive phase is needed for "probit" link. |
thin |
Thinning interval for monitored parameters. |
burnin |
Number of burn-in iterations ( |
taccept |
Ignored; no adaptive phase is needed for "probit" link. |
tuneadjust |
Ignored; no adaptive phase is needed for "probit" link. |
proppbeta |
Ignored; no adaptive phase is needed for "probit" link. |
propzp |
Ignored; no adaptive phase is needed for "probit" link. |
propsigmap |
Ignored; no adaptive phase is needed for "probit" link. |
propphibeta |
Ignored; no adaptive phase is needed for "probit" link. |
propzphi |
Ignored; no adaptive phase is needed for "probit" link. |
propsigmaphi |
Ignored; no adaptive phase is needed for "probit" link. |
maxnumbasis |
Maximum number of basis vectors to use when proposing latent history frequency updates. Default is |
pbeta0 |
Scaler or vector (of length k) specifying mean of pbeta ~ multivariateNormal(pbeta0, pSigma0) prior. If |
pSigma0 |
Scaler or k x k matrix specifying covariance matrix of pbeta ~ multivariateNormal(pbeta0, pSigma0) prior. If |
phibeta0 |
Scaler or vector (of length k) specifying mean of phibeta ~ multivariateNormal(phibeta0, phiSigma0) prior. If |
phiSigma0 |
Scaler or k x k matrix specifying covariance matrix of phibeta ~ multivariateNormal(phibeta0, phiSigma0) prior. If |
l0p |
Specifies "shape" parameter for [sigma2_zp] ~ invGamma(l0p,d0p) prior. Default is |
d0p |
Specifies "scale" parameter for [sigma2_zp] ~ invGamma(l0p,d0p) prior. Default is |
l0phi |
Specifies "shape" parameter for [sigma2_zphi] ~ invGamma(l0phi,d0phi) prior. Default is |
d0phi |
Specifies "scale" parameter for [sigma2_zphi] ~ invGamma(l0phi,d0phi) prior. Default is |
a0delta |
Scaler or vector (of length d) specifying the prior for the conditional (on detection) probability of type 1 (delta_1), type 2 (delta_2), and both type 1 and type 2 encounters (1-delta_1-delta_2). If |
a0alpha |
Specifies "shape1" parameter for [alpha] ~ Beta(a0alpha, b0alpha) prior. Only applicable when |
b0alpha |
Specifies "shape2" parameter for [alpha] ~ Beta(a0alpha, b0alpha) prior. Only applicable when |
a0psi |
Specifies "shape1" parameter for [psi] ~ Beta(a0psi,b0psi) prior. Default is |
b0psi |
Specifies "shape2" parameter for [psi] ~ Beta(a0psi,b0psi) prior. Default is |
initial.values |
Optional list of |
known |
Optional integer vector indicating whether the encounter history of an individual is known with certainty (i.e., the observed encounter history is the true encounter history). Encounter histories with at least one type 4 encounter are automatically assumed to be known, and |
link |
Link function for survival and capture probabilities. Only probit link is currently implemented. |
printlog |
Logical indicating whether to print the progress of chains and any errors to a log file in the working directory. Ignored when |
... |
Additional " |
The first time multimarkCJS
(or multimarkClosed
) is called, it will likely produce a firewall warning alerting users that R has requested the ability to accept incoming network connections. Incoming network connections are required to use parallel processing as implemented in multimarkCJS
. Note that setting parms="all"
is required for any multimarkCJS
model output to be used in multimodelCJS
.
A list containing the following:
mcmc |
Markov chain Monte Carlo object of class |
mod.p |
Model formula for detection probability (as specified by |
mod.phi |
Model formula for survival probability (as specified by |
mod.delta |
Formula always |
DM |
A list of design matrices for detection and survival probability respectively generated by |
initial.values |
A list containing the parameter and latent variable values at iteration |
mms |
An object of class |
Brett T. McClintock
Bonner, S. J., and Holmberg J. 2013. Mark-recapture with multiple, non-invasive marks. Biometrics 69: 766-775.
McClintock, B. T., Conn, P. B., Alonso, R. S., and Crooks, K. R. 2013. Integrated modeling of bilateral photo-identification data in mark-recapture analyses. Ecology 94: 1464-1471.
McClintock, B. T., Bailey, L. L., Dreher, B. P., and Link, W. A. 2014. Probit models for capture-recapture data subject to imperfect detection, individual heterogeneity and misidentification. The Annals of Applied Statistics 8: 2461-2484.
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Simulate open population data using defaults data <- simdataCJS() #Fit default open population model sim.dot <- multimarkCJS(data$Enc.Mat) #Posterior summary for monitored parameters summary(sim.dot$mcmc) plot(sim.dot$mcmc) #' #Fit ``age'' model with 2 age classes (e.g., juvenile and adult) for survival #using 'parameters' and 'right' arguments from RMark::make.design.data sim.age <- multimarkCJS(data$Enc.Mat,mod.phi=~age, parameters=list(Phi=list(age.bins=c(0,1,4))),right=FALSE) summary(getprobsCJS(sim.age))
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Simulate open population data using defaults data <- simdataCJS() #Fit default open population model sim.dot <- multimarkCJS(data$Enc.Mat) #Posterior summary for monitored parameters summary(sim.dot$mcmc) plot(sim.dot$mcmc) #' #Fit ``age'' model with 2 age classes (e.g., juvenile and adult) for survival #using 'parameters' and 'right' arguments from RMark::make.design.data sim.age <- multimarkCJS(data$Enc.Mat,mod.phi=~age, parameters=list(Phi=list(age.bins=c(0,1,4))),right=FALSE) summary(getprobsCJS(sim.age))
This function fits closed population abundance models for capture-mark-recapture data consisting of multiple non-invasive marks using Bayesian analysis methods. Markov chain Monte Carlo (MCMC) is used to draw samples from the joint posterior distribution.
multimarkClosed( Enc.Mat, data.type = "never", covs = data.frame(), mms = NULL, mod.p = ~1, mod.delta = ~type, parms = c("pbeta", "delta", "N"), nchains = 1, iter = 12000, adapt = 1000, bin = 50, thin = 1, burnin = 2000, taccept = 0.44, tuneadjust = 0.95, proppbeta = 0.1, propzp = 1, propsigmap = 1, npoints = 500, maxnumbasis = 1, a0delta = 1, a0alpha = 1, b0alpha = 1, a = 25, mu0 = 0, sigma2_mu0 = 1.75, a0psi = 1, b0psi = 1, initial.values = NULL, known = integer(), printlog = FALSE, ... )
multimarkClosed( Enc.Mat, data.type = "never", covs = data.frame(), mms = NULL, mod.p = ~1, mod.delta = ~type, parms = c("pbeta", "delta", "N"), nchains = 1, iter = 12000, adapt = 1000, bin = 50, thin = 1, burnin = 2000, taccept = 0.44, tuneadjust = 0.95, proppbeta = 0.1, propzp = 1, propsigmap = 1, npoints = 500, maxnumbasis = 1, a0delta = 1, a0alpha = 1, b0alpha = 1, a = 25, mu0 = 0, sigma2_mu0 = 1.75, a0psi = 1, b0psi = 1, initial.values = NULL, known = integer(), printlog = FALSE, ... )
Enc.Mat |
A matrix of observed encounter histories with rows corresponding to individuals and columns corresponding to sampling occasions (ignored unless |
data.type |
Specifies the encounter history data type. All data types include non-detections (type 0 encounter), type 1 encounter (e.g., left-side), and type 2 encounters (e.g., right-side). When both type 1 and type 2 encounters occur for the same individual within a sampling occasion, these can either be "non-simultaneous" (type 3 encounter) or "simultaneous" (type 4 encounter). Three data types are currently permitted:
|
covs |
A data frame of temporal covariates for detection probabilities (ignored unless |
mms |
An optional object of class |
mod.p |
Model formula for detection probability. For example, |
mod.delta |
Model formula for conditional probabilities of type 1 (delta_1) and type 2 (delta_2) encounters, given detection. Currently only |
parms |
A character vector giving the names of the parameters and latent variables to monitor. Possible parameters are logit-scale detection probability parameters (" |
nchains |
The number of parallel MCMC chains for the model. |
iter |
The number of MCMC iterations. |
adapt |
The number of iterations for proposal distribution adaptation. If |
bin |
Bin length for calculating acceptance rates during adaptive phase ( |
thin |
Thinning interval for monitored parameters. |
burnin |
Number of burn-in iterations ( |
taccept |
Target acceptance rate during adaptive phase ( |
tuneadjust |
Adjustment term during adaptive phase ( |
proppbeta |
Scaler or vector (of length k) specifying the initial standard deviation of the Normal(pbeta[j], proppbeta[j]) proposal distribution. If |
propzp |
Scaler or vector (of length M) specifying the initial standard deviation of the Normal(zp[i], propzp[i]) proposal distribution. If |
propsigmap |
Scaler specifying the initial Gamma(shape = 1/ |
npoints |
Number of Gauss-Hermite quadrature points to use for numerical integration. Accuracy increases with number of points, but so does computation time. |
maxnumbasis |
Maximum number of basis vectors to use when proposing latent history frequency updates. Default is |
a0delta |
Scaler or vector (of length d) specifying the prior for the conditional (on detection) probability of type 1 (delta_1), type 2 (delta_2), and both type 1 and type 2 encounters (1-delta_1-delta_2). If |
a0alpha |
Specifies "shape1" parameter for [alpha] ~ Beta(a0alpha, b0alpha) prior. Only applicable when |
b0alpha |
Specifies "shape2" parameter for [alpha] ~ Beta(a0alpha, b0alpha) prior. Only applicable when |
a |
Scale parameter for [sigma_z] ~ half-Cauchy(a) prior for the individual hetegeneity term sigma_zp = sqrt(sigma2_zp). Default is “uninformative” |
mu0 |
Scaler or vector (of length k) specifying mean of pbeta[j] ~ Normal(mu0[j], sigma2_mu0[j]) prior. If |
sigma2_mu0 |
Scaler or vector (of length k) specifying variance of pbeta[j] ~ Normal(mu0[j], sigma2_mu0[j]) prior. If |
a0psi |
Specifies "shape1" parameter for [psi] ~ Beta(a0psi,b0psi) prior. Default is |
b0psi |
Specifies "shape2" parameter for [psi] ~ Beta(a0psi,b0psi) prior. Default is |
initial.values |
Optional list of |
known |
Optional integer vector indicating whether the encounter history of an individual is known with certainty (i.e., the observed encounter history is the true encounter history). Encounter histories with at least one type 4 encounter are automatically assumed to be known, and |
printlog |
Logical indicating whether to print the progress of chains and any errors to a log file in the working directory. Ignored when |
... |
Additional " |
The first time multimarkClosed
(or multimarkCJS
) is called, it will likely produce a firewall warning alerting users that R has requested the ability to accept incoming network connections. Incoming network connections are required to use parallel processing as implemented in multimarkClosed
. Note that setting parms="all"
is required for any multimarkClosed
model output to be used in multimodelClosed
.
A list containing the following:
mcmc |
Markov chain Monte Carlo object of class |
mod.p |
Model formula for detection probability (as specified by |
mod.delta |
Model formula for conditional probability of type 1 or type 2 encounter, given detection (as specified by |
DM |
A list of design matrices for detection probability generated for model |
initial.values |
A list containing the parameter and latent variable values at iteration |
mms |
An object of class |
Brett T. McClintock
Bonner, S. J., and Holmberg J. 2013. Mark-recapture with multiple, non-invasive marks. Biometrics 69: 766-775.
McClintock, B. T., Conn, P. B., Alonso, R. S., and Crooks, K. R. 2013. Integrated modeling of bilateral photo-identification data in mark-recapture analyses. Ecology 94: 1464-1471.
McClintock, B. T., Bailey, L. L., Dreher, B. P., and Link, W. A. 2014. Probit models for capture-recapture data subject to imperfect detection, individual heterogeneity and misidentification. The Annals of Applied Statistics 8: 2461-2484.
bobcat
, processdata
, multimodelClosed
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Run single chain using the default model for bobcat data bobcat.dot<-multimarkClosed(bobcat) #Posterior summary for monitored parameters summary(bobcat.dot$mcmc) plot(bobcat.dot$mcmc)
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Run single chain using the default model for bobcat data bobcat.dot<-multimarkClosed(bobcat) #Posterior summary for monitored parameters summary(bobcat.dot$mcmc) plot(bobcat.dot$mcmc)
This function fits spatially-explicit population abundance models for capture-mark-recapture data consisting of multiple non-invasive marks using Bayesian analysis methods. Markov chain Monte Carlo (MCMC) is used to draw samples from the joint posterior distribution.
multimarkClosedSCR( Enc.Mat, trapCoords, studyArea = NULL, buffer = NULL, ncells = 1024, data.type = "never", covs = data.frame(), mms = NULL, mod.p = ~1, mod.delta = ~type, detection = "half-normal", parms = c("pbeta", "delta", "N"), nchains = 1, iter = 12000, adapt = 1000, bin = 50, thin = 1, burnin = 2000, taccept = 0.44, tuneadjust = 0.95, proppbeta = 0.1, propsigma = 1, propcenter = NULL, maxnumbasis = 1, a0delta = 1, a0alpha = 1, b0alpha = 1, sigma_bounds = NULL, mu0 = 0, sigma2_mu0 = 1.75, a0psi = 1, b0psi = 1, initial.values = NULL, known = integer(), scalemax = 10, printlog = FALSE, ... )
multimarkClosedSCR( Enc.Mat, trapCoords, studyArea = NULL, buffer = NULL, ncells = 1024, data.type = "never", covs = data.frame(), mms = NULL, mod.p = ~1, mod.delta = ~type, detection = "half-normal", parms = c("pbeta", "delta", "N"), nchains = 1, iter = 12000, adapt = 1000, bin = 50, thin = 1, burnin = 2000, taccept = 0.44, tuneadjust = 0.95, proppbeta = 0.1, propsigma = 1, propcenter = NULL, maxnumbasis = 1, a0delta = 1, a0alpha = 1, b0alpha = 1, sigma_bounds = NULL, mu0 = 0, sigma2_mu0 = 1.75, a0psi = 1, b0psi = 1, initial.values = NULL, known = integer(), scalemax = 10, printlog = FALSE, ... )
Enc.Mat |
A matrix containing the observed encounter histories with rows corresponding to individuals and ( |
trapCoords |
A matrix of dimension |
studyArea |
is a 3-column matrix containing the coordinates for the centroids of a contiguous grid of cells that define the study area and available habitat. Each row corresponds to a grid cell. The first 2 columns (“x” and “y”) indicate the Cartesian x- and y-coordinate for the centroid of each grid cell, and the third column (“avail”) indicates whether the cell is available habitat (=1) or not (=0). All cells must be square and have the same resolution. If |
buffer |
A scaler in same units as |
ncells |
The number of grid cells in the study area when |
data.type |
Specifies the encounter history data type. All data types include non-detections (type 0 encounter), type 1 encounter (e.g., left-side), and type 2 encounters (e.g., right-side). When both type 1 and type 2 encounters occur for the same individual within a sampling occasion, these can either be "non-simultaneous" (type 3 encounter) or "simultaneous" (type 4 encounter). Three data types are currently permitted:
|
covs |
A data frame of time- and/or trap-dependent covariates for detection probabilities (ignored unless |
mms |
An optional object of class |
mod.p |
Model formula for detection probability as a function of distance from activity centers. For example, |
mod.delta |
Model formula for conditional probabilities of type 1 (delta_1) and type 2 (delta_2) encounters, given detection. Currently only |
detection |
Model for detection probability as a function of distance from activity centers . Must be " |
parms |
A character vector giving the names of the parameters and latent variables to monitor. Possible parameters are cloglog-scale detection probability parameters (" |
nchains |
The number of parallel MCMC chains for the model. |
iter |
The number of MCMC iterations. |
adapt |
The number of iterations for proposal distribution adaptation. If |
bin |
Bin length for calculating acceptance rates during adaptive phase ( |
thin |
Thinning interval for monitored parameters. |
burnin |
Number of burn-in iterations ( |
taccept |
Target acceptance rate during adaptive phase ( |
tuneadjust |
Adjustment term during adaptive phase ( |
proppbeta |
Scaler or vector (of length k) specifying the initial standard deviation of the Normal(pbeta[j], proppbeta[j]) proposal distribution. If |
propsigma |
Scaler specifying the initial Gamma(shape = 1/ |
propcenter |
Scaler specifying the neighborhood distance when proposing updates to activity centers. When |
maxnumbasis |
Maximum number of basis vectors to use when proposing latent history frequency updates. Default is |
a0delta |
Scaler or vector (of length d) specifying the prior for the conditional (on detection) probability of type 1 (delta_1), type 2 (delta_2), and both type 1 and type 2 encounters (1-delta_1-delta_2). If |
a0alpha |
Specifies "shape1" parameter for [alpha] ~ Beta(a0alpha, b0alpha) prior. Only applicable when |
b0alpha |
Specifies "shape2" parameter for [alpha] ~ Beta(a0alpha, b0alpha) prior. Only applicable when |
sigma_bounds |
Positive vector of length 2 for the lower and upper bounds for the [sigma_scr] ~ Uniform(sigma_bounds[1], sigma_bounds[2]) (or [sqrt(lambda)] when |
mu0 |
Scaler or vector (of length k) specifying mean of pbeta[j] ~ Normal(mu0[j], sigma2_mu0[j]) prior. If |
sigma2_mu0 |
Scaler or vector (of length k) specifying variance of pbeta[j] ~ Normal(mu0[j], sigma2_mu0[j]) prior. If |
a0psi |
Specifies "shape1" parameter for [psi] ~ Beta(a0psi,b0psi) prior. Default is |
b0psi |
Specifies "shape2" parameter for [psi] ~ Beta(a0psi,b0psi) prior. Default is |
initial.values |
Optional list of |
known |
Optional integer vector indicating whether the encounter history of an individual is known with certainty (i.e., the observed encounter history is the true encounter history). Encounter histories with at least one type 4 encounter are automatically assumed to be known, and |
scalemax |
Upper bound for internal re-scaling of grid cell centroid coordinates. Default is |
printlog |
Logical indicating whether to print the progress of chains and any errors to a log file in the working directory. Ignored when |
... |
Additional " |
The first time multimarkSCRClosed
is called, it will likely produce a firewall warning alerting users that R has requested the ability to accept incoming network connections. Incoming network connections are required to use parallel processing as implemented in multimarkClosed
. Note that setting parms="all"
is required for any multimarkClosed
model output to be used in multimodelClosed
.
A list containing the following:
mcmc |
Markov chain Monte Carlo object of class |
mod.p |
Model formula for detection probability (as specified by |
mod.delta |
Model formula for conditional probability of type 1 or type 2 encounter, given detection (as specified by |
mod.det |
Model formula for detection function (as specified by |
DM |
A list of design matrices for detection probability generated for model |
initial.values |
A list containing the parameter and latent variable values at iteration |
mms |
An object of class |
Brett T. McClintock
Bonner, S. J., and Holmberg J. 2013. Mark-recapture with multiple, non-invasive marks. Biometrics 69: 766-775.
Gopalaswamy, A.M., Royle, J.A., Hines, J.E., Singh, P., Jathanna, D., Kumar, N. and Karanth, K.U. 2012. Program SPACECAP: software for estimating animal density using spatially explicit capture-recapture models. Methods in Ecology and Evolution 3:1067-1072.
King, R., McClintock, B. T., Kidney, D., and Borchers, D. L. 2016. Capture-recapture abundance estimation using a semi-complete data likelihood approach. The Annals of Applied Statistics 10: 264-285
McClintock, B. T., Conn, P. B., Alonso, R. S., and Crooks, K. R. 2013. Integrated modeling of bilateral photo-identification data in mark-recapture analyses. Ecology 94: 1464-1471.
McClintock, B. T., Bailey, L. L., Dreher, B. P., and Link, W. A. 2014. Probit models for capture-recapture data subject to imperfect detection, individual heterogeneity and misidentification. The Annals of Applied Statistics 8: 2461-2484.
Royle, J.A., Karanth, K.U., Gopalaswamy, A.M. and Kumar, N.S. 2009. Bayesian inference in camera trapping studies for a class of spatial capture-recapture models. Ecology 90: 3233-3244.
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Generate object of class "multimarkSCRsetup" from simulated data sim.data<-simdataClosedSCR() Enc.Mat <- sim.data$Enc.Mat trapCoords <- sim.data$spatialInputs$trapCoords studyArea <- sim.data$spatialInputs$studyArea #Run single chain using the default model for simulated data example.dot<-multimarkClosedSCR(Enc.Mat,trapCoords,studyArea) #Posterior summary for monitored parameters summary(example.dot$mcmc) plot(example.dot$mcmc)
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Generate object of class "multimarkSCRsetup" from simulated data sim.data<-simdataClosedSCR() Enc.Mat <- sim.data$Enc.Mat trapCoords <- sim.data$spatialInputs$trapCoords studyArea <- sim.data$spatialInputs$studyArea #Run single chain using the default model for simulated data example.dot<-multimarkClosedSCR(Enc.Mat,trapCoords,studyArea) #Posterior summary for monitored parameters summary(example.dot$mcmc) plot(example.dot$mcmc)
"multimarkSCRsetup"
A class of spatial 'mulitmark' model inputs
Enc.Mat
Object of class "matrix"
. The observed encounter histories (with rows corresponding to individuals and columns corresponding to sampling occasions).
data.type
Object of class "character"
. The encounter history data type ("never", "sometimes", or "always").
vAll.hists
Object of class "integer"
. An ordered vector containing all possible encounter histories in sequence.
Aprime
Object of class "sparseMatrix"
. Transpose of the A matrix mapping latent encounter histories to observed histories.
indBasis
Object of class "numeric"
.An ordered vector of the indices of the three encounter histories updated by each basis vector.
ncolbasis
Object of class "integer"
. The number of needed basis vectors.
knownx
Object of class "integer"
. Frequencies of known encounter histories.
C
Object of class "integer"
. Sampling occasion of first capture for each encounter history.
L
Object of class "integer"
. Sampling occasion of last capture for each encounter history.
naivex
Object of class "integer"
. “Naive” latent history frequencies assuming a one-to-one mapping with Enc.Mat
.
covs
Object of class "data.frame"
. Temporal covariates for detection probability (the number of rows in the data frame must equal the number of sampling occasions).
spatialInputs
Object of class "list"
. List is of length 4 containing trapCoords
and studyArea
after re-scaling coordinates based on maxscale
, as well as the original (not re-scaled) grid cell resolution (origCellRes
) and re-scaling range (Srange
).
Objects can be created by calls of the form processdata(Enc.Mat, ...)
or new("multimarkSCRsetup", ...)
.
No methods defined with class "multimarkSCRsetup".
Brett T. McClintock
showClass("multimarkSCRsetup")
showClass("multimarkSCRsetup")
"multimarksetup"
A class of 'mulitmark' model inputs
Enc.Mat
Object of class "matrix"
. The observed encounter histories (with rows corresponding to individuals and columns corresponding to sampling occasions).
data.type
Object of class "character"
. The encounter history data type ("never", "sometimes", or "always").
vAll.hists
Object of class "integer"
. An ordered vector containing all possible encounter histories in sequence.
Aprime
Object of class "sparseMatrix"
. Transpose of the A matrix mapping latent encounter histories to observed histories.
indBasis
Object of class "numeric"
.An ordered vector of the indices of the three encounter histories updated by each basis vector.
ncolbasis
Object of class "integer"
. The number of needed basis vectors.
knownx
Object of class "integer"
. Frequencies of known encounter histories.
C
Object of class "integer"
. Sampling occasion of first capture for each encounter history.
L
Object of class "integer"
. Sampling occasion of last capture for each encounter history.
naivex
Object of class "integer"
. “Naive” latent history frequencies assuming a one-to-one mapping with Enc.Mat
.
covs
Object of class "data.frame"
. Temporal covariates for detection probability (the number of rows in the data frame must equal the number of sampling occasions).
Objects can be created by calls of the form processdata(Enc.Mat, ...)
or new("multimarksetup", ...)
.
No methods defined with class "multimarksetup".
Brett T. McClintock
showClass("multimarksetup")
showClass("multimarksetup")
This function performs Bayesian multimodel inference for a set of 'multimark' open population survival (i.e., Cormack-Jolly-Seber) models using the reversible jump Markov chain Monte Carlo (RJMCMC) algorithm proposed by Barker & Link (2013).
multimodelCJS( modlist, modprior = rep(1/length(modlist), length(modlist)), monparms = "phi", miter = NULL, mburnin = 0, mthin = 1, M1 = NULL, pbetapropsd = 1, zppropsd = NULL, phibetapropsd = 1, zphipropsd = NULL, sigppropshape = 1, sigppropscale = 0.01, sigphipropshape = 1, sigphipropscale = 0.01, printlog = FALSE )
multimodelCJS( modlist, modprior = rep(1/length(modlist), length(modlist)), monparms = "phi", miter = NULL, mburnin = 0, mthin = 1, M1 = NULL, pbetapropsd = 1, zppropsd = NULL, phibetapropsd = 1, zphipropsd = NULL, sigppropshape = 1, sigppropscale = 0.01, sigphipropshape = 1, sigphipropscale = 0.01, printlog = FALSE )
modlist |
A list of individual model output lists returned by |
modprior |
Vector of length |
monparms |
Parameters to monitor. Only parameters common to all models can be monitored (e.g., " |
miter |
The number of RJMCMC iterations per chain. If |
mburnin |
Number of burn-in iterations ( |
mthin |
Thinning interval for monitored parameters. |
M1 |
Integer vector indicating the initial model for each chain, where |
pbetapropsd |
Scaler specifying the standard deviation of the Normal(0, pbetapropsd) proposal distribution for " |
zppropsd |
Scaler specifying the standard deviation of the Normal(0, zppropsd) proposal distribution for " |
phibetapropsd |
Scaler specifying the standard deviation of the Normal(0, phibetapropsd) proposal distribution for " |
zphipropsd |
Scaler specifying the standard deviation of the Normal(0, zphipropsd) proposal distribution for " |
sigppropshape |
Scaler specifying the shape parameter of the invGamma(shape = sigppropshape, scale = sigppropscale) proposal distribution for " |
sigppropscale |
Scaler specifying the scale parameter of the invGamma(shape = sigppropshape, scale = sigppropscale) proposal distribution for " |
sigphipropshape |
Scaler specifying the shape parameter of the invGamma(shape = sigphipropshape, scale = sigphipropscale) proposal distribution for " |
sigphipropscale |
Scaler specifying the scale parameter of the invGamma(shape = sigphipropshape, scale = sigphipropscale) proposal distribution for " |
printlog |
Logical indicating whether to print the progress of chains and any errors to a log file in the working directory. Ignored when |
Note that setting parms="all"
is required when fitting individual multimarkCJS
models to be included in modlist
.
A list containing the following:
rjmcmc |
Reversible jump Markov chain Monte Carlo object of class |
pos.prob |
A list of calculated posterior model probabilities for each chain, including the overall posterior model probabilities across all chains. |
Brett T. McClintock
Barker, R. J. and Link. W. A. 2013. Bayesian multimodel inference by RJMCMC: a Gibbs sampling approach. The American Statistician 67: 150-156.
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Generate object of class "multimarksetup" from simulated data data_type = "always" noccas <- 5 phibetaTime <- seq(2,0,length=noccas-1) # declining trend in survival data <- simdataCJS(noccas=5,phibeta=phibetaTime,data.type=data_type) setup <- processdata(data$Enc.Mat,data.type=data_type) #Run single chain using the default model. Note parms="all". sim.pdot.phidot <- multimarkCJS(mms=setup,parms="all",iter=1000,adapt=500,burnin=500) #Run single chain with temporal trend for phi. Note parms="all". sim.pdot.phiTime <- multimarkCJS(mms=setup,mod.phi=~Time,parms="all",iter=1000,adapt=500,burnin=500) #Perform RJMCMC using defaults modlist <- list(mod1=sim.pdot.phidot,mod2=sim.pdot.phiTime) sim.M <- multimodelCJS(modlist=modlist) #Posterior model probabilities sim.M$pos.prob #multimodel posterior summary for survival (display first cohort only) summary(sim.M$rjmcmc[,paste0("phi[1,",1:(noccas-1),"]")])
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Generate object of class "multimarksetup" from simulated data data_type = "always" noccas <- 5 phibetaTime <- seq(2,0,length=noccas-1) # declining trend in survival data <- simdataCJS(noccas=5,phibeta=phibetaTime,data.type=data_type) setup <- processdata(data$Enc.Mat,data.type=data_type) #Run single chain using the default model. Note parms="all". sim.pdot.phidot <- multimarkCJS(mms=setup,parms="all",iter=1000,adapt=500,burnin=500) #Run single chain with temporal trend for phi. Note parms="all". sim.pdot.phiTime <- multimarkCJS(mms=setup,mod.phi=~Time,parms="all",iter=1000,adapt=500,burnin=500) #Perform RJMCMC using defaults modlist <- list(mod1=sim.pdot.phidot,mod2=sim.pdot.phiTime) sim.M <- multimodelCJS(modlist=modlist) #Posterior model probabilities sim.M$pos.prob #multimodel posterior summary for survival (display first cohort only) summary(sim.M$rjmcmc[,paste0("phi[1,",1:(noccas-1),"]")])
This function performs Bayesian multimodel inference for a set of 'multimark' closed population abundance models using the reversible jump Markov chain Monte Carlo (RJMCMC) algorithm proposed by Barker & Link (2013).
multimodelClosed( modlist, modprior = rep(1/length(modlist), length(modlist)), monparms = "N", miter = NULL, mburnin = 0, mthin = 1, M1 = NULL, pbetapropsd = 1, zppropsd = NULL, sigppropshape = 6, sigppropscale = 4, printlog = FALSE )
multimodelClosed( modlist, modprior = rep(1/length(modlist), length(modlist)), monparms = "N", miter = NULL, mburnin = 0, mthin = 1, M1 = NULL, pbetapropsd = 1, zppropsd = NULL, sigppropshape = 6, sigppropscale = 4, printlog = FALSE )
modlist |
A list of individual model output lists returned by |
modprior |
Vector of length |
monparms |
Parameters to monitor. Only parameters common to all models can be monitored (e.g., " |
miter |
The number of RJMCMC iterations per chain. If |
mburnin |
Number of burn-in iterations ( |
mthin |
Thinning interval for monitored parameters. |
M1 |
Integer vector indicating the initial model for each chain, where |
pbetapropsd |
Scaler specifying the standard deviation of the Normal(0, pbetapropsd) proposal distribution for " |
zppropsd |
Scaler specifying the standard deviation of the Normal(0, zppropsd) proposal distribution for " |
sigppropshape |
Scaler specifying the shape parameter of the invGamma(shape = sigppropshape, scale = sigppropscale) proposal distribution for |
sigppropscale |
Scaler specifying the scale parameter of the invGamma(shape = sigppropshape, scale = sigppropscale) proposal distribution for |
printlog |
Logical indicating whether to print the progress of chains and any errors to a log file in the working directory. Ignored when |
Note that setting parms="all"
is required when fitting individual multimarkClosed
or markClosed
models to be included in modlist
.
A list containing the following:
rjmcmc |
Reversible jump Markov chain Monte Carlo object of class |
pos.prob |
A list of calculated posterior model probabilities for each chain, including the overall posterior model probabilities across all chains. |
Brett T. McClintock
Barker, R. J. and Link. W. A. 2013. Bayesian multimodel inference by RJMCMC: a Gibbs sampling approach. The American Statistician 67: 150-156.
multimarkClosed
, markClosed
, processdata
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Generate object of class "multimarksetup" setup <- processdata(bobcat) #Run single chain using the default model for bobcat data. Note parms="all". bobcat.dot <- multimarkClosed(mms=setup,parms="all",iter=1000,adapt=500,burnin=500) #Run single chain for bobcat data with time effects. Note parms="all". bobcat.time <- multimarkClosed(mms=setup,mod.p=~time,parms="all",iter=1000,adapt=500,burnin=500) #Perform RJMCMC using defaults modlist <- list(mod1=bobcat.dot,mod2=bobcat.time) bobcat.M <- multimodelClosed(modlist=modlist,monparms=c("N","p")) #Posterior model probabilities bobcat.M$pos.prob #multimodel posterior summary for abundance summary(bobcat.M$rjmcmc[,"N"])
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Generate object of class "multimarksetup" setup <- processdata(bobcat) #Run single chain using the default model for bobcat data. Note parms="all". bobcat.dot <- multimarkClosed(mms=setup,parms="all",iter=1000,adapt=500,burnin=500) #Run single chain for bobcat data with time effects. Note parms="all". bobcat.time <- multimarkClosed(mms=setup,mod.p=~time,parms="all",iter=1000,adapt=500,burnin=500) #Perform RJMCMC using defaults modlist <- list(mod1=bobcat.dot,mod2=bobcat.time) bobcat.M <- multimodelClosed(modlist=modlist,monparms=c("N","p")) #Posterior model probabilities bobcat.M$pos.prob #multimodel posterior summary for abundance summary(bobcat.M$rjmcmc[,"N"])
This function performs Bayesian multimodel inference for a set of 'multimark' spatial population abundance models using the reversible jump Markov chain Monte Carlo (RJMCMC) algorithm proposed by Barker & Link (2013).
multimodelClosedSCR( modlist, modprior = rep(1/length(modlist), length(modlist)), monparms = "N", miter = NULL, mburnin = 0, mthin = 1, M1 = NULL, pbetapropsd = 1, sigpropmean = 0.8, sigpropsd = 0.4, printlog = FALSE )
multimodelClosedSCR( modlist, modprior = rep(1/length(modlist), length(modlist)), monparms = "N", miter = NULL, mburnin = 0, mthin = 1, M1 = NULL, pbetapropsd = 1, sigpropmean = 0.8, sigpropsd = 0.4, printlog = FALSE )
modlist |
A list of individual model output lists returned by |
modprior |
Vector of length |
monparms |
Parameters to monitor. Only parameters common to all models can be monitored (e.g., " |
miter |
The number of RJMCMC iterations per chain. If |
mburnin |
Number of burn-in iterations ( |
mthin |
Thinning interval for monitored parameters. |
M1 |
Integer vector indicating the initial model for each chain, where |
pbetapropsd |
Scaler specifying the standard deviation of the Normal(0, pbetapropsd) proposal distribution for " |
sigpropmean |
Scaler specifying the mean of the inverse Gamma proposal distribution for |
sigpropsd |
Scaler specifying the standard deviation of the inverse Gamma proposal distribution for |
printlog |
Logical indicating whether to print the progress of chains and any errors to a log file in the working directory. Ignored when |
Note that setting parms="all"
is required when fitting individual multimarkClosedSCR
or markClosedSCR
models to be included in modlist
.
A list containing the following:
rjmcmc |
Reversible jump Markov chain Monte Carlo object of class |
pos.prob |
A list of calculated posterior model probabilities for each chain, including the overall posterior model probabilities across all chains. |
Brett T. McClintock
Barker, R. J. and Link. W. A. 2013. Bayesian multimodel inference by RJMCMC: a Gibbs sampling approach. The American Statistician 67: 150-156.
multimarkClosedSCR
, processdataSCR
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Generate object of class "multimarkSCRsetup" sim.data<-simdataClosedSCR() Enc.Mat<-sim.data$Enc.Mat trapCoords<-sim.data$spatialInputs$trapCoords studyArea<-sim.data$spatialInputs$studyArea setup<-processdataSCR(Enc.Mat,trapCoords,studyArea) #Run single chain using the default model for simulated data. Note parms="all". example.dot <- multimarkClosedSCR(mms=setup,parms="all",iter=1000,adapt=500,burnin=500) #Run single chain for simulated data with behavior effects. Note parms="all". example.c <- multimarkClosedSCR(mms=setup,mod.p=~c,parms="all",iter=1000,adapt=500,burnin=500) #Perform RJMCMC using defaults modlist <- list(mod1=example.dot,mod2=example.c) example.M <- multimodelClosedSCR(modlist=modlist,monparms=c("N","D","sigma2_scr")) #Posterior model probabilities example.M$pos.prob #multimodel posterior summary for abundance and density summary(example.M$rjmcmc[,c("N","D")])
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Generate object of class "multimarkSCRsetup" sim.data<-simdataClosedSCR() Enc.Mat<-sim.data$Enc.Mat trapCoords<-sim.data$spatialInputs$trapCoords studyArea<-sim.data$spatialInputs$studyArea setup<-processdataSCR(Enc.Mat,trapCoords,studyArea) #Run single chain using the default model for simulated data. Note parms="all". example.dot <- multimarkClosedSCR(mms=setup,parms="all",iter=1000,adapt=500,burnin=500) #Run single chain for simulated data with behavior effects. Note parms="all". example.c <- multimarkClosedSCR(mms=setup,mod.p=~c,parms="all",iter=1000,adapt=500,burnin=500) #Perform RJMCMC using defaults modlist <- list(mod1=example.dot,mod2=example.c) example.M <- multimodelClosedSCR(modlist=modlist,monparms=c("N","D","sigma2_scr")) #Posterior model probabilities example.M$pos.prob #multimodel posterior summary for abundance and density summary(example.M$rjmcmc[,c("N","D")])
This function plots the study area grid, available habitat, and trap coordinates for spatial capture-recapture studies. Activity centers and capture locations can also be plotted.
plotSpatialData( mms = NULL, trapCoords, studyArea, centers = NULL, trapLines = FALSE )
plotSpatialData( mms = NULL, trapCoords, studyArea, centers = NULL, trapLines = FALSE )
mms |
An optional object of class |
trapCoords |
A matrix of dimension |
studyArea |
A 3-column matrix defining the study area and available habitat. Each row corresponds to a grid cell. The first 2 columns indicate the Cartesian x- and y-coordinate for the centroid of each grid cell, and the third column indicates whether the cell is available habitat (=1) or not (=0). All cells must have the same resolution. Ignored unless |
centers |
An optional vector indicating the grid cell (i.e., the row of |
trapLines |
Logical indicating whether to draw lines from activity centers to respective traps at which each individual was captured. Default is |
Brett T. McClintock
#Plot the tiger example data plotSpatialData(trapCoords=tiger$trapCoords,studyArea=tiger$studyArea)
#Plot the tiger example data plotSpatialData(trapCoords=tiger$trapCoords,studyArea=tiger$studyArea)
This function generates an object of class multimarksetup
that is required to fit ‘multimark’ models.
processdata( Enc.Mat, data.type = "never", covs = data.frame(), known = integer() )
processdata( Enc.Mat, data.type = "never", covs = data.frame(), known = integer() )
Enc.Mat |
A matrix of observed encounter histories with rows corresponding to individuals and columns corresponding to sampling occasions (ignored unless |
data.type |
Specifies the encounter history data type. All data types include non-detections (type 0 encounter), type 1 encounter (e.g., left-side), and type 2 encounters (e.g., right-side). When both type 1 and type 2 encounters occur for the same individual within a sampling occasion, these can either be "non-simultaneous" (type 3 encounter) or "simultaneous" (type 4 encounter). Three data types are currently permitted:
|
covs |
A data frame of temporal covariates for detection probabilities (ignored unless |
known |
Optional integer vector indicating whether the encounter history of an individual is known with certainty (i.e., the observed encounter history is the true encounter history). Encounter histories with at least one type 4 encounter are automatically assumed to be known, and |
An object of class multimarksetup
.
Brett T. McClintock
Bonner, S. J., and Holmberg J. 2013. Mark-recapture with multiple, non-invasive marks. Biometrics 69: 766-775.
McClintock, B. T., Conn, P. B., Alonso, R. S., and Crooks, K. R. 2013. Integrated modeling of bilateral photo-identification data in mark-recapture analyses. Ecology 94: 1464-1471.
multimarksetup-class
, multimarkClosed
, bobcat
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Generate object of class "multimarksetup" setup <- processdata(bobcat) #Run single chain using the default model for bobcat data bobcat.dot<-multimarkClosed(mms=setup) #Run single chain for bobcat data with temporal effects (i.e., mod.p=~time) bobcat.time <- multimarkClosed(mms=setup,mod.p=~time)
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Generate object of class "multimarksetup" setup <- processdata(bobcat) #Run single chain using the default model for bobcat data bobcat.dot<-multimarkClosed(mms=setup) #Run single chain for bobcat data with temporal effects (i.e., mod.p=~time) bobcat.time <- multimarkClosed(mms=setup,mod.p=~time)
This function generates an object of class multimarkSCRsetup
that is required to fit spatial ‘multimark’ models.
processdataSCR( Enc.Mat, trapCoords, studyArea = NULL, buffer = NULL, ncells = NULL, data.type = "never", covs = data.frame(), known = integer(), scalemax = 10 )
processdataSCR( Enc.Mat, trapCoords, studyArea = NULL, buffer = NULL, ncells = NULL, data.type = "never", covs = data.frame(), known = integer(), scalemax = 10 )
Enc.Mat |
A matrix containing the observed encounter histories with rows corresponding to individuals and ( |
trapCoords |
A matrix of dimension |
studyArea |
is a 3-column matrix containing the coordinates for the centroids of a contiguous grid of cells that define the study area and available habitat. Each row corresponds to a grid cell. The first 2 columns indicate the Cartesian x- and y-coordinate for the centroid of each grid cell, and the third column indicates whether the cell is available habitat (=1) or not (=0). All cells must be square and have the same resolution. If |
buffer |
A scaler in same units as |
ncells |
The number of grid cells in the study area when |
data.type |
Specifies the encounter history data type. All data types include non-detections (type 0 encounter), type 1 encounter (e.g., left-side), and type 2 encounters (e.g., right-side). When both type 1 and type 2 encounters occur for the same individual within a sampling occasion, these can either be "non-simultaneous" (type 3 encounter) or "simultaneous" (type 4 encounter). Three data types are currently permitted:
|
covs |
A data frame of time- and/or trap-dependent covariates for detection probabilities (ignored unless |
known |
Optional integer vector indicating whether the encounter history of an individual is known with certainty (i.e., the observed encounter history is the true encounter history). Encounter histories with at least one type 4 encounter are automatically assumed to be known, and |
scalemax |
Upper bound for internal re-scaling of grid cell centroid coordinates. Default is |
An object of class multimarkSCRsetup
.
Brett T. McClintock
Bonner, S. J., and Holmberg J. 2013. Mark-recapture with multiple, non-invasive marks. Biometrics 69: 766-775.
Gopalaswamy, A.M., Royle, J.A., Hines, J.E., Singh, P., Jathanna, D., Kumar, N. and Karanth, K.U. 2012. Program SPACECAP: software for estimating animal density using spatially explicit capture-recapture models. Methods in Ecology and Evolution 3:1067-1072.
McClintock, B. T., Conn, P. B., Alonso, R. S., and Crooks, K. R. 2013. Integrated modeling of bilateral photo-identification data in mark-recapture analyses. Ecology 94: 1464-1471.
Royle, J.A., Karanth, K.U., Gopalaswamy, A.M. and Kumar, N.S. 2009. Bayesian inference in camera trapping studies for a class of spatial capture-recapture models. Ecology 90: 3233-3244.
multimarkSCRsetup-class
, multimarkClosedSCR
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Generate object of class "multimarksetup" from simulated data sim.data<-simdataClosedSCR() Enc.Mat <- sim.data$Enc.Mat trapCoords <- sim.data$spatialInputs$trapCoords studyArea <- sim.data$spatialInputs$studyArea setup <- processdataSCR(Enc.Mat,trapCoords,studyArea) #Run single chain using the default model for simulated data example.dot<-multimarkClosedSCR(mms=setup)
# This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin #Generate object of class "multimarksetup" from simulated data sim.data<-simdataClosedSCR() Enc.Mat <- sim.data$Enc.Mat trapCoords <- sim.data$spatialInputs$trapCoords studyArea <- sim.data$spatialInputs$studyArea setup <- processdataSCR(Enc.Mat,trapCoords,studyArea) #Run single chain using the default model for simulated data example.dot<-multimarkClosedSCR(mms=setup)
This function generates encounter histories from simulated open population capture-mark-recapture data consisting of multiple non-invasive marks.
simdataCJS( N = 100, noccas = 5, pbeta = -0.25, sigma2_zp = 0, phibeta = 1.6, sigma2_zphi = 0, delta_1 = 0.4, delta_2 = 0.4, alpha = 0.5, data.type = "never", link = "probit" )
simdataCJS( N = 100, noccas = 5, pbeta = -0.25, sigma2_zp = 0, phibeta = 1.6, sigma2_zphi = 0, delta_1 = 0.4, delta_2 = 0.4, alpha = 0.5, data.type = "never", link = "probit" )
N |
Number of individuals. |
noccas |
Number of sampling occasions. |
pbeta |
Logit- or probit-scale intercept term(s) for capture probability (p). Must be a scaler or vector of length |
sigma2_zp |
Logit- or probit-scale individual heterogeneity variance term for capture probability (p). |
phibeta |
Logit- or probit-scale intercept term(s) for survival probability ( |
sigma2_zphi |
Logit- or probit-scale individual heterogeneity variance term for survival probability ( |
delta_1 |
Conditional probability of type 1 encounter, given detection. |
delta_2 |
Conditional probability of type 2 encounter, given detection. |
alpha |
Conditional probability of simultaneous type 1 and type 2 detection, given both types encountered. Only applies when |
data.type |
Specifies the encounter history data type. All data types include non-detections (type 0 encounter), type 1 encounter (e.g., left-side), and type 2 encounters (e.g., right-side). When both type 1 and type 2 encounters occur for the same individual within a sampling occasion, these can either be "non-simultaneous" (type 3 encounter) or "simultaneous" (type 4 encounter). Three data types are currently permitted:
|
link |
Link function for detection probability. Must be " |
A list containing the following:
Enc.Mat |
A matrix containing the observed encounter histories with rows corresponding to individuals and columns corresponding to sampling occasions. |
trueEnc.Mat |
A matrix containing the true (latent) encounter histories with rows corresponding to individuals and columns corresponding to sampling occasions. |
Brett T. McClintock
Bonner, S. J., and Holmberg J. 2013. Mark-recapture with multiple, non-invasive marks. Biometrics 69: 766-775.
McClintock, B. T., Conn, P. B., Alonso, R. S., and Crooks, K. R. 2013. Integrated modeling of bilateral photo-identification data in mark-recapture analyses. Ecology 94: 1464-1471.
#simulate data for data.type="sometimes" using defaults data<-simdataCJS(data.type="sometimes")
#simulate data for data.type="sometimes" using defaults data<-simdataCJS(data.type="sometimes")
This function generates encounter histories from simulated closed population capture-mark-recapture data consisting of multiple non-invasive marks.
simdataClosed( N = 100, noccas = 5, pbeta = -0.4, tau = 0, sigma2_zp = 0, delta_1 = 0.4, delta_2 = 0.4, alpha = 0.5, data.type = "never", link = "logit" )
simdataClosed( N = 100, noccas = 5, pbeta = -0.4, tau = 0, sigma2_zp = 0, delta_1 = 0.4, delta_2 = 0.4, alpha = 0.5, data.type = "never", link = "logit" )
N |
True population size or abundance. |
noccas |
The number of sampling occasions. |
pbeta |
Logit- or probit-scale intercept term(s) for capture probability (p). Must be a scaler or vector of length |
tau |
Additive logit- or probit-scale behavioral effect term for recapture probability (c). |
sigma2_zp |
Logit- or probit-scale individual heterogeneity variance term. |
delta_1 |
Conditional probability of type 1 encounter, given detection. |
delta_2 |
Conditional probability of type 2 encounter, given detection. |
alpha |
Conditional probability of simultaneous type 1 and type 2 detection, given both types encountered. Only applies when |
data.type |
Specifies the encounter history data type. All data types include non-detections (type 0 encounter), type 1 encounter (e.g., left-side), and type 2 encounters (e.g., right-side). When both type 1 and type 2 encounters occur for the same individual within a sampling occasion, these can either be "non-simultaneous" (type 3 encounter) or "simultaneous" (type 4 encounter). Three data types are currently permitted:
|
link |
Link function for detection probability. Must be " |
A list containing the following:
Enc.Mat |
A matrix containing the observed encounter histories with rows corresponding to individuals and columns corresponding to sampling occasions. |
trueEnc.Mat |
A matrix containing the true (latent) encounter histories with rows corresponding to individuals and columns corresponding to sampling occasions. |
Brett T. McClintock
Bonner, S. J., and Holmberg J. 2013. Mark-recapture with multiple, non-invasive marks. Biometrics 69: 766-775.
McClintock, B. T., Conn, P. B., Alonso, R. S., and Crooks, K. R. 2013. Integrated modeling of bilateral photo-identification data in mark-recapture analyses. Ecology 94: 1464-1471.
#simulate data for data.type="sometimes" using defaults data<-simdataClosed(data.type="sometimes")
#simulate data for data.type="sometimes" using defaults data<-simdataClosed(data.type="sometimes")
This function generates encounter histories from spatially-explicit capture-mark-recapture data consisting of multiple non-invasive marks.
simdataClosedSCR( N = 30, ntraps = 9, noccas = 5, pbeta = 0.25, tau = 0, sigma2_scr = 0.75, lambda = 0.75, delta_1 = 0.4, delta_2 = 0.4, alpha = 0.5, data.type = "never", detection = "half-normal", spatialInputs = NULL, buffer = 3 * sqrt(sigma2_scr), ncells = 1024, scalemax = 10, plot = TRUE )
simdataClosedSCR( N = 30, ntraps = 9, noccas = 5, pbeta = 0.25, tau = 0, sigma2_scr = 0.75, lambda = 0.75, delta_1 = 0.4, delta_2 = 0.4, alpha = 0.5, data.type = "never", detection = "half-normal", spatialInputs = NULL, buffer = 3 * sqrt(sigma2_scr), ncells = 1024, scalemax = 10, plot = TRUE )
N |
True population size or abundance. |
ntraps |
The number of traps. If |
noccas |
Scaler indicating the number of sampling occasions per trap. |
pbeta |
Complementary loglog-scale intercept term for detection probability (p). Must be a scaler or vector of length |
tau |
Additive complementary loglog-scale behavioral effect term for recapture probability (c). |
sigma2_scr |
Complementary loglog-scale term for effect of distance in the “half-normal” detection function. Ignored unless |
lambda |
Complementary loglog-scale term for effect of distance in the “exponential” detection function. Ignored unless |
delta_1 |
Conditional probability of type 1 encounter, given detection. |
delta_2 |
Conditional probability of type 2 encounter, given detection. |
alpha |
Conditional probability of simultaneous type 1 and type 2 detection, given both types encountered. Only applies when |
data.type |
Specifies the encounter history data type. All data types include non-detections (type 0 encounter), type 1 encounter (e.g., left-side), and type 2 encounters (e.g., right-side). When both type 1 and type 2 encounters occur for the same individual within a sampling occasion, these can either be "non-simultaneous" (type 3 encounter) or "simultaneous" (type 4 encounter). Three data types are currently permitted:
|
detection |
Model for detection probability as a function of distance from activity centers. Must be " |
spatialInputs |
A list of length 3 composed of objects named
If |
buffer |
A scaler indicating the buffer around the bounding box of |
ncells |
The number of grid cells in the study area when |
scalemax |
Upper bound for grid cell centroid x- and y-coordinates. Default is |
plot |
Logical indicating whether to plot the simulated trap coordinates, study area, and activity centers using |
Please be very careful when specifying your own spatialInputs
; multimarkClosedSCR
and markClosedSCR
do little to verify that these make sense during model fitting.
A list containing the following:
Enc.Mat |
Matrix containing the observed encounter histories with rows corresponding to individuals and ( |
trueEnc.Mat |
Matrix containing the true (latent) encounter histories with rows corresponding to individuals and ( |
spatialInputs |
List of length 2 with objects named
|
centers |
|
Brett T. McClintock
Bonner, S. J., and Holmberg J. 2013. Mark-recapture with multiple, non-invasive marks. Biometrics 69: 766-775.
McClintock, B. T., Conn, P. B., Alonso, R. S., and Crooks, K. R. 2013. Integrated modeling of bilateral photo-identification data in mark-recapture analyses. Ecology 94: 1464-1471.
Royle, J.A., Karanth, K.U., Gopalaswamy, A.M. and Kumar, N.S. 2009. Bayesian inference in camera trapping studies for a class of spatial capture-recapture models. Ecology 90: 3233-3244.
processdataSCR
, multimarkClosedSCR
, markClosedSCR
#simulate data for data.type="sometimes" using defaults data<-simdataClosedSCR(data.type="sometimes")
#simulate data for data.type="sometimes" using defaults data<-simdataClosedSCR(data.type="sometimes")
Example tiger data for multimark
package.
These spatial capture-recapture data with a single mark type are summarized in a list of length 3 containing the following objects:
Enc.Mat
is a 44 x (noccas*ntraps) matrix containing observed encounter histories for 44 tigers across noccas=48
sampling occasions and ntraps=120
traps.
trapCoords
is a matrix of dimension ntraps
x (2 + noccas
) indicating the Cartesian coordinates and operating occasions for the traps, where rows correspond to trap, the first column the x-coordinate, and the second column the y-coordinate. The last noccas
columns indicate whether or not the trap was operating on each of the occasions, where ‘1’ indciates the trap was operating and ‘0’ indicates the trap was not operating.
studyArea
is a 3-column matrix containing the coordinates for the centroids of the contiguous grid of cells that define the study area and available habitat. Each row corresponds to a grid cell. The first 2 columns indicate the Cartesian x- and y-coordinate for the centroid of each grid cell, and the third column indicates whether the cell is available habitat (=1) or not (=0). The grid cells are 0.336 km^2 resolution.
These data were obtained from the R package SPACECAP
and modified by projecting onto a regular rectangular grid consisting of square grid cells (as is required by the spatial capture-recapture models in multimark
).
We thank Ullas Karanth, Wildlife Conservation Society, for providing the tiger data for use as an example with this package.
Gopalaswamy, A.M., Royle, J.A., Hines, J.E., Singh, P., Jathanna, D., Kumar, N. and Karanth, K.U. 2012. Program SPACECAP: software for estimating animal density using spatially explicit capture-recapture models. Methods in Ecology and Evolution 3:1067-1072.
Royle, J.A., Karanth, K.U., Gopalaswamy, A.M. and Kumar, N.S. 2009. Bayesian inference in camera trapping studies for a class of spatial capture-recapture models. Ecology 90: 3233-3244.
data(tiger) #plot the traps and available habitat within the study area plotSpatialData(trapCoords=tiger$trapCoords,studyArea=tiger$studyArea) # This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin # Fit spatial model to tiger data Enc.Mat<-tiger$Enc.Mat trapCoords<-tiger$trapCoords studyArea<-tiger$studyArea tiger.dot<-markClosedSCR(Enc.Mat,trapCoords,studyArea,iter=100,adapt=50,burnin=50) summary(tiger.dot$mcmc)
data(tiger) #plot the traps and available habitat within the study area plotSpatialData(trapCoords=tiger$trapCoords,studyArea=tiger$studyArea) # This example is excluded from testing to reduce package check time # Example uses unrealistically low values for nchain, iter, and burnin # Fit spatial model to tiger data Enc.Mat<-tiger$Enc.Mat trapCoords<-tiger$trapCoords studyArea<-tiger$studyArea tiger.dot<-markClosedSCR(Enc.Mat,trapCoords,studyArea,iter=100,adapt=50,burnin=50) summary(tiger.dot$mcmc)