Multivariate Point Process Package  0.1
Classes | Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
LowRankHawkesProcess Class Reference

LowRankHawkesProcess implements the standard multivariate Hawkes process. More...

#include "include/LowRankHawkesProcess.h"

Inheritance diagram for LowRankHawkesProcess:
IProcess

Classes

struct  OPTION
 Options used to configure the fitting of the Low-rank Hawkes process. More...
 

Public Types

enum  RegCoef { LAMBDA0, LAMBDA }
 Regularization coefficients. More...
 

Public Member Functions

 LowRankHawkesProcess (const unsigned &num_rows, const unsigned &num_cols, const Eigen::VectorXd &beta)
 The constructor. More...
 
void fit (const std::vector< Sequence > &data, const OPTION &options)
 Maximum likelihood estimation for the model parameters. More...
 
void fit (const std::vector< Sequence > &data, const LowRankHawkesProcess::OPTION &options, const Eigen::MatrixXd &TrueLambda0, const Eigen::MatrixXd &TrueAlpha)
 Maximum likelihood estimation for the model parameters. More...
 
virtual void NegLoglikelihood (double &objvalue, Eigen::VectorXd &gradient)
 Negative log-likelihood of the sequences from the observed user-item pairs. More...
 
virtual void Gradient (const unsigned &k, Eigen::VectorXd &gradient)
 Returns the gradient w.r.t. the model parameters on the k-th sequence. More...
 
virtual double Intensity (const double &t, const Sequence &data, Eigen::VectorXd &intensity_dim)
 Returns the summation \(\sum_{d=1}^D\lambda^*_d(t)\) of the intensity value \(\lambda^*_d(t)\) of each dimension in a given sequence data at the time t. More...
 
virtual double IntensityUpperBound (const double &t, const double &L, const Sequence &data, Eigen::VectorXd &intensity_upper_dim)
 Returns the upper bound of the summation of the intensity value on each dimension from time t to t + L given the history of past events in sequence data. Let \({\lambda_d^*(t)}\) be the conditional intensity function on the d-th dimension where \(d=1\dotso D\), and num_dims_ = D. This function returns

\begin{align} \lambda_0^D(t) \geq \sum_{d=1}^D\sup_{\tau\in[t, t + \tau(t)]}\lambda^*_d(\tau), \end{align}

where the returned value \(\lambda_0^D(t)\) will be used for Ogata's Thinning algorithm. More...

 
virtual double IntensityIntegral (const double &lower, const double &upper, const Sequence &data)
 Returns the integral of the intensity function \(\int_{a}^b\lambda^*(\tau)d\tau\) where \(a = lower\) and \(b = upper\). More...
 
virtual double PredictNextEventTime (const Sequence &data, const unsigned &num_simulations)
 Predict the next event timing by the expectation \(\int_{t_n}^\infty tf^*(t)dt\). Currently, we use the sample average by simulations to approximate the expectation since the conditional density \(f^*(t)\) normally does not have an analytic form. More...
 
- Public Member Functions inherited from IProcess
 IProcess (const unsigned &n, const unsigned &num_dims)
 The constructor. More...
 
const Eigen::VectorXd & GetParameters ()
 Return the column vector of model parameters. More...
 
unsigned GetNumDims ()
 Return the number of dimensions in the process. More...
 
void SetParameters (const Eigen::VectorXd &v)
 Set the model parameters. More...
 
void PlotIntensityFunction (const Sequence &data)
 Plots the intensity functions based on the given sequence. It plots the intensity function and the associated event points up of each dimension in the same figure. More...
 
void PlotIntensityFunction (const Sequence &data, const unsigned &dim_id)
 Plots the intensity function and the associated event points of the dimension dim_id. More...
 

Protected Member Functions

void Initialize (const std::vector< Sequence > &data)
 Initialize the temporal features event_intensity_features_ and integral_intensity_features_. More...
 
unsigned Vec2Ind (const unsigned &i, const unsigned &j)
 Helper function returns the respective index in a column vector of the given (row, column) pair. More...
 
void Ind2Vec (const unsigned &ind, unsigned &i, unsigned &j)
 Helper function returns the respective (row, column) pair given the index in the corresponding column vector. More...
 
- Protected Member Functions inherited from IProcess
void InitializeDimension (const std::vector< Sequence > &data)
 

Protected Attributes

Eigen::VectorXd beta_
 A column vector of size \(|\Omega|\) where each component is the decaying rate of the respective exponential triggering kernel for each observed user-item pair \((u,i)\in\Omega\). More...
 
