slope 6.5.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#include <cassert>
17
18namespace slope {
19
32{
33private:
34 std::vector<SlopeFit> fits;
35
36public:
40 SlopePath() = default;
41
48 SlopePath(const std::vector<SlopeFit>& fits)
49 : fits{ fits }
50 {
51 }
52
59 const SlopeFit& operator()(const size_t step) const
60 {
61 assert(step < fits.size());
62
63 return fits[step];
64 }
65
75 std::vector<Eigen::VectorXd> getIntercepts(
76 const bool original_scale = true) const
77 {
78 std::vector<Eigen::VectorXd> intercepts;
79
80 for (const auto& fit : fits) {
81 intercepts.emplace_back(fit.getIntercepts(original_scale));
82 }
83
84 return intercepts;
85 }
86
97 std::vector<Eigen::SparseMatrix<double>> getCoefs(
98 const bool original_scale = true) const
99 {
100 std::vector<Eigen::SparseMatrix<double>> coefs;
101
102 for (const auto& fit : fits) {
103 coefs.emplace_back(fit.getCoefs(original_scale));
104 }
105
106 return coefs;
107 }
108
114 std::vector<Clusters> getClusters() const
115 {
116 std::vector<Clusters> clusters;
117
118 for (const auto& fit : fits) {
119 clusters.emplace_back(fit.getClusters());
120 }
121
122 return clusters;
123 }
124
128 Eigen::ArrayXd getAlpha() const
129 {
130
131 Eigen::ArrayXd alphas(fits.size());
132
133 for (size_t i = 0; i < fits.size(); i++) {
134 alphas(i) = fits[i].getAlpha();
135 }
136
137 return alphas;
138 }
139
143 const Eigen::ArrayXd& getLambda() const { return fits.front().getLambda(); }
144
148 std::vector<double> getDeviance() const
149 {
150 std::vector<double> deviances;
151
152 for (const auto& fit : fits) {
153 deviances.emplace_back(fit.getDeviance());
154 }
155
156 return deviances;
157 }
158
162 double getNullDeviance() const { return fits.front().getNullDeviance(); }
163
167 std::vector<std::vector<double>> getPrimals() const
168 {
169 std::vector<std::vector<double>> primals;
170
171 for (const auto& fit : fits) {
172 primals.emplace_back(fit.getPrimals());
173 }
174
175 return primals;
176 }
177
181 std::vector<std::vector<double>> getDuals() const
182 {
183 std::vector<std::vector<double>> duals;
184
185 for (const auto& fit : fits) {
186 duals.emplace_back(fit.getDuals());
187 }
188
189 return duals;
190 }
191
195 std::vector<std::vector<double>> getTime() const
196 {
197 std::vector<std::vector<double>> time;
198
199 for (const auto& fit : fits) {
200 time.emplace_back(fit.getTime());
201 }
202
203 return time;
204 }
205
209 std::vector<int> getPasses() const
210 {
211 std::vector<int> passes;
212
213 for (const auto& fit : fits) {
214 passes.emplace_back(fit.getPasses());
215 }
216
217 return passes;
218 }
219
225 const std::vector<double> getDevianceRatios() const
226 {
227 std::vector<double> ratios;
228
229 for (const auto& fit : fits) {
230 ratios.emplace_back(fit.getDevianceRatio());
231 }
232
233 return ratios;
234 }
235
241 std::vector<std::vector<double>> getGaps() const
242 {
243 std::vector<std::vector<double>> gaps;
244
245 for (const auto& fit : fits) {
246 gaps.emplace_back(fit.getGaps());
247 }
248
249 return gaps;
250 }
251
256 std::size_t size() const { return fits.size(); }
257};
258
259} // namespace slope
A class representing the results of SLOPE (Sorted L1 Penalized Estimation) fitting.
Definition slope_fit.h:27
Container class for SLOPE regression solution paths.
Definition slope_path.h:32
std::vector< Eigen::VectorXd > getIntercepts(const bool original_scale=true) const
Returns the vector of intercept terms for each solution in the path.
Definition slope_path.h:75
std::vector< double > getDeviance() const
Gets the deviance values for each solution.
Definition slope_path.h:148
const SlopeFit & operator()(const size_t step) const
Get one of the coefficients along the path.
Definition slope_path.h:59
std::vector< std::vector< double > > getGaps() const
Computes the duality gaps (primal - dual objectives) for each solution.
Definition slope_path.h:241
std::vector< std::vector< double > > getDuals() const
Gets the dual objective values during optimization.
Definition slope_path.h:181
std::vector< int > getPasses() const
Gets the number of iterations for each solution.
Definition slope_path.h:209
std::vector< std::vector< double > > getPrimals() const
Gets the primal objective values during optimization.
Definition slope_path.h:167
std::vector< std::vector< double > > getTime() const
Gets the computation times for each solution.
Definition slope_path.h:195
const std::vector< double > getDevianceRatios() const
Computes the deviance ratio (explained deviance) for each solution.
Definition slope_path.h:225
SlopePath()=default
Constructs an empty SlopePath object.
const Eigen::ArrayXd & getLambda() const
Gets the lambda (regularization) weights.
Definition slope_path.h:143
Eigen::ArrayXd getAlpha() const
Gets the alpha parameter sequence.
Definition slope_path.h:128
std::size_t size() const
Gets the number of solutions in the path.
Definition slope_path.h:256
std::vector< Clusters > getClusters() const
Returns the clusters for each solution in the path.
Definition slope_path.h:114
std::vector< Eigen::SparseMatrix< double > > getCoefs(const bool original_scale=true) const
Returns the vector of coefficient matrices for each solution in the path.
Definition slope_path.h:97
SlopePath(const std::vector< SlopeFit > &fits)
Constructs a SlopePath object containing SLOPE regression solution path data.
Definition slope_path.h:48
double getNullDeviance() const
Gets the null model deviance.
Definition slope_path.h:162
Namespace containing SLOPE regression implementation.
Definition clusters.h:11
SLOPE (Sorted L-One Penalized Estimation) fitting results.