Index

Description

MatrixLM.MlmType
Mlm(B::Array{Float64,2}, varB::Array{Float64,2}, sigma::Array{Float64,2},    
    data::RawData, weights, targetType, lambda::Float64)

Type for storing the results of an mlm model fit.

source
MatrixLM.PredictorsType
Predictors(X::AbstractArray{Float64,2}, Z::AbstractArray{Float64,2},
           hasXIntercept::Bool, hasZIntercept::Bool)

Type for storing predictor (covariate) matrices. Also stores boolean variables hasXIntercept and hasZIntercept (if they are not supplied, they default to false).

source
MatrixLM.RawDataType
RawData(response::Response, predictors::Predictors)

Type for storing response and predictor matrices

Also stores dimensions of matrices as n, m, p, and q.

  • n : number of rows of X = number of rows of Y
  • m : number of rows of Z = number of columns of Y
  • p : number of columns of X
  • q : number of columns of Z

The constructor will compute n, m, p, and q based on the response and predictor matrices and assert that they are consistent.

source
MatrixLM.add_interceptMethod
add_intercept(A::AbstractArray{Float64,2})

Insert an intercept column (column of ones) at the beginning of a 2d array.

Arguments

  • A::AbstractArray{Float64,2}: 2d array of floats

Value

Returns A with an intercept column

source
MatrixLM.calc_coeffsMethod
calc_coeffs(X::AbstractArray{Float64,2}, Y::AbstractArray{Float64,2}, 
            Z::AbstractArray{Float64,2}, XTX::AbstractArray{Float64,2}, 
            ZTZ::AbstractArray{Float64,2})

Calculates the the coefficient estimates

Arguments

  • X::AbstractArray{Float64,2}: The row covariates, with all categorical variables coded in appropriate contrasts
  • Y::AbstractArray{Float64,2}: The multivariate response
  • Z::AbstractArray{Float64,2}: The column covariates, with all categorical variables coded in appropriate contrasts
  • XTX::AbstractArray{Float64,2}: X*transpose(X) product as a 2d array of floats
  • ZTZ::AbstractArray{Float64,2}: Z*transpose(Z) product as a 2d array of floats

Value

2d array of floats

source
MatrixLM.calc_preds!Method
calc_preds!(preds::AbstractArray{Float64,2}, 
            X::AbstractArray{Float64,2}, 
            Z::AbstractArray{Float64,2}, 
            B::AbstractArray{Float64,2})

Predict values in place

Arguments

  • preds::AbstractArray{Float64,2}: The predicted values, to be updated in place
  • X::AbstractArray{Float64,2}: The row covariates, standardized as necessary
  • Z::AbstractArray{Float64,2}: The column covariates, standardized as necessary
  • B::AbstractArray{Float64,2}: Coefficient estimates

Value

None; updates predicted values in place.

source
MatrixLM.calc_predsMethod
calc_preds(X::AbstractArray{Float64,2}, 
           Z::AbstractArray{Float64,2}, 
           B::AbstractArray{Float64,2})

Predict values

Arguments

  • X::AbstractArray{Float64,2}: The row covariates, standardized as necessary
  • Z::AbstractArray{Float64,2}: The column covariates, standardized as necessary
  • B::AbstractArray{Float64,2}: Coefficient estimates

Value

2d array of floats

source
MatrixLM.calc_resid!Method
calc_resid!(resid::AbstractArray{Float64,2}, 
                 X::AbstractArray{Float64,2}, 
                 Y::AbstractArray{Float64,2}, 
                 Z::AbstractArray{Float64,2}, 
                 B::AbstractArray{Float64,2})

Calculate residuals in place

Arguments

  • resid: 2d array of floats consisting of the residuals, to be updated in place
  • X: 2d array of floats consisting of the row covariates, standardized as necessary
  • Y: 2d array of floats consisting of the multivariate response observations, standardized as necessary
  • Z: 2d array of floats consisting of the column covariates, standardized as necessary
  • B: 2d array of floats consisting of coefficient estimates

Value

None; updates residuals in place.

source
MatrixLM.calc_residMethod
calc_resid(X::AbstractArray{Float64,2}, Y::AbstractArray{Float64,2}, 
       Z::AbstractArray{Float64,2}, B::AbstractArray{Float64,2})

Calculate residuals

Arguments

  • X::AbstractArray{Float64,2}: The row covariates, standardized as necessary
  • Y::AbstractArray{Float64,2}: The multivariate response observations, standardized as necessary
  • Z::AbstractArray{Float64,2}: The column covariates, standardized as necessary
  • B::AbstractArray{Float64,2}: The coefficient estimates

