slope 0.29.0
Loading...
Searching...
No Matches
slope::Poisson Class Reference

The Poisson class represents a Poisson regression loss function. More...

#include <poisson.h>

Inheritance diagram for slope::Poisson:
Collaboration diagram for slope::Poisson:

Public Member Functions

double loss (const Eigen::MatrixXd &eta, const Eigen::MatrixXd &y) override
 Calculates the negative log-likelihood loss for the Poisson regression.
 
double dual (const Eigen::MatrixXd &theta, const Eigen::MatrixXd &y, const Eigen::VectorXd &w) override
 Calculates the Fenchel conjugate (dual) of the Poisson loss.
 
Eigen::MatrixXd residual (const Eigen::MatrixXd &eta, const Eigen::MatrixXd &y) override
 Calculates the residual (negative gradient) for the Poisson regression.
 
void updateWeightsAndWorkingResponse (Eigen::VectorXd &w, Eigen::VectorXd &z, const Eigen::VectorXd &eta, const Eigen::VectorXd &y) override
 Updates the weights and working response for IRLS (Iteratively Reweighted Least Squares).
 
Eigen::MatrixXd preprocessResponse (const Eigen::MatrixXd &y) override
 Preprocesses the response for the Poisson model.
 
void updateIntercept (Eigen::VectorXd &beta0, const Eigen::MatrixXd &eta, const Eigen::MatrixXd &y) override
 Updates the intercept with a gradient descent update. Unlike the Quadratic and Logistic cases, the Poisson regression intercept update is not quite as simple since the gradient is not Lipschitz continuous. Instead we use a backtracking line search here.
 
Eigen::MatrixXd link (const Eigen::MatrixXd &mu) override
 The link function.
 
Eigen::MatrixXd inverseLink (const Eigen::MatrixXd &eta) override
 The inverse link function, also known as the mean function.
 
Eigen::MatrixXd predict (const Eigen::MatrixXd &eta) override
 Return predicted response, that is 0 or 1 depending on the predicted probabilities.
 
- Public Member Functions inherited from slope::Loss
virtual ~Loss ()=default
 Destructor for the Loss class.
 
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 loglikelihood of the null (intercept-only) model.
 

Additional Inherited Members

- Protected Member Functions inherited from slope::Loss
 Loss (double lipschitz_constant)
 Constructs an loss function with specified Lipschitz constant.
 

Detailed Description

The Poisson class represents a Poisson regression loss function.

The Poisson regression loss function is used for modeling count data. It assumes the response variable follows a Poisson distribution. The log-likelihood for a single observation is:

\[ \ell(y_i|\eta_i) = y_i\eta_i - e^{\eta_i} - \log(y_i!) \]

where \(\eta_i\) is the linear predictor and \(y_i\) is the observed count.

Definition at line 24 of file poisson.h.

Constructor & Destructor Documentation

◆ Poisson()

slope::Poisson::Poisson ( )
inlineexplicit

Definition at line 28 of file poisson.h.

Member Function Documentation

◆ dual()

double slope::Poisson::dual ( const Eigen::MatrixXd &  theta,
const Eigen::MatrixXd &  y,
const Eigen::VectorXd &  w 
)
overridevirtual

Calculates the Fenchel conjugate (dual) of the Poisson loss.

Parameters
thetaThe dual variables
yThe observed counts vector
wThe weights vector
Returns
The dual objective value

Implements slope::Loss.

Definition at line 13 of file poisson.cpp.

◆ inverseLink()

Eigen::MatrixXd slope::Poisson::inverseLink ( const Eigen::MatrixXd &  eta)
overridevirtual

The inverse link function, also known as the mean function.

Parameters
etaLinear predictor
Returns
\( \exp(\eta) \)

Implements slope::Loss.

Definition at line 69 of file poisson.cpp.

◆ link()

Eigen::MatrixXd slope::Poisson::link ( const Eigen::MatrixXd &  mu)
overridevirtual

The link function.

Parameters
muMean.
Returns
\( \log(\mu) \)

Implements slope::Loss.

Definition at line 63 of file poisson.cpp.

◆ loss()

double slope::Poisson::loss ( const Eigen::MatrixXd &  eta,
const Eigen::MatrixXd &  y 
)
overridevirtual

Calculates the negative log-likelihood loss for the Poisson regression.

Parameters
etaThe linear predictor vector \(\eta\)
yThe observed counts vector \(y\)
Returns
The negative log-likelihood: \(-\sum_i(y_i\eta_i - e^{\eta_i})\) (ignoring constant terms)

Implements slope::Loss.

Definition at line 7 of file poisson.cpp.

◆ predict()

Eigen::MatrixXd slope::Poisson::predict ( const Eigen::MatrixXd &  eta)
overridevirtual

Return predicted response, that is 0 or 1 depending on the predicted probabilities.

Parameters
etaThe linear predictor
Returns
The predicted response, which is the same as the inverse link in the case of the Poisson loss.

Implements slope::Loss.

Definition at line 75 of file poisson.cpp.

◆ preprocessResponse()

Eigen::MatrixXd slope::Poisson::preprocessResponse ( const Eigen::MatrixXd &  y)
overridevirtual

Preprocesses the response for the Poisson model.

Checks if the response is non-negative and throws an error otherwise.

Parameters
yPredicted values vector (n x 1)
Returns
Vector of residuals (n x 1)

Implements slope::Loss.

Definition at line 41 of file poisson.cpp.

◆ residual()

Eigen::MatrixXd slope::Poisson::residual ( const Eigen::MatrixXd &  eta,
const Eigen::MatrixXd &  y 
)
overridevirtual

Calculates the residual (negative gradient) for the Poisson regression.

Parameters
etaThe linear predictor vector \(\eta\)
yThe observed counts vector \(y\)
Returns
The residual vector: \(e^{\eta} - y\)

Implements slope::Loss.

Definition at line 25 of file poisson.cpp.

◆ updateIntercept()

void slope::Poisson::updateIntercept ( Eigen::VectorXd &  beta0,
const Eigen::MatrixXd &  eta,
const Eigen::MatrixXd &  y 
)
overridevirtual

Updates the intercept with a gradient descent update. Unlike the Quadratic and Logistic cases, the Poisson regression intercept update is not quite as simple since the gradient is not Lipschitz continuous. Instead we use a backtracking line search here.

Parameters
beta0The current intercept
etaThe current linear predictor
yThe observed counts vector

Reimplemented from slope::Loss.

Definition at line 51 of file poisson.cpp.

◆ updateWeightsAndWorkingResponse()

void slope::Poisson::updateWeightsAndWorkingResponse ( Eigen::VectorXd &  w,
Eigen::VectorXd &  z,
const Eigen::VectorXd &  eta,
const Eigen::VectorXd &  y 
)
overridevirtual

Updates the weights and working response for IRLS (Iteratively Reweighted Least Squares).

Parameters
wThe weights vector: \(w = e^{\eta}\)
zThe working response vector: \(z = \eta + (y - e^{\eta})/e^{\eta}\)
etaThe current linear predictor
yThe observed counts vector

Implements slope::Loss.

Definition at line 31 of file poisson.cpp.


The documentation for this class was generated from the following files: