Package 'degreenet'

Title: Models for Skewed Count Distributions Relevant to Networks
Description: Likelihood-based inference for skewed count distributions, typically of degrees used in network modeling. "degreenet" is a part of the "statnet" suite of packages for network analysis. See Jones and Handcock <doi:10.1098/rspb.2003.2369>.
Authors: Mark S. Handcock [aut, cre, cph]
Maintainer: Mark S. Handcock <[email protected]>
License: GPL-3 + file LICENSE
Version: 1.3-6
Built: 2024-10-26 03:39:37 UTC
Source: https://github.com/cran/degreenet

Help Index


Models for Skewed Count Distributions Relevant to Networks

Description

degreenet is a collection of functions to fit, diagnose, and simulate from distributions for skewed count data. The coverage of distributions is very selective, focusing on those that have been proposed to model the degree distribution on networks. For the rationale for this choice, see the papers in the references section below. For a list of functions type: help(package='degreenet')

For a complete list of the functions, use library(help="degreenet") or read the rest of the manual. For a simple demonstration, use demo(packages="degreenet").

The degreenet package is part of the statnet suite of packages. The suite was developed to facilitate the statistical analysis of network data.

When publishing results obtained using this package alone see the citation in citation(package="degreenet"). The citation for the original paper to use this package is Handcock and Jones (2003) and it should be cited for the theoretical development.

If you use other packages in the statnet suite, please cite it as:

Mark S. Handcock, David R. Hunter, Carter T. Butts, Steven M. Goodreau, and Martina Morris. 2003 statnet: Software tools for the Statistical Modeling of Network Data
https://statnet.org. For complete citation information, use
citation(package="statnet").

All programs derived from this or other statnet packages must cite them appropriately.

Details

See the Handcock and Jones (2003) reference (and the papers it cites and is cited by) for more information on the methodology.

Recent advances in the statistical modeling of random networks have had an impact on the empirical study of social networks. Statistical exponential family models (Strauss and Ikeda 1990) are a generalization of the Markov random network models introduced by Frank and Strauss (1986). These models recognize the complex dependencies within relational data structures. To date, the use of stochastic network models for networks has been limited by three interrelated factors: the complexity of realistic models, the lack of simulation tools for inference and validation, and a poor understanding of the inferential properties of nontrivial models.

This package relies on the network package which allows networks to be represented in R. The statnet suite of packages allows maximum likelihood estimates of exponential random network models to be calculated using Markov Chain Monte Carlo, as well as a broad range of statistical analysis of networks, such as tools for plotting networks, simulating networks and assessing model goodness-of-fit.

For detailed information on how to download and install the software, go to the statnet website: https://statnet.org. A tutorial, support newsgroup, references and links to further resources are provided there.

Author(s)

Mark S. Handcock [email protected]

Maintainer: Mark S. Handcock [email protected]

References

Frank, O., and Strauss, D.(1986). Markov graphs. Journal of the American Statistical Association, 81, 832-842.

Jones, J. H. and Handcock, M. S. (2003). An assessment of preferential attachment as a mechanism for human sexual network formation, Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

Handcock, M. S., Hunter, D. R., Butts, C. T., Goodreau, S. M., and Morris, M. (2003), statnet: Software tools for the Statistical Modeling of Network Data.,
URL https://statnet.org

Strauss, D., and Ikeda, M.(1990). Pseudolikelihood estimation for social networks. Journal of the American Statistical Association, 85, 204-212.


Conway Maxwell Poisson Modeling of Discrete Data

Description

Functions to Estimate the Conway Maxwell Poisson Discrete Probability Distribution via maximum likelihood.

Usage

acmpmle(x, cutoff = 1, cutabove = 1000, guess=c(7,3),
    method="BFGS", conc=FALSE, hellinger=FALSE, hessian=TRUE)

Arguments

x

A vector of counts (one per observation).

cutoff

Calculate estimates conditional on exceeding this value.

cutabove

Calculate estimates conditional on not exceeding this value.

guess

Initial estimate at the MLE.

method

Method of optimization. See "optim" for details.

conc

Calculate the concentration index of the distribution?

hellinger

Minimize Hellinger distance of the parametric model from the data instead of maximizing the likelihood.

hessian

