slope 0.28.0
Loading...
Searching...
No Matches
slope_fit.h
Go to the documentation of this file.
1
6#pragma once
7
8#include "losses/setup_loss.h"
9#include "normalize.h"
10#include "utils.h"
11#include <Eigen/Dense>
12#include <Eigen/SparseCore>
13#include <memory>
14
15namespace slope {
16
26{
27private:
28 Eigen::VectorXd intercepts;
29 Eigen::SparseMatrix<double> coefs;
30 std::vector<std::vector<int>>
31 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
46public:
48 SlopeFit() = default;
49
68 SlopeFit(const Eigen::VectorXd& intercepts,
69 const Eigen::SparseMatrix<double>& coefs,
70 const std::vector<std::vector<int>>& clusters,
71 const double alpha,
72 const Eigen::ArrayXd& lambda,
73 const double deviance,
74 const double null_deviance,
75 const std::vector<double>& primals,
76 const std::vector<double>& duals,
77 const std::vector<double>& time,
78 const int passes,
79 const std::string& centering_type,
80 const std::string& scaling_type)
81 : intercepts{ intercepts }
82 , coefs{ coefs }
83 , clusters{ clusters }
84 , alpha{ alpha }
85 , lambda{ lambda }
86 , deviance{ deviance }
87 , null_deviance{ null_deviance }
88 , primals{ primals }
89 , duals{ duals }
90 , time{ time }
91 , passes{ passes }
92 , centering_type{ centering_type }
93 , scaling_type{ scaling_type }
94 {
95 }
96
100 const Eigen::VectorXd& getIntercepts() const { return intercepts; }
101
105 const Eigen::SparseMatrix<double>& getCoefs() const { return coefs; }
106
110 const std::vector<std::vector<int>>& getClusters() const { return clusters; }
111
115 const Eigen::ArrayXd& getLambda() const { return lambda; }
116
120 double getAlpha() const { return alpha; }
121
125 double getDeviance() const { return deviance; }
126
130 const std::string& getLossType() const { return loss_type; }
131
135 double getNullDeviance() const { return null_deviance; }
136
140 const std::vector<double>& getPrimals() const { return primals; }
141
145 const std::vector<double>& getDuals() const { return duals; }
146
150 const std::vector<double>& getTime() const { return time; }
151
155 int getPasses() const { return passes; }
156
163 double getDevianceRatio() const { return 1.0 - deviance / null_deviance; }
164
171 std::vector<double> getGaps() const
172 {
173 std::vector<double> gaps(primals.size());
174 for (size_t i = 0; i < primals.size(); i++) {
175 gaps[i] = primals[i] - duals[i];
176 }
177 return gaps;
178 }
179
188 template<typename T>
189 Eigen::MatrixXd predict(T& x, const std::string& type = "response") const
190 {
191 validateOption(type, { "response", "linear" }, "type");
192
193 int n = x.rows();
194 int m = coefs.cols();
195
196 std::unique_ptr<Loss> loss = setupLoss(this->loss_type);
197
198 Eigen::VectorXd x_centers;
199 Eigen::VectorXd x_scales;
200 normalize(
201 x, x_centers, x_scales, this->centering_type, this->scaling_type, false);
202
203 Eigen::MatrixXd eta = Eigen::MatrixXd::Zero(n, m);
204
205 for (int k = 0; k < m; ++k) {
206 for (typename Eigen::SparseMatrix<double>::InnerIterator it(coefs, k); it;
207 ++it) {
208 int j = it.row();
209 eta.col(k) += x.col(j) * (it.value() / x_scales(j));
210 eta.col(k).array() -= it.value() * x_centers(j) / x_scales(j);
211 }
212
213 eta.col(k).array() += intercepts(k);
214 }
215
216 if (type == "linear") {
217 return eta;
218 }
219
220 // Return predictions
221 return loss->predict(eta);
222 }
223};
224
225} // namespace slope
A class representing the results of SLOPE (Sorted L1 Penalized Estimation) fitting.
Definition slope_fit.h:26
const std::vector< std::vector< int > > & getClusters() const
Gets the clusters.
Definition slope_fit.h:110
const Eigen::SparseMatrix< double > & getCoefs() const
Gets the sparse coefficient matrix for this fit.
Definition slope_fit.h:105
SlopeFit(const Eigen::VectorXd &intercepts, const Eigen::SparseMatrix< double > &coefs, const std::vector< std::vector< int > > &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)
Construct a new Slope Fit object.
Definition slope_fit.h:68
const Eigen::ArrayXd & getLambda() const
Gets the lambda (regularization) parameter used.
Definition slope_fit.h:115
const Eigen::VectorXd & getIntercepts() const
Gets the intercept terms for this SLOPE fit.
Definition slope_fit.h:100
double getAlpha() const
Gets the alpha (mixing) parameter used.
Definition slope_fit.h:120
const std::vector< double > & getPrimals() const
Gets the sequence of primal objective values during optimization.
Definition slope_fit.h:140
SlopeFit()=default
Default constructor.
const std::vector< double > & getDuals() const
Gets the sequence of dual objective values during optimization.
Definition slope_fit.h:145
const std::vector< double > & getTime() const
Gets the sequence of computation times during optimization.
Definition slope_fit.h:150
std::vector< double > getGaps() const
Calculate the duality gaps during optimization.
Definition slope_fit.h:171
Eigen::MatrixXd predict(T &x, const std::string &type="response") const
Predict the response for a given input matrix.
Definition slope_fit.h:189
const std::string & getLossType() const
Gets the model's loss type.
Definition slope_fit.h:130
double getDevianceRatio() const
Calculate the deviance ratio (1 - deviance/null_deviance)
Definition slope_fit.h:163
double getNullDeviance() const
Gets the null model deviance.
Definition slope_fit.h:135
double getDeviance() const
Gets the model deviance.
Definition slope_fit.h:125
int getPasses() const
Gets the total number of optimization iterations.
Definition slope_fit.h:155
Namespace containing SLOPE regression implementation.
Definition clusters.cpp:5
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.
Definition utils.cpp:7
JitNormalization normalize(Eigen::MatrixXd &x, Eigen::VectorXd &x_centers, Eigen::VectorXd &x_scales, const std::string &centering_type, const std::string &scaling_type, const bool modify_x)
Definition normalize.cpp:6
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.