Title: | Ridge Group Sparse Optimization Problem for Estimation of a Meta Model Based on Reproducing Kernel Hilbert Spaces |
---|---|
Description: | Estimates the Hoeffding decomposition of an unknown function by solving ridge group sparse optimization problem based on reproducing kernel Hilbert spaces, and approximates its sensitivity indices (see Kamari, H., Huet, S. and Taupin, M.-L. (2019) <arXiv:1905.13695>). |
Authors: | Halaleh Kamari |
Maintainer: | Halaleh Kamari <[email protected]> |
License: | GPL (>= 2) |
Version: | 1.1 |
Built: | 2024-11-03 06:44:41 UTC |
Source: | CRAN |
Fits a meta model to an unknown model by solving the ridge group sparse (or group lasso) optimization problem based on the reproducing kernel Hilbert spaces (RKHS), for the Gaussian regression model :
where variables are independent and uniformly distributed on
and are independent of
's.
We define the ridge group sparse criteria by :
and the group lasso criteria is obtained by setting in the criteria above. We set
to be the group lasso penalty parameter.
For each pair of the penalty parameters in the ridge group sparse criteria, one meta model, called RKHS meta model, is calculated. The RKHS meta model is an additive model with at most vMax groups. It satisfies the properties of the Hoeffding decomposition, and its terms estimate the terms in the Hoeffding decomposition of the function
.
These estimators are evaluated using a testing dataset. That is, the prediction error is calculated for each RKHS meta model and the one with the minimum prediction error is the "best" estimator for the true model . It provides a function that estimates the empirical sensitivity indices of the "best" RKHS meta model as an approximation of the true sensitivity indices.
Details.
Halaleh Kamari.
Maintainer: [email protected]
Kamari, H., Huet, S. and Taupin, M.-L. (2019) RKHSMetaMod : An R package to estimate the Hoeffding decomposition of an unknown function by solving RKHS Ridge Group Sparse optimization problem. <arXiv:1905.13695>
d <- 3 n <- 50 library(lhs) X <- maximinLHS(n, d) c <- c(0.2,0.6,0.8) F <- 1;for (a in 1:d) F <- F*(abs(4*X[,a]-2)+c[a])/(1+c[a]) epsilon <- rnorm(n,0,1);sigma <- 0.2 Y <- F + sigma*epsilon Dmax <- 3 kernel <- "matern" frc <- c(10,100) gamma <- c(.5,.01,.001,0) result <- RKHSMetMod(Y,X,kernel,Dmax,gamma,frc,FALSE)
d <- 3 n <- 50 library(lhs) X <- maximinLHS(n, d) c <- c(0.2,0.6,0.8) F <- 1;for (a in 1:d) F <- F*(abs(4*X[,a]-2)+c[a])/(1+c[a]) epsilon <- rnorm(n,0,1);sigma <- 0.2 Y <- F + sigma*epsilon Dmax <- 3 kernel <- "matern" frc <- c(10,100) gamma <- c(.5,.01,.001,0) result <- RKHSMetMod(Y,X,kernel,Dmax,gamma,frc,FALSE)
Calculates the Gram matrices for
vMax, and returns their associated eigenvalues and eigenvectors. The calculated Gram matrices may be not positive definite. The option "correction" of this function allows to replace the matrices
that are not positive definite by their "nearest positive definite" matrices.
calc_Kv(X, kernel, Dmax, correction, verbose, tol)
calc_Kv(X, kernel, Dmax, correction, verbose, tol)
X |
Matrix of observations with |
kernel |
Character, the type of the reproducing kernel: matern |
Dmax |
Integer, between |
correction |
Logical, if TRUE, the program makes the correction to the matrices |
verbose |
Logical, if TRUE, the group |
tol |
Scalar, used if correction is TRUE. For each matrix |
Let be the eigenvalues associated with matrix
. Set
and
. The eigenvalues of
that is not positive definite are replaced by
epsilon, with espilon
tol. The value of tol depends on the type of the kernel and it is chosen small.
List of two components "names.Grp" and "kv":
names.Grp |
Vector of size vMax, indicates the name of groups included in the meta model. |
kv |
List of vMax components with the same names as the vector names.Grp. Each element of the list is a list of two components "Evalues" and "Q": |
Evalues |
Vector of size |
Q |
Matrix with |
Note.
Halaleh Kamari
Kamari, H., Huet, S. and Taupin, M.-L. (2019) RKHSMetaMod : An R package to estimate the Hoeffding decomposition of an unknown function by solving RKHS Ridge Group Sparse optimization problem. <arXiv:1905.13695>
d <- 3 n <- 50 library(lhs) X <- maximinLHS(n, d) c <- c(0.2,0.6,0.8) F <- 1;for (a in 1:d) F <- F*(abs(4*X[,a]-2)+c[a])/(1+c[a]) epsilon <- rnorm(n,0,1);sigma <- 0.2 Y <- F + sigma*epsilon Dmax <- 3 kernel <- "matern" Kv <- calc_Kv(X, kernel, Dmax) names <- Kv$names.Grp Eigen.val1 <- Kv$kv$v1.$Evalues Eigen.vec1 <- Kv$kv$v1.$Q
d <- 3 n <- 50 library(lhs) X <- maximinLHS(n, d) c <- c(0.2,0.6,0.8) F <- 1;for (a in 1:d) F <- F*(abs(4*X[,a]-2)+c[a])/(1+c[a]) epsilon <- rnorm(n,0,1);sigma <- 0.2 Y <- F + sigma*epsilon Dmax <- 3 kernel <- "matern" Kv <- calc_Kv(X, kernel, Dmax) names <- Kv$names.Grp Eigen.val1 <- Kv$kv$v1.$Evalues Eigen.vec1 <- Kv$kv$v1.$Q
Fits a solution of the group lasso problem based on RKHS, with active groups in the obtained solution for the Gaussian regression model. It determines
, for which the number of active groups in the solution of the RKHS group lasso problem is equal to
, and returns the RKHS meta model associated with
.
grplasso_q(Y, Kv, q, rat, Num)
grplasso_q(Y, Kv, q, rat, Num)
Y |
Vector of response observations of size |
Kv |
List of eigenvalues and eigenvectors of positive definite Gram matrices |
q |
Integer, the number of active groups in the obtained solution. |
rat |
Positive scalar, used to restrict the minimum value of |
Num |
Integer, used to restrict the number of different values of the penalty parameter |
Input Kv should contain the eigenvalues and eigenvectors of positive definite Gram matrices . It is necessary to set input "correction" in the function
calc_Kv
equal to "TRUE".
List of components: "mus", "qs", "mu", "res":
mus |
Vector, values of the evaluated penalty parameters |
qs |
Vector, number of active groups associated with each value of |
mu |
Scalar, value of |
res |
An RKHS Group Lasso object: |
intercept |
Scalar, estimated value of intercept. |
teta |
Matrix with vMax rows and |
fit.v |
Matrix with |
fitted |
Vector of size |
Norm.H |
Vector of size vMax, estimated values of the penalty norm. |
supp |
Vector of active groups. |
Nsupp |
Vector of the names of the active groups. |
SCR |
Scalar, equals to |
crit |
Scalar, indicates the value of the penalized criteria. |
MaxIter |
Integer, number of iterations until convergence is reached. |
convergence |
TRUE or FALSE. Indicates whether the algorithm has converged or not. |
RelDiffCrit |
Scalar, value of the first convergence criteria at the last iteration, |
RelDiffPar |
Scalar, value of the second convergence criteria at the last iteration, |
Note.
Halaleh Kamari
Kamari, H., Huet, S. and Taupin, M.-L. (2019) RKHSMetaMod : An R package to estimate the Hoeffding decomposition of an unknown function by solving RKHS Ridge Group Sparse optimization problem. <arXiv:1905.13695>
d <- 3 n <- 50 library(lhs) X <- maximinLHS(n, d) c <- c(0.2,0.6,0.8) F <- 1;for (a in 1:d) F <- F*(abs(4*X[,a]-2)+c[a])/(1+c[a]) epsilon <- rnorm(n,0,1);sigma <- 0.2 Y <- F + sigma*epsilon Dmax <- 3 kernel <- "matern" Kv <- calc_Kv(X, kernel, Dmax, TRUE, TRUE) result <- grplasso_q(Y,Kv,5,100 ,Num=10) result$mu result$res$Nsupp
d <- 3 n <- 50 library(lhs) X <- maximinLHS(n, d) c <- c(0.2,0.6,0.8) F <- 1;for (a in 1:d) F <- F*(abs(4*X[,a]-2)+c[a])/(1+c[a]) epsilon <- rnorm(n,0,1);sigma <- 0.2 Y <- F + sigma*epsilon Dmax <- 3 kernel <- "matern" Kv <- calc_Kv(X, kernel, Dmax, TRUE, TRUE) result <- grplasso_q(Y,Kv,5,100 ,Num=10) result$mu result$res$Nsupp
Calculates the value of the penalty parameter in the RKHS group lasso problem when the first penalized parameter group enters the model.
mu_max(Y, matZ)
mu_max(Y, matZ)
Y |
Vector of response observations of size |
matZ |
List of vMax components. Each component includes the eigenvalues and eigenvectors of the positive definite Gram matrices |
Details.
An object of type numeric is returned.
Note.
Halaleh Kamari
Kamari, H., Huet, S. and Taupin, M.-L. (2019) RKHSMetaMod : An R package to estimate the Hoeffding decomposition of an unknown function by solving RKHS Ridge Group Sparse optimization problem. <arXiv:1905.13695>
Meier, L. Van de Geer, S. and Buhlmann, P. (2008) The group LASSO for logistic regression. Journal of the Royal Statistical Society Series B. 70. 53-71. 10.1111/j.1467-9868.2007.00627.x.
d <- 3 n <- 50 library(lhs) X <- maximinLHS(n, d) c <- c(0.2,0.6,0.8) F <- 1;for (a in 1:d) F <- F*(abs(4*X[,a]-2)+c[a])/(1+c[a]) epsilon <- rnorm(n,0,1);sigma <- 0.2 Y <- F + sigma*epsilon Dmax <- 3 kernel <- "matern" Kv <- calc_Kv(X, kernel, Dmax, TRUE,TRUE) matZ <- Kv$kv mumax <- mu_max(Y, matZ) mumax
d <- 3 n <- 50 library(lhs) X <- maximinLHS(n, d) c <- c(0.2,0.6,0.8) F <- 1;for (a in 1:d) F <- F*(abs(4*X[,a]-2)+c[a])/(1+c[a]) epsilon <- rnorm(n,0,1);sigma <- 0.2 Y <- F + sigma*epsilon Dmax <- 3 kernel <- "matern" Kv <- calc_Kv(X, kernel, Dmax, TRUE,TRUE) matZ <- Kv$kv mumax <- mu_max(Y, matZ) mumax
Fits the solution of the RKHS ridge group sparse optimization problem for the Gaussian regression model.
pen_MetMod(Y, Kv, gamma, mu, resg, gama_v, mu_v, maxIter, verbose, calcStwo)
pen_MetMod(Y, Kv, gamma, mu, resg, gama_v, mu_v, maxIter, verbose, calcStwo)
Y |
Vector of response observations of size |
Kv |
List, includes the eigenvalues and eigenvectors of the positive definite Gram matrices |
gamma |
Vector of positive scalars. Values of the penalty parameter |
mu |
Vector of positive scalars. Values of the penalty parameter |
resg |
List of initial parameters, includes the |
gama_v |
Scalar zero or vector of vMax positive scalars, considered as weights for the Ridge penalty. Set to zero, to consider no weights, i.e. all weights equal to |
mu_v |
Scalar zero or a vector with vMax scalars, considered as weigths of Sparse Group penalty. Set to zero, to consider no weights, i.e. all weights equal to |
maxIter |
Integer, shows the maximum number of loops through initial active groups at the first step and maximum number of loops through all groups at the second step. Set as |
verbose |
Logical, if TRUE, for each pair of penalty parameters |
calcStwo |
Logical, if TRUE, the program does a second step after convergence: the algorithm is done over all groups by taking the estimated parameters at the first step as initial values. Set as FALSE by default. |
Input Kv should contain the eigenvalues and eigenvectors of positive definite Gram matrices . It is necessary to set input "correction" in the function
calc_Kv
equal to "TRUE".
List of l components, with l equals to the number of pairs of the penalty parameters . Each component of the list is a list of
components "mu", "gamma" and "Meta-Model":
mu |
Positive scalar, an element of the input vector mu associated with the estimated Meta-Model. |
gamma |
Positive scalar, an element of the input vector gamma associated with the estimated Meta-Model. |
Meta-Model |
Estimated meta model associated with penalty parameters mu and gamma. List of |
intercept |
Scalar, estimated value of intercept. |
teta |
Matrix with vMax rows and |
fit.v |
Matrix with |
fitted |
Vector of size |
Norm.n |
Vector of size vMax, estimated values for the Ridge penalty norm. |
Norm.H |
Vector of size vMax, estimated values for the Group Sparse penalty norm. |
supp |
Vector of active groups. |
Nsupp |
Vector of the names of the active groups. |
SCR |
Scalar equals to |
crit |
Scalar indicates the value of the penalized criteria. |
gamma.v |
Vector of size vMax, coefficients of the Ridge penalty norm, |
mu.v |
Vector of size vMax, coefficients of the Group Sparse penalty norm, |
iter |
List of three components if calcStwo |
convergence |
TRUE or FALSE. Indicates whether the algorithm has converged or not. |
RelDiffCrit |
List of two components if calcStwo |
RelDiffPar |
List of two components if calcStwo |
Note.
Halaleh Kamari
Huet, S. and Taupin, M. L. (2017) Metamodel construction for sensitivity analysis. ESAIM: Procs 60, 27-69.
Kamari, H., Huet, S. and Taupin, M.-L. (2019) RKHSMetaMod : An R package to estimate the Hoeffding decomposition of an unknown function by solving RKHS Ridge Group Sparse optimization problem. <arXiv:1905.13695>
d <- 3 n <- 50 library(lhs) X <- maximinLHS(n, d) c <- c(0.2,0.6,0.8) F <- 1;for (a in 1:d) F <- F*(abs(4*X[,a]-2)+c[a])/(1+c[a]) epsilon <- rnorm(n,0,1);sigma <- 0.2 Y <- F + sigma*epsilon Dmax <- 3 kernel <- "matern" Kv <- calc_Kv(X, kernel, Dmax, TRUE,TRUE, tol = 1e-08) vMax <- length(Kv$names.Grp) matZ <- Kv$kv mumax <- mu_max(Y, matZ) mug1 <- mumax/10 mug2 <- mumax/100 gr1 <- RKHSgrplasso(Y,Kv, mug1) gr2 <- RKHSgrplasso(Y,Kv, mug2) gamma <- c(.5,.01,.001) #rescaling the penalty parameter mu <- c(mug1/sqrt(n),mug2/sqrt(n)) resg <- list(gr1,gr2) res <- pen_MetMod(Y,Kv,gamma,mu,resg,0,0) l <- length(res) for(i in 1:l){print(res[[i]]$mu)} for(i in 1:l){print(res[[i]]$gamma)} for(i in 1:l){print(res[[i]]$`Meta-Model`$Nsupp)} gama_v <- rep(1,vMax) mu_v <- rep(1,vMax) res.w <- pen_MetMod(Y,Kv,gamma,mu,resg,gama_v,mu_v) for(i in 1:l){print(res.w[[i]]$`Meta-Model`$Nsupp)}
d <- 3 n <- 50 library(lhs) X <- maximinLHS(n, d) c <- c(0.2,0.6,0.8) F <- 1;for (a in 1:d) F <- F*(abs(4*X[,a]-2)+c[a])/(1+c[a]) epsilon <- rnorm(n,0,1);sigma <- 0.2 Y <- F + sigma*epsilon Dmax <- 3 kernel <- "matern" Kv <- calc_Kv(X, kernel, Dmax, TRUE,TRUE, tol = 1e-08) vMax <- length(Kv$names.Grp) matZ <- Kv$kv mumax <- mu_max(Y, matZ) mug1 <- mumax/10 mug2 <- mumax/100 gr1 <- RKHSgrplasso(Y,Kv, mug1) gr2 <- RKHSgrplasso(Y,Kv, mug2) gamma <- c(.5,.01,.001) #rescaling the penalty parameter mu <- c(mug1/sqrt(n),mug2/sqrt(n)) resg <- list(gr1,gr2) res <- pen_MetMod(Y,Kv,gamma,mu,resg,0,0) l <- length(res) for(i in 1:l){print(res[[i]]$mu)} for(i in 1:l){print(res[[i]]$gamma)} for(i in 1:l){print(res[[i]]$`Meta-Model`$Nsupp)} gama_v <- rep(1,vMax) mu_v <- rep(1,vMax) res.w <- pen_MetMod(Y,Kv,gamma,mu,resg,gama_v,mu_v) for(i in 1:l){print(res.w[[i]]$`Meta-Model`$Nsupp)}
Computes the prediction error by considering a testing dataset.
PredErr(X, XT, YT, mu, gamma, res, kernel, Dmax)
PredErr(X, XT, YT, mu, gamma, res, kernel, Dmax)
X |
Matrix of observations with |
XT |
Matrix of observations of the testing dataset with |
YT |
Vector of response observations of testing dataset of size |
mu |
Vector of positive scalars. Values of the Group Sparse penalty parameter in decreasing order. See function |
gamma |
Vector of positive scalars. Values of the Ridge penalty parameter in decreasing order. See function |
res |
List, includes a squence of estimated meta models for the learning dataset, using RKHS Ridge Group Sparse or RKHS Group Lasso algorithm, associated with the penalty parameters mu and gamma. It should have the same format as the output of one of the functions: |
kernel |
Character, shows the type of the reproducing kernel: matern, brownian, gaussian, linear, quad. The same kernel should be chosen as the one used for the learning dataset. See function |
Dmax |
Integer between |
Details.
Matrix of the prediction errors is returned. Each element of the matrix is the obtained prediction error associated with one RKHS meta model in "res".
Note.
Halaleh Kamari
Kamari, H., Huet, S. and Taupin, M.-L. (2019) RKHSMetaMod : An R package to estimate the Hoeffding decomposition of an unknown function by solving RKHS Ridge Group Sparse optimization problem. <arXiv:1905.13695>
calc_Kv
, pen_MetMod
, RKHSMetMod
, RKHSMetMod_qmax
d <- 3 n <- 50 nT <- 50 library(lhs) X <- maximinLHS(n, d) XT <- maximinLHS(nT, d) c <- c(0.2,0.6,0.8) F <- 1;for (a in 1:d) F <- F*(abs(4*X[,a]-2)+c[a])/(1+c[a]) FT <- 1;for (a in 1:d) FT <- FT*(abs(4*XT[,a]-2)+c[a])/(1+c[a]) sigma <- 0.2 epsilon <- rnorm(n,0,1);Y <- F + sigma*epsilon epsilonT <- rnorm(nT,0,1);YT <- FT + sigma*epsilonT Dmax <- 3 kernel <- "matern" frc <- c(10,100) gamma <- c(.5,.01,.001) res <- RKHSMetMod(Y,X,kernel,Dmax,gamma,frc,FALSE) mu <- vector() l <- length(gamma) for(i in 1:length(frc)){mu[i]=res[[(i-1)*l+1]]$mu} error <- PredErr(X,XT, YT,mu,gamma, res, kernel,Dmax) error
d <- 3 n <- 50 nT <- 50 library(lhs) X <- maximinLHS(n, d) XT <- maximinLHS(nT, d) c <- c(0.2,0.6,0.8) F <- 1;for (a in 1:d) F <- F*(abs(4*X[,a]-2)+c[a])/(1+c[a]) FT <- 1;for (a in 1:d) FT <- FT*(abs(4*XT[,a]-2)+c[a])/(1+c[a]) sigma <- 0.2 epsilon <- rnorm(n,0,1);Y <- F + sigma*epsilon epsilonT <- rnorm(nT,0,1);YT <- FT + sigma*epsilonT Dmax <- 3 kernel <- "matern" frc <- c(10,100) gamma <- c(.5,.01,.001) res <- RKHSMetMod(Y,X,kernel,Dmax,gamma,frc,FALSE) mu <- vector() l <- length(gamma) for(i in 1:length(frc)){mu[i]=res[[(i-1)*l+1]]$mu} error <- PredErr(X,XT, YT,mu,gamma, res, kernel,Dmax) error
Fits the solution of an RKHS group lasso problem for the Gaussian regression model.
RKHSgrplasso(Y, Kv, mu, maxIter, verbose)
RKHSgrplasso(Y, Kv, mu, maxIter, verbose)
Y |
Vector of response observations of size |
Kv |
List, includes the eigenvalues and eigenvectors of the positive definite Gram matrices |
mu |
Positive scalar, value of the penalty parameter |
maxIter |
Integer, shows the maximum number of loops through all groups. Set as |
verbose |
Logical, if TRUE, prints: the number of current iteration, active groups and convergence criterias. Set as FALSE by default. |
Input Kv should contain the eigenvalues and eigenvectors of positive definite Gram matrices . It is necessary to set input correction in the function
calc_Kv
equal to "TRUE".
Estimated RKHS meta model, list with components:
intercept |
Scalar, estimated value of intercept. |
teta |
Matrix with vMax rows and |
fit.v |
Matrix with |
fitted |
Vector of size |
Norm.H |
Vector of size vMax, estimated values of the penalty norm. |
supp |
Vector of active groups. |
Nsupp |
Vector of the names of the active groups. |
SCR |
Scalar, equals to |
crit |
Scalar, indicates the value of penalized criteria. |
MaxIter |
Integer, number of iterations until convergence is reached. |
convergence |
TRUE or FALSE. Indicates whether the algorithm has converged or not. |
RelDiffCrit |
Scalar, value of the first convergence criteria at the last iteration, |
RelDiffPar |
Scalar, value of the second convergence criteria at the last iteration, |
Note.
Halaleh Kamari
Kamari, H., Huet, S. and Taupin, M.-L. (2019) RKHSMetaMod : An R package to estimate the Hoeffding decomposition of an unknown function by solving RKHS Ridge Group Sparse optimization problem. <arXiv:1905.13695>
Meier, L. Van de Geer, S. and Buhlmann, P. (2008) The group LASSO for logistic regression. Journal of the Royal Statistical Society Series B. 70. 53-71. 10.1111/j.1467-9868.2007.00627.x.
d <- 3 n <- 50 library(lhs) X <- maximinLHS(n, d) c <- c(0.2,0.6,0.8) F <- 1;for (a in 1:d) F <- F*(abs(4*X[,a]-2)+c[a])/(1+c[a]) epsilon <- rnorm(n,0,1);sigma <- 0.2 Y <- F + sigma*epsilon Dmax <- 3 kernel <- "matern" Kv <- calc_Kv(X, kernel, Dmax, TRUE, TRUE) matZ <- Kv$kv mumax <- mu_max(Y, matZ) mug <- mumax/10 gr <- RKHSgrplasso(Y,Kv, mug , 1000, FALSE) gr$Nsupp
d <- 3 n <- 50 library(lhs) X <- maximinLHS(n, d) c <- c(0.2,0.6,0.8) F <- 1;for (a in 1:d) F <- F*(abs(4*X[,a]-2)+c[a])/(1+c[a]) epsilon <- rnorm(n,0,1);sigma <- 0.2 Y <- F + sigma*epsilon Dmax <- 3 kernel <- "matern" Kv <- calc_Kv(X, kernel, Dmax, TRUE, TRUE) matZ <- Kv$kv mumax <- mu_max(Y, matZ) mug <- mumax/10 gr <- RKHSgrplasso(Y,Kv, mug , 1000, FALSE) gr$Nsupp
Calculates the Gram matrices for a chosen reproducing kernel and fits the solution of an RKHS ridge group sparse or an RKHS group lasso problem for each pair of penalty parameters
, for the Gaussian regression model.
RKHSMetMod(Y, X, kernel, Dmax, gamma, frc, verbose)
RKHSMetMod(Y, X, kernel, Dmax, gamma, frc, verbose)
Y |
Vector of response observations of size |
X |
Matrix of observations with |
kernel |
Character, indicates the type of the reproducing kernel: matern |
Dmax |
Integer, between |
gamma |
Vector of non negative scalars, values of the penalty parameter |
frc |
Vector of positive scalars. Each element of the vector sets a value to the penalty parameter |
verbose |
Logical, if TRUE, prints: the group |
Details.
List of l components, with l equals to the number of pairs of the penalty parameters . Each component of the list is a list of
components "mu", "gamma" and "Meta-Model":
mu |
Positive scalar, penalty parameter |
gamma |
Positive scalar, an element of the input vector gamma associated with the estimated Meta-Model. |
Meta-Model |
An RKHS Ridge Group Sparse or RKHS Group Lasso object associated with the penalty parameters mu and gamma: |
intercept |
Scalar, estimated value of intercept. |
teta |
Matrix with vMax rows and |
fit.v |
Matrix with |
fitted |
Vector of size |
Norm.n |
Vector of size vMax, estimated values for the Ridge penalty norm. |
Norm.H |
Vector of size vMax, estimated values of the Sparse Group penalty norm. |
supp |
Vector of active groups. |
Nsupp |
Vector of the names of the active groups. |
SCR |
Scalar, equals to |
crit |
Scalar, indicates the value of the penalized criteria. |
gamma.v |
Vector of size vMax, coefficients of the Ridge penalty norm, |
mu.v |
Vector of size vMax, coefficients of the Group Sparse penalty norm, |
iter |
List of two components: maxIter, and the number of iterations until the convergence is achieved. |
convergence |
TRUE or FALSE. Indicates whether the algorithm has converged or not. |
RelDiffCrit |
Scalar, value of the first convergence criteria at the last iteration, |
RelDiffPar |
Scalar, value of the second convergence criteria at the last iteration, |
For the case the outputs "mu"
and "Meta-Model" is the same as the one returned by the function
RKHSgrplasso
.
Halaleh Kamari
Kamari, H., Huet, S. and Taupin, M.-L. (2019) RKHSMetaMod : An R package to estimate the Hoeffding decomposition of an unknown function by solving RKHS Ridge Group Sparse optimization problem. <arXiv:1905.13695>
calc_Kv
, mu_max
, RKHSgrplasso
, pen_MetMod
d <- 3 n <- 50 library(lhs) X <- maximinLHS(n, d) c <- c(0.2,0.6,0.8) F <- 1;for (a in 1:d) F <- F*(abs(4*X[,a]-2)+c[a])/(1+c[a]) epsilon <- rnorm(n,0,1);sigma <- 0.2 Y <- F + sigma*epsilon Dmax <- 3 kernel <- "matern" frc <- c(10,100) gamma <- c(.5,.01,.001,0) result <- RKHSMetMod(Y,X,kernel,Dmax,gamma,frc,FALSE) l <- length(result) for(i in 1:l){print(result[[i]]$mu)} for(i in 1:l){print(result[[i]]$gamma)} for(i in 1:l){print(result[[i]]$`Meta-Model`$Nsupp)}
d <- 3 n <- 50 library(lhs) X <- maximinLHS(n, d) c <- c(0.2,0.6,0.8) F <- 1;for (a in 1:d) F <- F*(abs(4*X[,a]-2)+c[a])/(1+c[a]) epsilon <- rnorm(n,0,1);sigma <- 0.2 Y <- F + sigma*epsilon Dmax <- 3 kernel <- "matern" frc <- c(10,100) gamma <- c(.5,.01,.001,0) result <- RKHSMetMod(Y,X,kernel,Dmax,gamma,frc,FALSE) l <- length(result) for(i in 1:l){print(result[[i]]$mu)} for(i in 1:l){print(result[[i]]$gamma)} for(i in 1:l){print(result[[i]]$`Meta-Model`$Nsupp)}
Calculates the Gram matrices for a chosen kernel, determines
, note
, for which the number of active groups in the RKHS group lasso solution is equal to qmax, and fits a solution of an RKHS ridge group sparse or an RKHS group lasso problem for each pair of penalty parameters
, in the Gaussian regression model.
RKHSMetMod_qmax(Y, X, kernel, Dmax, gamma, qmax, rat, Num, verbose)
RKHSMetMod_qmax(Y, X, kernel, Dmax, gamma, qmax, rat, Num, verbose)
Y |
Vector of response observations of size |
X |
Matrix of observations with |
kernel |
Character, indicates the type of the reproducing kernel: matern |
Dmax |
Integer, between |
gamma |
Vector of non negative scalars, values of the penalty parameter |
qmax |
Integer, shows the maximum number of active groups in the obtained solution. |
rat |
Positive scalar, to restrict the minimum value of |
Num |
Integer, it is used to restrict the number of different values of the penalty parameter |
verbose |
Logical, if TRUE, prints: the group |
Details.
List of three components "mus", "qs", and "MetaModel":
mus |
Vector, values of the evaluated penalty parameters |
qs |
Vector, number of active groups associated with each element in mus. |
MetaModel |
List with the same length as the vector gamma. Each component of the list is a list of |
mu |
Scalar, the value |
gamma |
Positive scalar, element of the input vector gamma associated with the estimated Meta-Model. |
Meta-Model |
An RKHS Ridge Group Sparse or RKHS Group Lasso object associated with the penalty parameters mu and gamma: |
intercept |
Scalar, estimated value of intercept. |
teta |
Matrix with vMax rows and |
fit.v |
Matrix with |
fitted |
Vector of size |
Norm.n |
Vector of size vMax, estimated values for the Ridge penalty norm. |
Norm.H |
Vector of size vMax, estimated values of the Sparse Group penalty norm. |
supp |
Vector of active groups. |
Nsupp |
Vector of the names of the active groups. |
SCR |
Scalar, equals to |
crit |
Scalar, indicates the value of penalized criteria. |
gamma.v |
Vector, coefficients of the Ridge penalty norm, |
mu.v |
Vector, coefficients of the Group Sparse penalty norm, |
iter |
List of two components: maxIter, and the number of iterations until the convergence is achieved. |
convergence |
TRUE or FALSE. Indicates whether the algorithm has converged or not. |
RelDiffCrit |
Scalar, value of the first convergence criteria at the last iteration, |
RelDiffPar |
Scalar, value of the second convergence criteria at the last iteration, |
For the case the outputs "mu"
and "Meta-Model" is the same as the one returned by the function
RKHSgrplasso
.
Halaleh Kamari
Kamari, H., Huet, S. and Taupin, M.-L. (2019) RKHSMetaMod : An R package to estimate the Hoeffding decomposition of an unknown function by solving RKHS Ridge Group Sparse optimization problem. <arXiv:1905.13695>
calc_Kv
, mu_max
, RKHSgrplasso
, pen_MetMod
, grplasso_q
d <- 3 n <- 50 library(lhs) X <- maximinLHS(n, d) c <- c(0.2,0.6,0.8) F <- 1;for (a in 1:d) F <- F*(abs(4*X[,a]-2)+c[a])/(1+c[a]) epsilon <- rnorm(n,0,1);sigma <- 0.2 Y <- F + sigma*epsilon Dmax <- 3 kernel <- "matern" gamma <- c(.5,.01,.001,0) Num <- 10 rat <- 100 qmax <- 4 result <- RKHSMetMod_qmax(Y, X, kernel, Dmax, gamma, qmax, rat, Num,FALSE) names(result) result$mus result$qs l <- length(gamma) for(i in 1:l){print(result$MetaModel[[i]]$mu)} for(i in 1:l){print(result$MetaModel[[i]]$gamma)} for(i in 1:l){print(result$MetaModel[[i]]$`Meta-Model`$Nsupp)}
d <- 3 n <- 50 library(lhs) X <- maximinLHS(n, d) c <- c(0.2,0.6,0.8) F <- 1;for (a in 1:d) F <- F*(abs(4*X[,a]-2)+c[a])/(1+c[a]) epsilon <- rnorm(n,0,1);sigma <- 0.2 Y <- F + sigma*epsilon Dmax <- 3 kernel <- "matern" gamma <- c(.5,.01,.001,0) Num <- 10 rat <- 100 qmax <- 4 result <- RKHSMetMod_qmax(Y, X, kernel, Dmax, gamma, qmax, rat, Num,FALSE) names(result) result$mus result$qs l <- length(gamma) for(i in 1:l){print(result$MetaModel[[i]]$mu)} for(i in 1:l){print(result$MetaModel[[i]]$gamma)} for(i in 1:l){print(result$MetaModel[[i]]$`Meta-Model`$Nsupp)}
Calculates the empirical sensitivity indices.
SI_emp(res,ErrPred)
SI_emp(res,ErrPred)
res |
List, includes a sequence of estimated meta models, the solutions of the RKHS Ridge Group Sparse or RKHS Group Lasso problems. It should has the same format as the output of one of the functions: |
ErrPred |
Matrix or NULL. If matrix, each element of the matrix is the obtained prediction error associated with one RKHS meta model in "res". It should have the same format as the output of the function |
Details.
If input ErrPred"NULL", Vector of the empirical sensitivity incdices for the meta model with the minimum Prediction error is returned. If ErrPred
"NULL", a list of the vectors is returned. Each vector is the obtained sensitivity indices associated with one meta model in "res".
Note.
Halaleh Kamari
Kamari, H., Huet, S. and Taupin, M.-L. (2019) RKHSMetaMod : An R package to estimate the Hoeffding decomposition of an unknown function by solving RKHS Ridge Group Sparse optimization problem. <arXiv:1905.13695>
PredErr
, pen_MetMod
, RKHSMetMod
, RKHSMetMod_qmax
d <- 3 n <- 50;nT <- 50 library(lhs) X <- maximinLHS(n, d);XT <- maximinLHS(nT, d) c <- c(0.2,0.6,0.8) F <- 1;for (a in 1:d) F <- F*(abs(4*X[,a]-2)+c[a])/(1+c[a]) FT <- 1;for (a in 1:d) FT <- FT*(abs(4*XT[,a]-2)+c[a])/(1+c[a]) sigma <- 0.2 epsilon <- rnorm(n,0,1);Y <- F + sigma*epsilon epsilonT <- rnorm(nT,0,1);YT <- FT + sigma*epsilonT Dmax <- 3 kernel <- "matern" frc <- c(10) gamma=c(.5,.01,.001) res <- RKHSMetMod(Y,X,kernel,Dmax,gamma,frc,FALSE) mu <- vector() l <- length(gamma) for(i in 1:length(frc)){mu[i]=res[[(i-1)*l+1]]$mu} error <- PredErr(X,XT, YT,mu,gamma, res, kernel,Dmax) SI.minErr <- SI_emp(res, error) SI.minErr SI <- SI_emp(res, NULL) SI
d <- 3 n <- 50;nT <- 50 library(lhs) X <- maximinLHS(n, d);XT <- maximinLHS(nT, d) c <- c(0.2,0.6,0.8) F <- 1;for (a in 1:d) F <- F*(abs(4*X[,a]-2)+c[a])/(1+c[a]) FT <- 1;for (a in 1:d) FT <- FT*(abs(4*XT[,a]-2)+c[a])/(1+c[a]) sigma <- 0.2 epsilon <- rnorm(n,0,1);Y <- F + sigma*epsilon epsilonT <- rnorm(nT,0,1);YT <- FT + sigma*epsilonT Dmax <- 3 kernel <- "matern" frc <- c(10) gamma=c(.5,.01,.001) res <- RKHSMetMod(Y,X,kernel,Dmax,gamma,frc,FALSE) mu <- vector() l <- length(gamma) for(i in 1:length(frc)){mu[i]=res[[(i-1)*l+1]]$mu} error <- PredErr(X,XT, YT,mu,gamma, res, kernel,Dmax) SI.minErr <- SI_emp(res, error) SI.minErr SI <- SI_emp(res, NULL) SI