Many network outcomes are counts, such as the number of emails between individuals or the number of conflicts between countries. One way to model such outcomes is to treat them as ordinal and fit an ordinal probit model. This can be done with the ame
command using the family=ord
option. However, the MCMC algorithm for this model can take a long time to run, and can mix poorly if the data do not contain much information. In such cases, it may be preferable to use a more parametric approach, such as a Poisson regression model.
In this tutorial, I illustrate how to fit an overdispersed Poisson model to some network data with a count outcome, using the functions provided by the amen
package. The model we will fit is the following:
zi,j=β⊤xi,j+ai+bj+u⊤ivj+ϵi,jyi,j|zi,j∼Poisson(ezi,j) where Var(aibi)=Σab, Var(uivi)=Σuv, and Var(ϵi,jϵj,i)=Σϵ=σ2(1ρρ1). We can think of the parameter σ2 here as an overdispersion parameter, in the sense that if σ2 were zero, then yi,j would be Poisson with mean exp(β⊤xi,j+ai+bj+u⊤ivj). On the other hand, if σ2 is large then the variance of yi,j will be larger than exp(β⊤xi,j+ai+bj+u⊤ivj).
We will construct a Gibbs sampler to approximate the joint posterior distribution of all unknown quantities, which includes
the global parameters β,Σab,Σuv,σ2,ρ;
the latent nodal effects {(ai,bi,ui,vi):i=1,…,n};
the unobserved latent dyadic variables zi,j.
We will fit a rank-1 AMEN model to the dataset on sheep dominance encounters that is included with the amen
package:
library(amen)
$dom
Y<-sheep$age
age<-sheep
1:8,1:8] Y[
## V1 V2 V3 V4 V5 V6 V7 V8
## [1,] NA 0 0 0 0 0 0 1
## [2,] 0 NA 0 0 5 2 1 0
## [3,] 0 0 NA 0 7 4 0 0
## [4,] 0 0 8 NA 0 0 0 0
## [5,] 0 0 0 0 NA 1 0 0
## [6,] 0 0 0 7 0 NA 0 0
## [7,] 0 0 1 0 0 0 NA 0
## [8,] 0 1 1 4 5 3 0 NA
1:8] age[
## [1] 9 8 9 9 7 8 9 9
Here, yi,j is the number of times sheep i “dominated” sheep j over some time period, and xi is the age of sheep i in years.
boxplot( apply(Y,1,mean,na.rm=TRUE) ~ age )
boxplot( apply(Y,2,mean,na.rm=TRUE) ~ age )
Roughly speaking, older sheep dominate younger sheep.
Let’s fit a model for yi,j with main and interaction effects of age. Let xi=agei−¯age be the centered age variable, and let xi,j=(xi,xj,xixj).
The corresponding design matrix can be constructed as follows:
$age - mean(sheep$age)
x<-sheep
design_array(Xrow=x,Xcol=x,Xdyad=outer(x,x)) X<-
Now we need some starting values for some of the unknown quantities in the model:
log(Y+1) ; diag(Z)<-0
Z<-
diag(2)
Sab<-
1 ; Suv<-diag(2*R) ; U<-V<-matrix(0,nrow(Y),R)
R<-
1 ; rho<-0 s2<-
Now if we observed Z and assumed it followed the AME model above, we could just use the following MCMC algorithm to approximate the joint posterior distribution of the AME model parameters:
for(s in 1:10)
{ # update beta, a and b
rbeta_ab_fc(Z,Sab,rho,X,s2,offset=U%*%t(V))
tmp<-$beta ; a<-tmp$a ; b<-tmp$b
beta<-tmp
# update UV
rUV_fc(Z,U,V,Suv,rho,offset=Xbeta(X,beta)+outer(a,b,"+"))
tmp<-$U ; V<-tmp$V
U<-tmp
# update Sab
rSab_fc(a,b)
Sab<-
# update Suv
rSuv_fc(U,V)
Suv<-
# update s2
rs2_fc(Z,rho,offset=Xbeta(X,beta)+outer(a,b,"+")+U%*%t(V))
s2<-
# update rho
rrho_mh(Z,rho,s2,offset=Xbeta(X,beta)+outer(a,b,"+")+U%*%t(V))
rho<- }
However, in our model Z is the matrix of unknown log expectations of the elements of Y. To account for the fact that Z is unknown, we need to integrate/simulate over the possible values of Z in the Markov chain. This can be done by adding a Metropolis update for Z to the above MCMC algorithm. We do this as follows:
Simulate a candidate value for each zi,j: z∗i,j∼N(zi,j,cσ2), where c is a tuning parameter.
Replace the pair (zi,j,zj,i) with the proposed values (z∗i,j,z∗j,i) if logp(yi,j|z∗i,j)p(yj,i|z∗j,i)p(z∗i,j,z∗j,i)p(yi,j|zi,j)p(yj,i|zj,i)p(zi,j,zj,i)>logh where h∼U(0,1). Here p(zi,j,zj,i) is the joint density of (zi,j,zj,i) conditional on β,ai,bi,ui,vi,aj,bj,uj,vj and ρ and σ2. This joint density is not exotic: it is simply a bivariate normal density for (zi,j,zj,i) where the mean vector is (β⊤xi,j+ai+bj+u⊤ivj,β⊤xj,i+aj+bi+u⊤jvi) and covariance matrix Σϵ.
Calculation of the log numerator (and denominator) in the above ratio can be done with the amen
function ldZgbme
. The letters stand for “log density of Z from a gbme model”. You can search the web for more on gbme models. In R, updating Z using this procedure looks like the following:
# propose candidate Z
+matrix(rnorm(nrow(Y)^2),nrow(Y),nrow(Y))*sqrt(s2)
Zp<-Z
# compute acceptance ratio
Xbeta(X,beta)+outer(a,b,"+")+U%*%t(V)
EZ<-ldZgbme(Zp,Y,function(y,z){ dpois(y,exp(z),log=TRUE) },EZ,rho,s2) -
lr<- ldZgbme(Z, Y,function(y,z){ dpois(y,exp(z),log=TRUE) },EZ,rho,s2)
# simulate symmetric matrix of (log) uniform rvs
matrix(log(runif(nrow(Y)^2)),nrow(Y),nrow(Y))
lh<-lower.tri(lh)]<-t(lh)[lower.tri(lh)]
lh[
# update dyads for which lr>lh, and keep track of acceptances
>lh]<-Zp[lr>lh] Z[lr
Note that the ldZgbme
function takes as an argument a function that specifies the log density (or mass) function for each yi,j given zi,j. So for our Poisson model, this function is function(y,z){ dpois(y,exp(z),log=TRUE) }
. You can fit other types of models by specifying different log density functions.
Let’s implement this:
#### Starting values
log(Y+1) ; diag(Z)<-0
Z<-
diag(2)
Sab<-
1 ; Suv<-diag(2*R) ; U<-V<-matrix(0,nrow(Y),R)
R<-
1 ; rho<-0
s2<-
#### Parameter values to be saved
NULL
BETA<-VE<-LL<-matrix(0,nrow(Y),nrow(Y))
ACC<-
#### MCMC
set.seed(1)
for(s in 1:10000)
{ ## Update AMEN parameters
# update beta, a and b
rbeta_ab_fc(Z,Sab,rho,X,s2,offset=U%*%t(V))
tmp<-$beta ; a<-tmp$a ; b<-tmp$b
beta<-tmp
# update UV
rUV_fc(Z,U,V,Suv,rho,offset=Xbeta(X,beta)+outer(a,b,"+"))
tmp<-$U ; V<-tmp$V
U<-tmp
# update Sab
rSab_fc(a,b)
Sab<-
# update Suv
rSuv_fc(U,V)
Suv<-
# update s2
rs2_fc(Z,rho,offset=Xbeta(X,beta)+outer(a,b,"+")+U%*%t(V))
s2<-
# update rho
rrho_mh(Z,rho,s2,offset=Xbeta(X,beta)+outer(a,b,"+")+U%*%t(V))
rho<-
## Update Z
# propose candidate Z
+matrix(rnorm(nrow(Y)^2),nrow(Y),nrow(Y))*sqrt(s2)
Zp<-Z
# compute acceptance ratio
Xbeta(X,beta)+outer(a,b,"+")+U%*%t(V)
EZ<-ldZgbme(Zp,Y,function(y,z){ dpois(y,exp(z),log=TRUE) },EZ,rho,s2) -
lr<- ldZgbme(Z, Y,function(y,z){ dpois(y,exp(z),log=TRUE) },EZ,rho,s2)
# simulate symmetric matrix of (log) uniform rvs
matrix(log(runif(nrow(Y)^2)),nrow(Y),nrow(Y))
lh<-lower.tri(lh)]<-t(lh)[lower.tri(lh)]
lh[
# update dyads for which lr>lh, and keep track of acceptances
>lh]<-Zp[lr>lh]
Z[lr>lh]<-ACC[lr>lh]+1
ACC[lr
## Output
if(s%%25==0)
{ cat(s,range(ACC/s),"\n")
rbind(BETA,beta)
BETA<-rbind(VE,c(s2,rho))
VE<-
par(mfrow=c(1,3),mar=c(3,3,1,1),mgp=c(1.75,.75,0))
matplot(BETA,type="l")
c(LL,sum(dpois(Y,exp(Z),log=TRUE),na.rm=TRUE)) ; plot(LL,type="l")
LL<-hist(ACC/s)
} }
Let’s interpret the parameter estimates:
apply(BETA,2,mean)
## [1] -1.62872581 0.19926245 -0.33723549 0.05615933
apply(BETA,2,sd)
## [1] 0.31035583 0.08344497 0.05881178 0.01047963
So the older a sheep is, the more likely it is do dominate others, and the younger a sheep is, the more likely it is to be dominated. Also, the interaction term suggests some homophily on age in terms of number of dominance encounters.
apply(VE,2,mean)
## [1] 0.9692351 -0.7259700
apply(VE,2,sd)
## [1] 0.1199476 0.1731239
There is a strong negative dyadic correlation for these data, reflecting that if i dominates j frequently, j dominates i less frequently than expected based on the age effects and random effects.