slope 0.29.0
Loading...
Searching...
No Matches
logistic.cpp
1#include "logistic.h"
2#include "../constants.h"
3#include "../math.h"
4
5namespace slope {
6
7double
8Logistic::loss(const Eigen::MatrixXd& eta, const Eigen::MatrixXd& y)
9{
10 double loss =
11 eta.array().exp().log1p().sum() - y.reshaped().dot(eta.reshaped());
12 return loss / y.rows();
13}
14
15double
16Logistic::dual(const Eigen::MatrixXd& theta,
17 const Eigen::MatrixXd& y,
18 const Eigen::VectorXd&)
19{
20 using Eigen::log;
21
22 // Clamp probabilities to [p_min, 1-p_min]
23 Eigen::ArrayXd pr =
24 (theta + y).array().min(constants::P_MAX).max(constants::P_MIN);
25
26 return -(pr * log(pr) + (1.0 - pr) * log(1.0 - pr)).mean();
27}
28
29Eigen::MatrixXd
30Logistic::residual(const Eigen::MatrixXd& eta, const Eigen::MatrixXd& y)
31{
32 return 1.0 / (1.0 + (-eta).array().exp()) - y.array();
33}
34
35Eigen::MatrixXd
36Logistic::preprocessResponse(const Eigen::MatrixXd& y)
37{
38 // Check if the response is in {0, 1} and convert it otherwise
39 Eigen::MatrixXd y_clamped = y.array().min(1.0).max(0.0);
40
41 // Throw an error if the response is not binary
42 if ((y_clamped.array() != 0.0 && y_clamped.array() != 1.0).any()) {
43 throw std::invalid_argument("Response must be binary");
44 }
45
46 return y_clamped;
47}
48
49void
51 Eigen::VectorXd& z,
52 const Eigen::VectorXd& eta,
53 const Eigen::VectorXd& y)
54{
55 Eigen::ArrayXd pr = (1.0 / (1.0 + (-eta.array()).exp()))
57 .max(constants::P_MIN);
58 w = pr * (1.0 - pr);
59 z = eta.array() + (y.array() - pr) / w.array();
60}
61
62Eigen::MatrixXd
63Logistic::link(const Eigen::MatrixXd& mu)
64{
65 return mu.unaryExpr([](const double& x) {
66 return logit(std::clamp(x, constants::P_MIN, constants::P_MAX));
67 });
68}
69
70Eigen::MatrixXd
71Logistic::inverseLink(const Eigen::MatrixXd& eta)
72{
73 return eta.unaryExpr(
74 [](const double& x) { return 1.0 / (1.0 + std::exp(-x)); });
75}
76
77Eigen::MatrixXd
78Logistic::predict(const Eigen::MatrixXd& eta)
79{
80 Eigen::MatrixXd prob = inverseLink(eta);
81 return prob.unaryExpr([](double pr) { return pr > 0.5 ? 1.0 : 0.0; });
82}
83
84} // namespace slope
Eigen::MatrixXd preprocessResponse(const Eigen::MatrixXd &y)
Preprocesses the response for the quadratic model.
Definition logistic.cpp:36
double dual(const Eigen::MatrixXd &theta, const Eigen::MatrixXd &y, const Eigen::VectorXd &w)
Calculates the dual for the logistic loss function.
Definition logistic.cpp:16
Eigen::MatrixXd inverseLink(const Eigen::MatrixXd &eta)
The inverse link function, also known as the mean function.
Definition logistic.cpp:71
void updateWeightsAndWorkingResponse(Eigen::VectorXd &w, Eigen::VectorXd &z, const Eigen::VectorXd &eta, const Eigen::VectorXd &y)
Updates the weights and working response for the logistic loss function.
Definition logistic.cpp:50
Eigen::MatrixXd residual(const Eigen::MatrixXd &eta, const Eigen::MatrixXd &y)
Calculates the residual for the logistic loss function.
Definition logistic.cpp:30
double loss(const Eigen::MatrixXd &eta, const Eigen::MatrixXd &y)
Calculates the loss for the logistic loss function.
Definition logistic.cpp:8
Eigen::MatrixXd predict(const Eigen::MatrixXd &eta)
Return predicted response, that is 0 or 1 depending on the predicted probabilities.
Definition logistic.cpp:78
Eigen::MatrixXd link(const Eigen::MatrixXd &mu)
The link function.
Definition logistic.cpp:63
Logistic loss function implementation for SLOPE algorithm.
constexpr double P_MIN
Minimum allowed probability value to avoid numerical underflow.
Definition constants.h:30
constexpr double P_MAX
Maximum allowed probability value to avoid numerical issues near 1.0.
Definition constants.h:33
Namespace containing SLOPE regression implementation.
Definition clusters.cpp:5
T logit(const T &x)
Definition math.h:83