Eigen::VectorXd event_intensity_features_
 Temporal features associated with each event. More...
 
Eigen::VectorXd integral_intensity_features_
 Temporal features associated with each user-item pair. More...
 
Eigen::SparseMatrix< double > pair_event_map_
 A sparse bitmap mapping matrix. More...
 
Eigen::VectorXd observation_window_T_
 observation_window_T_ stores the observation window for each sequence. More...
 
Eigen::VectorXi observed_idx_
 Stores the index of each observed user-item pair in the total \(m\) by \(n\) grid. More...
 
unsigned num_rows_
 Total number of rows(or users). More...
 
unsigned num_cols_
 Total number of columns(or items). More...
 
LowRankHawkesProcess::OPTION options_
 A configuration object which saves the optimization options. More...
 
- Protected Attributes inherited from IProcess
Eigen::VectorXd parameters_
 A column vector represents all model parameters of the process. More...
 
unsigned num_dims_
 The total number of dimensions of the process. More...
 
std::vector< std::vector< std::vector< double > > > all_timestamp_per_dimension_
 all_timestamp_per_dimension_ is a 3-d array where all_timestamp_per_dimension_[c][n][i] records the i-th event on the n-th dimension in the c-th sequence. More...
 

Detailed Description

LowRankHawkesProcess implements the standard multivariate Hawkes process.

The Low-rank Hawkes Process is an $mn$-dimensional Hawkes process arranged in an \(m\)-by- \(n\) grid. For instance, we have \(m\) users and \(n\) items. The conditional intensity function of the entry \((u,i)\) between user \(u\) and item \(i\) is \(\lambda^{u,i}(t) = \boldsymbol{\Lambda}_0(u,i) + \boldsymbol{A}(u,i)\sum_{t^{u,i}_k < t}\exp(-\beta_{u,i}(t - t^{u,i}_k))\) where \(\boldsymbol{\Lambda}_0\) and \(\boldsymbol{A}\) are the low-rank base intensity matrix and the self-exciting matrix, respectively.

Definition at line 20 of file LowRankHawkesProcess.h.

Member Enumeration Documentation

Regularization coefficients.

Enumerator
LAMBDA0 

Regularization coefficient for \(\|\boldsymbol{\Lambda}_0\|\)

LAMBDA 

Regularization coefficient for \(\|\mathbf{A}\|\)

Definition at line 101 of file LowRankHawkesProcess.h.

Constructor & Destructor Documentation

LowRankHawkesProcess::LowRankHawkesProcess ( const unsigned &  num_rows,
const unsigned &  num_cols,
const Eigen::VectorXd &  beta 
)
inline

The constructor.

Parameters
[in]num_rowsthe total number of rows(or users).
[in]num_dimsthe total number of columns(or items).
[in]betaa column vector of size \(|\Omega|\) where each component is the decaying rate of the respective exponential triggering kernel for each observed user-item pair \((u,i)\in\Omega\).

Definition at line 159 of file LowRankHawkesProcess.h.

Member Function Documentation

void LowRankHawkesProcess::fit ( const std::vector< Sequence > &  data,
const OPTION options 
)

Maximum likelihood estimation for the model parameters.

Parameters
[in]datavectors of observed sequences.
[in]optionsdata structure containing different configurations for the optimization algorithm and the respective regularizations.

Definition at line 135 of file LowRankHawkesProcess.cc.

void LowRankHawkesProcess::fit ( const std::vector< Sequence > &  data,
const LowRankHawkesProcess::OPTION options,
const Eigen::MatrixXd &  TrueLambda0,
const Eigen::MatrixXd &  TrueAlpha 
)

Maximum likelihood estimation for the model parameters.

Parameters
[in]datavectors of observed sequences.
[in]optionsdata structure containing different configurations for the optimization algorithm and the respective regularizations.
[in]TrueLambda0ground truth \(\boldsymbol{\Lambda}_0\).
[in]TrueAlphaground truth \(\boldsymbol{A}\).

Definition at line 147 of file LowRankHawkesProcess.cc.

void LowRankHawkesProcess::Gradient ( const unsigned &  k,
Eigen::VectorXd &  gradient 
)
virtual

Returns the gradient w.r.t. the model parameters on the k-th sequence.

Parameters
[in]ksequence index.
[out]gradientthe gradient vector w.r.t. the model parameters.

Implements IProcess.

Definition at line 130 of file LowRankHawkesProcess.cc.

