slope 0.29.0
Loading...
Searching...
No Matches
normalize.h
Go to the documentation of this file.
1
6#pragma once
7
8#include "jit_normalization.h"
9#include "math.h"
10#include <Eigen/SparseCore>
11
12namespace slope {
13
31template<typename T>
32void
33computeCenters(Eigen::VectorXd& x_centers, const T& x, const std::string& type)
34{
35 int p = x.cols();
36
37 if (type == "manual") {
38 if (x_centers.size() != p) {
39 throw std::invalid_argument("Invalid dimensions in centers");
40 }
41
42 if (!x_centers.allFinite()) {
43 throw std::invalid_argument("Centers must be finite");
44 }
45
46 } else if (type == "mean") {
47 x_centers = means(x);
48 } else if (type == "min") {
49 x_centers = mins(x);
50 } else if (type != "none") {
51 throw std::invalid_argument("Invalid centering type");
52 }
53}
54
72template<typename T>
73void
74computeScales(Eigen::VectorXd& x_scales, const T& x, const std::string& type)
75{
76 int p = x.cols();
77
78 if (type == "manual") {
79 if (x_scales.size() != p) {
80 throw std::invalid_argument("Invalid dimensions in scales");
81 }
82 if (!x_scales.allFinite()) {
83 throw std::invalid_argument("Scales must be finite");
84 }
85 } else if (type == "sd") {
86 x_scales = stdDevs(x);
87 } else if (type == "l1") {
88 x_scales = l1Norms(x);
89 } else if (type == "l2") {
90 x_scales = l2Norms(x);
91 } else if (type == "max_abs") {
92 x_scales = maxAbs(x);
93 } else if (type == "range") {
94 x_scales = ranges(x);
95 } else if (type != "none") {
96 throw std::invalid_argument("Invalid scaling type");
97 }
98}
99
122normalize(Eigen::MatrixXd& x,
123 Eigen::VectorXd& x_centers,
124 Eigen::VectorXd& x_scales,
125 const std::string& centering_type,
126 const std::string& scaling_type,
127 const bool modify_x);
128
153normalize(Eigen::SparseMatrix<double>& x,
154 Eigen::VectorXd& x_centers,
155 Eigen::VectorXd& x_scales,
156 const std::string& centering_type,
157 const std::string& scaling_type,
158 const bool modify_x);
159
180std::tuple<Eigen::VectorXd, Eigen::MatrixXd>
181rescaleCoefficients(const Eigen::VectorXd& beta0,
182 const Eigen::VectorXd& beta,
183 const Eigen::VectorXd& x_centers,
184 const Eigen::VectorXd& x_scales,
185 const bool intercept);
186
187} // namespace slope
Enums to control predictor standardization behavior.
Mathematical support functions for the slope package.
Namespace containing SLOPE regression implementation.
Definition clusters.cpp:5
Eigen::VectorXd l1Norms(const T &x)
Computes the L1 (Manhattan) norms for each column of a matrix.
Definition math.h:439
Eigen::VectorXd maxAbs(const Eigen::SparseMatrix< double > &x)
Computes the maximum absolute value for each column of a matrix.
Definition math.cpp:126
Eigen::VectorXd means(const Eigen::SparseMatrix< double > &x)
Computes the arithmetic mean for each column of a sparse matrix.
Definition math.cpp:76
JitNormalization
Enums to control predictor standardization behavior.
Eigen::VectorXd stdDevs(const Eigen::SparseMatrix< double > &x)
Computes the standard deviation for each column of a matrix.
Definition math.cpp:180
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
Eigen::VectorXd ranges(const Eigen::SparseMatrix< double > &x)
Computes the range (max - min) for each column of a matrix.
Definition math.cpp:97
void computeScales(Eigen::VectorXd &x_scales, const T &x, const std::string &type)
Definition normalize.h:74
Eigen::VectorXd l2Norms(const Eigen::SparseMatrix< double > &x)
Computes the L2 (Euclidean) norms for each column of a sparse matrix.
Definition math.cpp:56
Eigen::VectorXd mins(const Eigen::SparseMatrix< double > &x)
Computes the minimum value for each column of a sparse matrix.
Definition math.cpp:153
void computeCenters(Eigen::VectorXd &x_centers, const T &x, const std::string &type)
Definition normalize.h:33
std::tuple< Eigen::VectorXd, Eigen::MatrixXd > rescaleCoefficients(const Eigen::VectorXd &beta0, const Eigen::VectorXd &beta, const Eigen::VectorXd &x_centers, const Eigen::VectorXd &x_scales, const bool intercept)
Rescales the coefficients using the given parameters.
Definition normalize.cpp:94