slope 0.29.0
Loading...
Searching...
No Matches
diagnostics.h
Go to the documentation of this file.
1
6#pragma once
7
8#include "jit_normalization.h"
9#include "losses/loss.h"
10#include "sorted_l1_norm.h"
11#include <Eigen/Dense>
12#include <memory>
13#include <numeric>
14
15namespace slope {
16
39template<typename MatrixType>
40double
41computeDual(const Eigen::VectorXd& beta,
42 const Eigen::MatrixXd& residual,
43 const std::unique_ptr<Loss>& loss,
44 const SortedL1Norm& sl1_norm,
45 const Eigen::ArrayXd& lambda,
46 const MatrixType& x,
47 const Eigen::MatrixXd& y,
48 const Eigen::VectorXd& x_centers,
49 const Eigen::VectorXd& x_scales,
50 const JitNormalization& jit_normalization,
51 const bool intercept)
52{
53 int n = x.rows();
54 int pm = beta.size();
55
56 Eigen::VectorXd gradient(pm);
57
58 std::vector<int> full_set(pm);
59 std::iota(full_set.begin(), full_set.end(), 0);
60
61 updateGradient(gradient,
62 x,
63 residual,
64 full_set,
65 x_centers,
66 x_scales,
67 Eigen::VectorXd::Ones(n),
68 jit_normalization);
69
70 Eigen::MatrixXd theta = residual;
71
72 // First compute gradient with potential offset for intercept case
73 Eigen::VectorXd dual_gradient = gradient;
74
75 // TODO: Can we avoid this copy? Maybe revert offset afterwards or,
76 // alternatively, solve intercept until convergence and then no longer
77 // need the offset at all.
78 if (intercept) {
79 Eigen::VectorXd theta_mean = theta.colwise().mean();
80 theta.rowwise() -= theta_mean.transpose();
81
82 offsetGradient(dual_gradient,
83 x,
84 theta_mean,
85 full_set,
86 x_centers,
87 x_scales,
88 jit_normalization);
89 }
90
91 // Common scaling operation
92 double dual_norm = sl1_norm.dualNorm(dual_gradient, lambda);
93 theta.array() /= std::max(1.0, dual_norm);
94
95 double dual = loss->dual(theta, y, Eigen::VectorXd::Ones(n));
96
97 return dual;
98}
99
100} // namespace slope
Class representing the Sorted L1 Norm.
double dualNorm(const Eigen::VectorXd &a, const Eigen::ArrayXd &lambda) const
Computes the dual norm of a vector.
Enums to control predictor standardization behavior.
The declartion of the Objctive class and its subclasses, which represent the data-fitting part of the...
Namespace containing SLOPE regression implementation.
Definition clusters.cpp:5
JitNormalization
Enums to control predictor standardization behavior.
double computeDual(const Eigen::VectorXd &beta, const Eigen::MatrixXd &residual, const std::unique_ptr< Loss > &loss, const SortedL1Norm &sl1_norm, const Eigen::ArrayXd &lambda, const MatrixType &x, const Eigen::MatrixXd &y, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, const JitNormalization &jit_normalization, const bool intercept)
Computes the dual objective function value for SLOPE optimization.
Definition diagnostics.h:41
void updateGradient(Eigen::VectorXd &gradient, const T &x, const Eigen::MatrixXd &residual, const std::vector< int > &active_set, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, const Eigen::VectorXd &w, const JitNormalization jit_normalization)
Definition math.h:223
void offsetGradient(Eigen::VectorXd &gradient, const T &x, const Eigen::VectorXd &offset, const std::vector< int > &active_set, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, const JitNormalization jit_normalization)
Definition math.h:295
The declaration of the SortedL1Norm class.