slope 6.0.1
Loading...
Searching...
No Matches
loss.h
Go to the documentation of this file.
1
7#pragma once
8
9#include <Eigen/Core>
10
11namespace slope {
12
20class Loss
21{
22public:
26 virtual ~Loss() = default;
27
38 virtual double loss(const Eigen::MatrixXd& eta, const Eigen::MatrixXd& y) = 0;
39
51 virtual double dual(const Eigen::MatrixXd& theta,
52 const Eigen::MatrixXd& y,
53 const Eigen::VectorXd& w) = 0;
54
67 Eigen::MatrixXd residual(const Eigen::MatrixXd& eta,
68 const Eigen::MatrixXd& y);
69 // TODO: This function could be removed since it is just
70 // g^{-1}(eta) - y for all families
71
81 virtual Eigen::MatrixXd hessianDiagonal(const Eigen::MatrixXd& eta) = 0;
82
94 virtual void updateWeightsAndWorkingResponse(Eigen::MatrixXd& w,
95 Eigen::MatrixXd& z,
96 const Eigen::MatrixXd& eta,
97 const Eigen::MatrixXd& y);
98
104 virtual Eigen::MatrixXd preprocessResponse(const Eigen::MatrixXd& y) = 0;
105
114 virtual void updateIntercept(Eigen::VectorXd& beta0,
115 const Eigen::MatrixXd& eta,
116 const Eigen::MatrixXd& y)
117 {
118 Eigen::MatrixXd residual = this->residual(eta, y);
119 beta0 -= residual.colwise().mean() / this->lipschitz_constant;
120 };
121 // TODO: This function is not currently used anywhere and could
122 // probably be removed.
123
130 virtual Eigen::MatrixXd link(const Eigen::MatrixXd& mu) = 0;
131
138 virtual Eigen::MatrixXd inverseLink(const Eigen::MatrixXd& eta) = 0;
139
145 virtual Eigen::MatrixXd predict(const Eigen::MatrixXd& eta) = 0;
146
154 virtual double deviance(const Eigen::MatrixXd& eta, const Eigen::MatrixXd& y)
155 {
156 return 2.0 * (loss(eta, y) - loss(link(y), y));
157 }
158
159protected:
167 explicit Loss(double lipschitz_constant)
168 : lipschitz_constant(lipschitz_constant)
169 {
170 }
171
172private:
173 const double lipschitz_constant;
174};
175
176} // namespace slope
Loss function interface.
Definition loss.h:21
Loss(double lipschitz_constant)
Constructs an loss function with specified Lipschitz constant.
Definition loss.h:167
virtual double loss(const Eigen::MatrixXd &eta, const Eigen::MatrixXd &y)=0
Calculates the loss function.
virtual Eigen::MatrixXd preprocessResponse(const Eigen::MatrixXd &y)=0
Preprocess response.
virtual double deviance(const Eigen::MatrixXd &eta, const Eigen::MatrixXd &y)
Computes deviance, which is 2 times the difference between the loglikelihood of the model and the log...
Definition loss.h:154
virtual void updateIntercept(Eigen::VectorXd &beta0, const Eigen::MatrixXd &eta, const Eigen::MatrixXd &y)
Updates the intercept with a gradient descent update.
Definition loss.h:114
Eigen::MatrixXd residual(const Eigen::MatrixXd &eta, const Eigen::MatrixXd &y)
Calculates the generalized residual.
virtual void updateWeightsAndWorkingResponse(Eigen::MatrixXd &w, Eigen::MatrixXd &z, const Eigen::MatrixXd &eta, const Eigen::MatrixXd &y)
Updates weights and working response.
virtual Eigen::MatrixXd inverseLink(const Eigen::MatrixXd &eta)=0
The inverse link function, also known as the mean function.
virtual Eigen::MatrixXd predict(const Eigen::MatrixXd &eta)=0
Return predicted response.
virtual double dual(const Eigen::MatrixXd &theta, const Eigen::MatrixXd &y, const Eigen::VectorXd &w)=0
Calculates the dual loss.
virtual Eigen::MatrixXd hessianDiagonal(const Eigen::MatrixXd &eta)=0
Calculates the hessian diagonal of the loss function.
virtual ~Loss()=default
Destructor for the Loss class.
virtual Eigen::MatrixXd link(const Eigen::MatrixXd &mu)=0
The link function.
Namespace containing SLOPE regression implementation.
Definition clusters.h:11