void LowRankHawkesProcess::Ind2Vec ( const unsigned &  ind,
unsigned &  i,
unsigned &  j 
)
protected

Helper function returns the respective (row, column) pair given the index in the corresponding column vector.

Parameters
[in]indindex in the long column vector.
[out]ireturned row index.
[out]jreturned column index.

Definition at line 21 of file LowRankHawkesProcess.cc.

void LowRankHawkesProcess::Initialize ( const std::vector< Sequence > &  data)
protected

Initialize the temporal features event_intensity_features_ and integral_intensity_features_.

Parameters
[in]datacollection of input sequences induced from each observed user-item pair.

Definition at line 28 of file LowRankHawkesProcess.cc.

double LowRankHawkesProcess::Intensity ( const double &  t,
const Sequence data,
Eigen::VectorXd &  intensity_dim 
)
virtual

Returns the summation \(\sum_{d=1}^D\lambda^*_d(t)\) of the intensity value \(\lambda^*_d(t)\) of each dimension in a given sequence data at the time t.

Parameters
[in]tthe given time.
[in]datathe given sequence of the past events until time t.
[out]intensity_dima column vector of size num_dims_ where each component stores the intensity value of the respetive dimension at time t given the past sequence in data.
Returns
the summation of the intensity value from each dimension.

Implements IProcess.

Definition at line 89 of file LowRankHawkesProcess.cc.

double LowRankHawkesProcess::IntensityIntegral ( const double &  lower,
const double &  upper,
const Sequence data 
)
virtual

Returns the integral of the intensity function \(\int_{a}^b\lambda^*(\tau)d\tau\) where \(a = lower\) and \(b = upper\).

Parameters
[in]lowerstarting point of the integral.
[in]upperending point of the integral.
[in]datasequence of past events.
Returns
\(\int_{a}^b\lambda^*(\tau)d\tau\) where \(a = lower\) and \(b = upper\).

Implements IProcess.

Definition at line 164 of file LowRankHawkesProcess.cc.

double LowRankHawkesProcess::IntensityUpperBound ( const double &  t,
const double &  L,
const Sequence data,
Eigen::VectorXd &  intensity_upper_dim 
)
virtual

Returns the upper bound of the summation of the intensity value on each dimension from time t to t + L given the history of past events in sequence data. Let \({\lambda_d^*(t)}\) be the conditional intensity function on the d-th dimension where \(d=1\dotso D\), and num_dims_ = D. This function returns

\begin{align} \lambda_0^D(t) \geq \sum_{d=1}^D\sup_{\tau\in[t, t + \tau(t)]}\lambda^*_d(\tau), \end{align}

where the returned value \(\lambda_0^D(t)\) will be used for Ogata's Thinning algorithm.

Parameters
tthe starting time.
Lthe duration.
datathe given sequence of the past events until time t.
intensity_upper_dima column vector of size num_dims_ storing the upper bound of the intensity function on each dimension from time t to t + L.
Returns
the summation of the upper-bound of each intensity function from the respetive dimension within the interval [t, t + L].

Implements IProcess.

Definition at line 94 of file LowRankHawkesProcess.cc.

void LowRankHawkesProcess::NegLoglikelihood ( double &  objvalue,
Eigen::VectorXd &  gradient 
)
virtual

Negative log-likelihood of the sequences from the observed user-item pairs.

\begin{align} - \frac{1}{|\Omega|}\sum_{\mathcal{T}^{u,i}}\bigg\{\sum_{t^{u,i}_j\in\mathcal{T}^{u,i}}\log\bigg(\boldsymbol{\Lambda}_0(u,i) + \boldsymbol{A}(u,i)\sum_{t^{u,i}_k < t^{u,i}_j}\exp(-\beta_{u,i}(t^{u,i}_j - t^{u,i}_k))\bigg) - T^{u,i}\boldsymbol{\Lambda}_0(u,i) - \boldsymbol{A}(u,i)\sum_{t^{u,i}_k < T^{u,i}}\int_{t^{u,i}_k}^{T^{u,i}}\lambda^{u,i}(t)dt\bigg\}, \end{align}

where \(\mathcal{T}^{u,i}\) is the sequence of events induced from the pair \((u,i)\).

Parameters
[out]objvaluenegative loglikelihood.
[out]gradientgradient of the parameters.

Implements IProcess.

Definition at line 99 of file LowRankHawkesProcess.cc.

double LowRankHawkesProcess::PredictNextEventTime ( const Sequence data,
const unsigned &  num_simulations 
)
virtual

