slope 0.28.0
Loading...
Searching...
No Matches
slope_path.h
Go to the documentation of this file.
1
11#pragma once
12
13#include "slope_fit.h"
14#include <Eigen/Dense>
15#include <Eigen/SparseCore>
16
17namespace slope {
18
31{
32private:
33 std::vector<SlopeFit> fits;
34
35public:
39 SlopePath() = default;
40
47 SlopePath(const std::vector<SlopeFit>& fits)
48 : fits{ fits }
49 {
50 }
51
58 const SlopeFit& operator()(const size_t step)
59 {
60 assert(step < fits.size());
61
62 return fits[step];
63 }
64
74 std::vector<Eigen::VectorXd> getIntercepts() const
75 {
76 std::vector<Eigen::VectorXd> intercepts;
77
78 for (const auto& fit : fits) {
79 intercepts.emplace_back(fit.getIntercepts());
80 }
81
82 return intercepts;
83 }
84
95 std::vector<Eigen::SparseMatrix<double>> getCoefs() const
96 {
97 std::vector<Eigen::SparseMatrix<double>> coefs;
98
99 for (const auto& fit : fits) {
100 coefs.emplace_back(fit.getCoefs());
101 }
102
103 return coefs;
104 }
105
115 std::vector<std::vector<std::vector<int>>> getClusters() const
116 {
117 std::vector<std::vector<std::vector<int>>> clusters;
118
119 for (const auto& fit : fits) {
120 clusters.emplace_back(fit.getClusters());
121 }
122
123 return clusters;
124 }
125
135 Eigen::SparseMatrix<double> getCoefs(const std::size_t i) const
136 {
137 assert(i < fits.size() && "Index out of bounds");
138
139 return fits[i].getCoefs();
140 }
141
145 Eigen::ArrayXd getAlpha() const
146 {
147
148 Eigen::ArrayXd alphas(fits.size());
149
150 for (size_t i = 0; i < fits.size(); i++) {
151 alphas(i) = fits[i].getAlpha();
152 }
153
154 return alphas;
155 }
156
160 const Eigen::ArrayXd& getLambda() const { return fits.front().getLambda(); }
161
165 std::vector<double> getDeviance() const
166 {
167 std::vector<double> deviances;
168
169 for (const auto& fit : fits) {
170 deviances.emplace_back(fit.getDeviance());
171 }
172
173 return deviances;
174 }
175
179 double getNullDeviance() const { return fits.front().getNullDeviance(); }
180
184 std::vector<std::vector<double>> getPrimals() const
185 {
186 std::vector<std::vector<double>> primals;
187
188 for (const auto& fit : fits) {
189 primals.emplace_back(fit.getPrimals());
190 }
191
192 return primals;
193 }
194
198 std::vector<std::vector<double>> getDuals() const
199 {
200 std::vector<std::vector<double>> duals;
201
202 for (const auto& fit : fits) {
203 duals.emplace_back(fit.getDuals());
204 }
205
206 return duals;
207 }
208
212 std::vector<std::vector<double>> getTime() const
213 {
214 std::vector<std::vector<double>> time;
215
216 for (const auto& fit : fits) {
217 time.emplace_back(fit.getTime());
218 }
219
220 return time;
221 }
222
226 std::vector<int> getPasses() const
227 {
228 std::vector<int> passes;
229
230 for (const auto& fit : fits) {
231 passes.emplace_back(fit.getPasses());
232 }
233
234 return passes;
235 }
236
242 const std::vector<double> getDevianceRatios() const
243 {
244 std::vector<double> ratios;
245
246 for (const auto& fit : fits) {
247 ratios.emplace_back(fit.getDevianceRatio());
248 }
249
250 return ratios;
251 }
252
258 std::vector<std::vector<double>> getGaps() const
259 {
260 std::vector<std::vector<double>> gaps;
261
262 for (const auto& fit : fits) {
263 gaps.emplace_back(fit.getGaps());
264 }
265
266 return gaps;
267 }
268};
269
270} // namespace slope
A class representing the results of SLOPE (Sorted L1 Penalized Estimation) fitting.
Definition slope_fit.h:26
Container class for SLOPE regression solution paths.
Definition slope_path.h:31
std::vector< double > getDeviance() const
Gets the deviance values for each solution.
Definition slope_path.h:165
std::vector< std::vector< double > > getGaps() const
Computes the duality gaps (primal - dual objectives) for each solution.
Definition slope_path.h:258
const SlopeFit & operator()(const size_t step)
Get one of the coefficients along the path.
Definition slope_path.h:58
std::vector< std::vector< double > > getDuals() const
Gets the dual objective values during optimization.
Definition slope_path.h:198
std::vector< int > getPasses() const
Gets the number of iterations for each solution.
Definition slope_path.h:226
std::vector< std::vector< double > > getPrimals() const
Gets the primal objective values during optimization.
Definition slope_path.h:184
std::vector< std::vector< double > > getTime() const
Gets the computation times for each solution.
Definition slope_path.h:212
const std::vector< double > getDevianceRatios() const
Computes the deviance ratio (explained deviance) for each solution.
Definition slope_path.h:242
SlopePath()=default
Constructs an empty SlopePath object.
std::vector< Eigen::VectorXd > getIntercepts() const
Returns the vector of intercept terms for each solution in the path.
Definition slope_path.h:74
const Eigen::ArrayXd & getLambda() const
Gets the lambda (regularization) weights.
Definition slope_path.h:160
Eigen::SparseMatrix< double > getCoefs(const std::size_t i) const
Gets the sparse coefficient matrix for a specific solution in the path.
Definition slope_path.h:135
Eigen::ArrayXd getAlpha() const
Gets the alpha parameter sequence.
Definition slope_path.h:145
std::vector< Eigen::SparseMatrix< double > > getCoefs() const
Returns the vector of coefficient matrices for each solution in the path.
Definition slope_path.h:95
SlopePath(const std::vector< SlopeFit > &fits)
Constructs a SlopePath object containing SLOPE regression solution path data.
Definition slope_path.h:47
double getNullDeviance() const
Gets the null model deviance.
Definition slope_path.h:179
std::vector< std::vector< std::vector< int > > > getClusters() const
Returns the clusters for each solution in the path.
Definition slope_path.h:115
Namespace containing SLOPE regression implementation.
Definition clusters.cpp:5
SLOPE (Sorted L-One Penalized Estimation) fitting results.