The jSDM_binomial_probit function performs a Binomial probit regression in a Bayesian framework. The function calls a Gibbs sampler written in 'C++' code which uses conjugate priors to estimate the conditional posterior distribution of model's parameters.
burnin: The number of burnin iterations for the sampler.
mcmc: The number of Gibbs iterations for the sampler. Total number of Gibbs iterations is equal to burnin+mcmc. burnin+mcmc must be divisible by 10 and superior or equal to 100 so that the progress bar can be displayed.
thin: The thinning interval used in the simulation. The number of mcmc iterations must be divisible by this value.
presence_data: A matrix nsitexnspecies indicating the presence by a 1 (or the absence by a 0) of each species on each site.
site_formula: A one-sided formula of the form '~x1+...+xp' specifying the explanatory variables for the suitability process of the model, used to form the design matrix X of size nsitexnp.
trait_data: A data frame containing the species traits which can be included as part of the model. Default to NULL to fit a model without species traits.
trait_formula: A one-sided formula of the form '~ t1 + ... + tk + x1:t1 + ... + xp:tk' specifying the interactions between the environmental variables and the species traits to be considered in the model, used to form the trait design matrix Tr of size nspeciesxnt and to set to 0 the γ parameters corresponding to interactions not taken into account according to the formula. Default to NULL to fit a model with all possible interactions between species traits found in trait_data and environmental variables defined by site_formula.
site_data: A data frame containing the model's explanatory variables by site.
n_latent: An integer which specifies the number of latent variables to generate. Defaults to 0.
site_effect: A string indicating whether row effects are included as fixed effects ("fixed"), as random effects ("random"), or not included ("none") in the model. If fixed effects, then for parameter identifiability the first row effect is set to zero, which analogous to acting as a reference level when dummy variables are used. If random effects, they are drawn from a normal distribution with mean zero and unknown variance, analogous to a random intercept in mixed models. Defaults to "none".
lambda_start: Starting values for λ parameters corresponding to the latent variables for each species must be either a scalar or a nlatentxnspecies upper triangular matrix with strictly positive values on the diagonal, ignored if n_latent=0. If lambda_start takes a scalar value, then that value will serve for all of the λ parameters except those concerned by the constraints explained above.
W_start: Starting values for latent variables must be either a scalar or a nsitexnlatent matrix, ignored if n_latent=0. If W_start takes a scalar value, then that value will serve for all of the Wil with i=1,...,nsite and l=1,...,nlatent.
beta_start: Starting values for β parameters of the suitability process for each species must be either a scalar or a npxnspecies matrix. If beta_start takes a scalar value, then that value will serve for all of the β parameters.
alpha_start: Starting values for random site effect parameters must be either a scalar or a nsite-length vector, ignored if site_effect="none". If alpha_start takes a scalar value, then that value will serve for all of the α parameters.
gamma_start: Starting values for γ parameters that represent the influence of species-specific traits on species' responses β, gamma_start must be either a scalar, a vector of length nt, np or nt.np or a ntxnp matrix. If gamma_start takes a scalar value, then that value will serve for all of the γ parameters. If gamma_start is a vector of length nt or nt.np the resulting ntxnp matrix is filled by column with specified values, if a np-length vector is given, the matrix is filled by row.
V_alpha: Starting value for variance of random site effect if site_effect="random" or constant variance of the Gaussian prior distribution for the fixed site effect if site_effect="fixed". Must be a strictly positive scalar, ignored if site_effect="none".
mu_beta: Means of the Normal priors for the β parameters of the suitability process. mu_beta must be either a scalar or a np-length vector. If mu_beta takes a scalar value, then that value will serve as the prior mean for all of the β parameters. The default value is set to 0 for an uninformative prior, ignored if trait_data is specified.
V_beta: Variances of the Normal priors for the β parameters of the suitability process. V_beta must be either a scalar or a npxnp symmetric positive semi-definite square matrix. If V_beta takes a scalar value, then that value will serve as the prior variance for all of the β parameters, so the variance covariance matrix used in this case is diagonal with the specified value on the diagonal. The default variance is large and set to 10 for an uninformative flat prior.
mu_lambda: Means of the Normal priors for the λ parameters corresponding to the latent variables. mu_lambda must be either a scalar or a nlatent-length vector. If mu_lambda takes a scalar value, then that value will serve as the prior mean for all of the λ parameters. The default value is set to 0 for an uninformative prior.
V_lambda: Variances of the Normal priors for the λ parameters corresponding to the latent variables. V_lambda must be either a scalar or a nlatentxnlatent symmetric positive semi-definite square matrix. If V_lambda takes a scalar value, then that value will serve as the prior variance for all of λ parameters, so the variance covariance matrix used in this case is diagonal with the specified value on the diagonal. The default variance is large and set to 10 for an uninformative flat prior.
mu_gamma: Means of the Normal priors for the γ parameters. mu_gamma must be either a scalar, a vector of length nt, np or nt.np or a ntxnp matrix. If mu_gamma takes a scalar value, then that value will serve as the prior mean for all of the γ parameters. If mu_gamma is a vector of length nt or nt.np the resulting ntxnp matrix is filled by column with specified values, if a np-length vector is given, the matrix is filled by row. The default value is set to 0 for an uninformative prior, ignored if trait_data=NULL.
V_gamma: Variances of the Normal priors for the γ parameters. V_gamma must be either a scalar, a vector of length nt, np or nt.np or a ntxnp positive matrix. If V_gamma takes a scalar value, then that value will serve as the prior variance for all of the γ parameters. If V_gamma is a vector of length nt or nt.np the resulting ntxnp matrix is filled by column with specified values, if a np-length vector is given, the matrix is filled by row. The default variance is large and set to 10 for an uninformative flat prior, ignored if trait_data=NULL.
shape_Valpha: Shape parameter of the Inverse-Gamma prior for the random site effect variance V_alpha, ignored if site_effect="none" or site_effect="fixed". Must be a strictly positive scalar. Default to 0.5 for weak informative prior.
rate_Valpha: Rate parameter of the Inverse-Gamma prior for the random site effect variance V_alpha, ignored if site_effect="none" or site_effect="fixed"
Must be a strictly positive scalar. Default to 0.0005 for weak informative prior.
seed: The seed for the random number generator. Default to 1234.
verbose: A switch (0,1) which determines whether or not the progress of the sampler is printed to the screen. Default is 1: a progress bar is printed, indicating the step (in %) reached by the Gibbs sampler.
Returns
An object of class "jSDM" acting like a list including : - mcmc.alpha: An mcmc object that contains the posterior samples for site effects α, not returned if site_effect="none".
mcmc.V_alpha: An mcmc object that contains the posterior samples for variance of random site effect, not returned if site_effect="none" or site_effect="fixed".
mcmc.latent: A list by latent variable of mcmc objects that contains the posterior samples for latent variables Wl with l=1,...,nlatent, not returned if n_latent=0.
mcmc.sp: A list by species of mcmc objects that contains the posterior samples for species effects βj and λj if n_latent>0 with j=1,...,nspecies.
mcmc.gamma: A list by covariates of mcmc objects that contains the posterior samples for γp parameters with p=1,...,np if trait_data is specified.
mcmc.Deviance: The posterior sample of the deviance (D) is also provided, with D defined as : D=−2log(∏ijP(yij∣βj,λj,αi,Wi)).
Z_latent: Predictive posterior mean of the latent variable Z.
probit_theta_latent: Predictive posterior mean of the probability to each species to be present on each site, transformed by probit link function.
theta_latent: Predictive posterior mean of the probability to each species to be present on each site.
model_spec: Various attributes of the model fitted, including the response and model matrix used, distributional assumptions as link function, family and number of latent variables, hyperparameters used in the Bayesian estimation and mcmc, burnin and thin.
The mcmc. objects can be summarized by functions provided by the coda package.
Details
We model an ecological process where the presence or absence of species j on site i is explained by habitat suitability.
Ecological process:
yij∼Bernoulli(θij)yijBernoulli(θij),
where
if n_latent=0 and site_effect="none"
probit (θij)=Xiβj
if n_latent>0 and site_effect="none"
probit (θij)=Xiβj+Wiλj
if n_latent=0 and site_effect="random"
probit (θij)=Xiβj+αi and αiN(0,Vα)
if n_latent>0 and site_effect="fixed"
probit (θij)=Xiβj+Wiλj+αi
if n_latent=0 and site_effect="fixed"
probit (θij)=Xiβj+αi
if n_latent>0 and site_effect="random"
probit (θij)=Xiβj+Wiλj+αi and αiN(0,Vα)
In the absence of data on species traits (trait_data=NULL), the effect of species j: βj; follows the same a priori Gaussian distribution such that βjNnp(μβ,Vβ), for each species.
If species traits data are provided, the effect of species j: βj; follows an a priori Gaussian distribution such that βjNnp(μ\betaj,Vβ), where muβ.jp=Σ(k=1,...,nt)tjk.γkp, takes different values for each species.
We assume that γ.kpN(μγ.kp,Vγ.kp) as prior distribution.
We define the matrix γ=(γkp)fork=1,...,ntandp=1,...,np such as :
x_0
x_1
...
x_p
...
x_np
__________
________
________
________
________
________
t_0 |
gamma_(0,0)
γ(0,1)
...
γ(0,p)
...
γ(0,np)
{ effect of
|
intercept
environmental
|
variables
t_1 |
γ(1,0)
γ(1,1)
...
γ(1,p)
...
γ(1,np)
... |
...
...
...
...
...
...
t_k |
γ(k,0)
γ(k,1)
...
γ(k,p)
...
γ(k,np)
... |
...
...
...
...
...
...
t_nt |
γ(nt,0)
γ(nt,1)
...
γ(nt,p)
...
γ(nt,np)
average
trait effect
interaction
traits
environment
Examples
#======================================# jSDM_binomial_probit()# Example with simulated data#====================================#=================#== Load librarieslibrary(jSDM)#==================#' #== Data simulation#= Number of sitesnsite <-150#= Set seed for repeatabilityseed <-1234set.seed(seed)#= Number of speciesnsp<-20#= Number of latent variablesn_latent <-2#= Ecological process (suitability)x1 <- rnorm(nsite,0,1)x2 <- rnorm(nsite,0,1)X <- cbind(rep(1,nsite),x1,x2)np <- ncol(X)#= Latent variables WW <- matrix(rnorm(nsite*n_latent,0,1), nsite, n_latent)#= Fixed species effect beta beta.target <- t(matrix(runif(nsp*np,-2,2), byrow=TRUE, nrow=nsp))#= Factor loading lambda lambda.target <- matrix(0, n_latent, nsp)mat <- t(matrix(runif(nsp*n_latent,-2,2), byrow=TRUE, nrow=nsp))lambda.target[upper.tri(mat, diag=TRUE)]<- mat[upper.tri(mat, diag=TRUE)]diag(lambda.target)<- runif(n_latent,0,2)#= Variance of random site effect V_alpha.target <-0.5#= Random site effect alphaalpha.target <- rnorm(nsite,0, sqrt(V_alpha.target))# Simulation of response data with probit linkprobit_theta <- X%*%beta.target + W%*%lambda.target + alpha.target
theta <- pnorm(probit_theta)e <- matrix(rnorm(nsp*nsite,0,1),nsite,nsp)# Latent variable Z Z_true <- probit_theta + e
# Presence-absence matrix YY <- matrix (NA, nsite,nsp)for(i in1:nsite){for(j in1:nsp){if( Z_true[i,j]>0){Y[i,j]<-1}else{Y[i,j]<-0}}}#==================================#== Site-occupancy model# Increase number of iterations (burnin and mcmc) to get convergencemod<-jSDM_binomial_probit(# Iteration burnin=200, mcmc=200, thin=1,# Response variable presence_data=Y,# Explanatory variables site_formula=~x1+x2, site_data = X, n_latent=2, site_effect="random",# Starting values alpha_start=0, beta_start=0, lambda_start=0, W_start=0, V_alpha=1,# Priors shape_Valpha=0.5, rate_Valpha=0.0005, mu_beta=0, V_beta=1, mu_lambda=0, V_lambda=1, seed=1234, verbose=1)# ===================================================# Result analysis# ===================================================#==========#== Outputsoldpar <- par(no.readonly =TRUE)#= Parameter estimates## beta_jmean_beta <- matrix(0,nsp,ncol(X))pdf(file=file.path(tempdir(),"Posteriors_beta_jSDM_probit.pdf"))par(mfrow=c(ncol(X),2))for(j in1:nsp){ mean_beta[j,]<- apply(mod$mcmc.sp[[j]][,1:ncol(X)],2, mean)for(p in1:ncol(X)){ coda::traceplot(mod$mcmc.sp[[j]][,p]) coda::densplot(mod$mcmc.sp[[j]][,p], main = paste(colnames(mod$mcmc.sp[[j]])[p],", species : ",j)) abline(v=beta.target[p,j],col='red')}}dev.off()## lambda_jmean_lambda <- matrix(0,nsp,n_latent)pdf(file=file.path(tempdir(),"Posteriors_lambda_jSDM_probit.pdf"))par(mfrow=c(n_latent*2,2))for(j in1:nsp){ mean_lambda[j,]<- apply(mod$mcmc.sp[[j]][,(ncol(X)+1):(ncol(X)+n_latent)],2, mean)for(l in1:n_latent){ coda::traceplot(mod$mcmc.sp[[j]][,ncol(X)+l]) coda::densplot(mod$mcmc.sp[[j]][,ncol(X)+l], main=paste(colnames(mod$mcmc.sp[[j]])[ncol(X)+l],", species : ",j)) abline(v=lambda.target[l,j],col='red')}}dev.off()# Species effects beta and factor loadings lambdapar(mfrow=c(1,2))plot(t(beta.target), mean_beta, main="species effect beta", xlab ="obs", ylab ="fitted")abline(a=0,b=1,col='red')plot(t(lambda.target), mean_lambda, main="factor loadings lambda", xlab ="obs", ylab ="fitted")abline(a=0,b=1,col='red')## W latent variablespar(mfrow=c(1,2))for(l in1:n_latent){ plot(W[,l], summary(mod$mcmc.latent[[paste0("lv_",l)]])[[1]][,"Mean"], main = paste0("Latent variable W_", l), xlab ="obs", ylab ="fitted") abline(a=0,b=1,col='red')}## alphapar(mfrow=c(1,3))plot(alpha.target, summary(mod$mcmc.alpha)[[1]][,"Mean"], xlab ="obs", ylab ="fitted", main="site effect alpha")abline(a=0,b=1,col='red')## Valphacoda::traceplot(mod$mcmc.V_alpha)coda::densplot(mod$mcmc.V_alpha)abline(v=V_alpha.target,col='red')## Deviancesummary(mod$mcmc.Deviance)plot(mod$mcmc.Deviance)#= Predictions## Zpar(mfrow=c(1,2))plot(Z_true,mod$Z_latent, main="Z_latent", xlab="obs", ylab="fitted")abline(a=0,b=1,col='red')## probit_thetaplot(probit_theta,mod$probit_theta_latent, main="probit(theta)",xlab="obs",ylab="fitted")abline(a=0,b=1,col='red')## probabilities thetapar(mfrow=c(1,1))plot(theta,mod$theta_latent, main="theta",xlab="obs",ylab="fitted")abline(a=0,b=1,col='red')par(oldpar)
References
Chib, S. and Greenberg, E. (1998) Analysis of multivariate probit models. Biometrika, 85, 347-361.
Warton, D. I.; Blanchet, F. G.; O'Hara, R. B.; O'Hara, R. B.; Ovaskainen, O.; Taskinen, S.; Walker, S. C. and Hui, F. K. C. (2015) So Many Variables: Joint Modeling in Community Ecology. Trends in Ecology & Evolution, 30, 766-779.
Ovaskainen, O., Tikhonov, G., Norberg, A., Blanchet, F. G., Duan, L., Dunson, D., Roslin, T. and Abrego, N. (2017) How to make more out of community data? A conceptual framework and its implementation as models and software. Ecology Letters, 20, 561-576.