Calculate the hessian of the information matrix (for use with calculating the standard errors.

Value

theta

vector of MLE of the parameters.

asycov

asymptotic covariance matrix.

asycor

asymptotic correlation matrix.

se

vector of standard errors for the MLE.

conc

The value of the concentration index (if calculated).

Note

See the papers on https://handcock.github.io/?q=Holland for details.

Based on the C code in the package compoisson written by Jeffrey Dunn (2008).

References

compoisson: Conway-Maxwell-Poisson Distribution, Jeffrey Dunn, 2008, R package version 0.3

See Also

ayulemle, awarmle, simcmp

Examples

# Simulate a Conway Maxwell Poisson distribution over 100
# observations with mean of 7 and variance of 3
# This leads to a mean of 1

set.seed(1)
s4 <- simcmp(n=100, v=c(7,3))
table(s4)

#
# Calculate the MLE and an asymptotic confidence
# interval for the parameters
#

acmpmle(s4)

Discrete version of q-Exponential Modeling of Discrete Data

Description

Functions to Estimate the Discrete version of q-Exponential Probability Distribution via maximum likelihood.

Usage

adqemle(x, cutoff = 1, cutabove = 1000, guess = c(3.5,1),
    method = "BFGS", conc = FALSE, hellinger = FALSE, hessian=TRUE)

Arguments

x

A vector of counts (one per observation).

cutoff

Calculate estimates conditional on exceeding this value.

cutabove

Calculate estimates conditional on not exceeding this value.

guess

Initial estimate at the MLE.

conc

Calculate the concentration index of the distribution?

method

Method of optimization. See "optim" for details.

hellinger

Minimize Hellinger distance of the parametric model from the data instead of maximizing the likelihood.

hessian

Calculate the hessian of the information matrix (for use with calculating the standard errors.

Value

theta

vector of MLE of the parameters.

asycov

asymptotic covariance matrix.

asycor

asymptotic correlation matrix.

se

vector of standard errors for the MLE.

conc

The value of the concentration index (if calculated).

Note

See the papers on https://handcock.github.io/?q=Holland for details

References

Jones, J. H. and Handcock, M. S. "An assessment of preferential attachment as a mechanism for human sexual network formation," Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

See Also

ayulemle, adqemle, simdqe

Examples

# Simulate a Discrete version of q-Exponential distribution over 100
# observations with a PDF exponent of 3.5 and a 
# sigma scale of 1

set.seed(1)
s4 <- simdqe(n=100, v=c(3.5,1))
table(s4)

#
# Calculate the MLE and an asymptotic confidence
# interval for the parameters
#

s4est <- adqemle(s4)
s4est

# Calculate the MLE and an asymptotic confidence
# interval for rho under the Yule model
#

s4yuleest <- ayulemle(s4)
s4yuleest

#
# Compare the AICC and BIC for the two models
#

lldqeall(v=s4est$theta,x=s4)
llyuleall(v=s4yuleest$theta,x=s4)

Poisson Lognormal Modeling of Discrete Data

Description

Functions to Estimate the Poisson Lognormal Discrete Probability Distribution via maximum likelihood.

Usage

aplnmle(x, cutoff = 1, cutabove = 1000, guess = c(0.6,1.2),
    method = "BFGS", conc = FALSE, hellinger = FALSE, hessian=TRUE,logn=TRUE)

Arguments

x

A vector of counts (one per observation).

cutoff

Calculate estimates conditional on exceeding this value.

cutabove

Calculate estimates conditional on not exceeding this value.

guess

Initial estimate at the MLE.

method

Method of optimization. See "optim" for details.

conc

Calculate the concentration index of the distribution?

hellinger

Minimize Hellinger distance of the parametric model from the data instead of maximizing the likelihood.

hessian

Calculate the hessian of the information matrix (for use with calculating the standard errors.

logn

Use logn parametrization, that is, mean and variance on the observation scale.

Value

theta

vector of MLE of the parameters.

asycov

asymptotic covariance matrix.

asycor

asymptotic correlation matrix.

se

vector of standard errors for the MLE.

conc

The value of the concentration index (if calculated).

Note

See the papers on https://handcock.github.io/?q=Holland for details

References

Jones, J. H. and Handcock, M. S. "An assessment of preferential attachment as a mechanism for human sexual network formation," Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

See Also

ayulemle, awarmle, simpln

Examples

# Simulate a Poisson Lognormal distribution over 100
# observations with lognormal mean of -1 and lognormal variance of 1
# This leads to a mean of 1

set.seed(1)
s4 <- simpln(n=100, v=c(-1,1))
table(s4)

#
# Calculate the MLE and an asymptotic confidence
# interval for the parameters
#

s4est <- aplnmle(s4)
s4est

# Calculate the MLE and an asymptotic confidence
# interval for rho under the Yule model
#

s4yuleest <- ayulemle(s4)
s4yuleest

# Calculate the MLE and an asymptotic confidence
# interval for rho under the Waring model
#

s4warest <- awarmle(s4)
s4warest

#
# Compare the AICC and BIC for the three models
#

llplnall(v=s4est$theta,x=s4)
llyuleall(v=s4yuleest$theta,x=s4)
llwarall(v=s4warest$theta,x=s4)

Waring Modeling of Discrete Data

Description

Functions to Estimate the Waring Discrete Probability Distribution via maximum likelihood.

Usage

awarmle(x, cutoff = 1, cutabove = 1000, guess = c(3.5,0.1),
    method = "BFGS", conc = FALSE, hellinger = FALSE, hessian=TRUE)

Arguments

x

A vector of counts (one per observation).

cutoff

Calculate estimates conditional on exceeding this value.

cutabove

Calculate estimates conditional on not exceeding this value.

guess

Initial estimate at the MLE.

conc

Calculate the concentration index of the distribution?

method

Method of optimization. See "optim" for details.

hellinger

Minimize Hellinger distance of the parametric model from the data instead of maximizing the likelihood.

hessian

Calculate the hessian of the information matrix (for use with calculating the standard errors.

Value

theta

vector of MLE of the parameters.

asycov

asymptotic covariance matrix.

asycor

asymptotic correlation matrix.

se

vector of standard errors for the MLE.

conc

The value of the concentration index (if calculated).

Note

See the papers on https://handcock.github.io/?q=Holland for details

References

Jones, J. H. and Handcock, M. S. "An assessment of preferential attachment as a mechanism for human sexual network formation," Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

See Also

ayulemle, awarmle, simwar

Examples

# Simulate a Waring distribution over 100
# observations with a PDf exponent of 3.5 and a 
# probability of including a new actor of 0.1

set.seed(1)
s4 <- simwar(n=100, v=c(3.5,0.1))
table(s4)

#
# Calculate the MLE and an asymptotic confidence
# interval for the parameters
#

s4est <- awarmle(s4)
s4est

# Calculate the MLE and an asymptotic confidence
# interval for rho under the Yule model
#

s4yuleest <- ayulemle(s4)
s4yuleest

#
# Compare the AICC and BIC for the two models
#

llwarall(v=s4est$theta,x=s4)
llyuleall(v=s4yuleest$theta,x=s4)

Yule Distribution Modeling of Discrete Data

Description

Functions to Estimate the Yule Discrete Probability Distribution via maximum likelihood.

Usage

ayulemle(x, cutoff = 1, cutabove = 1000, guess = 3.5, conc = FALSE,
     method = "BFGS", hellinger = FALSE, hessian = TRUE, weights = rep(1, length(x)))

Arguments

x

A vector of counts (one per observation).

cutoff

Calculate estimates conditional on exceeding this value.

cutabove

Calculate estimates conditional on not exceeding this value.

guess

Initial estimate at the MLE.

conc

Calculate the concentration index of the distribution?

method

Method of optimization. See "optim" for details.

hellinger

Minimize Hellinger distance of the parametric model from the data instead of maximizing the likelihood.

hessian

Calculate the hessian of the information matrix (for use with calculating the standard errors.

weights

sample weights on the observed counts.

Value

theta

vector of MLE of the parameters.

asycov

asymptotic covariance matrix.

asycor

asymptotic correlation matrix.

se

vector of standard errors for the MLE.

conc

The value of the concentration index (if calculated).

Note

See the papers on https://handcock.github.io/?q=Holland for details

References

Jones, J. H. and Handcock, M. S. "An assessment of preferential attachment as a mechanism for human sexual network formation," Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

See Also

ayulemle, awarmle, simyule

Examples

# Simulate a Yule distribution over 100
# observations with PDf exponent of 3.5

set.seed(1)
s4 <- simyule(n=100, rho=3.5)
table(s4)

#
# Calculate the MLE and an asymptotic confidence
# interval for the parameters
#

s4est <- ayulemle(s4)
s4est

#
# Compute the AICC and BIC for the model
#

llyuleall(v=s4est$theta,x=s4)

Calculate Bootstrap Estimates and Confidence Intervals for the Discrete Pareto Distribution

Description

Uses the parametric bootstrap to estimate the bias and confidence interval of the MLE of the Discrete Pareto Distribution.

Usage

bsdp(x, cutoff=1, m=200, np=1, alpha=0.95)
bootstrapdp(x,cutoff=1,cutabove=1000,
                          m=200,alpha=0.95,guess=3.31,hellinger=FALSE,
                          mle.meth="adpmle")

Arguments

x

A vector of counts (one per observation).

cutoff

Calculate estimates conditional on exceeding this value.

m

Number of bootstrap samples to draw.

np

Number of parameters in the model (1 by default).

alpha

Type I error for the confidence interval.

hellinger

Minimize Hellinger distance of the parametric model from the data instead of maximizing the likelihood.

cutabove

Calculate estimates conditional on not exceeding this value.

guess

Initial estimate at the MLE.

mle.meth

Method to use to compute the MLE.

Value

dist

matrix of sample CDFs, one per row.

obsmle

The Discrete Pareto MLE of the PDF exponent.

bsmles

Vector of bootstrap MLE.

quantiles

Quantiles of the bootstrap MLEs.

pvalue

p-value of the Anderson-Darling statistics relative to the bootstrap MLEs.

obsmands

Observed Anderson-Darling Statistic.

meanmles

Mean of the bootstrap MLEs.

guess

Initial estimate at the MLE.

mle.meth

Method to use to compute the MLE.

Note

See the papers on https://handcock.github.io/?q=Holland for details

References

Jones, J. H. and Handcock, M. S. "An assessment of preferential attachment as a mechanism for human sexual network formation," Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

See Also

anbmle, simdp, lldp

Examples

## Not run: 
# Now, simulate a Discrete Pareto distribution over 100
# observations with expected count 1 and probability of another
# of 0.2

set.seed(1)
s4 <- simdp(n=100, v=3.31)
table(s4)

#
# Calculate the MLE and an asymptotic confidence
# interval for the parameter.
#

s4est <- adpmle(s4)
s4est

#
# Use the bootstrap to compute a confidence interval rather than using the 
# asymptotic confidence interval for the parameter.
#

bsdp(s4, m=20)

## End(Not run)

Calculate Bootstrap Estimates and Confidence Intervals for the Negative Binomial Distribution

Description

Uses the parametric bootstrap to estimate the bias and confidence interval of the MLE of the Negative Binomial Distribution.

Usage

bsnb(x, cutoff=1, m=200, np=2, alpha=0.95, hellinger=FALSE)
bootstrapnb(x,cutoff=1,cutabove=1000,
                          m=200,alpha=0.95,guess=c(5, 0.2),
                          file="none")

Arguments

x

A vector of counts (one per observation).

cutoff

Calculate estimates conditional on exceeding this value.

m

Number of bootstrap samples to draw.

np

Number of parameters in the model (1 by default).

alpha

Type I error for the confidence interval.

hellinger

Minimize Hellinger distance of the parametric model from the data instead of maximizing the likelihood.

cutabove

Calculate estimates conditional on not exceeding this value.

guess

Guess at the parameter value.

file

Name of the file to store the results. By default do not save the results.

Value

dist

matrix of sample CDFs, one per row.

obsmle

The Negative Binomial MLE of the PDF exponent.

bsmles

Vector of bootstrap MLE.

quantiles

Quantiles of the bootstrap MLEs.

pvalue

p-value of the Anderson-Darling statistics relative to the bootstrap MLEs.

obsmands

Observed Anderson-Darling Statistic.

meanmles

Mean of the bootstrap MLEs.

guess

Initial estimate at the MLE.

mle.meth

Method to use to compute the MLE.

Note

See the papers on https://handcock.github.io/?q=Holland for details

References

Jones, J. H. and Handcock, M. S. "An assessment of preferential attachment as a mechanism for human sexual network formation," Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

See Also

anbmle, simnb, llnb

Examples

# Now, simulate a Negative Binomial distribution over 100
# observations with expected count 1 and probability of another
# of 0.2

set.seed(1)
s4 <- simnb(n=100, v=c(5,0.2))
table(s4)

#
# Calculate the MLE and an asymptotic confidence
# interval for the parameter.
#

s4est <- anbmle(s4)
s4est

#
# Use the bootstrap to compute a confidence interval rather than using the 
# asymptotic confidence interval for the parameter.
#

bsnb(s4, m=20)

Calculate Bootstrap Estimates and Confidence Intervals for the Poisson Lognormal Distribution

Description

Uses the parametric bootstrap to estimate the bias and confidence interval of the MLE of the Poisson Lognormal Distribution.

Usage

bspln(x, cutoff=1, m=200, np=2, alpha=0.95, v=NULL,
                   hellinger=FALSE)
bootstrappln(x,cutoff=1,cutabove=1000,
                          m=200,alpha=0.95,guess=c(0.6,1.2), file = "none")

Arguments

x

A vector of counts (one per observation).

cutoff

Calculate estimates conditional on exceeding this value.

m

Number of bootstrap samples to draw.

np

Number of parameters in the model (1 by default).

alpha

Type I error for the confidence interval.

v

Parameter value to use for the bootstrap distribution. By default it is the MLE of the data.

hellinger

Minimize Hellinger distance of the parametric model from the data instead of maximizing the likelihood.

cutabove

Calculate estimates conditional on not exceeding this value.

guess

Initial estimate at the MLE.

file

Name of the file to store the results. By default do not save the results.

Value

dist

matrix of sample CDFs, one per row.

obsmle

The Poisson Lognormal MLE of the PDF exponent.

bsmles

Vector of bootstrap MLE.

quantiles

Quantiles of the bootstrap MLEs.

pvalue

p-value of the Anderson-Darling statistics relative to the bootstrap MLEs.

obsmands

Observed Anderson-Darling Statistic.

meanmles

Mean of the bootstrap MLEs.

Note

See the papers on https://handcock.github.io/?q=Holland for details

References

Jones, J. H. and Handcock, M. S. "An assessment of preferential attachment as a mechanism for human sexual network formation," Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

See Also

anbmle, simpln, llpln

Examples

# Now, simulate a Poisson Lognormal distribution over 100
# observations with expected count 1 and probability of another
# of 0.2

set.seed(1)
s4 <- simpln(n=100, v=c(5,0.2))
table(s4)

#
# Calculate the MLE and an asymptotic confidence
# interval for the parameter.
#

s4est <- aplnmle(s4)
s4est

#
# Use the bootstrap to compute a confidence interval rather than using the 
# asymptotic confidence interval for the parameter.
#

bspln(s4, m=5)

Calculate Bootstrap Estimates and Confidence Intervals for the Waring Distribution

Description

Uses the parametric bootstrap to estimate the bias and confidence interval of the MLE of the Waring Distribution.

Usage

bswar(x, cutoff=1, m=200, np=2, alpha=0.95, v=NULL,
                   hellinger=FALSE)
bootstrapwar(x,cutoff=1,cutabove=1000,
             m=200,alpha=0.95,guess=c(3.31, 0.1),file="none",
             conc = FALSE)

Arguments

x

A vector of counts (one per observation).

cutoff

Calculate estimates conditional on exceeding this value.

m

Number of bootstrap samples to draw.

np

Number of parameters in the model (1 by default).

alpha

Type I error for the confidence interval.

v

Parameter value to use for the bootstrap distribution. By default it is the MLE of the data.

hellinger

Minimize Hellinger distance of the parametric model from the data instead of maximizing the likelihood.

cutabove

Calculate estimates conditional on not exceeding this value.

guess

Guess at the parameter value.

file

Name of the file to store the results. By default do not save the results.

conc

Calculate the concentration index of the distribution?

Value

dist

matrix of sample CDFs, one per row.

obsmle

The Waring MLE of the PDF exponent.

bsmles

Vector of bootstrap MLE.

quantiles

Quantiles of the bootstrap MLEs.

pvalue

p-value of the Anderson-Darling statistics relative to the bootstrap MLEs.

obsmands

Observed Anderson-Darling Statistic.

meanmles

Mean of the bootstrap MLEs.

guess

Initial estimate at the MLE.

mle.meth

Method to use to compute the MLE.

Note

See the papers on https://handcock.github.io/?q=Holland for details

References

Jones, J. H. and Handcock, M. S. "An assessment of preferential attachment as a mechanism for human sexual network formation," Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

See Also

anbmle, simwar, llwar

Examples

# Now, simulate a Waring distribution over 100
# observations with expected count 1 and probability of another
# of 0.2

set.seed(1)
s4 <- simwar(n=100, v=c(5,0.2))
table(s4)

#
# Calculate the MLE and an asymptotic confidence
# interval for the parameter.
#

s4est <- awarmle(s4)
s4est

#
# Use the bootstrap to compute a confidence interval rather than using the 
# asymptotic confidence interval for the parameter.
#

bswar(s4, m=20)

Calculate Bootstrap Estimates and Confidence Intervals for the Yule Distribution

Description

Uses the parametric bootstrap to estimate the bias and confidence interval of the MLE of the Yule Distribution.

Usage

bsyule(x, cutoff=1, m=200, np=1, alpha=0.95, v=NULL,
                   hellinger=FALSE, cutabove=1000)
bootstrapyule(x,cutoff=1,cutabove=1000,
                          m=200,alpha=0.95,guess=3.31,hellinger=FALSE,
                          mle.meth="ayulemle")

Arguments

x

A vector of counts (one per observation).

cutoff

Calculate estimates conditional on exceeding this value.

m

Number of bootstrap samples to draw.

np

Number of parameters in the model (1 by default).

alpha

Type I error for the confidence interval.

v

Parameter value to use for the bootstrap distribution. By default it is the MLE of the data.

hellinger

Minimize Hellinger distance of the parametric model from the data instead of maximizing the likelihood.

cutabove

Calculate estimates conditional on not exceeding this value.

guess

Initial estimate at the MLE.

mle.meth

Method to use to compute the MLE.

Value

dist

matrix of sample CDFs, one per row.

obsmle

The Yule MLE of the PDF exponent.

bsmles

Vector of bootstrap MLE.

quantiles

Quantiles of the bootstrap MLEs.

pvalue

p-value of the Anderson-Darling statistics relative to the bootstrap MLEs.

obsmands

Observed Anderson-Darling Statistic.

meanmles

Mean of the bootstrap MLEs.

Note

See the papers on https://handcock.github.io/?q=Holland for details

References

Jones, J. H. and Handcock, M. S. "An assessment of preferential attachment as a mechanism for human sexual network formation," Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

See Also

ayulemle, simyule, llyule

Examples

# Now, simulate a Yule distribution over 100
# observations with rho=4.0

set.seed(1)
s4 <- simyule(n=100, rho=4)
table(s4)

#
# Calculate the MLE and an asymptotic confidence
# interval for rho
#

s4est <- ayulemle(s4)
s4est

#
# Use the bootstrap to compute a confidence interval rather than using the 
# asymptotic confidence interval for rho.
#

bsyule(s4, m=20)

Models for Count Distributions

Description

Functions to Estimate Parametric Discrete Probability Distributions via maximum likelihood Based on categorical response

Usage

gyulemle(x, cutoff = 1, cutabove = 1000, guess = 3.5, conc = FALSE, 
    method = "BFGS", hellinger = FALSE, hessian=TRUE)

Arguments

x

A vector of categories for counts (one per observation). The values of x and the categories are: 0=0, 1=1, 2=2, 3=3, 4=4, 5=5-10, 6=11-20, 7=21-100, 8=>100

cutoff

Calculate estimates conditional on exceeding this value.

cutabove

Calculate estimates conditional on not exceeding this value.

guess

Initial estimate at the MLE.

conc

Calculate the concentration index of the distribution?

method

Method of optimization. See "optim" for details.

hellinger

Minimize Hellinger distance of the parametric model from the data instead of maximizing the likelihood.

hessian

Calculate the hessian of the information matrix (for use with calculating the standard errors.

Value

result

vector of parameter estimates - lower 95% confidence value, upper 95% confidence value, the PDF MLE, the asymptotic standard error, and the number of data values >=cutoff and <=cutabove.

theta

The Yule MLE of the PDF exponent.

value

The maximized value of the function.

conc

The value of the concentration index (if calculated).

Note

See the papers on https://handcock.github.io/?q=Holland for details

References

Jones, J. H. and Handcock, M. S. "An assessment of preferential attachment as a mechanism for human sexual network formation," Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

Examples

#
# Simulate a Yule distribution over 100
# observations with rho=4.0
#

set.seed(1)
s4 <- simyule(n=100, rho=4)
table(s4)

#
# Recode it as categorical
#

s4[s4 >  4 & s4 < 11] <- 5
s4[s4 > 100] <- 8
s4[s4 >  20] <- 7
s4[s4 >  10] <- 6

#
# Calculate the MLE and an asymptotic confidence
# interval for rho
#

s4est <- gyulemle(s4)
s4est

#
# Calculate the MLE and an asymptotic confidence
# interval for rho under the Waring model (i.e., rho=4, p=2/3)
#

s4warest <- gwarmle(s4)
s4warest

#
# Compare the AICC and BIC for the two models
#

llgyuleall(v=s4est$theta,x=s4)
llgwarall(v=s4warest$theta,x=s4)

Calculate the Conditional log-likelihood for Count Distributions

Description

Functions to Estimate the Conditional Log-likelihood for Discrete Probability Distributions. The likelihood is calcualted condition on the count being at least the cutoff value and less than or equal to the cutabove value.

Usage

llgyule(v, x, cutoff=1,cutabove=1000,xr=1:10000,hellinger=FALSE)

Arguments

v

A vector of parameters for the Yule (a 1-vector - the scaling exponent).

x

A vector of categories for counts (one per observation). The values of x and the categories are: 0=0, 1=1, 2=2, 3=3, 4=4, 5=5-10, 6=11-20, 7=21-100, 8=>100

cutoff

Calculate estimates conditional on exceeding this value.

cutabove

Calculate estimates conditional on not exceeding this value.

xr

range of count values to use to approximate the set of all realistic counts.

hellinger

Calculate the Hellinger distance of the parametric model from the data instead of the log-likelihood?

Value

the log-likelihood for the data x at parameter value v (or the Hellinder distance if hellinger=TRUE).

Note

See the papers on https://handcock.github.io/?q=Holland for details

References

Jones, J. H. and Handcock, M. S. "An assessment of preferential attachment as a mechanism for human sexual network formation," Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

See Also

gyulemle, llgyuleall, dyule

Examples

#
# Simulate a Yule distribution over 100
# observations with rho=4.0
#

set.seed(1)
s4 <- simyule(n=100, rho=4)
table(s4)

#
# Recode it as categorical
#

s4[s4 >  4 & s4 < 11] <- 5
s4[s4 > 100] <- 8
s4[s4 >  20] <- 7
s4[s4 >  10] <- 6

#
# Calculate the MLE and an asymptotic confidence
# interval for rho
#

s4est <- gyulemle(s4)
s4est

#
# Calculate the MLE and an asymptotic confidence
# interval for rho under the Waring model (i.e., rho=4, p=2/3)
#

s4warest <- gwarmle(s4)
s4warest

#
# Compare the log-likelihoods for the two models
#

llgyule(v=s4est$theta,x=s4)
llgwar(v=s4warest$theta,x=s4)

Calculate the log-likelihood for Count Distributions

Description

Functions to Estimate the Log-likelihood for Discrete Probability Distributions Based on Categorical Response.

Usage

llgyuleall(v, x, cutoff = 2, cutabove = 1000,  np=1)

Arguments

v

A vector of parameters for the Yule (a 1-vector - the scaling exponent).

x

A vector of categories for counts (one per observation). The values of x and the categories are: 0=0, 1=1, 2=2, 3=3, 4=4, 5=5-10, 6=11-20, 7=21-100, 8=>100

cutoff

Calculate estimates conditional on exceeding this value.

cutabove

Calculate estimates conditional on not exceeding this value.

np

wnumber of parameters in the model. For the Yule this is 1.

Value

the log-likelihood for the data x at parameter value v.

Note

See the papers on https://handcock.github.io/?q=Holland for details

References

Jones, J. H. and Handcock, M. S. "An assessment of preferential attachment as a mechanism for human sexual network formation," Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

See Also

gyulemle, llgyule, dyule, llgwarall

Examples

#
# Simulate a Yule distribution over 100
# observations with rho=4.0
#

set.seed(1)
s4 <- simyule(n=100, rho=4)
table(s4)

#
# Recode it as categorical
#

s4[s4 >  4 & s4 < 11] <- 5
s4[s4 > 100] <- 8
s4[s4 >  20] <- 7
s4[s4 >  10] <- 6


#
# Calculate the MLE and an asymptotic confidence
# interval for rho
#

s4est <- gyulemle(s4)
s4est

# Calculate the MLE and an asymptotic confidence
# interval for rho under the Waring model (i.e., rho=4, p=2/3)
#

s4warest <- gwarmle(s4)
s4warest

#
# Compare the AICC and BIC for the two models
#

llgyuleall(v=s4est$theta,x=s4)
llgwarall(v=s4warest$theta,x=s4)

Calculate the Conditional log-likelihood for the Poisson Lognormal Distributions

Description

Compute the Conditional Log-likelihood for the Poisson Lognormal Discrete Probability Distribution. The likelihood is calculated conditionl on the count being at least the cutoff value and less than or equal to the cutabove value.

Usage

llpln(v, x, cutoff=1,cutabove=1000,xr=1:10000,hellinger=FALSE,logn = TRUE)

Arguments

v

A vector of parameters for the Yule (a 1-vector - the scaling exponent).

x

A vector of counts (one per observation).

cutoff

Calculate estimates conditional on exceeding this value.

cutabove

Calculate estimates conditional on not exceeding this value.

xr

range of count values to use to approximate the set of all realistic counts.

hellinger

Calculate the Hellinger distance of the parametric model from the data instead of the log-likelihood?

logn

Use logn parametrization, that is, mean and variance on the observation scale.

Value

the log-likelihood for the data x at parameter value v (or the Hellinder distance if hellinger=TRUE).

Note

See the papers on https://handcock.github.io/?q=Holland for details

References

Jones, J. H. and Handcock, M. S. "An assessment of preferential attachment as a mechanism for human sexual network formation," Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

See Also

aplnmle, llplnall, dpln

Examples

# Simulate a Poisson Lognormal distribution over 100
# observations with lognormal mean -1 and logormal standard deviation 1.

set.seed(1)
s4 <- simpln(n=100, v=c(-1,1))
table(s4)

#
# Calculate the MLE and an asymptotic confidence
# interval for rho
#

s4est <- aplnmle(s4)
s4est

#
# Calculate the MLE and an asymptotic confidence
# interval for rho under the Waring model
#

s4warest <- awarmle(s4)
s4warest

#
# Compare the log-likelihoods for the two models
#

llpln(v=s4est$theta,x=s4)
llwar(v=s4warest$theta,x=s4)

Calculate the Conditional log-likelihood for Count Distributions

Description

Functions to Estimate the Conditional Log-likelihood for Discrete Probability Distributions. The likelihood is calcualted condition on the count being at least the cutoff value and less than or equal to the cutabove value.

Usage

llyule(v, x, cutoff=1,cutabove=1000, xr=1:10000 ,hellinger=FALSE,
       weights = rep(1, length(x)))

Arguments

v

A vector of parameters for the Yule (a 1-vector - the scaling exponent).

x

A vector of counts (one per observation).

cutoff

Calculate estimates conditional on exceeding this value.

cutabove

Calculate estimates conditional on not exceeding this value.

xr

range of count values to use to approximate the set of all realistic counts.

hellinger

Calculate the Hellinger distance of the parametric model from the data instead of the log-likelihood?

weights

sample weights on the observed counts.

Value

the log-likelihood for the data x at parameter value v (or the Hellinder distance if hellinger=TRUE).

Note

See the papers on https://handcock.github.io/?q=Holland for details

References

Jones, J. H. and Handcock, M. S. "An assessment of preferential attachment as a mechanism for human sexual network formation," Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

See Also

ayulemle, llyuleall, dyule

Examples

# Simulate a Yule distribution over 100
# observations with rho=4.0

set.seed(1)
s4 <- simyule(n=100, rho=4)
table(s4)

#
# Calculate the MLE and an asymptotic confidence
# interval for rho
#

s4est <- ayulemle(s4)
s4est

#
# Calculate the MLE and an asymptotic confidence
# interval for rho under the Waring model (i.e., rho=4, p=2/3)
#

s4warest <- awarmle(s4)
s4warest

#
# Compare the log-likelihoods for the two models
#

llyule(v=s4est$theta,x=s4)
llwar(v=s4warest$theta,x=s4)

Calculate the log-likelihood for Count Distributions

Description

Functions to Estimate the Log-likelihood for Discrete Probability Distributions.

Usage

llyuleall(v, x, cutoff = 2, cutabove = 1000,  np=1)

Arguments

v

A vector of parameters for the Yule (a 1-vector - the scaling exponent).

x

A vector of counts (one per observation).

cutoff

Calculate estimates conditional on exceeding this value.

cutabove

Calculate estimates conditional on not exceeding this value.

np

wnumber of parameters in the model. For the Yule this is 1.

Value

the log-likelihood for the data x at parameter value v.

Note

See the papers on https://handcock.github.io/?q=Holland for details

References

Jones, J. H. and Handcock, M. S. "An assessment of preferential attachment as a mechanism for human sexual network formation," Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

See Also

ayulemle, llyule, dyule, llwarall

Examples

# Simulate a Yule distribution over 100
# observations with rho=4.0

set.seed(1)
s4 <- simyule(n=100, rho=4)
table(s4)

#
# Calculate the MLE and an asymptotic confidence
# interval for rho
#

s4est <- ayulemle(s4)
s4est

# Calculate the MLE and an asymptotic confidence
# interval for rho under the Waring model (i.e., rho=4, p=2/3)
#

s4warest <- awarmle(s4)
s4warest

#
# Compare the AICC and BIC for the two models
#

llyuleall(v=s4est$theta,x=s4)
llwarall(v=s4warest$theta,x=s4)

Generate a undirected network with a given sequence of degrees

Description

Generate a undirected network where the degree of each actor is specified. The degree is the number of actors the actor is tied to. This returns a network object and requires the igraph package.

Usage

reedmolloy(deg, maxit=10, verbose=TRUE)

Arguments

deg

vector of counts where element ii is the degree of actor ii. Its sum should be even.

maxit

integer; maximum number of jitterings of the degree sequence to find a valid network.

verbose

Print out details of the progress of the algorithm.

Value

The network is returned as a network object.

Note

See the papers on https://handcock.github.io/?q=Holland for details

References

Jones, J. H. and Handcock, M. S. "An assessment of preferential attachment as a mechanism for human sexual network formation," Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

See Also

ayulemle, dyule

Examples

# Now, simulate a Poisson Lognormal distribution over 100
# observations with mean = -1 and s.d. = 1.

set.seed(2)
s4 <- simpln(n=100, v=c(-1,1))
table(s4)
#
simr <- reedmolloy(s4)
simr

Rounded Poisson Lognormal Modeling of Discrete Data

Description

Functions to Estimate the Rounded Poisson Lognormal Discrete Probability Distribution via maximum likelihood.

Usage

rplnmle(x, cutoff = 1, cutabove = 1000, guess = c(0.6,1.2),
    method = "BFGS", conc = FALSE, hellinger = FALSE, hessian=TRUE)

Arguments

x

A vector of counts (one per observation).

cutoff

Calculate estimates conditional on exceeding this value.

cutabove

Calculate estimates conditional on not exceeding this value.

guess

Initial estimate at the MLE.

conc

Calculate the concentration index of the distribution?

method

Method of optimization. See "optim" for details.

hellinger

Minimize Hellinger distance of the parametric model from the data instead of maximizing the likelihood.

hessian

Calculate the hessian of the information matrix (for use with calculating the standard errors.

Value

theta

vector of MLE of the parameters.

asycov

asymptotic covariance matrix.

asycor

asymptotic correlation matrix.

se

vector of standard errors for the MLE.

conc

The value of the concentration index (if calculated).

Note

See the papers on https://handcock.github.io/?q=Holland for details

References

Jones, J. H. and Handcock, M. S. "An assessment of preferential attachment as a mechanism for human sexual network formation," Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

See Also

aplnmle

Examples

# Simulate a Poisson Lognormal distribution over 100
# observations with lognormal mean of -1 and lognormal variance of 1
# This leads to a mean of 1

set.seed(1)
s4 <- simpln(n=100, v=c(-1,1))
table(s4)

#
# Calculate the MLE and an asymptotic confidence
# interval for the parameters
#

s4est <- rplnmle(s4)
s4est

Generate a (non-random) network from a Yule Distribution

Description

Generate a network with a given number of actors having a degree distribution draw from a Yule distribution. The resultant network is not random - that is, is not a random draw from all such networks.

Usage

ryule(n=20,rho=2.5, maxdeg=n-1, maxit=10, verbose=FALSE)

Arguments

n

Number of actors in the network.

rho

PDF exponent of the Yule distribution.

maxdeg

Maximum degree to sample (using truncation of the distribution). If this is greater then n-1 then n-1 is used.

maxit

integer; maximum number of resamplings of the degree sequence to find a valid network.

verbose

Print out details of the progress of the algorithm.

Value

If the network package is available, the network is returned as a network object. If not a sociomatrix is returned.

Note

See the papers on https://handcock.github.io/?q=Holland for details

References

Jones, J. H. and Handcock, M. S. "An assessment of preferential attachment as a mechanism for human sexual network formation," Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

See Also

ayulemle, dyule, reedmolloy

Examples

# Now, simulate a Yule network of 30
# actors with rho=4.0
ryule(n=30, rho=4)

Simulate from a Conway Maxwell Poisson Distribution

Description

Functions to generate random samples from a Conway Maxwell Poisson Probability Distribution

Usage

simcmp(n=100, v=c(7,2.6), maxdeg=10000)

Arguments

n

number of samples to draw.

v

Conway Maxwell Poisson parameters: lognormal mean and lognormal s.d.

maxdeg

Maximum degree to sample (using truncation of the distribution).

Value

vector of random draws or samples.

Note

See the papers on https://handcock.github.io/?q=Holland for details

References

Jones, J. H. and Handcock, M. S. "An assessment of preferential attachment as a mechanism for human sexual network formation," Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

See Also

acmpmle, dcmp

Examples

# Now, simulate a Conway Maxwell Poisson distribution over 100
# observations with lognormal mean -1 and lognormal standard deviation 1.

set.seed(1)
s4 <- simcmp(n=100, v=c(7,3))
table(s4)

Simulate from a Discrete Pareto Distribution

Description

Functions to generate random samples from a Discrete Pareto Probability Distribution

Usage

simdp(n=100, v=3.5, maxdeg=10000)

Arguments

n

number of samples to draw.

v

Discrete Pareto parameters: PDF exponent.

maxdeg

Maximum degree to sample (using truncation of the distribution).

Value

vector of random draws or samples.

Note

See the papers on https://handcock.github.io/?q=Holland for details

References

Jones, J. H. and Handcock, M. S. "An assessment of preferential attachment as a mechanism for human sexual network formation," Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

See Also

adpmle, ddp

Examples

## Not run: 
# Now, simulate a Discrete Pareto distribution over 100
# observations with lognormal mean -1 and lognormal standard deviation 1.

set.seed(1)
s4 <- simdp(n=100, v=3.5)
table(s4)

## End(Not run)

Simulate from a Negative Binomial Distribution

Description

Functions to generate random samples from a Negative Binomial Probability Distribution

Usage

simnb(n=100, v=c(5,0.2), maxdeg=10000)

Arguments

n

number of samples to draw.

v

Negative Binomial parameters: expected count and probability of another.

maxdeg

Maximum degree to sample (using truncation of the distribution).

Value

vector of random draws or samples.

Note

See the papers on https://handcock.github.io/?q=Holland for details

References

Jones, J. H. and Handcock, M. S. "An assessment of preferential attachment as a mechanism for human sexual network formation," Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

See Also

anbmle, dnb

Examples

# Now, simulate a Negative Binomial distribution over 100
# observations with lognormal mean -1 and lognormal standard deviation 1.

set.seed(1)
s4 <- simnb(n=100, v=c(5,0.2))
table(s4)

Simulate from a Poisson Lognormal Distribution

Description

Functions to generate random samples from a Poisson Lognormal Probability Distribution

Usage

simpln(n=100, v=c(0.6,1.2), maxdeg=10000, cutoff=1)

Arguments

n

number of samples to draw.

v

Poisson Lognormal parameters: lognormal mean and lognormal s.d.

maxdeg

Maximum degree to sample (using truncation of the distribution).

cutoff

Calculate estimates conditional on exceeding this value.

Value

vector of random draws or samples.

Note

See the papers on https://handcock.github.io/?q=Holland for details

References

Jones, J. H. and Handcock, M. S. "An assessment of preferential attachment as a mechanism for human sexual network formation," Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

See Also

aplnmle, dpln

Examples

# Now, simulate a Poisson Lognormal distribution over 100
# observations with lognormal mean -1 and lognormal standard deviation 1.

set.seed(1)
s4 <- simpln(n=100, v=c(-1,1))
table(s4)

Simulate from a Waring Distribution

Description

Functions to generate random samples from a Waring Probability Distribution

Usage

simwar(n=100, v=c(3.5, 0.1), maxdeg=10000)

Arguments

n

number of samples to draw.

v

Waring parameters: scaling exponent and probability of a new actor.

maxdeg

Maximum degree to sample (using truncation of the distribution).

Value

vector of random draws or samples.

Note

See the papers on https://handcock.github.io/?q=Holland for details

References

Jones, J. H. and Handcock, M. S. "An assessment of preferential attachment as a mechanism for human sexual network formation," Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

See Also

awarmle, dwar

Examples

# Now, simulate a Waring distribution over 100
# observations with Waring with exponent 3.5 and probability of a new
# actor 0.1.

set.seed(1)
s4 <- simwar(n=100, v=c(3.5, 0.1))
table(s4)

Simulate from a Yule Distribution

Description

Functions to generate random samples from a Yule Probability Distribution

Usage

simyule(n=100, rho=4, maxdeg=10000)

Arguments

n

number of samples to draw.

rho

Yule PDF exponent.

maxdeg

Maximum degree to sample (using truncation of the distribution).

Value

vector of random draws or samples.

Note

See the papers on https://handcock.github.io/?q=Holland for details

References

Jones, J. H. and Handcock, M. S. "An assessment of preferential attachment as a mechanism for human sexual network formation," Proceedings of the Royal Society, B, 2003, 270, 1123-1128.

See Also

ayulemle, dyule

Examples

# Now, simulate a Yule distribution over 100
# observations with rho=4.0

set.seed(1)
s4 <- simyule(n=100, rho=4)
table(s4)

Number of sex partners in the last 12 months for men and women in Sweden

Description

This is a data set used in Jones and Handcock (2002) and

The data are counts of the numbers of sex partners for men and women in the last twelve months. The data from the 1996 “Sex in Sweden" survey based on a nationwide probability sample and financed by the Swedish National Board of Health.

Usage

data(sweden)

Source

We thanks Dr. Bo Lewin, Professor of Sociology, Uppsala University and head of the research team responsible for the “Sex in Sweden" study for providing the Swedish data used in this study. This research supported by Grant 7R01DA012831-02 from NIDA and Grant 1R01HD041877 from NICHD.

References

Lewin, B. (1996). Sex in Sweden, Stockholm: National Institute of Public Health.

Handcock, Mark S. and Jones, James Holland (2004), “Likelihood-Based Inference for Stochastic Models of Sexual Network Formation" Theoretical Population Biology, doi:10.1016/j.tpb.2003.09.006.

Jones, James Holland and Handcock, Mark S. (2003), Nature, 423, 6940, 605-606.

Handcock, Mark S. and Jones, James Holland (2003), “An assessment of preferential attachment as a mechanism for human sexual network formation" Proceedings of the Royal Society, B., 270, 1123-1128.

See Also

ayulemle