slope 6.0.1
Loading...
Searching...
No Matches
slope_fit.h
Go to the documentation of this file.
1
6#pragma once
7
8#include "clusters.h"
9#include "losses/setup_loss.h"
10#include "normalize.h"
11#include "utils.h"
12#include <Eigen/Dense>
13#include <Eigen/SparseCore>
14#include <memory>
15
16namespace slope {
17
27{
28private:
29 Eigen::VectorXd intercepts;
30 Eigen::SparseMatrix<double> coefs;
31 Clusters clusters;
32 double alpha;
33 Eigen::ArrayXd lambda;
34 double deviance;
35 double null_deviance;
36 std::vector<double>
37 primals;
38 std::vector<double>
39 duals;
40 std::vector<double> time;
41 int passes;
42 std::string centering_type;
43 std::string scaling_type;
44 std::string loss_type;
45 bool has_intercept;
46 Eigen::VectorXd x_centers;
47 Eigen::VectorXd x_scales;
48
49public:
50 SlopeFit() = default;
51
73 SlopeFit(const Eigen::VectorXd& intercepts,
74 const Eigen::SparseMatrix<double>& coefs,
75 const Clusters& clusters,
76 const double alpha,
77 const Eigen::ArrayXd& lambda,
78 const double deviance,
79 const double null_deviance,
80 const std::vector<double>& primals,
81 const std::vector<double>& duals,
82 const std::vector<double>& time,
83 const int passes,
84 const std::string& centering_type,
85 const std::string& scaling_type,
86 const bool has_intercept,
87 const Eigen::VectorXd& x_centers,
88 const Eigen::VectorXd& x_scales)
89 : intercepts{ intercepts }
90 , coefs{ coefs }
91 , clusters{ clusters }
92 , alpha{ alpha }
93 , lambda{ lambda }
94 , deviance{ deviance }
95 , null_deviance{ null_deviance }
96 , primals{ primals }
97 , duals{ duals }
98 , time{ time }
99 , passes{ passes }
100 , centering_type{ centering_type }
101 , scaling_type{ scaling_type }
102 , has_intercept{ has_intercept }
103 , x_centers{ x_centers }
104 , x_scales{ x_scales }
105 {
106 }
107
113 Eigen::VectorXd getIntercepts(const bool original_scale = true) const
114 {
115 // TODO: Scale intercepts independently of coefficients
116 if (original_scale) {
117 auto [beta0_out, beta_out] = rescaleCoefficients(intercepts,
118 coefs,
119 this->x_centers,
120 this->x_scales,
121 this->has_intercept);
122 return beta0_out;
123 }
124
125 return intercepts;
126 }
127
133 Eigen::SparseMatrix<double> getCoefs(const bool original_scale = true) const
134 {
135 if (original_scale) {
136 // TODO: Scale coefficients independently of intercepts
137 auto [beta0_out, beta_out] = rescaleCoefficients(intercepts,
138 coefs,
139 this->x_centers,
140 this->x_scales,
141 this->has_intercept);
142 return beta_out.sparseView();
143 }
144
145 return coefs;
146 }
147
151 const Clusters& getClusters() const { return clusters; }
152
156 const Eigen::ArrayXd& getLambda() const { return lambda; }
157
161 double getAlpha() const { return alpha; }
162
166 double getDeviance() const { return deviance; }
167
171 const std::string& getLossType() const { return loss_type; }
172
176 double getNullDeviance() const { return null_deviance; }
177
181 const std::vector<double>& getPrimals() const { return primals; }
182
186 const std::vector<double>& getDuals() const { return duals; }
187
191 const std::vector<double>& getTime() const { return time; }
192
196 int getPasses() const { return passes; }
197
204 double getDevianceRatio() const { return 1.0 - deviance / null_deviance; }
205
212 std::vector<double> getGaps() const
213 {
214 std::vector<double> gaps(primals.size());
215 for (size_t i = 0; i < primals.size(); i++) {
216 gaps[i] = primals[i] - duals[i];
217 }
218 return gaps;
219 }
220
226 bool hasIntercept() const { return has_intercept; }
227
236 template<typename T>
237 Eigen::MatrixXd predict(Eigen::EigenBase<T>& x,
238 const std::string& type = "response") const
239 {
240 validateOption(type, { "response", "linear" }, "type");
241
242 Eigen::MatrixXd eta = x.derived() * getCoefs();
243
244 if (has_intercept) {
245 eta.rowwise() += getIntercepts().transpose();
246 }
247
248 if (type == "linear") {
249 return eta;
250 }
251
252 // Return predictions
253 std::unique_ptr<Loss> loss = setupLoss(this->loss_type);
254
255 return loss->predict(eta);
256 }
257};
258
259} // namespace slope
Representation of the clusters in SLOPE.
Definition clusters.h:18
A class representing the results of SLOPE (Sorted L1 Penalized Estimation) fitting.
Definition slope_fit.h:27
const Clusters & getClusters() const
Gets the clusters.
Definition slope_fit.h:151
bool hasIntercept() const
Checks if model has intercept.
Definition slope_fit.h:226
Eigen::VectorXd getIntercepts(const bool original_scale=true) const
Gets the intercept terms for this SLOPE fit.
Definition slope_fit.h:113
const Eigen::ArrayXd & getLambda() const
Gets the lambda (regularization) parameter used.
Definition slope_fit.h:156
SlopeFit(const Eigen::VectorXd &intercepts, const Eigen::SparseMatrix< double > &coefs, const Clusters &clusters, const double alpha, const Eigen::ArrayXd &lambda, const double deviance, const double null_deviance, const std::vector< double > &primals, const std::vector< double > &duals, const std::vector< double > &time, const int passes, const std::string &centering_type, const std::string &scaling_type, const bool has_intercept, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales)
Construct a new Slope Fit object.
Definition slope_fit.h:73
double getAlpha() const
Gets the alpha (mixing) parameter used.
Definition slope_fit.h:161
Eigen::SparseMatrix< double > getCoefs(const bool original_scale=true) const
Gets the sparse coefficient matrix for this fit.
Definition slope_fit.h:133
const std::vector< double > & getPrimals() const
Gets the sequence of primal objective values during optimization.
Definition slope_fit.h:181
const std::vector< double > & getDuals() const
Gets the sequence of dual objective values during optimization.
Definition slope_fit.h:186
const std::vector< double > & getTime() const
Gets the sequence of computation times during optimization.
Definition slope_fit.h:191
Eigen::MatrixXd predict(Eigen::EigenBase< T > &x, const std::string &type="response") const
Predict the response for a given input matrix.
Definition slope_fit.h:237
std::vector< double > getGaps() const
Calculate the duality gaps during optimization.
Definition slope_fit.h:212
const std::string & getLossType() const
Gets the model's loss type.
Definition slope_fit.h:171
double getDevianceRatio() const
Calculate the deviance ratio (1 - deviance/null_deviance)
Definition slope_fit.h:204
double getNullDeviance() const
Gets the null model deviance.
Definition slope_fit.h:176
double getDeviance() const
Gets the model deviance.
Definition slope_fit.h:166
int getPasses() const
Gets the total number of optimization iterations.
Definition slope_fit.h:196
The declaration of the Clusters class.
Namespace containing SLOPE regression implementation.
Definition clusters.h:11
std::unique_ptr< Loss > setupLoss(const std::string &loss)
Factory function to create the appropriate loss function based on the distribution family.
void validateOption(const std::string &value, const std::set< std::string > &valid_options, const std::string &parameter_name)
Validates if a given value exists in a set of valid options.
std::tuple< Eigen::VectorXd, Eigen::MatrixXd > rescaleCoefficients(const Eigen::VectorXd &beta0, const Eigen::SparseMatrix< double > &beta, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, const bool intercept)
Rescales the coefficients using the given parameters.
Functions to normalize the design matrix and rescale coefficients in case the design was normalized.
Factory function to create the appropriate loss function based on.
Various utility functions.