Value

2d array of floats

source
MatrixLM.calc_sigmaMethod
calc_sigma(resid::AbstractArray{Float64,2}, targetType::AbstractString)

Estimates variance of errors and the shrinkage coefficient, with variance shrinkage.

Arguments

  • resid::AbstractArray{Float64,2}: 2d array of floats consisting of the residuals
  • targetType::AbstractString: Indicating the target type toward which to shrink the variance. Acceptable inputs are "A", "B", "C", and "D".
    • "A": Target is identity matrix
    • "B": Target is diagonal matrix with constant diagonal
    • "C": Target is has same diagonal element, and same off-diagonal element
    • "D": Target is diagonal matrix with unequal entries

Value

Tuple

  • sigma: 2d array of floats; shrunk estimated variance of errors
  • lambda: floating scalar; estimated shrinkage coefficient (0 = no shrinkage, 1 = complete shrinkage)
source
MatrixLM.calc_sigmaMethod
calc_sigma(resid::AbstractArray{Float64,2}, targetType::Nothing)

Estimates variance of errors and the shrinkage coefficient, without variance shrinkage.

Arguments

  • resid::AbstractArray{Float64,2}: 2d array of floats consisting of the residuals
  • targetType : nothing

Value

Tuple

  • sigma: 2d array of floats; estimated variance of errors
  • lambda: 0.0

Some notes

Since this version of calc_sigma does not implement variance shrinkage, the shrinkage coefficient lambda is 0.

source
MatrixLM.calc_varMethod
calc_var(X::AbstractArray{Float64,2}, Z::AbstractArray{Float64,2},
         XTX::AbstractArray{Float64,2}, ZTZ::AbstractArray{Float64,2}, 
         sigma::AbstractArray{Float64,2})

Calculate the variance (diagonal of the covariance matrix) of the coefficient estimates.

Arguments

  • X::AbstractArray{Float64,2}: The row covariates, with all categorical variables coded in appropriate contrasts
  • Z::AbstractArray{Float64,2}: The column covariates, with all categorical variables coded in appropriate contrasts
  • XTX::AbstractArray{Float64,2}: X*transpose(X) product as a 2d array of floats
  • ZTZ::AbstractArray{Float64,2}: Z*transpose(Z) product as a 2d array of floats
  • sigma::AbstractArray{Float64,2}: 2d array of floats consisting of the estimated sigma

Value

2d array of floats

source
MatrixLM.centerMethod
center(A::AbstractArray{Float64,2})

Centers columns of a 2d array

Arguments

  • A::AbstractArray{Float64,2}: 2d array of floats

Value

2d array of floats

source
MatrixLM.coefMethod
coef(MLM::Mlm)

Extracts coefficients from Mlm object

Arguments

  • MLM::Mlm: Mlm object

Value

2d array of floats

source
MatrixLM.contrFunction
contr(df::DataFrames.DataFrame, cVars::AbstractArray{Symbol,1}, 
      cTypes::AbstractArray{String,1}=repeat(["treat"], inner=length(cVars)), 
      trtRefs::AbstractArray= repeat([nothing], inner=length(cVars)))

Converts categorical variables in a DataFrame to specified contrast types. All other variables are left as-is.

Arguments

  • df::DataFrames.DataFrame: DataFrame of variables
  • cVar::Symbol: symbol for the categorical variable in df to be converted
  • cTypes::AbstractArray{String,1}: 1d array of character strings of the same length as cVars, indicating the types of contrasts to use. Defaults to treatment contrasts ("treat") for all variables in cVars. Other options include "sum" for sum contrasts, "noint" for treatment contrasts with no intercept, and "sumnoint" for sum contrasts with no intercept. For "treat" cTypes, you can also specify the level to use as the reference treatment using trtRefs.
  • trtRefs::AbstractArray: optional 1d array of character strings of the same length as cVars, specifying the level to use as the references for treatment contrasts. Defaults to nothing for all variables in cVars.

Value

DataFrame with same variables as the original DataFrame, but categorical variables converted to dummy contrasts.

Some notes

If cVars consists of only an empty Symbol, i.e. cVars=[Symbol()], this will signal to the function that no contrasts should be created. The original DataFrame will be returned.

source
MatrixLM.cov_estMethod
cov_est(resid::AbstractArray{Float64,2})

Estimates error variance and its variance/covariance

Arguments

  • resid::AbstractArray{Float64,2}: 2d array of floats consisting of the residuals

Value

Tuple

  • est: 2d array of floats; estimate
  • varest: 2d array of floats; variance/covariance estimate

