edit

Public Documentation

Documentation for MixedModels.jl's public interface.

MixedModels

# MixedModels.GeneralizedLinearMixedModelType.

GeneralizedLinearMixedModel

Generalized linear mixed-effects model representation

Members:

  • LMM: a LinearMixedModel - the local approximation to the GLMM.
  • β: the fixed-effects vector
  • β₀: similar to β. User in the PIRLS algorithm if step-halving is needed.
  • θ: covariance parameter vector
  • b: similar to u, equivalent to broadcast!(*, b, LMM.Λ, u)
  • u: a vector of matrices of random effects
  • u₀: similar to u. Used in the PIRLS algorithm if step-halving is needed.
  • resp: a GlmResp object
  • η: the linear predictor
  • wt: vector of prior case weights, a value of T[] indicates equal weights.

source

# MixedModels.LinearMixedModelType.

LinearMixedModel

Linear mixed-effects model representation

Members

  • formula: the formula for the model
  • fixefnames: names of the fixed effects (for displaying coefficients)
  • wttrms: a length nt vector of weighted model matrices. The last two elements are X and y.
  • trms: a vector of unweighted model matrices. If isempty(sqrtwts) the same object as wttrms
  • Λ: a length nt - 2 vector of lower triangular matrices
  • sqrtwts: the Diagonal matrix of the square roots of the case weights. Allowed to be size 0
  • A: an nt × nt symmetric matrix of matrices representing hcat(Z,X,y)'hcat(Z,X,y)
  • R: a nt × nt matrix of matrices - the upper Cholesky factor of Λ'AΛ+I
  • opt: an OptSummary object

source

# MixedModels.OptSummaryType.

OptSummary

Summary of an NLopt optimization

Members

  • initial: a copy of the initial parameter values in the optimization
  • lowerbd: lower bounds on the parameter values
  • ftol_rel: as in NLopt
  • ftol_abs: as in NLopt
  • xtol_rel: as in NLopt
  • xtol_abs: as in NLopt
  • final: a copy of the final parameter values from the optimization
  • fmin: the final value of the objective
  • feval: the number of function evaluations
  • optimizer: the name of the optimizer used, as a Symbol
  • returnvalue: the return value, as a Symbol

source

# MixedModels.ReMatType.

ReMat

A representation of the model matrix for a random-effects term

source

# MixedModels.ScalarReMatType.

ScalarReMat

The representation of the model matrix for a scalar random-effects term

Members

  • f: the grouping factor as a CategoricalVector
  • z: the raw random-effects model matrix as a Vector
  • fnm: the name of the grouping factor as a Symbol
  • cnms: a Vector of column names

source

# MixedModels.VarCorrType.

VarCorr

An encapsulation of information on the fitted random-effects variance-covariance matrices.

Members

  • Λ: the vector of lower triangular matrices from the MixedModel
  • fnms: a Vector{ASCIIString} of grouping factor names
  • cnms: a Vector{Vector{ASCIIString}} of column names
  • s: the estimate of σ, the standard deviation of the per-observation noise

The main purpose of defining this type is to isolate the logic in the show method.

source

# MixedModels.VectorReMatType.

VectorReMat

The representation of the model matrix for a vector-valued random-effects term

Members

  • f: the grouping factor as a CategoricalVector
  • z: the transposed raw random-effects model matrix
  • fnm: the name of the grouping factor as a Symbol
  • cnms: a Vector of column names (row names after transposition) of z

source

# Base.LinAlg.condMethod.

cond(m::MixedModel)

Returns the vector of the condition numbers of the blocks of m.Λ

source

# Base.stdMethod.

std{T}(m::MixedModel{T})

The estimated standard deviations of the variance components as a Vector{Vector{T}}.

source

# MixedModels.LaplaceDevianceMethod.

LaplaceDeviance{T}(m::GeneralizedLinearMixedModel{T})::T

Return the Laplace approximation to the deviance of m.

If the distribution D does not have a scale parameter the Laplace approximation is defined as the squared length of the conditional modes, u, plus the determinant of Λ'Z'ZΛ + 1, plus the sum of the squared deviance residuals.

source

# MixedModels.bootstrap!Method.

bootstrap!{T}(r::Matrix{T}, m::LinearMixedModel{T}, f!::Function;
    β=fixef(m), σ=sdest(m), θ=getθ(m))

Overwrite columns of r with the results of applying the mutating extractor f! to parametric bootstrap replications of model m.

The signature of f! should be f!{T}(v::AbstractVector{T}, m::LinearMixedModel{T})

Named Arguments

β::Vector{T}, σ::T, and θ::Vector{T} are the values of the parameters in m for simulation of the responses.

source

# MixedModels.bootstrapMethod.

bootstrap{T}(N, m::LinearMixedModel{T},
    β::Vector{T}=fixef(m), σ::T=sdest(m), θ::Vector{T}=getθ(m))

Perform N parametric bootstrap replication fits of m, returning a data frame with column names :obj, the objective function at convergence, , the estimated standard deviation of the residuals, βᵢ, i = 1,...,p, the fixed-effects coefficients, and θᵢ, i = 1,...,k the covariance parameters.

Named Arguments

β::Vector{T}, σ::T, and θ::Vector{T} are the values of the parameters in m for simulation of the responses.

source

# MixedModels.cfactor!Method.

cfactor!(A::AbstractMatrix)

