Title: | Missingness in Multi-Task Regression with Network Estimation |
---|---|
Description: | Efficient procedures for fitting conditional graphical lasso models that link a set of predictor variables to a set of response variables (or tasks), even when the response data may contain missing values. 'missoNet' simultaneously estimates the predictor coefficients for all tasks by leveraging information from one another, in order to provide more accurate predictions in comparison to modeling them individually. Additionally, 'missoNet' estimates the response network structure influenced by conditioning predictor variables using a L1-regularized conditional Gaussian graphical model. Unlike most penalized multi-task regression methods (e.g., MRCE), 'missoNet' is capable of obtaining estimates even when the response data is corrupted by missing values. The method automatically enjoys the theoretical and computational benefits of convexity, and returns solutions that are comparable to the estimates obtained without missingness. |
Authors: | Yixiao Zeng [aut, cre, cph], Celia Greenwood [ths, aut], Archer Yang [ths, aut] |
Maintainer: | Yixiao Zeng <[email protected]> |
License: | GPL-2 |
Version: | 1.2.0 |
Built: | 2024-11-08 06:45:23 UTC |
Source: | CRAN |
Package: | missoNet |
Type: | Package |
Version: | 1.0.0 |
Date: | 2022-10-01 |
License: | GPL-2 |
Fit a series of ‘missoNet
’ models with user-supplied regularization parameter pairs for the lasso penalties, {(,
)}.
Perform k-fold cross-validation for ‘missoNet
’ over a grid of (auto-computed) regularization parameter pairs.
S3 method for plotting the cross-validation errors from a fitted 'cv.missoNet'
object.
S3 method for making predictions of response values from a fitted 'cv.missoNet'
object.
Quickly generate synthetic data for simulation studies.
This function performs k-fold cross-validation for ‘missoNet
’. The regularization path is computed
for all possible combinations of values given in the two regularization parameter sequences, namely and
.
‘
cv.missoNet
’ will select the most suitable model among all cross-validated fits along the path.
See the details of ‘missoNet
’ for the model definition.
To help users, the ‘cv.missoNet
’ function is designed to automatically determine the likely ranges
of the regularization parameters over which the cross-validation searches.
cv.missoNet( X, Y, kfold = 5, rho = NULL, lambda.Beta = NULL, lambda.Theta = NULL, lamBeta.min.ratio = NULL, lamTheta.min.ratio = NULL, n.lamBeta = NULL, n.lamTheta = NULL, lamBeta.scale.factor = 1, lamTheta.scale.factor = 1, Beta.maxit = 1000, Beta.thr = 1e-04, eta = 0.8, Theta.maxit = 1000, Theta.thr = 1e-04, eps = 1e-08, penalize.diagonal = TRUE, diag.penalty.factor = NULL, standardize = TRUE, standardize.response = TRUE, fit.1se = FALSE, fit.relax = FALSE, permute = TRUE, with.seed = NULL, parallel = FALSE, cl = NULL, verbose = 1 )
cv.missoNet( X, Y, kfold = 5, rho = NULL, lambda.Beta = NULL, lambda.Theta = NULL, lamBeta.min.ratio = NULL, lamTheta.min.ratio = NULL, n.lamBeta = NULL, n.lamTheta = NULL, lamBeta.scale.factor = 1, lamTheta.scale.factor = 1, Beta.maxit = 1000, Beta.thr = 1e-04, eta = 0.8, Theta.maxit = 1000, Theta.thr = 1e-04, eps = 1e-08, penalize.diagonal = TRUE, diag.penalty.factor = NULL, standardize = TRUE, standardize.response = TRUE, fit.1se = FALSE, fit.relax = FALSE, permute = TRUE, with.seed = NULL, parallel = FALSE, cl = NULL, verbose = 1 )
X |
Numeric predictor matrix ( |
Y |
Numeric response matrix ( |
kfold |
Number of folds for cross-validation – the default is |
rho |
(Optional) A scalar or a numeric vector of length |
lambda.Beta |
(Optional) Numeric vector: a user-supplied sequence of non-negative values for { |
lambda.Theta |
(Optional) Numeric vector: a user-supplied sequence of non-negative values for { |
lamBeta.min.ratio |
The smallest value of |
lamTheta.min.ratio |
The smallest value of |
n.lamBeta |
The number of |
n.lamTheta |
The number of |
lamBeta.scale.factor |
A positive multiplication factor for scaling the entire |
lamTheta.scale.factor |
A positive multiplication factor for scaling the entire |
Beta.maxit |
The maximum number of iterations of the fast iterative shrinkage-thresholding algorithm (FISTA) for updating |
Beta.thr |
The convergence threshold of the FISTA algorithm for updating |
eta |
The backtracking line search shrinkage factor; the default is |
Theta.maxit |
The maximum number of iterations of the ‘ |
Theta.thr |
The convergence threshold of the ‘ |
eps |
A numeric tolerance level for the L1 projection of the empirical covariance matrix; the default is |
penalize.diagonal |
Logical: should the diagonal elements of |
diag.penalty.factor |
Numeric: a separate penalty multiplication factor for the diagonal elements of |
standardize |
Logical: should the columns of |
standardize.response |
Logical: should the columns of |
fit.1se |
Logical: the default is |
fit.relax |
Logical: the default is |
permute |
Logical: should the subject indices for the cross-validation be permuted? The default is |
with.seed |
A random number seed for the permutation. |
parallel |
Logical: the default is |
cl |
A cluster object created by ‘ |
verbose |
Value of |
The ‘cv.missoNet
’ function fits ‘missoNet
’ models ('kfold'
*
'n.lamBeta'
*
'n.lamTheta'
)
times in the whole cross-validation process. That is, for the th-fold (
) computation, the models are fitted at each of
the all (
'n.lamBeta'
*
'n.lamTheta'
) possible combinations of the regularization parameters (,
), with the
th
fold of the training data omitted. The errors are accumulated, and the averaged errors as well as the standard deviations are computed over all folds.
Note that the results of ‘
cv.missoNet
’ are random, since the samples are randomly split into k-folds. Users can eliminate this randomness
by setting 'permute = FALSE'
, or by explicitly assigning a seed to the permutation of samples.
A user-supplied sequence for {} and/or {
} is permitted,
otherwise the program computes an appropriate range of values for the regularization parameters using other control arguments.
Note that ‘
cv.missoNet
’ standardizes 'X'
and 'Y'
to have unit variances before computing its
sequences (and then unstandardizes the resulting coefficients); if you wish to reproduce/compare results with those of other softwares,
it is best to supply pre-standardized
'X'
and 'Y'
. If the algorithm is running slowly, track its progress with 'verbose = 2'
.
In most cases, choosing a sparser grid for the regularization parameters (e.g. smaller 'n.lamBeta'
and/or 'n.lamTheta'
) or setting 'Beta.thr = 1.0E-3'
(or even bigger)
allows the algorithm to make faster progress.
After cross-validation, the regression coefficient matrix and the precision matrix
can be
estimated at three special
pairs, by reapplying ‘
missoNet
’ to the entire training dataset:
"lambda.min
" := (), at which the minimum mean cross-validated error is achieved;
"lambda.1se.Beta
" := (), where
is the largest
at which the mean cross-validated error is within one standard error of the minimum;
"lambda.1se.Theta
" := (), where
is the largest
at which the mean cross-validated error is within one standard error of the minimum.
The corresponding estimates, along with the values, are stored in three separate lists inside the returned object:
'est.min'
, 'est.1se.B'
and 'est.1se.Tht'
('est.1se.B'
and 'est.1se.Tht'
are 'NULL'
unless the argument 'fit.1se = TRUE'
when calling ‘cv.missoNet
’).
The ‘cv.missoNet
’ function returns an R object of S3 class 'cv.missoNet'
for which there are a set of accessory functions available.
The plotting function ‘plot.cv.missoNet
’ can be used to graphically identify the optimal pair of the regularization parameters,
and the prediction function ‘predict.cv.missoNet
’ can be used to make predictions of response values given new input 'X'
.
See the vignette for examples.
This function returns a 'cv.missoNet'
object containing a named 'list'
with all the ingredients of the cross-validated fit:
est.min |
A
|
est.1se.B |
A |
est.1se.Tht |
A |
rho |
A vector of length |
fold.index |
The subject indices identifying which fold each observation is in. |
lambda.Beta.vec |
A flattened vector of length ( |
lambda.Theta.vec |
A flattened vector of length ( |
cvm |
A flattened vector of length ( |
cvup |
Upper cross-validated errors. |
cvlo |
Lower cross-validated errors. |
penalize.diagonal |
Logical: whether the diagonal elements of |
diag.penalty.factor |
The additional penalty multiplication factor for the diagonal elements of |
Yixiao Zeng [email protected], Celia M.T. Greenwood and Archer Yi Yang.
## Simulate a dataset with response values missing completely at random (MCAR), ## the overall missing rate is around 10%. set.seed(123) # reproducibility sim.dat <- generateData(n = 300, p = 50, q = 20, rho = 0.1, missing.type = "MCAR") tr <- 1:240 # training set indices tst <- 241:300 # test set indices X.tr <- sim.dat$X[tr, ] # predictor matrix Y.tr <- sim.dat$Z[tr, ] # corrupted response matrix ## Perform a five-fold cross-validation WITH specified 'lambda.Beta' and 'lambda.Theta'. ## 'standardize' and 'standardize.response' are 'TRUE' by default. lamB.vec <- 10^(seq(from = 0, to = -1, length.out = 5)) lamTht.vec <- 10^(seq(from = 0, to = -1, length.out = 5)) cvfit <- cv.missoNet(X = X.tr, Y = Y.tr, kfold = 5, lambda.Beta = lamB.vec, lambda.Theta = lamTht.vec) ## Perform a five-fold cross-validation WITHOUT specified 'lambda.Beta' and 'lambda.Theta'. ## In this case, a grid of 'lambda.Beta' and 'lambda.Theta' values in a (hopefully) reasonable ## range will be computed and used by the program. ## ## < This simplest command should be a good start for most analyses. > cvfit <- cv.missoNet(X = X.tr, Y = Y.tr, kfold = 5) ## Alternatively, compute the cross-validation folds in parallel on a cluster with 2 cores. ## ## 'fit.1se = TRUE' tells the program to make additional estimations of the parameters at the ## largest value of 'lambda.Beta' / 'lambda.Theta' that gives the most regularized model such ## that the cross-validated error is within one standard error of the minimum. cl <- parallel::makeCluster(min(parallel::detectCores()-1, 2)) cvfit <- cv.missoNet(X = X.tr, Y = Y.tr, kfold = 5, fit.1se = TRUE, parallel = TRUE, cl = cl, permute = TRUE, with.seed = 486) # permute with seed for reproducibility parallel::stopCluster(cl) ## Use PRE-STANDARDIZED training data if you wish to compare the results with other softwares. X.tr.std <- scale(X.tr, center = TRUE, scale = TRUE) Y.tr.std <- scale(Y.tr, center = TRUE, scale = TRUE) cvfit.std <- cv.missoNet(X = X.tr.std, Y = Y.tr.std, kfold = 5, standardize = FALSE, standardize.response = FALSE) ## Plot the (standardized) mean cross-validated errors in a heatmap. plot(cvfit, type = "cv.heatmap") ## Plot the (standardized) mean cross-validated errors in a 3D scatterplot. plot(cvfit, type = "cv.scatter", plt.surf = TRUE) ## Extract the estimates at "lambda.min". Beta.hat1 <- cvfit$est.min$Beta Theta.hat1 <- cvfit$est.min$Theta ## Extract the estimates at "lambda.1se.Beta" (if 'fit.1se' = TRUE). Beta.hat2 <- cvfit$est.1se.B$Beta Theta.hat2 <- cvfit$est.1se.B$Theta ## Extract the estimates at "lambda.1se.Theta" (if 'fit.1se' = TRUE). Beta.hat3 <- cvfit$est.1se.Tht$Beta Theta.hat3 <- cvfit$est.1se.Tht$Theta ## Make predictions of response values on the test set. newy1 <- predict(cvfit, newx = sim.dat$X[tst, ], s = "lambda.min") newy2 <- predict(cvfit, newx = sim.dat$X[tst, ], s = "lambda.1se.Beta") # 'fit.1se' = TRUE newy3 <- predict(cvfit, newx = sim.dat$X[tst, ], s = "lambda.1se.Theta") # 'fit.1se' = TRUE
## Simulate a dataset with response values missing completely at random (MCAR), ## the overall missing rate is around 10%. set.seed(123) # reproducibility sim.dat <- generateData(n = 300, p = 50, q = 20, rho = 0.1, missing.type = "MCAR") tr <- 1:240 # training set indices tst <- 241:300 # test set indices X.tr <- sim.dat$X[tr, ] # predictor matrix Y.tr <- sim.dat$Z[tr, ] # corrupted response matrix ## Perform a five-fold cross-validation WITH specified 'lambda.Beta' and 'lambda.Theta'. ## 'standardize' and 'standardize.response' are 'TRUE' by default. lamB.vec <- 10^(seq(from = 0, to = -1, length.out = 5)) lamTht.vec <- 10^(seq(from = 0, to = -1, length.out = 5)) cvfit <- cv.missoNet(X = X.tr, Y = Y.tr, kfold = 5, lambda.Beta = lamB.vec, lambda.Theta = lamTht.vec) ## Perform a five-fold cross-validation WITHOUT specified 'lambda.Beta' and 'lambda.Theta'. ## In this case, a grid of 'lambda.Beta' and 'lambda.Theta' values in a (hopefully) reasonable ## range will be computed and used by the program. ## ## < This simplest command should be a good start for most analyses. > cvfit <- cv.missoNet(X = X.tr, Y = Y.tr, kfold = 5) ## Alternatively, compute the cross-validation folds in parallel on a cluster with 2 cores. ## ## 'fit.1se = TRUE' tells the program to make additional estimations of the parameters at the ## largest value of 'lambda.Beta' / 'lambda.Theta' that gives the most regularized model such ## that the cross-validated error is within one standard error of the minimum. cl <- parallel::makeCluster(min(parallel::detectCores()-1, 2)) cvfit <- cv.missoNet(X = X.tr, Y = Y.tr, kfold = 5, fit.1se = TRUE, parallel = TRUE, cl = cl, permute = TRUE, with.seed = 486) # permute with seed for reproducibility parallel::stopCluster(cl) ## Use PRE-STANDARDIZED training data if you wish to compare the results with other softwares. X.tr.std <- scale(X.tr, center = TRUE, scale = TRUE) Y.tr.std <- scale(Y.tr, center = TRUE, scale = TRUE) cvfit.std <- cv.missoNet(X = X.tr.std, Y = Y.tr.std, kfold = 5, standardize = FALSE, standardize.response = FALSE) ## Plot the (standardized) mean cross-validated errors in a heatmap. plot(cvfit, type = "cv.heatmap") ## Plot the (standardized) mean cross-validated errors in a 3D scatterplot. plot(cvfit, type = "cv.scatter", plt.surf = TRUE) ## Extract the estimates at "lambda.min". Beta.hat1 <- cvfit$est.min$Beta Theta.hat1 <- cvfit$est.min$Theta ## Extract the estimates at "lambda.1se.Beta" (if 'fit.1se' = TRUE). Beta.hat2 <- cvfit$est.1se.B$Beta Theta.hat2 <- cvfit$est.1se.B$Theta ## Extract the estimates at "lambda.1se.Theta" (if 'fit.1se' = TRUE). Beta.hat3 <- cvfit$est.1se.Tht$Beta Theta.hat3 <- cvfit$est.1se.Tht$Theta ## Make predictions of response values on the test set. newy1 <- predict(cvfit, newx = sim.dat$X[tst, ], s = "lambda.min") newy2 <- predict(cvfit, newx = sim.dat$X[tst, ], s = "lambda.1se.Beta") # 'fit.1se' = TRUE newy3 <- predict(cvfit, newx = sim.dat$X[tst, ], s = "lambda.1se.Theta") # 'fit.1se' = TRUE
The ‘generateData
’ function is used to readily produce synthetic data with randomly/systematically-missing values from a conditional Gaussian graphical model.
This function supports three types of missing mechanisms that can be specified by users – missing completely at random (MCAR), missing at random (MAR) and
missing not at random (MNAR).
generateData( X = NULL, Beta = NULL, E = NULL, Theta = NULL, Sigma.X = NULL, n, p, q, rho, missing.type = "MCAR", Beta.row.sparsity = 0.2, Beta.elm.sparsity = 0.2, with.seed = NULL )
generateData( X = NULL, Beta = NULL, E = NULL, Theta = NULL, Sigma.X = NULL, n, p, q, rho, missing.type = "MCAR", Beta.row.sparsity = 0.2, Beta.elm.sparsity = 0.2, with.seed = NULL )
X |
(Optional) a user-supplied predictor matrix ( |
Beta |
(Optional) a user-supplied regression coefficient matrix |
E |
(Optional) a user-supplied error matrix ( |
Theta |
(Optional) a user-supplied positive definite precision (inverse covariance) matrix |
Sigma.X |
(Optional) A user-supplied positive definite covariance matrix |
n |
Sample size. |
p |
The dimensionality of the predictors. |
q |
The dimensionality of the responses. |
rho |
A scalar or a numeric vector of length |
missing.type |
Character string: can be " |
Beta.row.sparsity |
A Bernoulli parameter between 0 and 1 controlling the approximate proportion of rows where at least one element could be nonzero in |
Beta.elm.sparsity |
A Bernoulli parameter between 0 and 1 controlling the approximate proportion of nonzero elements among the rows of |
with.seed |
A random number seed for the generative process. |
The dataset is simulated through the following steps:
If 'X = NULL'
(default), the function ‘MASS::mvrnorm(n, mean = rep(0, p), sigma = Sigma.X)
’ is used to simulate 'n'
samples from a 'p'
-variate Gaussian distribution for generating a predictor matrix 'X'
;
If 'Beta = NULL'
(default), the function ‘stats::rnorm(p*q, 0, 1)
’ is used to fill an empty () dimensional matrix
'Beta'
, of which the row sparsity and element sparsity are later controlled by the auxiliary arguments 'Beta.row.sparsity'
and 'Beta.elm.sparsity'
, respectively;
If 'E = NULL'
(default), the function ‘MASS::mvrnorm(n, mean = rep(0, q), sigma = solve(Theta))
’ is used to simulate 'n'
samples from a 'q'
-variate Gaussian distribution for generating an error matrix 'E'
;
A complete response matrix 'Y'
without missing values is then generated by the command 'Y = X %*% Beta + E'
;
To get a response matrix 'Z'
:= (
'Y'
) corrupted by missing data, the values in 'Y'
are partially replaced with 'NA'
s following the strategy specified by the arguments 'missing.type'
and 'rho'
.
To better illustrate the step 5 above, suppose for all i = 1,...,n
and j = 1,...,q
: 'Y[i, j]'
is replaced with 'NA'
if 'M[i, j] == 1'
, where 'M'
is an indicator matrix of missingness having the same dimension as 'Y'
.
The value of 'M[i, j]'
is partially controlled by the arguments 'missing.type'
and 'rho'
.
Below we sum up the three built-in missing mechanisms supported by the ‘generateData
’ function:
'missing.type'
== "MCAR
": 'Y[i, j] <- NA'
if 'M[i, j] == 1'
, where 'M[i, j] = rbinom(0, rho[j])'
;
'missing.type'
== "MAR
": 'Y[i, j] <- NA'
if 'M[i, j] == 1'
, where 'M[i, j] = rbinom(0, (rho[j] * c / (1 + exp(-(X %*% Beta)[i, j]))))'
, in which c
is a constant correcting the missing rate of the j
th column of 'Y'
to 'rho[j]'
;
'missing.type'
== "MNAR
": 'Y[i, j] <- NA'
if 'M[i, j] == 1'
, where 'M[i, j] = 1 * (Y[i, j] < Tj)'
, in which 'Tj = quantile(Y[ , j], rho[j])'
.
Of the aforementioned missing mechanisms, "MCAR
" is random, and the other two are systematic.
under "MCAR
", 'M[i, j]'
is not related to 'Y'
or to 'X'
;
under "MAR
", 'M[i, j]'
is related to 'X'
, but not related to 'Y'
after 'X'
is controlled;
under "MNAR
", 'M[i, j]'
is related to 'Y'
itself, even after 'X'
is controlled.
This function returns a 'list'
consisting of the following components:
X |
A simulated (or the user-supplied) predictor matrix ( |
Y |
A simulated response matrix without missing values ( |
Z |
A simulated response matrix with missing values coded as |
Beta |
The regression coefficient matrix |
Theta |
The precision matrix |
rho |
A vector of length |
missing.type |
Character string: the type of missing mechanism used to generate missing values in the response matrix. |
Yixiao Zeng [email protected], Celia M.T. Greenwood and Archer Yi Yang.
## Simulate a dataset with response values missing completely at random (MCAR), ## the overall missing rate is around 10%. sim.dat <- generateData(n = 300, p = 50, q = 20, rho = 0.1, missing.type = "MCAR") ## ------------------------------------------------------------------------------- ## Fit a missoNet model using the simulated dataset. X <- sim.dat$X # predictor matrix Y <- sim.dat$Z # corrupted response matrix fit <- missoNet(X = X, Y = Y, lambda.Beta = 0.1, lambda.Theta = 0.1) ## Simulate a dataset with response values missing at random (MAR), the approximate ## missing rate for each column of the response matrix is specified through a vector 'rho'. ## ## The row sparsity and element sparsity of the auto-generated 'Beta' could be ## adjusted correspondingly by using 'Beta.row.sparsity' and 'Beta.elm.sparsity'. n <- 300; p <- 50; q <- 20 rho <- runif(q, min = 0, max = 0.2) sim.dat <- generateData(n = n, p = p, q = q, rho = rho, missing.type = "MAR", Beta.row.sparsity = 0.3, Beta.elm.sparsity = 0.2) ## Simulate a dataset with response values missing not at random (MNAR), ## using the user-supplied 'Beta' and 'Theta'. n <- 300; p <- 50; q <- 20 Beta <- matrix(rnorm(p*q, 0, 1), p, q) # a nonsparse 'Beta' (p x q) Theta <- diag(q) # a diagonal 'Theta' (q x q) sim.dat <- generateData(Beta = Beta, Theta = Theta, n = n, p = p, q = q, rho = 0.1, missing.type = "MNAR") ## --------------------------------------------------------------------- ## Specifying just one of 'Beta' and 'Theta' is also allowed. sim.dat <- generateData(Theta = Theta, n = n, p = p, q = q, rho = 0.1, missing.type = "MNAR") ## User-supplied 'X', 'Beta' and 'E', in which case 'Y' is deterministic. n <- 300; p <- 50; q <- 20 X <- matrix(rnorm(n*p, 0, 1), n, p) Beta <- matrix(rnorm(p*q, 0, 1), p, q) E <- mvtnorm::rmvnorm(n, rep(0, q), sigma = diag(q)) sim.dat <- generateData(X = X, Beta = Beta, E = E, n = n, p = p, q = q, rho = 0.1, missing.type = "MCAR")
## Simulate a dataset with response values missing completely at random (MCAR), ## the overall missing rate is around 10%. sim.dat <- generateData(n = 300, p = 50, q = 20, rho = 0.1, missing.type = "MCAR") ## ------------------------------------------------------------------------------- ## Fit a missoNet model using the simulated dataset. X <- sim.dat$X # predictor matrix Y <- sim.dat$Z # corrupted response matrix fit <- missoNet(X = X, Y = Y, lambda.Beta = 0.1, lambda.Theta = 0.1) ## Simulate a dataset with response values missing at random (MAR), the approximate ## missing rate for each column of the response matrix is specified through a vector 'rho'. ## ## The row sparsity and element sparsity of the auto-generated 'Beta' could be ## adjusted correspondingly by using 'Beta.row.sparsity' and 'Beta.elm.sparsity'. n <- 300; p <- 50; q <- 20 rho <- runif(q, min = 0, max = 0.2) sim.dat <- generateData(n = n, p = p, q = q, rho = rho, missing.type = "MAR", Beta.row.sparsity = 0.3, Beta.elm.sparsity = 0.2) ## Simulate a dataset with response values missing not at random (MNAR), ## using the user-supplied 'Beta' and 'Theta'. n <- 300; p <- 50; q <- 20 Beta <- matrix(rnorm(p*q, 0, 1), p, q) # a nonsparse 'Beta' (p x q) Theta <- diag(q) # a diagonal 'Theta' (q x q) sim.dat <- generateData(Beta = Beta, Theta = Theta, n = n, p = p, q = q, rho = 0.1, missing.type = "MNAR") ## --------------------------------------------------------------------- ## Specifying just one of 'Beta' and 'Theta' is also allowed. sim.dat <- generateData(Theta = Theta, n = n, p = p, q = q, rho = 0.1, missing.type = "MNAR") ## User-supplied 'X', 'Beta' and 'E', in which case 'Y' is deterministic. n <- 300; p <- 50; q <- 20 X <- matrix(rnorm(n*p, 0, 1), n, p) Beta <- matrix(rnorm(p*q, 0, 1), p, q) E <- mvtnorm::rmvnorm(n, rep(0, q), sigma = diag(q)) sim.dat <- generateData(X = X, Beta = Beta, E = E, n = n, p = p, q = q, rho = 0.1, missing.type = "MCAR")
This function fits the conditional graphical lasso models to datasets with missing response values.
‘missoNet
’ computes the regularization path for the lasso penalties sequentially along the
bivariate regularization parameter sequence provided by the user.
missoNet( X, Y, lambda.Beta, lambda.Theta, rho = NULL, Beta.maxit = 10000, Beta.thr = 1e-08, eta = 0.8, Theta.maxit = 10000, Theta.thr = 1e-08, eps = 1e-08, penalize.diagonal = TRUE, diag.penalty.factor = NULL, standardize = TRUE, standardize.response = TRUE, fit.relax = FALSE, parallel = FALSE, cl = NULL, verbose = 1 )
missoNet( X, Y, lambda.Beta, lambda.Theta, rho = NULL, Beta.maxit = 10000, Beta.thr = 1e-08, eta = 0.8, Theta.maxit = 10000, Theta.thr = 1e-08, eps = 1e-08, penalize.diagonal = TRUE, diag.penalty.factor = NULL, standardize = TRUE, standardize.response = TRUE, fit.relax = FALSE, parallel = FALSE, cl = NULL, verbose = 1 )
X |
Numeric predictor matrix ( |
Y |
Numeric response matrix ( |
lambda.Beta |
A scalar or a numeric vector: a user-supplied sequence of non-negative value(s) for { |
lambda.Theta |
A scalar or a numeric vector: a user-supplied sequence of non-negative value(s) for { |
rho |
(Optional) A scalar or a numeric vector of length |
Beta.maxit |
The maximum number of iterations of the fast iterative shrinkage-thresholding algorithm (FISTA) for updating |
Beta.thr |
The convergence threshold of the FISTA algorithm for updating |
eta |
The backtracking line search shrinkage factor; the default is |
Theta.maxit |
The maximum number of iterations of the ‘ |
Theta.thr |
The convergence threshold of the ‘ |
eps |
A numeric tolerance level for the L1 projection of the empirical covariance matrix; the default is |
penalize.diagonal |
Logical: should the diagonal elements of |
diag.penalty.factor |
Numeric: a separate penalty multiplication factor for the diagonal elements of |
standardize |
Logical: should the columns of |
standardize.response |
Logical: should the columns of |
fit.relax |
Logical: the default is |
parallel |
Logical: the default is |
cl |
A cluster object created by ‘ |
verbose |
Value of |
‘missoNet
’ is the main model-fitting function which is specifically proposed to fit the conditional
graphical lasso models / penalized multi-task Gaussian regressions to (corrupted) datasets with response values missing at random (MAR).
To facilitate the interpretation of the model, let's temporarily assume that there are no missing values
in the data used to fit the model. Suppose we have observations of both a
-variate predictor
and a
-variate response
, for the
th sample (
),
‘
missoNet
’ assumes the model
where and
are one
realization of the
responses and the
predictors, respectively.
is an error vector drawn from a multivariate Gaussian distribution.
The regression coefficient matrix that mapping predictors to responses and
the precision (inverse covariance) matrix
that revealing the
responses' conditional dependencies are the parameters to be estimated by solving a penalized MLE problem
where
The response matrix has
th row (
),
and the predictor matrix
has
th row (
).
The intercept
is canceled out because of centering of the data matrices
and
.
denotes the indicator function for whether penalizing the diagonal elements of
or not.
When
, a global minimizer of the objective function defined above does not exist without the diagonal penalization.
Missingness in real data is inevitable. In this instance, the estimates based only on complete cases are likely to be biased,
and the objective function is likely to no longer be a biconvex optimization problem. In addition, many algorithms cannot be directly employed since they
require complete datasets as inputs. ‘missoNet
’ aims to handle the specific situation where the response matrix contains values that
are missing at random (MAR. Please refer to the vignette or other resources for more information about the differences between MAR, missing completely at
random (MCAR) and missing not at random (MNAR)). As it should be, ‘
missoNet
’ is also applicable to datasets with MCAR response values or without any missing values.
The method provides a unified framework for automatically solving a convex modification of the multi-task learning problem defined above,
using corrupted datasets. Moreover, ‘missoNet
’ enjoys the theoretical and computational benefits of convexity and returns
solutions that are comparable/close to the clean conditional graphical lasso estimates. Please refer to the original manuscript (coming soon) for more details of our method.
This function returns a 'list'
consisting of the following components:
est.list |
A named
|
rho |
A vector of length |
penalize.diagonal |
Logical: whether the diagonal elements of |
diag.penalty.factor |
The additional penalty multiplication factor for the diagonal elements of |
Yixiao Zeng [email protected], Celia M.T. Greenwood and Archer Yi Yang.
## Simulate a dataset with response values missing completely at random (MCAR), ## the overall missing rate is around 10%. set.seed(123) # reproducibility sim.dat <- generateData(n = 300, p = 50, q = 20, rho = 0.1, missing.type = "MCAR") tr <- 1:240 # training set indices tst <- 241:300 # test set indices X.tr <- sim.dat$X[tr, ] # predictor matrix Y.tr <- sim.dat$Z[tr, ] # corrupted response matrix ## Fit one missoNet model with two scalars for 'lambda.Beta' and 'lambda.Theta'. fit1 <- missoNet(X = X.tr, Y = Y.tr, lambda.Beta = 0.1, lambda.Theta = 0.2) ## Fit a series of missoNet models with the lambda pairs := (lambda.Beta, lambda.Theta) ## sequentially extracted from the 'lambda.Beta' and 'lambda.Theta' vectors, note that the ## two vectors must have the same length. lamB.vec <- 10^(seq(from = 0, to = -1, length.out = 5)) lamTht.vec <- rep(0.1, 5) fit2 <- missoNet(X = X.tr, Y = Y.tr, lambda.Beta = lamB.vec, lambda.Theta = lamTht.vec) ## Parallelization on a cluster with two cores. cl <- parallel::makeCluster(2) fit2 <- missoNet(X = X.tr, Y = Y.tr, lambda.Beta = lamB.vec, lambda.Theta = lamTht.vec, parallel = TRUE, cl = cl) parallel::stopCluster(cl) ## Extract the estimates at ('lamB.vec[1]', 'lamTht.vec[1]'). ## The estimates at the subsequent lambda pairs could be accessed in the same way. Beta.hat <- fit2$est.list[[1]]$Beta Theta.hat <- fit2$est.list[[1]]$Theta lambda.Beta <- fit2$est.list[[1]]$lambda.Beta # equal to 'lamB.vec[1]' lambda.Theta <- fit2$est.list[[1]]$lambda.Theta # equal to 'lamTht.vec[1]' ## Fit a series of missoNet models using PRE-STANDARDIZED training data ## if you wish to compare the results with other softwares. X.tr.std <- scale(X.tr, center = TRUE, scale = TRUE) Y.tr.std <- scale(Y.tr, center = TRUE, scale = TRUE) fit3 <- missoNet(X = X.tr.std, Y = Y.tr.std, lambda.Beta = lamB.vec, lambda.Theta = lamTht.vec, standardize = FALSE, standardize.response = FALSE)
## Simulate a dataset with response values missing completely at random (MCAR), ## the overall missing rate is around 10%. set.seed(123) # reproducibility sim.dat <- generateData(n = 300, p = 50, q = 20, rho = 0.1, missing.type = "MCAR") tr <- 1:240 # training set indices tst <- 241:300 # test set indices X.tr <- sim.dat$X[tr, ] # predictor matrix Y.tr <- sim.dat$Z[tr, ] # corrupted response matrix ## Fit one missoNet model with two scalars for 'lambda.Beta' and 'lambda.Theta'. fit1 <- missoNet(X = X.tr, Y = Y.tr, lambda.Beta = 0.1, lambda.Theta = 0.2) ## Fit a series of missoNet models with the lambda pairs := (lambda.Beta, lambda.Theta) ## sequentially extracted from the 'lambda.Beta' and 'lambda.Theta' vectors, note that the ## two vectors must have the same length. lamB.vec <- 10^(seq(from = 0, to = -1, length.out = 5)) lamTht.vec <- rep(0.1, 5) fit2 <- missoNet(X = X.tr, Y = Y.tr, lambda.Beta = lamB.vec, lambda.Theta = lamTht.vec) ## Parallelization on a cluster with two cores. cl <- parallel::makeCluster(2) fit2 <- missoNet(X = X.tr, Y = Y.tr, lambda.Beta = lamB.vec, lambda.Theta = lamTht.vec, parallel = TRUE, cl = cl) parallel::stopCluster(cl) ## Extract the estimates at ('lamB.vec[1]', 'lamTht.vec[1]'). ## The estimates at the subsequent lambda pairs could be accessed in the same way. Beta.hat <- fit2$est.list[[1]]$Beta Theta.hat <- fit2$est.list[[1]]$Theta lambda.Beta <- fit2$est.list[[1]]$lambda.Beta # equal to 'lamB.vec[1]' lambda.Theta <- fit2$est.list[[1]]$lambda.Theta # equal to 'lamTht.vec[1]' ## Fit a series of missoNet models using PRE-STANDARDIZED training data ## if you wish to compare the results with other softwares. X.tr.std <- scale(X.tr, center = TRUE, scale = TRUE) Y.tr.std <- scale(Y.tr, center = TRUE, scale = TRUE) fit3 <- missoNet(X = X.tr.std, Y = Y.tr.std, lambda.Beta = lamB.vec, lambda.Theta = lamTht.vec, standardize = FALSE, standardize.response = FALSE)
S3 method for plotting the cross-validation error surface from a fitted 'cv.missoNet'
object.
## S3 method for class 'cv.missoNet' plot( x, type = c("cv.heatmap", "cv.scatter"), detailed.axes = TRUE, plt.surf = TRUE, ... )
## S3 method for class 'cv.missoNet' plot( x, type = c("cv.heatmap", "cv.scatter"), detailed.axes = TRUE, plt.surf = TRUE, ... )
x |
A fitted |
type |
A character string for the type of plot, can be either " |
detailed.axes |
Logical: whether the detailed axes should be plotted. The default is |
plt.surf |
Logical: whether to draw the error surface. The default is |
... |
Other graphical arguments used by ‘ComplexHeatmap::Heatmap’ ( |
The plot object.
Yixiao Zeng [email protected], Celia M.T. Greenwood and Archer Yi Yang.
## Simulate a dataset. set.seed(123) # reproducibility sim.dat <- generateData(n = 200, p = 10, q = 10, rho = 0.1, missing.type = "MCAR") ## Perform a five-fold cross-validation on the simulated dataset. cvfit <- cv.missoNet(X = sim.dat$X, Y = sim.dat$Z, kfold = 5, fit.1se = TRUE, permute = TRUE, with.seed = 486) ## Plot the (standardized) mean cross-validated errors in a heatmap. plot(cvfit, type = "cv.heatmap") ## Plot the (standardized) mean cross-validated errors in a 3D scatterplot. plot(cvfit, type = "cv.scatter", plt.surf = TRUE)
## Simulate a dataset. set.seed(123) # reproducibility sim.dat <- generateData(n = 200, p = 10, q = 10, rho = 0.1, missing.type = "MCAR") ## Perform a five-fold cross-validation on the simulated dataset. cvfit <- cv.missoNet(X = sim.dat$X, Y = sim.dat$Z, kfold = 5, fit.1se = TRUE, permute = TRUE, with.seed = 486) ## Plot the (standardized) mean cross-validated errors in a heatmap. plot(cvfit, type = "cv.heatmap") ## Plot the (standardized) mean cross-validated errors in a 3D scatterplot. plot(cvfit, type = "cv.scatter", plt.surf = TRUE)
S3 method for making predictions of response values from a fitted 'cv.missoNet'
object.
## S3 method for class 'cv.missoNet' predict(object, newx = NULL, s = "lambda.min", ...)
## S3 method for class 'cv.missoNet' predict(object, newx = NULL, s = "lambda.min", ...)
object |
A fitted |
newx |
A predictor matrix of new values at which predictions are to be made. The columns of |
s |
Character string, the regularization parameter pair |
... |
Not used. Other arguments for predicting. |
The matrix of predicted values: 'newy = mu_hat + newx %*% Beta_hat'
.
Yixiao Zeng [email protected], Celia M.T. Greenwood and Archer Yi Yang.
## Simulate a dataset. set.seed(123) # reproducibility sim.dat <- generateData(n = 300, p = 10, q = 10, rho = 0.1, missing.type = "MCAR") tr <- 1:240 # training set indices tst <- 241:300 # test set indices ## Perform a five-fold cross-validation on the training set. cvfit <- cv.missoNet(X = sim.dat$X[tr, ], Y = sim.dat$Z[tr, ], kfold = 5, fit.1se = TRUE, permute = TRUE, with.seed = 486) ## Make predictions of response values on the test set. newy1 <- predict(cvfit, newx = sim.dat$X[tst, ], s = "lambda.min") newy2 <- predict(cvfit, newx = sim.dat$X[tst, ], s = "lambda.1se.Beta") # 'fit.1se' = TRUE newy3 <- predict(cvfit, newx = sim.dat$X[tst, ], s = "lambda.1se.Theta") # 'fit.1se' = TRUE
## Simulate a dataset. set.seed(123) # reproducibility sim.dat <- generateData(n = 300, p = 10, q = 10, rho = 0.1, missing.type = "MCAR") tr <- 1:240 # training set indices tst <- 241:300 # test set indices ## Perform a five-fold cross-validation on the training set. cvfit <- cv.missoNet(X = sim.dat$X[tr, ], Y = sim.dat$Z[tr, ], kfold = 5, fit.1se = TRUE, permute = TRUE, with.seed = 486) ## Make predictions of response values on the test set. newy1 <- predict(cvfit, newx = sim.dat$X[tst, ], s = "lambda.min") newy2 <- predict(cvfit, newx = sim.dat$X[tst, ], s = "lambda.1se.Beta") # 'fit.1se' = TRUE newy3 <- predict(cvfit, newx = sim.dat$X[tst, ], s = "lambda.1se.Theta") # 'fit.1se' = TRUE