slope 0.29.0
Loading...
Searching...
No Matches
slope::Loss Class Referenceabstract

Loss function interface. More...

#include <loss.h>

Inheritance diagram for slope::Loss:

Public Member Functions

virtual ~Loss ()=default
 Destructor for the Loss class.
 
virtual double loss (const Eigen::MatrixXd &eta, const Eigen::MatrixXd &y)=0
 Calculates the loss function.
 
virtual double dual (const Eigen::MatrixXd &theta, const Eigen::MatrixXd &y, const Eigen::VectorXd &w)=0
 Calculates the dual loss.
 
virtual Eigen::MatrixXd residual (const Eigen::MatrixXd &eta, const Eigen::MatrixXd &y)=0
 Calculates the residual.
 
virtual void updateWeightsAndWorkingResponse (Eigen::VectorXd &w, Eigen::VectorXd &z, const Eigen::VectorXd &eta, const Eigen::VectorXd &y)=0
 Updates the weights and working response.
 
virtual Eigen::MatrixXd preprocessResponse (const Eigen::MatrixXd &y)=0
 Preprocess response.
 
virtual void updateIntercept (Eigen::VectorXd &beta0, const Eigen::MatrixXd &eta, const Eigen::MatrixXd &y)
 Updates the intercept with a gradient descent update. Also updates the linear predictor (but not the residual).
 
virtual Eigen::MatrixXd link (const Eigen::MatrixXd &mu)=0
 The link function.
 
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 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.
 

Protected Member Functions

 Loss (double lipschitz_constant)
 Constructs an loss function with specified Lipschitz constant.
 

Detailed Description

Loss function interface.

This class defines the interface for an loss function, which is used in optimization algorithms. The loss function calculates the loss, dual, residual, and updates the weights and working response.

Definition at line 20 of file loss.h.

Constructor & Destructor Documentation

◆ Loss()

slope::Loss::Loss ( double  lipschitz_constant)
inlineexplicitprotected

Constructs an loss function with specified Lipschitz constant.

Parameters
lipschitz_constantThe Lipschitz constant for the loss function

The Lipschitz constant is used to ensure convergence in gradient-based optimization by bounding the rate of change of the gradient.

Definition at line 150 of file loss.h.

Member Function Documentation

◆ deviance()

virtual double slope::Loss::deviance ( const Eigen::MatrixXd &  eta,
const Eigen::MatrixXd &  y 
)
inlinevirtual

Computes deviance, which is 2 times the difference between the loglikelihood of the model and the loglikelihood of the null (intercept-only) model.

Parameters
etaThe predicted values.
yThe true values.

Definition at line 137 of file loss.h.

◆ dual()

virtual double slope::Loss::dual ( const Eigen::MatrixXd &  theta,
const Eigen::MatrixXd &  y,
const Eigen::VectorXd &  w 
)
pure virtual

Calculates the dual loss.

This function calculates the dual function given the estimated parameters (theta) and the true values (y).

Parameters
thetaThe estimated parameters.
yThe true values.
wWeights.
Returns
The dual value.

Implemented in slope::Logistic, slope::Multinomial, slope::Quadratic, and slope::Poisson.

◆ inverseLink()

virtual Eigen::MatrixXd slope::Loss::inverseLink ( const Eigen::MatrixXd &  eta)
pure virtual

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

Parameters
etaMean
Returns
The result of applying the inverse link function, which depends on the loss function used.

Implemented in slope::Logistic, slope::Multinomial, slope::Quadratic, and slope::Poisson.

◆ link()

virtual Eigen::MatrixXd slope::Loss::link ( const Eigen::MatrixXd &  mu)
pure virtual

The link function.

Parameters
muMean.
Returns
The result of applying the link function, which depends on the loss function used.

Implemented in slope::Logistic, slope::Multinomial, slope::Quadratic, and slope::Poisson.

◆ loss()

virtual double slope::Loss::loss ( const Eigen::MatrixXd &  eta,
const Eigen::MatrixXd &  y 
)
pure virtual

Calculates the loss function.

This function calculates the loss function given the predicted values (eta) and the true values (y).

Parameters
etaThe predicted values.
yThe true values.
Returns
The loss value.

Implemented in slope::Logistic, slope::Multinomial, slope::Quadratic, and slope::Poisson.

◆ predict()

virtual Eigen::MatrixXd slope::Loss::predict ( const Eigen::MatrixXd &  eta)
pure virtual

Return predicted response.

Parameters
etaThe linear predictor
Returns
The predicted response.

Implemented in slope::Logistic, slope::Multinomial, slope::Quadratic, and slope::Poisson.

◆ preprocessResponse()

virtual Eigen::MatrixXd slope::Loss::preprocessResponse ( const Eigen::MatrixXd &  y)
pure virtual

Preprocess response.

Parameters
yThe response

Implemented in slope::Logistic, slope::Multinomial, slope::Quadratic, and slope::Poisson.

◆ residual()

virtual Eigen::MatrixXd slope::Loss::residual ( const Eigen::MatrixXd &  eta,
const Eigen::MatrixXd &  y 
)
pure virtual

Calculates the residual.

This function calculates the residual given the predicted values (eta) and the true values (y).

Parameters
etaThe predicted values.
yThe true values.
Returns
The residual vector.

Implemented in slope::Logistic, slope::Multinomial, slope::Quadratic, and slope::Poisson.

◆ updateIntercept()

virtual void slope::Loss::updateIntercept ( Eigen::VectorXd &  beta0,
const Eigen::MatrixXd &  eta,
const Eigen::MatrixXd &  y 
)
inlinevirtual

Updates the intercept with a gradient descent update. Also updates the linear predictor (but not the residual).

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

Reimplemented in slope::Poisson.

Definition at line 99 of file loss.h.

◆ updateWeightsAndWorkingResponse()

virtual void slope::Loss::updateWeightsAndWorkingResponse ( Eigen::VectorXd &  w,
Eigen::VectorXd &  z,
const Eigen::VectorXd &  eta,
const Eigen::VectorXd &  y 
)
pure virtual

Updates the weights and working response.

This function updates the weights and working response given the predicted values (eta) and the true values (y).

Parameters
wThe weights to be updated.
zThe working response to be updated.
etaThe predicted values.
yThe true values.

Implemented in slope::Logistic, slope::Multinomial, slope::Quadratic, and slope::Poisson.


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