2d array of floats

source
MatrixLM.design_matrixMethod
design_matrix(frml, df::DataFrame,cntrst::Dict{Symbol, AbstractContrasts})

Returns the design matrix based on the formula terms and the data source.

# Arguments 

- `frml`: formula terms generated by the macro `@mlmformula`
- `df`: dataframe containing the data source table
- `cntrst`: dictionnary describing encoding method for categorical or ordinal variables, based on `StatsModels.jl`

design_matrix(frml, df::DataFrame, cntrst::Vector)

Returns the design matrix based on the formula terms and the data source.

# Arguments 

- `frml`: formula terms generated by the macro `@mlmformula`
- `df`: dataframe containing the data source table
- `cntrst`: A vector containing tuples of variable names and corresponding  encoding function.

design_matrix(frml, df::DataFrame)

Returns the default design matrix based on the formula terms and the data source, where
all the categorical variables are dummy coded.

# Arguments 

- `frml`: formula terms generated by the macro `@mlmformula`
- `df`: dataframe containing the data source table
source
MatrixLM.design_matrix_namesMethod
design_matrix_names(frml, df,cntrst::Dict{Symbol, AbstractContrasts})

Returns the columns names of the design matrix based on the formula terms and the data source.

# Arguments 

- `frml`: formula terms generated by the macro `@mlmformula`
- `df`: dataframe containing the data source table
- `cntrst`: dictionnary describing encoding method for categorical or ordinal variables, based on `StatsModels.jl`

design_matrix_names(frml, df::DataFrame, cntrst::Vector)

Returns the design matrix based on the formula terms and the data source.

# Arguments 

- `frml`: formula terms generated by the macro `@mlmformula`
- `df`: dataframe containing the data source table
- `cntrst`: A vector containing tuples of variable names and corresponding  encoding function.

design_matrix_names(frml, df::DataFrame)

Returns the columns names of the design matrix based on the formula terms and the data source, where
all the categorical variables are dummy coded.

# Arguments 

- `frml`: formula terms generated by the macro `@mlmformula`
- `df`: dataframe containing the data source table
source
MatrixLM.diagonalMethod
diagonal(A::AbstractArray{Float64,2})

Get the diagonal of a 2d array of floats. This just calls the base diag function.

Arguments

  • A::AbstractArray{Float64,2}: 2d array of floats

Value

1d array of floats

source
MatrixLM.diagonalMethod
diagonal(A::AbstractArray{Float64,1})

Get the diagonal of a 1d array of floats. Behaves like an identity function (returns itself).

Arguments

  • A::AbstractArray{Float64,1}: 1d array of floats

Value

1d array of floats

Some notes

Originally intended for use when A is a 1 by 1 array, so may have unintended consequences for a 1d array of length > 1.

source
MatrixLM.diagonalMethod
diagonal(A::Float64)

Get the diagonal of a single scalar (float) value. Behaves like an identity function (returns itself).

Arguments

  • A::Float64: floating scalar

Value

Floating scalar

source
MatrixLM.fittedMethod
fitted(MLM::Mlm)

Calculate fitted values of an Mlm object

Arguments

  • MLM::Mlm: Mlm object

Value

Response object

source
MatrixLM.get_XMethod
get_X(data::RawData)

Extract X matrix from RawData object

Arguments

  • data::RawData: RawData object

Value

2d array

source
MatrixLM.get_YMethod
get_Y(data::RawData)

Extract Y matrix from RawData object

Arguments

  • data::RawData: RawData object

Value

2d array

source
MatrixLM.get_ZMethod
get_Z(data::RawData)

Extract Z matrix from RawData object

Arguments

  • data::RawData: RawData object

Value

2d array

source
MatrixLM.get_dummyMethod
get_dummy(df::DataFrames.DataFrame, cVar::Symbol, 
          cType::String, trtRef::Nothing)

Convert categorical variable to dummy indicators using specified contrast type. This covers all cases except for treatment contrasts with a specified reference level.

Arguments

  • df::DataFrames.DataFrame: DataFrame of variables
  • cVar::Symbol: symbol for the categorical variable in df to be converted
  • cType::String: character string indicating the type of contrast to use for cVar
  • trtRef::Nothing: nothing

Value

DataFrame of dummy variables for the specified categorical variable

source
MatrixLM.get_dummyMethod
get_dummy(df::DataFrames.DataFrame, cVar::Symbol, 
          cType::String, trtRef::String)

Convert categorical variables to for treatment contrasts with a specified reference level.

