Elastic net¶
The elastic net (EN) seeks an approximate solution to a (usually undetermined) system of equations which penalizes both the one and two norms of the solution:
One could perhaps motivate the elastic net by recognizing that an LQ decomposition of an underdetermined systems of equations finds an \(x\) which exactly satisfies \(A x = b\) and simultaneously minimizes \(\| x \|_2\). From this perspective, for underdetermined systems, the elastic net trades the exact satisfaction of \(A x = b\), which would be provided via an LQ decomposition, for the hope of a sparse \(x\) via the addition of a penalty on \(\| x \|_1\).
Real instances of the problem are expressable as an (affine) quadratic program in a manner which is nearly identical to that of basis pursuit denoising (BPDN): by splitting \(x\) into its positive and negative components, say \(x = u - v\), and introducing the residual \(r = b - A x\), one arrives at
By default, Elemental solves this equation using a Mehrotra Predictor-Corrector primal-dual Interior Point Method.
Python API¶
-
EN
(A, b, lambda1, lambda2[, ctrl=None])¶ - Parameters
A – dense or sparse, sequential or distributed matrix
b – dense right-hand side vector (with type compatible to
A
)lambda1 – penalty on the one-norm of the solution
lambda2 – penalty on the square of the two-norm of the solution
ctrl – (optional)
QPAffineCtrl
instance
- Return type
dense solution vector (with type matching that of
b
)
C++ API¶
-
void
EN
(const Matrix<Real> &A, const Matrix<Real> &b, Real lambda1, Real lambda2, Matrix<Real> &x, const qp::affine::Ctrl<Real> &ctrl = qp::affine::Ctrl<Real>())¶
-
void
EN
(const ElementalMatrix<Real> &A, const ElementalMatrix<Real> &b, Real lambda1, Real lambda2, ElementalMatrix<Real> &x, const qp::affine::Ctrl<Real> &ctrl = qp::affine::Ctrl<Real>())¶
-
void
EN
(const SparseMatrix<Real> &A, const Matrix<Real> &b, Real lambda1, Real lambda2, Matrix<Real> &x, const qp::affine::Ctrl<Real> &ctrl = qp::affine::Ctrl<Real>())¶
-
void
EN
(const DistSparseMatrix<Real> &A, const DistMultiVec<Real> &b, Real lambda1, Real lambda2, DistMultiVec<Real> &x, const qp::affine::Ctrl<Real> &ctrl = qp::affine::Ctrl<Real>())¶
C API¶
Single-precision¶
-
ElError
ElENDist_s
(ElConstDistMatrix_s A, ElConstDistMatrix_s b, float lambda1, float lambda2, ElDistMatrix_s x)¶
Double-precision¶
-
ElError
ElENDist_d
(ElConstDistMatrix_d A, ElConstDistMatrix_d b, double lambda1, double lambda2, ElDistMatrix_d x)¶
-
ElError
ElENSparse_d
(ElConstSparseMatrix_d A, ElConstMatrix_d b, double lambda1, double lambda2, ElMatrix_d x)¶
-
ElError
ElENDistSparse_d
(ElConstDistSparseMatrix_d A, ElConstDistMultiVec_d b, double lambda1, double lambda2, ElDistMultiVec_d x)¶
Expert interface¶
Single-precision¶
-
ElError
ElENX_s
(ElConstMatrix_s A, ElConstMatrix_s b, float lambda1, float lambda2, ElMatrix_s x, ElQPAffineCtrl_s ctrl)¶
-
ElError
ElENXDist_s
(ElConstDistMatrix_s A, ElConstDistMatrix_s b, float lambda1, float lambda2, ElDistMatrix_s x, ElQPAffineCtrl_s ctrl)¶
Double-precision¶
-
ElError
ElENX_d
(ElConstMatrix_d A, ElConstMatrix_d b, double lambda1, double lambda2, ElMatrix_d x, ElQPAffineCtrl_d ctrl)¶
-
ElError
ElENXDist_d
(ElConstDistMatrix_d A, ElConstDistMatrix_d b, double lambda1, double lambda2, ElDistMatrix_d x, ElQPAffineCtrl_d ctrl)¶