A slightly modified version of chol! from Base

Uses inject! (as opposed to copy!), downdate! (as opposed to syrk! or gemm!) and recursive calls to cfactor!.

Note: The cfactor! method for dense matrices calls LAPACK.potrf! directly to avoid errors being thrown when A is computationally singular

source

# MixedModels.condVarMethod.

condVar(m::MixedModel)

Return the conditional variances matrices of the random effects.

The random effects are returned by ranef as a vector of length k, where k is the number of random effects terms. The ith element is a matrix of size vᵢ × ℓᵢ where vᵢ is the size of the vector-valued random effects for each of the ℓᵢ levels of the grouping factor. Technically those values are the modes of the conditional distribution of the random effects given the observed data.

This function returns an array of k three dimensional arrays, where the ith array is of size vᵢ × vᵢ × ℓᵢ. These are the diagonal blocks from the conditional variance-covariance matrix,

s² Λ(Λ'Z'ZΛ + I)⁻¹Λ'

source

# MixedModels.fixefMethod.

fixef(m::MixedModel)

Returns the fixed-effects parameter vector estimate.

source

# MixedModels.getθMethod.

getθ(A::LowerTriangular{T, Matrix{T}})

Return a vector of the elements of the lower triangle of A (column-major ordering)

source

# MixedModels.glmmMethod.

glmm(f::Formula, fr::ModelFrame, d::Distribution[, l::GLM.Link])

Return a GeneralizedLinearMixedModel object.

The value is ready to be fit! but has not yet been fit.

source

# MixedModels.lmmMethod.

lmm(f::DataFrames.Formula, fr::DataFrames.DataFrame; weights = [])

Create a LinearMixedModel from f, which contains both fixed-effects terms and random effects, and fr.

The return value is ready to be fit! but has not yet been fit.

source

# MixedModels.lmmMethod.

lmm(m::MixedModel)

Extract the LinearMixedModel from a MixedModel.

If m is a LinearMixedModel return m. If m is a GeneralizedLinearMixedModel return m.LMM.

source

# MixedModels.lowerbdMethod.

lowerbd{T}(A::LowerTriangular{T,Matrix{T}})

Return the vector of lower bounds on the parameters, θ.

These are the elements in the lower triangle in column-major ordering. Diagonals have a lower bound of 0. Off-diagonals have a lower-bound of -Inf.

source

# MixedModels.lowerbdMethod.

lowerbd(m::LinearMixedModel)

Return the vector of lower bounds on the covariance parameter vector θ

source

# MixedModels.objectiveMethod.

objective(m::LinearMixedModel)

Return negative twice the log-likelihood of model m

source

# MixedModels.pirls!Function.

pirls!(m::GeneralizedLinearMixedModel)

Use Penalized Iteratively Reweighted Least Squares (PIRLS) to determine the conditional modes of the random effects.

When varyβ is true both u and β are optimized with PIRLS. Otherwise only u is optimized and β is held fixed.

Passing verbose = true provides verbose output of the iterations.

source

# MixedModels.pwrssMethod.

pwrss(m::LinearMixedModel)

The penalized residual sum-of-squares.

source

# MixedModels.ranefMethod.

ranef{T}(m::MixedModel{T}, uscale=false)

Returns, as a Vector{Matrix{T}}, the conditional modes of the random effects in model m.

If uscale is true the random effects are on the spherical (i.e. u) scale, otherwise on the original scale.

source

# MixedModels.refit!Method.

refit!{T}(m::LinearMixedModel{T}[, y::Vector{T}])

Refit the model m after installing response y.

If y is omitted the current response vector is used.

source

# MixedModels.rematMethod.

 remat(e::Expr, df::DataFrames.DataFrame)

A factory for ReMat objects.

e should be of the form :(e1 | e2) where e1 is a valid rhs of a Formula and pool(e2) can be evaluated within df. The result is a ScalarReMat or a VectorReMat, as appropriate.

source

# MixedModels.sdestMethod.

sdest(m::LinearMixedModel)

Return the estimate of σ, the standard deviation of the per-observation noise.

source

# MixedModels.setθ!Method.

setθ!{T}(m::LinearMixedModel{T}, v::Vector{T})

Install v as the θ parameters in m. Changes m.Λ only.

source

# MixedModels.simulate!Method.

simulate!(m::LinearMixedModel; β=fixef(m), σ=sdest(m), θ=getθ(m))

Overwrite the response (i.e. m.trms[end]) with a simulated response vector from model m.

source

# MixedModels.varestMethod.

varest(m::LinearMixedModel)

Returns the estimate of σ², the variance of the conditional distribution of Y given B.

source

# StatsBase.fit!Function.

fit!(m::LinearMixedModel[, verbose::Bool=false[, optimizer::Symbol=:LN_BOBYQA]])

Optimize the objective of a LinearMixedModel.

A value for optimizer should be the name of an NLopt derivative-free optimizer allowing for box constraints.

source

# StatsBase.fit!Method.

fit!(m::GeneralizedLinearMixedModel[, verbose = false, fast = false])

Optimize the objective function for m.

When fast is true a potentially much faster but slightly less accurate algorithm, in which pirls! optimizes both the random effects and the fixed-effects parameters, is used.

source

# StatsBase.vcovMethod.

 vcov(m::MixedModel)

Returns the estimated covariance matrix of the fixed-effects estimator.

source