Arguments

  • df::DataFrames.DataFrame: DataFrame of variables
  • cVar::Symbol: symbol for the categorical variable in df to be converted
  • cType::String: character string indicating the type of contrast to use for cVar
  • trtRef::String: character string specifying the level in cVar to use as the reference

Value

DataFrame of dummy variables for the specified categorical variable

source
MatrixLM.kron_diagMethod
kron_diag(A, B)

Compute the diagonal of the Kronecker product of arrays or scalars

Arguments

  • A: square 2d array of floats, a 1d array of floats, or a scalar
  • B: square 2d array of floats, a 1d array of floats, or a scalar

Value

2d array of floats

source
MatrixLM.mlmMethod
mlm(data::RawData; addXIntercept::Bool=true, addZIntercept::Bool=true, weights=nothing, targetType=nothing)

Matrix linear model using least squares method. Column weighted least squares and shrinkage of the variance of the errors are options.

Arguments

  • data::RawData: RawData object

Keyword arguments

  • addXIntercept::Bool : boolean flag indicating whether or not to include an X intercept (row main effects). Defaults to true.
  • addZIntercept::Bool : boolean flag indicating whether or not to include a Z intercept (column main effects). Defaults to true.
  • weights : 1d array of floats to use as column weights for Y, or nothing. If the former, must be the same length as the number of columns of Y. Defaults to nothing.
  • targetType : string indicating the target type toward which to shrink the error variance, or nothing. If the former, acceptable inputs are "A", "B", "C", and "D". Defaults to nothing.
    • "A": Target is identity matrix
    • "B": Target is diagonal matrix with constant diagonal
    • "C": Target is has same diagonal element, and same off-diagonal element
    • "D": Target is diagonal matrix with unequal entries

Value

An Mlm object

source
MatrixLM.mlm_fitMethod
mlm_fit(data::RawData, weights::Nothing, targetType)

Matrix linear model using least squares method. Optionally incorporates shrinkage of the variance of the errors.

Arguments

  • data::RawData: RawData object
  • weights::Nothing: nothing
  • targetType: string indicating the target type toward which to shrink the error variance, or nothing. If the former, acceptable inputs are "A", "B", "C", and "D".
    • "A": Target is identity matrix
    • "B": Target is diagonal matrix with constant diagonal
    • "C": Target is has same diagonal element, and same off-diagonal element
    • "D": Target is diagonal matrix with unequal entries

Value

An Mlm object

source
MatrixLM.mlm_fitMethod
mlm_fit(data::RawData, weights::Array{Float64,1}, targetType)

Matrix linear model using column weighted least squares method. Optionally incorporates shrinkage of the variance of the errors.

Arguments

  • data::RawData : RawData object
  • weights::Array{Float64,1} : 1d array of floats to use as column weights for Y. Must be the same length as the number of columns of Y.
  • targetType : string indicating the target type toward which to shrink the error variance, or nothing. If the former, acceptable inputs are "A", "B", "C", and "D".
    • "A": Target is identity matrix
    • "B": Target is diagonal matrix with constant diagonal
    • "C": Target is has same diagonal element, and same off-diagonal element
    • "D": Target is diagonal matrix with unequal entries

Value

An Mlm object

source
MatrixLM.mlm_permsFunction
mlm_perms(data::RawData, nPerms::Int64=1000; 
          permFun::Function=shuffle_rows, 
          addXIntercept::Bool=true, addZIntercept::Bool=true, 
          weights=nothing, targetType=nothing, isMainEff::Bool=false)

Obtains permutation p-values for MLM t-statistics.

Arguments

  • data::RawData: RawData object
  • nPerms::Int64=1000: Number of permutations. Defaults to 1000.

Keyword arguments

  • permFun::Function: function used to permute Y. Defaults to shuffle_rows (shuffles rows of Y).
  • addXIntercept::Bool=true: Boolean flag indicating whether or not to include an X intercept (row main effects). Defaults to true.
  • addZIntercept::Bool=true: Boolean flag indicating whether or not to include a Z intercept (column main effects). Defaults to true.
  • weights: 1d array of floats to use as column weights for Y, or nothing. If the former, must be the same length as the number of columns of Y. Defaults to nothing.
  • targetType: string indicating the target type toward which to shrink the error variance, or nothing. If the former, acceptable inputs are "A", "B", "C", and "D". Defaults to nothing.
    • "A": Target is identity matrix
    • "B": Target is diagonal matrix with constant diagonal
    • "C": Target is has same diagonal element, and same off-diagonal element
    • "D": Target is diagonal matrix with unequal entries
  • isMainEff::Bool: boolean flag indicating whether or not to include p-values for the main effects