Predict the next event timing by the expectation \(\int_{t_n}^\infty tf^*(t)dt\). Currently, we use the sample average by simulations to approximate the expectation since the conditional density \(f^*(t)\) normally does not have an analytic form.

Parameters
[in]datathe sequence of past events.
[in]num_simulationsnumber of simulations we use to calculate the sample average.
Returns
the prediction of the next event timing.

Implements IProcess.

Definition at line 159 of file LowRankHawkesProcess.cc.

unsigned LowRankHawkesProcess::Vec2Ind ( const unsigned &  i,
const unsigned &  j 
)
protected

Helper function returns the respective index in a column vector of the given (row, column) pair.

Parameters
[in]irow index.
[in]jcolumn index.
Returns
the index of the pair in a column vector.

Definition at line 16 of file LowRankHawkesProcess.cc.

Member Data Documentation

Eigen::VectorXd LowRankHawkesProcess::beta_
protected

A column vector of size \(|\Omega|\) where each component is the decaying rate of the respective exponential triggering kernel for each observed user-item pair \((u,i)\in\Omega\).

Definition at line 28 of file LowRankHawkesProcess.h.

Eigen::VectorXd LowRankHawkesProcess::event_intensity_features_
protected

Temporal features associated with each event.

Suppose we have observed a collection \(\Omega = \{(u,i)\}\) of user-item pairs. Each user-item pair \((u,i)\) induces a sequence of \(n_{u,i}\) events. The total number of events is thus \(N = \sum_{(u,i)\in\Omega}n_{u,i}\). event_intensity_features_ is a column vector of size \(N\), which can be regarded as the sequential concatenation of each sequence from the respective user-item pair. Each component of event_intensity_features_ is the summation of the influence from the past events to the current event in the respective sequence, that is, \(\sum_{t^{u,i}_k < t^{u,i}_j}\exp(-\beta_{u,i}(t^{u,i}_j - t^{u,i}_k))\).

Definition at line 35 of file LowRankHawkesProcess.h.

Eigen::VectorXd LowRankHawkesProcess::integral_intensity_features_
protected

Temporal features associated with each user-item pair.

Suppose we have observed a collection \(\Omega = \{(u,i)\}\) of user-item pairs. integral_intensity_features_ is a column vector of size \(|\Omega|\) where each component stores \(\sum_{t^{u,i}_k < T^{u,i}}\int_{t^{u,i}_k}^{T^{u,i}}\lambda^{u,i}(t)dt\) for the respective user-item pair \((u,i)\) where \(T^{u,i}\) is the corresponding observation window.

Definition at line 42 of file LowRankHawkesProcess.h.

unsigned LowRankHawkesProcess::num_cols_
protected

Total number of columns(or items).

Definition at line 71 of file LowRankHawkesProcess.h.

unsigned LowRankHawkesProcess::num_rows_
protected

Total number of rows(or users).

Definition at line 66 of file LowRankHawkesProcess.h.

Eigen::VectorXd LowRankHawkesProcess::observation_window_T_
protected

observation_window_T_ stores the observation window for each sequence.

For each sequence induced from the pair \((u,i)\), \(\text{observation_window_T_}(c) = T^{u,i}\).

Definition at line 56 of file LowRankHawkesProcess.h.

Eigen::VectorXi LowRankHawkesProcess::observed_idx_
protected

Stores the index of each observed user-item pair in the total \(m\) by \(n\) grid.

Definition at line 61 of file LowRankHawkesProcess.h.

LowRankHawkesProcess::OPTION LowRankHawkesProcess::options_
protected

A configuration object which saves the optimization options.

Definition at line 148 of file LowRankHawkesProcess.h.

Eigen::SparseMatrix<double> LowRankHawkesProcess::pair_event_map_
protected

A sparse bitmap mapping matrix.

pair_event_map_ is a sparse binary matrix of size \(|\Omega|\) by \(\sum_{(u,i)\in\Omega}n_{u,i}\). For each pair \((u,i)\), the \(n_{u,i}\) entries (out of the total \(\sum_{(u,i)\in\Omega}n_{u,i}\) columns) in the respective row are set to be one. Because each row is of size \(\sum_{(u,i)\in\Omega}n_{u,i}\) which is the total number of events induced from all observed user-item pairs, the non-zero entries in each row marks the location of the events belonging to the respective \((u,i)\) pair.

Definition at line 49 of file LowRankHawkesProcess.h.


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