Index
MatrixLMnet.Mlmnet_bic
MatrixLMnet.Mlmnet_cv
MatrixLMnet.admm!
MatrixLMnet.backtransform!
MatrixLMnet.calc_avg_mse
MatrixLMnet.calc_avg_prop_zero
MatrixLMnet.calc_bic
MatrixLMnet.calc_grad
MatrixLMnet.calc_grad!
MatrixLMnet.calc_mse
MatrixLMnet.calc_mse
MatrixLMnet.calc_prop_zero
MatrixLMnet.calc_prop_zero
MatrixLMnet.cd!
MatrixLMnet.cd_active!
MatrixLMnet.coef
MatrixLMnet.coef
MatrixLMnet.coef_3d
MatrixLMnet.criterion
MatrixLMnet.findnotin
MatrixLMnet.fista!
MatrixLMnet.fista_bt!
MatrixLMnet.fitted
MatrixLMnet.fitted
MatrixLMnet.get_func
MatrixLMnet.inner_update_cd!
MatrixLMnet.inner_update_cd!
MatrixLMnet.ista!
MatrixLMnet.lambda_min
MatrixLMnet.make_folds
MatrixLMnet.make_folds_conds
MatrixLMnet.minimize_rows
MatrixLMnet.mlmnet
MatrixLMnet.mlmnet_bic
MatrixLMnet.mlmnet_bic_summary
MatrixLMnet.mlmnet_cv
MatrixLMnet.mlmnet_cv_summary
MatrixLMnet.mlmnet_pathwise
MatrixLMnet.mlmnet_perms
MatrixLMnet.normalize!
MatrixLMnet.outer_update_fista_bt!
MatrixLMnet.predict
MatrixLMnet.predict
MatrixLMnet.println_verbose
MatrixLMnet.prox
MatrixLMnet.prox
MatrixLMnet.prox
MatrixLMnet.prox
MatrixLMnet.prox_mat
MatrixLMnet.prox_mat
MatrixLMnet.prox_mat
MatrixLMnet.resid
MatrixLMnet.resid
MatrixLMnet.update_admm!
MatrixLMnet.update_cd_active_cyclic!
MatrixLMnet.update_cd_active_random!
MatrixLMnet.update_cd_cyclic!
MatrixLMnet.update_cd_random!
MatrixLMnet.update_fista!
MatrixLMnet.update_fista!
MatrixLMnet.update_fista2!
MatrixLMnet.update_fista2!
MatrixLMnet.update_ista!
MatrixLMnet.update_ista!
MatrixLMnet.valid_reduce2
Description
MatrixLMnet.Mlmnet_bic
— TypeMlmnet_bic(MLMNet, lambdas, alphas, data)
Type for storing the results of running BIC validation for mlmnet
MatrixLMnet.Mlmnet_cv
— TypeMlmnet_cv(MLMNets::Array{Mlmnet,1} , lambdas::Array{Float64,1}, alphas::Array{Float64,1}, data::RawData, rowFolds::Array{Array,1} , colFolds::Array{Array,1} )
Type for storing the results of running cross-validation for mlmnet
MatrixLMnet.admm!
— Methodadmm!(X::AbstractArray{Float64,2}, Y::AbstractArray{Float64,2},
Z::AbstractArray{Float64,2}, lambda::Float64, alpha::Float64,
B::AbstractArray{Float64,2},
regXidx::AbstractArray{Int64,1},
regZidx::AbstractArray{Int64,1}, reg::BitArray{2}, norms,
Qx::AbstractArray{Float64,2}, Qz::AbstractArray{Float64,2},
U::AbstractArray{Float64,2}, L::AbstractArray{Float64,2};
isVerbose::Bool=true, stepsize::Float64=0.01,
rho::Float64=1.0, setRho::Bool=true,
thresh::Float64=10.0^(-7), maxiter::Int=10^10,
tau_incr::Float64=2.0, tau_decr::Float64=2.0, mu::Float64=10.0)
Performs ADMM.
Arguments
- X = 2d array of floats consisting of the row covariates, with all categorical variables coded in appropriate contrasts
- Y = 2d array of floats consisting of the multivariate response
- Z = 2d array of floats consisting of the column covariates, with all categorical variables coded in appropriate contrasts
- lambda = lambda penalty, a floating scalar
- alpha = parameter (ϵ[0, 1]) determining the mix of penalties between L1 and L2
- B = 2d array of floats consisting of starting coefficient estimates
- regXidx = 1d array of indices corresponding to regularized X covariates
- regZidx = 1d array of indices corresponding to regularized Z covariates
- reg = 2d array of bits, indicating whether or not to regularize each of the coefficients
- norms = 2d array of floats consisting of the norms corresponding to each coefficient or
nothing
- Qx = 2d array of floats consisting of the eigenvectors of X
- Qz = 2d array of floats consisting of the eigenvectors of Z
- U = 2d array of floats consisting of the transformed Y matrix
- L = 2d array of floats consisting of the kronecker product of the eigenvalues of X and Z
Keyword arguments
- isVerbose = boolean flag indicating whether or not to print messages. Defaults to
true
. - stepsize = float; step size for updates (irrelevant for ADMM). Defaults to
0.01
. - rho = float; parameter that controls ADMM tuning. Defaults to
1.0
. - setRho = boolean flag indicating whether the ADMM tuning parameter
rho
should be calculated. Defaults totrue
. - thresh = threshold at which the coefficients are considered to have converged, a floating scalar. Defaults to
10^(-7)
. - maxiter = maximum number of update iterations. Defaults to
10^10
. - tau_incr = float; parameter that controls the factor at which rho increases. Defaults to 2.0.
- tau_decr = float; parameter that controls the factor at which rho decreases. Defaults to 2.0.
- mu = float; parameter that controls the factor at which the primal and dual residuals should be within each other. Defaults to 10.0.
Value
None; updates coefficients in place
Some notes
Convergence is determined as when the log ratio of the current and previous criteria is less than the threshold thresh
.
rho
controls ADMM tuning and can be specified by the user.
MatrixLMnet.backtransform!
— Methodbacktransform!(B::AbstractArray{Float64,4},
addXIntercept::Bool, addZIntercept::Bool,
meansX::AbstractArray{Float64,2},
meansZ::AbstractArray{Float64,2},
normsX::AbstractArray{Float64,2},
normsZ::AbstractArray{Float64,2})
Back-transform coefficient estimates B in place if X and Z were standardized prior to the estimation– when not including intercept columns for either X or Z.
Arguments
- B = 4d array of coefficient estimates
- addXIntercept = boolean flag indicating whether or not to X has an intercept column
- addZIntercept = boolean flag indicating whether or not to Z has an intercept column
- meansX = 2d array of column means of X, obtained prior to standardizing X
- meansZ = 2d array of column means of Z, obtained prior to standardizing Z
- normsX = 2d array of column norms of X, obtained prior to standardizing X
- normsZ = 2d array of column norms of Z, obtained prior to standardizing Z
Value
None; back-transforms B in place.
Some notes
B is a 4d array in which each coefficient matrix is stored along the third and fourth dimension.
MatrixLMnet.calc_avg_mse
— Methodcalc_avg_mse(MLMNet_cv::Mlmnet_cv)
Calculates average test MSE across folds.
Arguments
- MLMNetcv = MLMNetcv object
Value
2d array of floats
MatrixLMnet.calc_avg_prop_zero
— Methodcalc_avg_prop_zero(MLMNet_cv::Mlmnet_cv)
Calculates average proportion of zero interaction coefficients across folds.
Arguments
- MLMNetcv = MLMNetcv object
Value
1d array of floats
MatrixLMnet.calc_bic
— Methodcalc_bic(MLMNet::Mlmnet)
Calculates BIC for each model according to the lambda-alpha pair parameter.
Arguments
- MLMNets = Mlmnet object resulting from
mlmnet()
function.
Value
2d array of floats with dimensions equal to the number of lambdas by the number of alphas.
MatrixLMnet.calc_grad!
— Methodcalc_grad!(grad::AbstractArray{Float64,2},
X::AbstractArray{Float64,2},
Z::AbstractArray{Float64,2},
resid::AbstractArray{Float64,2})
Calculate gradient in place
Arguments
- gradient = 2d array of floats consisting of the gradient, to be updated in place
- X = 2d array of floats consisting of the row covariates, standardized as necessary
- Z = 2d array of floats consisting of the column covariates, standardized as necessary
- resid = 2d array of floats consisting of the residuals
Value
None; updates gradient in place.
MatrixLMnet.calc_grad
— Methodcalc_grad!(Xi::AbstractArray{Float64,1}, Zj::AbstractArray{Float64,1},
resid::AbstractArray{Float64,2})
Calculate gradient at a single coefficient
Arguments
- Xi = 1d array of floats consisting of the row covariates for the coefficient, standardized as necessary
- Zj = 1d array of floats consisting of the column covariates for the coefficient, standardized as necessary
- resid = 2d array of floats consisting of the residuals
Value
A floating scalar
MatrixLMnet.calc_mse
— Method calc_mse(MLMNets::AbstractArray{Mlmnet,1}, data::RawData,
lambdas::AbstractArray{Float64,1},
alphas::AbstractArray{Float64,1},
rowFolds::Array{Array{Int64,1},1},
colFolds::Array{Array{Int64,1},1})
Calculates test MSE for each of the CV folds for each lambda.
Arguments
- MLMNets = 1d array of Mlmnet objects resulting from running cross validation
- data = RawData object used to generate MLMNets
- lambdas = 1d array of floats consisting of the total penalties in descending order. If they are not in descending order, they will be sorted.
- alphas = 1d array of floats consisting of the penalty ratio that determines the mix of penalties between L1 and L2
- rowFolds = 1d array of arrays containing booleans for the row folds
- colFolds = 1d array of arrays containing booleans for the column folds
Value
2d array of floats with dimensions equal to the number of lambdas by the number of folds.
MatrixLMnet.calc_mse
— Method calc_mse(MLMNet::Mlmnet)
Calculates test MSE for each pair of lambda-alpha.
Arguments
- MLMNet = Mlmnet object
Value
Matrix of floats with dimensions equal to the number of lambdas by the number of alphas.
MatrixLMnet.calc_prop_zero
— Methodcalc_prop_zero(MLMNets::AbstractArray{Mlmnet,1},
lambdas::AbstractArray{Float64,1},
alphas::AbstractArray{Float64,1};
dig::Int64=12)
Calculates proportion of zero interaction coefficients for each of the CV folds for each lambda.
Arguments
- MLMNets = 1d array of Mlmnet objects resulting from running cross validation
- lambdas = 1d array of floats consisting of lambda penalties used to generate MLMNets
Keyword arguments
- dig = integer; digits of precision for zero coefficients. Defaults to 12.
Value
2d array of floats with dimensions equal to the number of lambdas by the number of folds.
MatrixLMnet.calc_prop_zero
— Methodcalc_prop_zero(MLMNet::Mlmnet; dig::Int64=12)
Calculates proportion of zero interaction coefficients for each of the CV folds for each lambda.
Arguments
- MLMNet = Mlmnet object
Keyword arguments
- dig = integer; digits of precision for zero coefficients. Defaults to 12.
Value
Matrix of floats with dimensions equal to the number of lambdas by the number of alphas.
MatrixLMnet.cd!
— Methodcd!(X::AbstractArray{Float64,2}, Y::AbstractArray{Float64,2},
Z::AbstractArray{Float64,2}, lambda::Float64, alpha::Float64,
B::AbstractArray{Float64,2},
regXidx::AbstractArray{Int64,1},
regZidx::AbstractArray{Int64,1}, reg::BitArray{2}, norms;
isVerbose::Bool=true, stepsize::Float64=0.01,
isRandom::Bool=true, thresh::Float64=10.0^(-7),
maxiter::Int=10^10)
Performs coordinate descent using either random or cyclic updates. Does NOT take advantage of the active set; see cd_active!
.
Arguments
- X = 2d array of floats consisting of the row covariates, with all categorical variables coded in appropriate contrasts
- Y = 2d array of floats consisting of the multivariate response
- Z = 2d array of floats consisting of the column covariates, with all categorical variables coded in appropriate contrasts
- lambda = lambda penalty, a floating scalar
- alpha = parameter (ϵ[0, 1]) determining the mix of penalties between L1 and L2
- B = 2d array of floats consisting of starting coefficient estimates
- regXidx = 1d array of indices corresponding to regularized X covariates
- regZidx = 1d array of indices corresponding to regularized Z covariates
- reg = 2d array of bits, indicating whether or not to regularize each of the coefficients
- norms = 2d array of floats consisting of the norms corresponding to each coefficient or
nothing
Keyword arguments
- isVerbose = boolean flag indicating whether or not to print messages. Defaults to
true
. - stepsize = float; step size for updates (irrelevant for coordinate descent). Defaults to
0.01
. - isRandom = boolean flag indicating whether to use random or cyclic updates. Defaults to
true
. - thresh = threshold at which the coefficients are considered to have converged, a floating scalar. Defaults to
10^(-7)
. - maxiter = maximum number of update iterations. Defaults to
10^10
.
Value
None; updates coefficients in place
Some notes
Convergence is determined as when the log ratio of the current and previous criteria is less than the threshold thresh
.
MatrixLMnet.cd_active!
— Methodcd_active!(X::AbstractArray{Float64,2}, Y::AbstractArray{Float64,2},
Z::AbstractArray{Float64,2}, lambda::Float64, alpha::Float64,
B::AbstractArray{Float64,2},
regXidx::AbstractArray{Int64,1},
regZidx::AbstractArray{Int64,1}, reg::BitArray{2}, norms;
isVerbose::Bool=true, stepsize::Float64=0.01,
isRandom::Bool=true, thresh::Float64=10.0^(-7),
maxiter::Int=10^10)
Performs coordinate descent, taking advantage of the active set, using either random or cyclic updates.
Arguments
- X = 2d array of floats consisting of the row covariates, with all categorical variables coded in appropriate contrasts
- Y = 2d array of floats consisting of the multivariate response
- Z = 2d array of floats consisting of the column covariates, with all categorical variables coded in appropriate contrasts
- lambda = lambda penalty, a floating scalar
- alpha = parameter (ϵ[0, 1]) determining the mix of penalties between L1 and L2
- B = 2d array of floats consisting of starting coefficient estimates
- regXidx = 1d array of indices corresponding to regularized X covariates
- regZidx = 1d array of indices corresponding to regularized Z covariates
- reg = 2d array of bits, indicating whether or not to regularize each of the coefficients
- norms = 2d array of floats consisting of the norms corresponding to each coefficient or
nothing
Keyword arguments
- isVerbose = boolean flag indicating whether or not to print messages. Defaults to
true
. - stepsize = float; step size for updates (irrelevant for coordinate descent). Defaults to
0.01
. - isRandom = boolean flag indicating whether to use random or cyclic updates. Defaults to
true
. - thresh = threshold at which the coefficients are considered to have converged, a floating scalar. Defaults to
10^(-7)
. - maxiter = maximum number of update iterations. Defaults to
10^10
.
Value
None; updates coefficients in place
Some notes
Convergence is determined as when the log ratio of the current and previous criteria is less than the threshold thresh
.
MatrixLMnet.coef
— Methodcoef(MLMNet::Mlmnet, lambda::Float64, alpha::Float64)
Extract coefficients from Mlmnet object at a given lambda
Arguments
- MLMNet = Mlmnet object
- lambda = lambda penalty to use, a floating scalar
- alpha = alpha penalty to determine the mix of penalties between L1 and L2 a floating scalar
Value
2d array of coefficients
MatrixLMnet.coef
— Methodcoef(MLMNet::Mlmnet)
Extract all coefficients from Mlmnet object
Arguments
- MLMNet = Mlmnet object
Value
3d array of coefficients
MatrixLMnet.coef_3d
— Methodcoef_3d(MLMNet::Mlmnet)
Extract coefficients from Mlmnet object as a flattened 2d array
Arguments
- MLMNet = Mlmnet object
Value
2d array of flattened coefficients, where each column corresponds to a different lambda and alpha
MatrixLMnet.criterion
— Methodcriterion(B::AbstractArray{Float64,2},
resid::AbstractArray{Float64,2},
lambdaL1::Float64, lambdaL2::Float64, crit_denom::AbstractArray{Int64,1})
Calculate the criterion for the Elastic-net penalty
Arguments
- B = 2d array of floats consisting of regularized coefficient estimates
- resid = 2d array of floats consisting of the residuals
- lambdaL1 = l1 penalty, a floating scalar
- lambdaL2 = l2 penalty, a floating scalar
- crit_denom = 1d array of 2 integers, the denominators of the criterion
Value
A floating scalar
MatrixLMnet.findnotin
— Methodfindnotin(a::AbstractArray{Int64,1}, b::AbstractArray{Int64,1})
Returns elements of b
that are not present in a
.
Arguments
- a = 1d array of integers
- b = 1d array of integers
Value
1d array of integers
MatrixLMnet.fista!
— Methodfista!(X::AbstractArray{Float64,2}, Y::AbstractArray{Float64,2},
Z::AbstractArray{Float64,2},
lambda::Float64, alpha::Float64,
B::AbstractArray{Float64,2},
regXidx::AbstractArray{Int64,1},
regZidx::AbstractArray{Int64,1}, reg::BitArray{2}, norms;
isVerbose::Bool=true, stepsize::Float64=0.01,
thresh::Float64=10.0^(-7), maxiter::Int=10^10)
Performs the Elastic-net version FISTA with fixed step size.
Arguments
- X = 2d array of floats consisting of the row covariates, with all categorical variables coded in appropriate contrasts
- Y = 2d array of floats consisting of the multivariate response
- Z = 2d array of floats consisting of the column covariates, with all categorical variables coded in appropriate contrasts
- lambda = penalty parameter, a floating scalar
- alpha = parameter (ϵ[0, 1]) determining the mix of penalties between L1 and L2
- B = 2d array of floats consisting of starting coefficient estimates
- regXidx = 1d array of indices corresponding to regularized X covariates
- regZidx = 1d array of indices corresponding to regularized Z covariates
- reg = 2d array of bits, indicating whether or not to regularize each of the coefficients
- norms = 2d array of floats consisting of the norms corresponding to each coefficient or
nothing
Keyword arguments
- isVerbose = boolean flag indicating whether or not to print messages. Defaults to
true
. - stepsize = float; step size for updates. Defaults to
0.01
. - thresh = threshold at which the coefficients are considered to have converged, a floating scalar. Defaults to
10^(-7)
. - maxiter = maximum number of update iterations. Defaults to
10^10
.
Value
None; updates coefficients in place
Some notes
Convergence is determined as when the log ratio of the current and previous criteria is less than the threshold thresh
.
The default method for choosing the fixed step size for fista!
is to use the reciprocal of the product of the maximum eigenvalues of X*transpose(X)
and Z*transpose(Z)
. This is computed by the mlmnet
function when fista!
is passed into the fun
argument and setStepSize
is set to true
. If setStepSize
is set to false
, the value of the stepsize
argument will be used as the fixed step size. Note that obtaining the eigenvalues when X
and/or Z
are very large may exceed computational limitations.
MatrixLMnet.fista_bt!
— Methodfista_bt!(X::AbstractArray{Float64,2}, Y::AbstractArray{Float64,2},
Z::AbstractArray{Float64,2}, lambda::Float64, alpha::Float64,
B::AbstractArray{Float64,2},
regXidx::AbstractArray{Int64,1},
regZidx::AbstractArray{Int64,1}, reg::BitArray{2}, norms;
isVerbose::Bool=true, stepsize::Float64=0.01,
gamma::Float64=0.5, thresh::Float64=10.0^(-7),
maxiter::Int=10^10)
Performs FISTA with backtracking.
Arguments
- X = 2d array of floats consisting of the row covariates, with all categorical variables coded in appropriate contrasts
- Y = 2d array of floats consisting of the multivariate response
- Z = 2d array of floats consisting of the column covariates, with all categorical variables coded in appropriate contrasts
- lambda = lambda penalty, a floating scalar
- alpha = parameter (ϵ[0, 1]) determining the mix of penalties between L1 and L2
- B = 2d array of floats consisting of starting coefficient estimates
- regXidx = 1d array of indices corresponding to regularized X covariates
- regZidx = 1d array of indices corresponding to regularized Z covariates
- reg = 2d array of bits, indicating whether or not to regularize each of the coefficients
- norms = 2d array of floats consisting of the norms corresponding to each coefficient or
nothing
Keyword arguments
- isVerbose = boolean flag indicating whether or not to print messages. Defaults to
true
. - stepsize = float; step size for updates. Defaults to
0.01
. - gamma = float; multiplying factor for step size backtracking/line search. Defaults to
0.5
. - thresh = threshold at which the coefficients are considered to have converged, a floating scalar. Defaults to
10^(-7)
. - maxiter = maximum number of update iterations. Defaults to
10^10
.
Value
None; updates coefficients in place
Some notes
Convergence is determined as when the log ratio of the current and previous criteria is less than the threshold thresh
.
Specifying a good starting step size (stepsize
) and multiplying factor (gamma
) in the mlmnet
function when fista_bt!
is passed into the fun
argument can be difficult. Shrinking the step size too gradually can result in slow convergence. Doing so too quickly can cause the criterion to diverge. We have found that setting stepsize
to 0.01 often works well in practice; choice of gamma
appears to be less consequential.
MatrixLMnet.fitted
— Methodfitted(MLMNet::Mlmnet, lambda::Float64, alpha::Float64)
Calculate fitted values of an Mlmnet object, given a lambda
Arguments
- MLMNet = Mlmnet object
- lambda = lambda penalty to use, a floating scalar
- alpha = alpha penalty to determine the mix of penalties between L1 and L2 a floating scalar
Value
2d array of fitted values
MatrixLMnet.fitted
— Methodfitted(MLMNet::Mlmnet)
Calculate fitted values of an Mlmnet object
Arguments
- MLMNet = Mlmnet object
Value
3d array of fitted values
MatrixLMnet.get_func
— Methodget_func(method::String)
Return actual module function name according to method name according to a dictionnary.
Arguments
- method = String describing selected method. The method can be "ista",
"fista", "fista_bt" or "admm".
Value
A function
MatrixLMnet.inner_update_cd!
— Methodinner_update_cd!(i::Int64, j::Int64, B::AbstractArray{Float64,2},
resid::AbstractArray{Float64,2},
X::AbstractArray{Float64,2},
Z::AbstractArray{Float64,2},
norms::AbstractArray{Float64,2}, lambda::Float64,
reg::BitArray{2})
Updates a single coefficient estimate in place (to be called by update_cd_cyclic!
, update_cd_random!
, update_cd_active_cyclic!
, or update_cd_active_random!
) when X
and Z
are not standardized.
Arguments
- i = row index of the coefficient to update
- j = column index of the coefficient to update
- B = 2d array of floats consisting of coefficient estimates
- resid = 2d array of floats consisting of the residuals
- X = 2d array of floats consisting of the row covariates, with all categorical variables coded in appropriate contrasts
- Z = 2d array of floats consisting of the column covariates, with all categorical variables coded in appropriate contrasts
- norms = 2d array of floats consisting of the norms corresponding to each coefficient
- lambda = lambda penalty, a floating scalar
- reg = 2d array of bits, indicating whether or not to regularize each of the coefficients
Value
None; updates a coefficient in place
MatrixLMnet.inner_update_cd!
— Methodinner_update_cd!(i::Int64, j::Int64, B::AbstractArray{Float64,2},
resid::AbstractArray{Float64,2},
X::AbstractArray{Float64,2},
Z::AbstractArray{Float64,2},
norms::Nothing, lambda::Float64, reg::BitArray{2})
Updates a single coefficient estimate in place (to be called by update_cd_cyclic!
, update_cd_random!
, update_cd_active_cyclic!
, or update_cd_active_random!
) when X
and Z
are both standardized.
Arguments
- i = row index of the coefficient to update
- j = column index of the coefficient to update
- B = 2d array of floats consisting of coefficient estimates
- resid = 2d array of floats consisting of the residuals
- X = 2d array of floats consisting of the row covariates, with all categorical variables coded in appropriate contrasts
- Z = 2d array of floats consisting of the column covariates, with all categorical variables coded in appropriate contrasts
- norms =
nothing
- lambda = lambda penalty, a floating scalar
- reg = 2d array of bits, indicating whether or not to regularize each of the coefficients
Value
None; updates a coefficient in place
MatrixLMnet.ista!
— MethodistaNet!(X::AbstractArray{Float64,2}, Y::AbstractArray{Float64,2},
Z::AbstractArray{Float64,2}, lambda::Float64, alpha::Float64,
B::AbstractArray{Float64,2},
regXidx::AbstractArray{Int64,1},
regZidx::AbstractArray{Int64,1}, reg::BitArray{2}, norms;
isVerbose::Bool=true, stepsize::Float64=0.01,
thresh::Float64=10.0^(-7), maxiter::Int=10^10)
Performs the Elastic-net version ISTA with fixed step size.
Arguments
- X = 2d array of floats consisting of the row covariates, with all categorical variables coded in appropriate contrasts
- Y = 2d array of floats consisting of the multivariate response
- Z = 2d array of floats consisting of the column covariates, with all categorical variables coded in appropriate contrasts
- lambda = penalty parameter, a floating scalar
- alpha = parameter (ϵ[0, 1]) determining the mix of penalties between L1 and L2
- B = 2d array of floats consisting of starting coefficient estimates
- regXidx = 1d array of indices corresponding to regularized X covariates
- regZidx = 1d array of indices corresponding to regularized Z covariates
- reg = 2d array of bits, indicating whether or not to regularize each of the coefficients
- norms = 2d array of floats consisting of the norms corresponding to each coefficient or
nothing
Keyword arguments
- isVerbose = boolean flag indicating whether or not to print messages. Defaults to
true
. - stepsize = float; step size for updates. Defaults to
0.01
. - thresh = threshold at which the coefficients are considered to have converged, a floating scalar. Defaults to
10^(-7)
. - maxiter = maximum number of update iterations. Defaults to
10^10
.
Value
None; updates coefficients in place
Some notes
Convergence is determined as when the log ratio of the current and previous criteria is less than the threshold thresh
.
The default method for choosing the fixed step size for ista!
is to use the reciprocal of the product of the maximum eigenvalues of X*transpose(X)
and Z*transpose(Z)
. This is computed by the mlmnet
function when ista!
is passed into the fun
argument and setStepSize
is set to true
. If setStepSize
is set to false
, the value of the stepsize
argument will be used as the fixed step size. Note that obtaining the eigenvalues when X
and/or Z
are very large may exceed computational limitations.
MatrixLMnet.lambda_min
— Methodlambda_min(MLMNet_cv::Mlmnet_cv)
Returns summary information for lambdas corresponding to the minimum average test MSE across folds and the MSE one that is standard error greater.
Arguments
- MLMNetcv = MLMNetcv object
Value
DataFrame from mlmnetcvsummary restricted to the lambdas and alphas that correspond to the minimum average test MSE across folds and the MSE that is one standard error greater.
MatrixLMnet.make_folds
— Functionmake_folds(n::Int64, k::Int64=10, k2::Int64=k)
Generate k
non-overlapping folds.
Arguments
- n = Total number of observations to split into folds.
- k = Number of folds to create. Defaults to 10. If
k=1
, then all the data (along this dimension) will be included in each fold. - k2 = If
k=1
, then all the data (along this dimension) will be included in each fold.k2
specifies how many folds there are. Defaults tok
, which is kind of silly, but there needs to be a placeholder.
Value
1d array of length k
of arrays of indices
MatrixLMnet.make_folds_conds
— Functionmake_folds_conds(conds::AbstractArray{String,1},
k::Int64=10, prop::Float64=1/k)
Generate k
folds for a set of conditions, making sure each level of each condition is represented in each fold.
Arguments
- conds = 1d array of conditions (strings)
- k = Number of folds to create. Defaults to 10.
- prop = Proportion of each condition level's replicates to include in each fold. Defaults to 1/k. Each fold will contain at least one replicate of each condition level.
Value
1d array of length k
of arrays of indices
MatrixLMnet.minimize_rows
— Methodminimize_rows(indices::Vector{CartesianIndex{2}})
Processes a vector of CartesianIndex
objects representing positions in a 2D matrix and returns a new vector of CartesianIndex objects. Each element in the resulting vector should represent the smallest row index for each unique column index.
Arguments
- indices = 1d array of
CartesianIndex
objects representing positions in a 2D matrix
Value
1d array of CartesianIndex
objects representing the smallest row index for each unique column index.
Example
julia> input_indices = [CartesianIndex(1, 1), CartesianIndex(2, 1), CartesianIndex(3, 2), CartesianIndex(1, 2)]
julia> output_indices = minimize_rows(input_indices)
2-element Vector{CartesianIndex{2}}:
CartesianIndex(1, 2)
CartesianIndex(1, 1)
MatrixLMnet.mlmnet
— Methodmlmnet(data::RawData,
lambdas::AbstractArray{Float64,1}, alphas::AbstractArray{Float64,1};
method::String = "ista",
isNaive::Bool=false,
addXIntercept::Bool=true, addZIntercept::Bool=true,
toXReg::BitArray{1}=trues(data.p),
toZReg::BitArray{1}=trues(data.q),
toXInterceptReg::Bool=false, toZInterceptReg::Bool=false,
toNormalize::Bool=true, isVerbose::Bool=true,
stepsize::Float64=0.01, setStepsize::Bool=true,
funArgs...)
Centers and normalizes X and Z predictor matrices, calculates fixed step size, performs the supplied method on two descending lists of lambdas and alphas (each for L1 and L2) using ``warm starts'', and backtransforms resulting coefficients, as is deemed necessary by the user inputs.
Arguments
- data = RawData object
- lambdas = 1d array of floats consisting of the total penalties in descending order. If they are not in descending order, they will be sorted.
- alphas = 1d array of floats consisting of the penalty ratio that determines the mix of penalties between L1 and L2
Keyword arguments
- methods = function name that applies the Elastic-net penalty estimate method; default is
ista
, and the other methods arefista
,fista_bt
,admm
andcd
- isNaive = boolean flag indicating whether to solve the Naive or non-Naive Elastic-net problem
- addXIntercept = boolean flag indicating whether or not to include an
X
intercept (row main effects). Defaults totrue
. - addZIntercept = boolean flag indicating whether or not to include a
Z
intercept (column main effects). Defaults totrue
. - toXReg = 1d array of bit flags indicating whether or not to regularize each of the
X
(row) effects. Defaults to 2d array oftrue
s with length equal to the number ofX
effects (equivalent todata.p
). - toZReg = 1d array of bit flags indicating whether or not to regularize each of the
Z
(column) effects. Defaults to 2d array oftrue
s with length equal to the number ofZ
effects (equivalent todata.q
). - toXInterceptReg = boolean flag indicating whether or not to regularize the
X
intercept Defaults tofalse
. - toZInterceptReg = boolean flag indicating whether or not to regularize the
Z
intercept. Defaults tofalse
. - toNormalize = boolean flag indicating if the columns of
X
andZ
should be standardized (to mean 0, standard deviation 1). Defaults totrue
. - isVerbose = boolean flag indicating whether or not to print messages. Defaults to
true
. - stepsize = float; step size for updates (irrelevant for coordinate descent and when
setStepsize
is set totrue
forista!
andfista!
). Defaults to0.01
. - setStepsize = boolean flag indicating whether the fixed step size should be calculated (for
ista!
andfista!
). Defaults totrue
. - funArgs = variable keyword arguments to be passed into
fun
Value
An Mlmnet object
Some notes
The default method for choosing the fixed step size for fista!
or istaNet!
is to use the reciprocal of the product of the maximum eigenvalues of X*transpose(X)
and Z*transpose(Z)
. This is computed when fista!
or ista!
is passed into the fun
argument and setStepsize
is set to true
. If setStepsize
is set to false
, the value of the stepsize
argument will be used as the fixed step size. Note that obtaining the eigenvalues when X
and/or Z
are very large may exceed computational limitations.
Specifying a good starting step size (stepsize
) and multiplying factor (gamma
) when fista_bt!
is passed into the fun
argument can be difficult. Shrinking the step size too gradually can result in slow convergence. Doing so too quickly can cause the criterion to diverge. We have found that setting stepsize
to 0.01 often works well in practice; choice of gamma
appears to be less consequential.
MatrixLMnet.mlmnet_bic
— Methodmlmnet_bic(data::RawData,
lambdas::AbstractArray{Float64,1},
alphas::AbstractArray{Float64,1};
method::String="ista", isNaive::Bool=false,
addXIntercept::Bool=true, addZIntercept::Bool=true,
toXReg::BitArray{1}=trues(size(get_X(data), 2)),
toZReg::BitArray{1}=trues(size(get_Z(data), 2)),
toXInterceptReg::Bool=false, toZInterceptReg::Bool=false,
toNormalize::Bool=true, isVerbose::Bool=true,
stepsize::Float64=0.01, setStepsize::Bool=true,
dig::Int64=12, funArgs...)
Performs BIC validation for mlmnet
.
Arguments
- data = RawData object
- lambdas = 1d array of floats consisting of the total penalties in descending order. If they are not in descending order, they will be sorted.
- alphas = 1d array of floats consisting of the penalty ratio that determines the mix of penalties between L1 and L2
Keyword arguments
- methods = function name that applies the Elastic-net penalty estimate method; default is
ista
, and the other methods arefista
,fista_bt
,admm
andcd
- isNaive = boolean flag indicating whether to solve the Naive or non-Naive Elastic-net problem
- addXIntercept = boolean flag indicating whether or not to include an
X
intercept (row main effects). Defaults totrue
. - addZIntercept = boolean flag indicating whether or not to include a
Z
intercept (column main effects). Defaults totrue
. - toXReg = 1d array of bit flags indicating whether or not to regularize each of the
X
(row) effects. Defaults to 2d array oftrue
s with length equal to the number ofX
effects (equivalent todata.p
). - toZReg = 1d array of bit flags indicating whether or not to regularize each of the
Z
(column) effects. Defaults to 2d array oftrue
s with length equal to the number ofZ
effects (equivalent todata.q
). - toXInterceptReg = boolean flag indicating whether or not to regularize the
X
intercept Defaults tofalse
. - toZInterceptReg = boolean flag indicating whether or not to regularize the
Z
intercept. Defaults tofalse
. - toNormalize = boolean flag indicating if the columns of
X
andZ
should be standardized (to mean 0, standard deviation 1). Defaults totrue
. - isVerbose = boolean flag indicating whether or not to print messages. Defaults to
true
. - stepsize = float; step size for updates (irrelevant for coordinate descent and when
setStepsize
is set totrue
forista!
andfista!
). Defaults to0.01
. - setStepsize = boolean flag indicating whether the fixed step size should be calculated (for
ista!
andfista!
). Defaults totrue
. - dig = integer; digits of precision for zero coefficients. Defaults to 12.
- funArgs = variable keyword arguments to be passed into
fun
Value
An Mlmnet_bic object.
Some notes
This is the base mlmnet_bic
function that all other variants call.
MatrixLMnet.mlmnet_bic_summary
— Methodmlmnet_bic_summary(MLMNet_bic::Mlmnet_bic)
Summarizes results of BIC-validation by returning a table with:
- Lambdas used
- MSE across folds for each lambda
- Proportion of zero interaction coefficients across each pair of lambda and alpha
Arguments
- MLMNetbic = Mlmnetbic object
Value
DataFrame summarizing BIC, MSE, proportion of zero interactions across each pair of lambda and alpha.
MatrixLMnet.mlmnet_cv
— Methodmlmnet_cvmlmnet_cv(data::RawData,
lambdas::AbstractArray{Float64,1},
alphas::AbstractArray{Float64,1},
rowFolds::Array{Array{Int64,1},1},
colFolds::Array{Array{Int64,1},1};
method::String="ista", isNaive::Bool=false,
addXIntercept::Bool=true, addZIntercept::Bool=true,
toXReg::BitArray{1}=trues(size(get_X(data), 2)),
toZReg::BitArray{1}=trues(size(get_Z(data), 2)),
toXInterceptReg::Bool=false, toZInterceptReg::Bool=false,
toNormalize::Bool=true, isVerbose::Bool=true,
stepsize::Float64=0.01, setStepsize::Bool=true,
dig::Int64=12, funArgs...)
Performs cross-validation for mlmnet
using row and column folds from user input.
Arguments
- data = RawData object
- lambdas = 1d array of floats consisting of the total penalties in descending order. If they are not in descending order, they will be sorted.
- alphas = 1d array of floats consisting of the penalty ratio that determines the mix of penalties between L1 and L2
- rowFolds = 1d array of arrays (one array for each fold), each containing the indices for a row fold; must be same length as colFolds. Can be generated with a call to
make_folds
, which is based onKfold
from the MLBase package. - colFolds = 1d array of arrays (one array for each fold), each containing the indices for a column fold; must be same length as rowFolds. Can be generated with a call to
make_folds
, which is based onKfold
from the MLBase package
Keyword arguments
- methods = function name that applies the Elastic-net penalty estimate method; default is
ista
, and the other methods arefista
,fista_bt
,admm
andcd
- isNaive = boolean flag indicating whether to solve the Naive or non-Naive Elastic-net problem
- addXIntercept = boolean flag indicating whether or not to include an
X
intercept (row main effects). Defaults totrue
. - addZIntercept = boolean flag indicating whether or not to include a
Z
intercept (column main effects). Defaults totrue
. - toXReg = 1d array of bit flags indicating whether or not to regularize each of the
X
(row) effects. Defaults to 2d array oftrue
s with length equal to the number ofX
effects (equivalent todata.p
). - toZReg = 1d array of bit flags indicating whether or not to regularize each of the
Z
(column) effects. Defaults to 2d array oftrue
s with length equal to the number ofZ
effects (equivalent todata.q
). - toXInterceptReg = boolean flag indicating whether or not to regularize the
X
intercept Defaults tofalse
. - toZInterceptReg = boolean flag indicating whether or not to regularize the
Z
intercept. Defaults tofalse
. - toNormalize = boolean flag indicating if the columns of
X
andZ
should be standardized (to mean 0, standard deviation 1). Defaults totrue
. - isVerbose = boolean flag indicating whether or not to print messages. Defaults to
true
. - stepsize = float; step size for updates (irrelevant for coordinate descent and when
setStepsize
is set totrue
forista!
andfista!
). Defaults to0.01
. - setStepsize = boolean flag indicating whether the fixed step size should be calculated (for
ista!
andfista!
). Defaults totrue
. - dig = integer; digits of precision for zero coefficients. Defaults to 12.
- funArgs = variable keyword arguments to be passed into
fun
Value
An Mlmnet_cv object.
Some notes
This is the base mlmnet_cv
function that all other variants call. Folds are computed in parallel when possible.
MatrixLMnet.mlmnet_cv_summary
— Methodmlmnet_cv_summary(MLMNet_cv::Mlmnet_cv)
Summarizes results of cross-validation by returning a table with:
- Lambdas used
- Average MSE across folds for each lambda
- Average proportion of zero interaction coefficeints across folds for each lambda
Arguments
- MLMNetcv = MLMNetcv object
Value
DataFrame summarizing average MSE and proportion of zero interactions across folds for each lambda.
MatrixLMnet.mlmnet_pathwise
— Methodmlmnet_pathwise(fun::Function, X::AbstractArray{Float64,2},
Y::AbstractArray{Float64,2},
Z::AbstractArray{Float64,2},
lambdas::AbstractArray{Float64,1},
alphas::AbstractArray{Float64, 1},
regXidx::AbstractArray{Int64,1},
regZidx::AbstractArray{Int64,1},
reg::BitArray{2}, norms; isVerbose::Bool=true,
stepsize::Float64=0.01, funArgs...)
Performs the supplied method on two descending lists of lambdas (for l1 and l2) using ``warm starts''.
Arguments
- fun = function that applies the Elastic-net pentalty estimate method
- X = 2d array of floats consisting of the row covariates, with all categorical variables coded in appropriate contrasts
- Y = 2d array of floats consisting of the multivariate response
- Z = 2d array of floats consisting of the column covariates, with all categorical variables coded in appropriate contrasts
- lambdas = 1d array of floats consisting of the total penalties in descending order. If they are not in descending order, they will be sorted.
- alphas = 1d array of floats consisting of the penalty ratio that determines the mix of penalties between L1 and L2
- regXidx = 1d array of indices corresponding to regularized X covariates
- regZidx = 1d array of indices corresponding to regularized Z covariates
- reg = 2d array of bits, indicating whether or not to regularize each of the coefficients
- norms = 2d array of floats consisting of the norms corresponding to each coefficient or
nothing
Keyword arguments
- isVerbose = boolean flag indicating whether or not to print messages. Defaults to
true
. - stepsize = float; step size for updates. Defaults to
0.01
. - funArgs = variable keyword arguments to be passed into
fun
Value
A 4d array consisting of the coefficient estimates, with the different lambdas and alphas along the first and second dimensions respectively
Some notes
Assumes that all necessary standardizations have been performed on X, Y, and Z, including adding on intercepts. To be called by mlmnet
, which performs standardization and backtransforming.
MatrixLMnet.mlmnet_perms
— Methodmlmnet_perms(data::RawData,
lambdas::AbstractArray{Float64,1}, alphas::AbstractArray{Float64,1};
method::String = "ista", isNaive::Bool=false,
permFun::Function=shuffle_rows,
addXIntercept::Bool=true, addZIntercept::Bool=true,
toXReg::BitArray{1}=trues(data.p),
toZReg::BitArray{1}=trues(data.q),
toXInterceptReg::Bool=false,
toZInterceptReg::Bool=false,
toNormalize::Bool=true, isVerbose::Bool=true,
stepsize::Float64=0.01, setStepsize=true, funArgs...)
Permutes response matrix Y in RawData object and then calls the mlmnet core function.
Arguments
- data = RawData object
- lambdas = 1d array of floats consisting of the total penalties in descending order. If they are not in descending order, they will be sorted.
- alphas = 1d array of floats consisting of the penalty ratio that determines the mix of penalties between L1 and L2
Keyword arguments
- methods = function name that applies the Elastic-net penalty estimate method; default is
ista
, and the other methods arefista
,fista_bt
,admm
andcd
- isNaive = boolean flag indicating whether to solve the Naive or non-Naive Elastic-net problem
- permFun = function used to permute
Y
. Defaults toshuffle_rows
(shuffles rows ofY
). - addXIntercept = boolean flag indicating whether or not to include an
X
intercept (row main effects). Defaults totrue
. - addZIntercept = boolean flag indicating whether or not to include a
Z
intercept (column main effects). Defaults totrue
. - toXReg = 1d array of bit flags indicating whether or not to regularize each of the
X
(row) effects. Defaults to 2d array oftrue
s with length equal to the number ofX
effects (equivalent todata.p
). - toZReg = 1d array of bit flags indicating whether or not to regularize each of the
Z
(column) effects. Defaults to 2d array oftrue
s with length equal to the number ofZ
effects (equivalent todata.q
). - toXInterceptReg = boolean flag indicating whether or not to regularize the
X
intercept Defaults tofalse
. - toZInterceptReg = boolean flag indicating whether or not to regularize the
Z
intercept. Defaults tofalse
. - toNormalize = boolean flag indicating if the columns of
X
andZ
should be standardized (to mean 0, standard deviation 1). Defaults totrue
. - isVerbose = boolean flag indicating whether or not to print messages. Defaults to
true
. - setStepsize = boolean flag indicating whether the fixed step size should be calculated (for
ista!
andfista!
). Defaults totrue
. - stepsize = float; step size for updates (irrelevant for coordinate descent and when
setStepsize
is set totrue
forista!
andfista!
). Defaults to0.01
. - funArgs = variable keyword arguments to be passed into
fun
Value
An MLMnet object
Some notes
The default method for choosing the fixed step size for fista!
or ista!
is to use the reciprocal of the product of the maximum eigenvalues of X*transpose(X)
and Z*transpose(Z)
. This is computed when fista!
or ista!
is passed into the fun
argument and setStepsize
is set to true
. If setStepsize
is set to false
, the value of the stepsize
argument will be used as the fixed step size. Note that obtaining the eigenvalues when X
and/or Z
are very large may exceed computational limitations.
Specifying a good starting step size (stepsize
) and multiplying factor (gamma
) when fista_bt!
is passed into the fun
argument can be difficult. Shrinking the step size too gradually can result in slow convergence. Doing so too quickly can cause the criterion to diverge. We have found that setting stepsize
to 0.01 often works well in practice; choice of gamma
appears to be less consequential.
MatrixLMnet.normalize!
— Methodnormalize!(A::AbstractArray{Float64,2}, hasIntercept::Bool)
Centers and normalizes the columns of A in place
Arguments
- A = 2d array of floats
- hasIntercept = boolean flag indicating whether the first column of A is the intercept
Value
Centers and normalizes A in place and returns 2d arrays of the column means and L2 norms of A before standardization.
MatrixLMnet.outer_update_fista_bt!
— Methodouter_update_fista_bt!(B::AbstractArray{Float64,2},
B_prev::AbstractArray{Float64,2},
A::AbstractArray{Float64,2},
resid::AbstractArray{Float64,2},
resid_B::AbstractArray{Float64,2},
grad::AbstractArray{Float64,2},
X::AbstractArray{Float64,2},
Y::AbstractArray{Float64,2},
Z::AbstractArray{Float64,2},
norms, lambdaL1::Float64, lambdaL2::Float64,
reg::BitArray{2},
iter::AbstractArray{Int64,1},
stepsize::AbstractArray{Float64,1},
gamma::Float64)
Uses backtracking to update step size for FISTA.
Arguments
- B = 2d array of floats consisting of coefficient estimates
- B_prev = 2d array of floats consisting of coefficient estimates saved from the previous iteration
- A = 2d array of floats consisting of extrapolated coefficients
- resid = 2d array of floats consisting of the residuals calculated from the extrapolated coefficients
- resid_B = 2d array of floats consisting of the residuals calculated from the coefficient estimates
- grad = 2d array of floats consisting of the gradient
- X = 2d array of floats consisting of the row covariates, with all categorical variables coded in appropriate contrasts
- Y = 2d array of floats consisting of the multivariate response observations
- Z = 2d array of floats consisting of the column covariates, with all categorical variables coded in appropriate contrasts
- norms = 2d array of floats consisting of the norms corresponding to each coefficient or
nothing
- lambdaL1 = l1 penalty, a floating scalar
- lambdaL2 = l2 penalty, a floating scalar
- reg = 2d array of bits, indicating whether or not to regularize each of the coefficients
- iter = 1d array consisting of a single integer keeping track of how many iterations have been computed
- stepsize = 1d array consisting of a float; step size of updates
- gamma = float; multiplying factor for step size backtracking/line search
Value
None; updates coefficients and step size in place
MatrixLMnet.predict
— Functionpredict(MLMNet::Mlmnet, newPredictors::Predictors=MLMNet.data.predictors)
Calculate new predictions based on Mlmnet object
Arguments
- MLMNet = Mlmnet object
- newPredictors = Predictors object. Defaults to the data.predictors field in the MLM object used to fit the model.
Value
4d array of predicted values
MatrixLMnet.predict
— Functionpredict(MLMNet::Mlmnet, lambda::Float64, alpha::Float64,
newPredictors::Predictors=MLMNet.data.predictors)
Calculate new predictions based on Mlmnet object and given a lambda
Arguments
- MLMNet = Mlmnet object
- lambda = lambda penalty to use, a floating scalar
- alpha = alpha penalty to determine the mix of penalties between L1 and L2 a floating scalar
- newPredictors = Predictors object. Defaults to the data.predictors field in the MLM object used to fit the model.
Value
2d array of predicted values
MatrixLMnet.println_verbose
— Functionprintln_verbose(x, isVerbose::Bool=true)
Version of println that only prints when isVerbose flag is true
Arguments
- x = something that can be printed
- isVerbose = boolean flag indicating whether or not to print messages. Defaults to
true
.
Value
None; prints to console
MatrixLMnet.prox
— Methodprox(b::Float64, gradient::Float64, b2sign::Float64,
lambda::Float64, norm::Float64)
Proximal (soft-thresholding) operator when step size is 1
Arguments
- b = coefficient to update, a float
- gradient = gradient of b, a float
- b2sign = sign of b + stepsize*gradient, a float
- lambda = lambda penalty , a float
- norm = norm corresponding to b, a float
Value
A floating scalar
MatrixLMnet.prox
— Methodprox(b, gradient, b2sign, lambda, norm, stepsize)
Proximal (soft-thresholding) operator
Arguments
- b = coefficient to update, a float
- gradient = gradient of b, a float
- b2sign = sign of b + stepsize*gradient, a float
- lambda = lambda penalty , a float
- norm = norm corresponding to b, a float
- stepsize = step size to multiply updates, a float
Value
A floating scalar
MatrixLMnet.prox
— Methodprox(b::Float64, gradient::Float64, b2sign::Float64,
lambda::Float64, norm::Nothing, stepsize::Float64)
Proximal (soft-thresholding) operator when not incorporating the norms (norms=1)
Arguments
- b = coefficient to update, a float
- gradient = gradient of b, a float
- b2sign = sign of b + stepsize*gradient, a float
- lambda = lambda penalty , a float
- norm = Nothing
- stepsize = step size to multiply updates, a float
Value
A floating scalar
MatrixLMnet.prox
— Methodprox(b, gradient, b2sign, lambda, norm)
Proximal (soft-thresholding) operator when not incorporating the norms (norms=1) and step size is 1
Arguments
- b = coefficient to update, a float
- gradient = gradient of b, a float
- b2sign = sign of b + stepsize*gradient, a float
- lambda = lambda penalty , a float
- norm = Nothing
Value
A floating scalar
MatrixLMnet.prox_mat
— Methodprox_mat(b::AbstractArray{Float64,2}, b2sign::AbstractArray{Float64,2},
lambda::Float64, norm::AbstractArray{Float64,2}, stepsize::Float64)
Proximal (soft-thresholding) operator
Arguments
- b = coefficient to update, a float
- b2sign = sign of b + stepsize*gradient, a float
- lambda = lambda penalty , a float
- norm = norm corresponding to b, a float
- stepsize = step size to multiply updates, a float
Value
A floating scalar
MatrixLMnet.prox_mat
— Methodprox_mat(b::AbstractArray{Float64,2}, b2sign::AbstractArray{Float64,2},
lambda::Float64, norm::AbstractArray{Float64,2})
Proximal (soft-thresholding) operator
Arguments
- b = coefficient to update, a float
- b2sign = sign of b + stepsize*gradient, a float
- lambda = lambda penalty , a float
- norm = norm corresponding to b, a float
Value
A floating scalar
MatrixLMnet.prox_mat
— Methodprox_mat(b::Float64, gradient::Float64, b2sign::Float64,
lambda::Float64, norm::Nothing, stepsize::Float64)
Proximal (soft-thresholding) operator when not incorporating the norms (norms=1)
Arguments
- b = coefficient to update, a float
- gradient = gradient of b, a float
- b2sign = sign of b + stepsize*gradient, a float
- lambda = lambda penalty , a float
- norm = Nothing
- stepsize = step size to multiply updates, a float
Value
A floating scalar
MatrixLMnet.resid
— Functionresid(MLMNet::Mlmnet, newData::RawData=MLMNet.data)
Calculate residuals of an MLMNet object
Arguments
- MLMNet = Mlmnet object
- newData = RawData object. Defaults to the data field in the MLM object used to fit the model.
Value
3d array of residuals
MatrixLMnet.resid
— Functionresid(MLMNet::Mlmnet, lambda::Float64, alpha::Float64, newData::RawData=MLMNet.data)
Calculate residuals of an MLMNet object, given a lambda
Arguments
- MLMNet = Mlmnet object
- lambda = lambda penalty to use, a floating scalar
- alpha = alpha penalty to determine the mix of penalties between L1 and L2 a floating scalar
- newData = RawData object. Defaults to the data field in the MLM object used to fit the model.
Value
2d array of residuals
MatrixLMnet.update_admm!
— Methodupdate_admm!update_admm!(B::AbstractArray{Float64,2},
B0::AbstractArray{Float64,2},
B2::AbstractArray{Float64,2},
resid::AbstractArray{Float64,2},
X::AbstractArray{Float64,2},
Y::AbstractArray{Float64,2},
Z::AbstractArray{Float64,2},
Qx::AbstractArray{Float64,2},
Qz::AbstractArray{Float64,2},
U::AbstractArray{Float64,2},
L::AbstractArray{Float64,2},
lambdaL1::Float64, lambdaL2::Float64,
regXidx::AbstractArray{Int64,1},
regZidx::AbstractArray{Int64,1},
rho::AbstractArray{Float64,1},
r::AbstractArray{Float64,2},
s::AbstractArray{Float64,2},
tau_incr::Float64, tau_decr::Float64, mu::Float64)
Updates coefficient estimates in place for each ADMM iteration.
Arguments
- B = 2d array of floats consisting of coefficient estimates for L1 updates
- B0 = 2d array of floats consisting of coefficient estimates for L2 updates
- B2 = 2d array of floats consisting of coefficient estimates for dual updates
- resid = 2d array of floats consisting of the residuals
- X = 2d array of floats consisting of the row covariates, with all categorical variables coded in appropriate contrasts
- Y = 2d array of floats consisting of the multivariate response
- Z = 2d array of floats consisting of the column covariates, with all categorical variables coded in appropriate contrasts
- Qx = 2d array of floats consisting of the eigenvectors of X
- Qz = 2d array of floats consisting of the eigenvectors of Z
- U = 2d array of floats consisting of the transformed Y matrix
- L = 2d array of floats consisting of the kronecker product of the eigenvalues of X and Z
- lambdaL1 = l1 penalty, a floating scalar
- lambdaL2 = l2 penalty, a floating scalar
- regXidx = 1d array of indices corresponding to regularized X covariates
- regZidx = 1d array of indices corresponding to regularized Z covariates
- rho = float; parameter that controls ADMM tuning.
- r = 2d array of floats consisting of the primal residuals.
- s = 2d array of floats consisting of the dual residuals.
- tau_incr = float; parameter that controls the factor at which rho increases. Defaults to 2.0.
- tau_decr = float; parameter that controls the factor at which rho decreases. Defaults to 2.0.
- mu = float; parameter that controls the factor at which the primal and dual residuals should be within each other. Defaults to 10.0.
Value
None; updates coefficients in place
Some notes
Convergence is determined as when the log ratio of the current and previous criteria is less than the threshold thresh
.
rho
controls ADMM tuning and can be specified by the user.
MatrixLMnet.update_cd_active_cyclic!
— Methodupdate_cd_active_cyclic!(B::AbstractArray{Float64,2},
resid::AbstractArray{Float64,2},
X::AbstractArray{Float64,2},
Z::AbstractArray{Float64,2},
norms, lambda::Float64, reg::BitArray{2},
nonreg_idx::Tuple{AbstractArray{Int64,1},
AbstractArray{Int64,1}},
active_idx::Tuple{AbstractArray{Int64,1},
AbstractArray{Int64,1}})
Cyclically updates active set of coefficients in place for each coordinate descent iteration.
Arguments
- B = 2d array of floats consisting of coefficient estimates
- resid = 2d array of floats consisting of the residuals
- X = 2d array of floats consisting of the row covariates, with all categorical variables coded in appropriate contrasts
- Z = 2d array of floats consisting of the column covariates, with all categorical variables coded in appropriate contrasts
- norms = 2d array of floats consisting of the norms corresponding to each coefficient or
nothing
- lambda = lambda penalty, a floating scalar
- reg = 2d array of bits, indicating whether or not to regularize each of the coefficients
- nonreg_idx = tuple with the 2d indices of the non-regularized coefficients as two 1d arrays of integers
- active_idx = tuple with the 2d indices of the active coefficients as two 1d arrays of integers
Value
None; updates coefficients in place
Some notes
Given that you pass in the indices for the non-regularized and active (regularized) coefficients separately, this function can be further optimized so that you don't check for regularization when updating each coefficient with inner_update!
.
MatrixLMnet.update_cd_active_random!
— Methodupdate_cd_active_random!(B::AbstractArray{Float64,2},
resid::AbstractArray{Float64,2},
X::AbstractArray{Float64,2},
Z::AbstractArray{Float64,2},
norms, lambda::Float64, reg::BitArray{2},
nonreg_idx::Tuple{AbstractArray{Int64,1},
AbstractArray{Int64,1}},
active_idx::Tuple{AbstractArray{Int64,1},
AbstractArray{Int64,1}})
Randomly updates active set of coefficients in place for each coordinate descent iteration.
Arguments
- B = 2d array of floats consisting of coefficient estimates
- resid = 2d array of floats consisting of the residuals
- X = 2d array of floats consisting of the row covariates, with all categorical variables coded in appropriate contrasts
- Z = 2d array of floats consisting of the column covariates, with all categorical variables coded in appropriate contrasts
- norms = 2d array of floats consisting of the norms corresponding to each coefficient or
nothing
- lambda = lambda penalty, a floating scalar
- reg = 2d array of bits, indicating whether or not to regularize each of the coefficients
- nonreg_idx = tuple with the 2d indices of the non-regularized coefficients as two 1d arrays of integers
- active_idx = tuple with the 2d indices of the active coefficients as two 1d arrays of integers
Value
None; updates coefficients in place
MatrixLMnet.update_cd_cyclic!
— Methodupdate_cd_cyclic!(B::AbstractArray{Float64,2},
resid::AbstractArray{Float64,2},
X::AbstractArray{Float64,2},
Z::AbstractArray{Float64,2},
norms, lambda::Float64, reg::BitArray{2})
Cyclically updates coefficients in place for each coordinate descent iteration.
Arguments
- B = 2d array of floats consisting of coefficient estimates
- resid = 2d array of floats consisting of the residuals
- X = 2d array of floats consisting of the row covariates, with all categorical variables coded in appropriate contrasts
- Z = 2d array of floats consisting of the column covariates, with all categorical variables coded in appropriate contrasts
- norms = 2d array of floats consisting of the norms corresponding to each coefficient or
nothing
- lambda = lambda penalty, a floating scalar
- reg = 2d array of bits, indicating whether or not to regularize each of the coefficients
Value
None; updates coefficients in place
MatrixLMnet.update_cd_random!
— Methodupdate_cd_random!(B::AbstractArray{Float64,2},
resid::AbstractArray{Float64,2},
X::AbstractArray{Float64,2},
Z::AbstractArray{Float64,2},
norms, lambda::Float64, reg::BitArray{2})
Randomly updates coefficients in place for each coordinate descent iteration.
Arguments
- B = 2d array of floats consisting of coefficient estimates
- resid = 2d array of floats consisting of the residuals
- X = 2d array of floats consisting of the row covariates, with all categorical variables coded in appropriate contrasts
- Z = 2d array of floats consisting of the column covariates, with all categorical variables coded in appropriate contrasts
- norms = 2d array of floats consisting of the norms corresponding to each coefficient or
nothing
- lambda = lambda penalty, a floating scalar
- reg = 2d array of bits, indicating whether or not to regularize each of the coefficients
Value
None; updates coefficients in place
MatrixLMnet.update_fista!
— Methodupdate_fista!(B::AbstractArray{Float64,2},
B_prev::AbstractArray{Float64,2},
A::AbstractArray{Float64,2},
resid::AbstractArray{Float64,2},
resid_B::AbstractArray{Float64,2},
grad::AbstractArray{Float64,2},
X::AbstractArray{Float64,2},
Y::AbstractArray{Float64,2},
Z::AbstractArray{Float64,2},
norms::AbstractArray{Float64,2},
lambdaL1::Float64, lambdaL2::Float64,
reg::BitArray{2},
iter::AbstractArray{Int64,1},
stepsize::AbstractArray{Float64,1})
Updates coefficient estimates in place for each FISTA iteration when X
and Z
are not standardized.
Arguments
- B = 2d array of floats consisting of coefficient estimates
- B_prev = 2d array of floats consisting of coefficient estimates saved from the previous iteration
- A = 2d array of floats consisting of extrapolated coefficients
- resid = 2d array of floats consisting of the residuals calculated from the extrapolated coefficients
- resid_B = 2d array of floats consisting of the residuals calculated from the coefficient estimates
- grad = 2d array of floats consisting of the gradient
- X = 2d array of floats consisting of the row covariates, with all categorical variables coded in appropriate contrasts
- Y = 2d array of floats consisting of the multivariate response observations
- Z = 2d array of floats consisting of the column covariates, with all categorical variables coded in appropriate contrasts
- norms = 2d array of floats consisting of the norms corresponding to each coefficient
- lambdaL1 = l1 penalty, a floating scalar
- lambdaL2 = l2 penalty, a floating scalar
- reg = 2d array of bits, indicating whether or not to regularize each of the coefficients
- iter = 1d array consisting of a single integer keeping track of how many iterations have been computed
- stepsize = 1d array consisting of a float; step size of updates
Value
None; updates coefficients in place
MatrixLMnet.update_fista!
— Methodupdate_fista!(B::AbstractArray{Float64,2},
B_prev::AbstractArray{Float64,2},
A::AbstractArray{Float64,2},
resid::AbstractArray{Float64,2},
resid_B::AbstractArray{Float64,2},
grad::AbstractArray{Float64,2},
X::AbstractArray{Float64,2},
Y::AbstractArray{Float64,2},
Z::AbstractArray{Float64,2},
norms::Nothing, lambdaL1::Float64, lambdaL2::Float64,
reg::BitArray{2},
iter::AbstractArray{Int64,1},
stepsize::AbstractArray{Float64,1})
Updates coefficient estimates in place for each FISTA iteration based on the Elastic-net silution, when X
and Z
are both standardized.
Arguments
- B = 2d array of floats consisting of coefficient estimates
- B_prev = 2d array of floats consisting of coefficient estimates saved from the previous iteration
- A = 2d array of floats consisting of extrapolated coefficients
- resid = 2d array of floats consisting of the residuals calculated from the extrapolated coefficients
- resid_B = 2d array of floats consisting of the residuals calculated from the coefficient estimates
- grad = 2d array of floats consisting of the gradient
- X = 2d array of floats consisting of the row covariates, with all categorical variables coded in appropriate contrasts
- Y = 2d array of floats consisting of the multivariate response observations
- Z = 2d array of floats consisting of the column covariates, with all categorical variables coded in appropriate contrasts
- norms =
nothing
- lambdaL1 = l1 penalty, a floating scalar
- lambdaL2 = l2 penalty, a floating scalar
- reg = 2d array of bits, indicating whether or not to regularize each of the coefficients
- iter = 1d array consisting of a single integer keeping track of how many iterations have been computed
- stepsize = 1d array consisting of a float; step size of updates
Value
None; updates coefficients in place
MatrixLMnet.update_fista2!
— Methodupdate_fista2!(B::AbstractArray{Float64,2},
A::AbstractArray{Float64,2},
resid_B::AbstractArray{Float64,2},
grad::AbstractArray{Float64,2},
X::AbstractArray{Float64,2},
Y::AbstractArray{Float64,2},
Z::AbstractArray{Float64,2},
norms::AbstractArray{Float64,2},
lambdaL1::Float64, lambdaL2::Float64,
reg::BitArray{2},
stepsize::AbstractArray{Float64,1})
Updates coefficient estimates in place for each FISTA iteration when X
and Z
are not standardized, but without updating the extrapolated coefficients.
Arguments
- B = 2d array of floats consisting of coefficient estimates
- A = 2d array of floats consisting of extrapolated coefficients
- resid_B = 2d array of floats consisting of the residuals calculated from the coefficient estimates
- grad = 2d array of floats consisting of the gradient
- X = 2d array of floats consisting of the row covariates, with all categorical variables coded in appropriate contrasts
- Y = 2d array of floats consisting of the multivariate response observations
- Z = 2d array of floats consisting of the column covariates, with all categorical variables coded in appropriate contrasts
- norms = 2d array of floats consisting of the norms corresponding to each coefficient
- lambda = lambda penalty, a floating scalar
- reg = 2d array of bits, indicating whether or not to regularize each of the coefficients
- stepsize = 1d array consisting of a float; step size of updates
Value
None; updates coefficients in place
MatrixLMnet.update_fista2!
— Methodupdate_fista2!(B::AbstractArray{Float64,2},
A::AbstractArray{Float64,2},
resid_B::AbstractArray{Float64,2},
grad::AbstractArray{Float64,2},
X::AbstractArray{Float64,2},
Y::AbstractArray{Float64,2},
Z::AbstractArray{Float64,2},
norms::Nothing, lambdaL1::Float64, lambdaL2::Float64,
reg::BitArray{2},
stepsize::AbstractArray{Float64,1})
Updates coefficient estimates in place for each FISTA iteration when X
and Z
are both standardized, but without updating the extrapolated coefficients.
Arguments
- B = 2d array of floats consisting of coefficient estimates
- A = 2d array of floats consisting of extrapolated coefficients
- resid_B = 2d array of floats consisting of the residuals calculated from the coefficient estimates
- grad = 2d array of floats consisting of the gradient
- X = 2d array of floats consisting of the row covariates, with all categorical variables coded in appropriate contrasts
- Y = 2d array of floats consisting of the multivariate response observations
- Z = 2d array of floats consisting of the column covariates, with all categorical variables coded in appropriate contrasts
- norms =
nothing
- lambdaL1 = l1 penalty, a floating scalar
- lambdaL2 = l2 penalty, a floating scalar
- reg = 2d array of bits, indicating whether or not to regularize each of the coefficients
- stepsize = 1d array consisting of a float; step size of updates
Value
None; updates coefficients in place
MatrixLMnet.update_ista!
— Methodupdate_istaNet!(B::AbstractArray{Float64,2},
resid::AbstractArray{Float64,2},
grad::AbstractArray{Float64,2},
X::AbstractArray{Float64,2},
Y::AbstractArray{Float64,2},
Z::AbstractArray{Float64,2},
norms::AbstractArray{Float64,2},
lambdaL1::Float64, lambdaL2::Float64, reg::BitArray{2},
stepsize::AbstractArray{Float64,1})
Updates coefficient estimates in place for each ISTA iteration based on the Elastic-net solution, when X
and Z
are not standardized.
Arguments
- B = 2d array of floats consisting of coefficient estimates
- resid = 2d array of floats consisting of the residuals
- grad = 2d array of floats consisting of the gradient
- X = 2d array of floats consisting of the row covariates, with all categorical variables coded in appropriate contrasts
- Y = 2d array of floats consisting of the multivariate response observations
- Z = 2d array of floats consisting of the column covariates, with all categorical variables coded in appropriate contrasts
- norms = 2d array of floats consisting of the norms corresponding to each coefficient
- lambdaL1 = penalty, a floating scalar
- lambdaL2 = penalty, a floating scalar
- reg = 2d array of bits, indicating whether or not to regularize each of the coefficients
- stepsize = 1d array consisting of a float; step size of updates
Value
None; updates coefficients in place
MatrixLMnet.update_ista!
— Methodupdate_istaNet!(B::AbstractArray{Float64,2},
resid::AbstractArray{Float64,2},
grad::AbstractArray{Float64,2},
X::AbstractArray{Float64,2},
Y::AbstractArray{Float64,2},
Z::AbstractArray{Float64,2},
norms::Nothing, lambdaL1::Float64, lambdaL2::Float64,
reg::BitArray{2},
stepsize::AbstractArray{Float64,1})
Updates coefficient estimates in place for each ISTA iteration based on the Elastic-net solution, when X
and Z
are both standardized.
Arguments
- B = 2d array of floats consisting of coefficient estimates
- resid = 2d array of floats consisting of the residuals
- grad = 2d array of floats consisting of the gradient
- X = 2d array of floats consisting of the row covariates, with all categorical variables coded in appropriate contrasts
- Y = 2d array of floats consisting of the multivariate response observations
- Z = 2d array of floats consisting of the column covariates, with all categorical variables coded in appropriate contrasts
- norms =
nothing
- lambdaL1 = l1 penalty, a floating scalar
- lambdaL1 = l2 penalty, a floating scalar
- reg = 2d array of bits, indicating whether or not to regularize each of the coefficients
- stepsize = 1d array consisting of a float; step size of updates
Value
None; updates coefficients in place
MatrixLMnet.valid_reduce2
— Functionvalid_reduce2(A::Array{Float64,3}, fun::Function=mean)
Reduce a 2d matrix across its columns using a given function, but ignoring NaN, Inf, and -Inf.
Arguments
- A = 2d array of floats
- fun = function with which to reduce across the columns of A
Value
2d array of floats