Value

Tuple

  • tStats: 2d array of floats; t-statistics
  • pvals: 2d array of floats; permutation p-values

Some notes

Permutations are computed in parallel when possible.

source
MatrixLM.perm_pvalsFunction
perm_pvals(fun::Function, data::RawData, nPerms::Int64=1000; 
           permFun::Function=shuffle_rows, funArgs...)

Obtains permutation p-values.

Arguments

  • fun::Function: function that returns a test statistic
  • data::RawData: RawData object
  • nPerms::Int64: number of permutations. Defaults to 1000.

Keyword arguments

  • permFun::Function: function used to permute Y. Defaults to shuffle_rows (shuffles rows of Y).
  • funArgs: variable keyword arguments to be passed into fun

Value

Tuple

  • testStats: 2d array of floats; t-statistics
  • pvals: 2d array of floats; permutation p-values

Some notes

Permutations are computed in parallel when possible.

source
MatrixLM.predictFunction
predict(MLM::Mlm, newPredictors::Predictors=MLM.data.predictors)

Calculates new predictions based on Mlm object

Arguments

  • MLM::Mlm: Mlm object
  • newPredictors::Predictors: Predictors object. Defaults to the data.predictors field in the Mlm object used to fit the model.

Value

Response object

source
MatrixLM.remove_interceptMethod
remove_intercept(A::AbstractArray{Float64,2})

Remove the intercept column, assumed to be the first column of a 2d array.

Arguments

  • A::AbstractArray{Float64,2}: 2d array of floats

Value

Returns A without the intercept column

source
MatrixLM.residFunction
resid(MLM::Mlm, newData::RawData=MLM.data)

Calculates residuals of an Mlm object

Arguments

  • MLM::Mlm: Mlm object
  • newData::RawData: RawData object. Defaults to the data field in the Mlm object used to fit the model.

Value

2d array of floats

source
MatrixLM.shrink_sigmaMethod
shrink_sigma(resid::AbstractArray{Float64,2}, targetType::String)

Estimates variance of errors and the shrinkage coefficient

Arguments

  • resid::AbstractArray{Float64,2}: 2d array of floats consisting of the residuals
  • targetType::String: string indicating the target type toward which to shrink the variance. Acceptable inputs are "A", "B", "C", and "D".
    • "A": Target is identity matrix
    • "B": Target is diagonal matrix with constant diagonal
    • "C": Target is has same diagonal element, and same off-diagonal element
    • "D": Target is diagonal matrix with unequal entries

Value

Tuple

  • sigma: 2d array of floats; shrunk estimated variance of errors
  • lambda: floating scalar; estimated shrinkage coefficient (0 = no shrinkage, 1 = complete shrinkage)

Reference

Ledoit, O., & Wolf, M. (2003). Improved estimation of the covariance matrix of stock returns with an application to portfolio selection. Journal of empirical finance, 10(5), 603-621.

source
MatrixLM.shuffle_colsMethod
shuffle_cols(A::AbstractArray{Float64,2})

Shuffles columns of a 2d array

Arguments

  • A::AbstractArray{Float64,2}: 2d array of floats

Value

Returns A with columns shuffled

source
MatrixLM.shuffle_rowsMethod
shuffle_rows(A::AbstractArray{Float64,2})

Shuffles rows of a 2d array

Arguments

  • A::AbstractArray{Float64,2}: 2d array of floats

Value

Returns A with rows shuffled

source
MatrixLM.t_statFunction
t_stat(MLM::Mlm, isMainEff::Bool=false)

Calculates t-statistics of an Mlm object

Arguments

  • MLM::Mlm: Mlm object
  • isMainEff::Bool : boolean flag indicating whether or not to include t-statistics for the main effects

Value

2d array of floats

source
MatrixLM.@mlmformulaMacro
mlmformula(ex)

Capture and parse a formula expression for matrix linear model.

The `@mlmformula` domain-specific language serves the purpose of facilitating table-to-matrix transformations.
It is structured to be intuitive for users who have experience with other statistical software.
An elementary formula in this language consists of individual terms. These terms may either be symbols that reference
data columns or literal numbers `0` or `1`. They are combined by the operators `+`, `&`, and `*`. 
To ensure correct parsing of the formula, the `@mlmformula`` macro needs to be invoked within parentheses. 
This macro is built upon the `@formula` macro from the `StatsModels.jl` package.

# Example
```julia
julia> @mlmformula(1 + varA * VarB)
1
varA(unknown)
VarB(unknown)
varA(unknown) & VarB(unknown)
```
source