Multivariate Point Process Package
0.1
|
HawkesLearningTriggeringKernel implements the multivariate Hawkes process where the triggering kernel can be learned from the data. More...
#include "include/HawkesLearningTriggeringKernel.h"
Classes | |
struct | OPTION |
Options used to configure the fitting of the general Hawkes Process with learned triggering kernels. More... | |
Public Types | |
enum | Regularizer { L22, GROUP, NONE } |
Supported regularizations used to fit the standard Hawkes Process. More... | |
enum | RegCoef { LAMBDA0, LAMBDA } |
Regularization coefficients. More... | |
Public Member Functions | |
HawkesLearningTriggeringKernel (const unsigned &n, const unsigned &num_dims, const Eigen::VectorXd &tau, const Eigen::VectorXd &sigma) | |
The constructor. More... | |
HawkesLearningTriggeringKernel (const unsigned &n, const unsigned &num_dims, const Graph *graph, const Eigen::VectorXd &tau, const Eigen::VectorXd &sigma) | |
The constructor. More... | |
void | fit (const std::vector< Sequence > &data, const OPTION &options) |
Maximum likelihood estimation for the model parameters. More... | |
virtual void | NegLoglikelihood (double &objvalue, Eigen::VectorXd &gradient) |
Negative loglikelihood of general Hawkes process. 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... | |
void | PlotTriggeringKernel (const unsigned &dim_m, const unsigned &dim_n, const double &T, const double &delta) |
Visualize the learned triggering kernel between dimension dim_m and dimension dim_n. 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 arrayK and arrayG from the input sequences. More... | |
void | InitializeWithGraph (const std::vector< Sequence > &data) |
Initialize the temporal features arrayK and arrayG from the input sequences where the dependency structure among the dimensions is given. More... | |
void | InitializeConstants () |
Initialize all helper variables. More... | |
void | PostProcessing () |
Post process the learned dependency structure. More... | |
void | GetNegLoglikelihood (double &objvalue, Eigen::VectorXd &gradient) |
Compute the negative loglikelihood. More... | |
Protected Member Functions inherited from IProcess | |
void | InitializeDimension (const std::vector< Sequence > &data) |
Protected Attributes | |
const double | PI = 3.14159265358979323846 |
Eigen::VectorXd | tau_ |
A column vector of size num_rbfs_. Each component of tau_ is the location of the respective RBF basis function. More... | |
unsigned | num_rbfs_ |
Total number of RBF basis functions. More... | |
Eigen::VectorXd | sigma_ |
A column vector of size num_rbfs_. Each component of sigma_ is the bandwidth of the respective RBF basis function. More... | |
Eigen::VectorXd | sqrt2sigma_ |
Helper variable defined as \(\sqrt{2}\)sigma_. More... | |
Eigen::VectorXd | sqrt2PIsigma_ |
Helper variable defined as \(0.5\sqrt{2\pi}\)sigma_. More... | |
Eigen::VectorXd | erfctau_sigma_ |
Helper variable defined as \(\text{erfc}\bigg(\frac{tau\_}{\sqrt{2}\text{sigma_}}\bigg)\). More... | |
std::vector< std::vector< std::vector< Eigen::MatrixXd > > > | arrayK |
Temporal features associated with the intensity function. More... | |
std::vector< std::vector< Eigen::MatrixXd > > | arrayG |
Temporal features derived from the integral of the intensity. More... | |
Eigen::VectorXd | observation_window_T_ |
A column vector of length \(C\) which is the total number of sequences. Each component records the observation window in the respective sequence. More... | |
unsigned | num_sequences_ |
Total number of observed sequences. More... | |
OPTION | options_ |
A configuration object which saves the optimization options. More... | |
const Graph * | graph_ |
A graph object represents the dependency structure among the dimensions. 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... | |
HawkesLearningTriggeringKernel implements the multivariate Hawkes process where the triggering kernel can be learned from the data.
A multivariate Hawkes Process is a process where the occurrence of an event to a dimension will trigger more events on this dimension and other related dimensions in the near future. The intensity function of each dimension of the Hawkes process is defined as the following:
\begin{align} \lambda^n(t) = \lambda_0^n + \sum_{m=1}^D\sum_{t^m_j < t}\gamma_{mn}(t, t^m_j), \end{align}
where \(\lambda_0^n\geq 0\) is the base intensity, \(D\) is the number of dimensions, and \(\gamma_{mn}(t, t^m_j)\geq 0\) captures the extent to which an event on dimension \(m\) at the time \(t^m_j\) can trigger an event on dimension \(n\) at the time \(t\). We parametrize \(\gamma_{mn}(t, t^m_j) = \boldsymbol{k}_{mn}(t - t^m_j)^\top\boldsymbol{\alpha}_{mn}\). \(\boldsymbol{k}_{mn}(t)\) is a column vector where the \(l\)-th component is defined as \(\boldsymbol{k}_{mnl}(t)=\exp(-(\frac{\tau_l - t}{\sqrt{2}\sigma_l})^2)\), and \(\tau_l\) and \(\sigma_l\) are the parameters of the basis function.
Definition at line 28 of file HawkesLearningTriggeringKernel.h.
Regularization coefficients.
Enumerator | |
---|---|
LAMBDA0 |
Regularization coefficient for \(\|\boldsymbol{\lambda}_0\|\) |
LAMBDA |
Regularization coefficient for GROUP LASSO |
Definition at line 55 of file HawkesLearningTriggeringKernel.h.
Supported regularizations used to fit the standard Hawkes Process.
Enumerator | |
---|---|
L22 |
L22 norm \(\|\cdot\|_2^2\) |
GROUP |
GROUP LASSO \(\sum_{m}\|\boldsymbol{\alpha}_{mn}\|_2\) |
NONE |
No regularization |
Definition at line 36 of file HawkesLearningTriggeringKernel.h.
|
inline |
The constructor.
[in] | n | the number of parameters in total. |
[in] | num_dims | the number of dimensions in the process. |
[in] | tau | the location of each basis function. |
[in] | sigma | the bandwidth of each basis function. |
Definition at line 188 of file HawkesLearningTriggeringKernel.h.
|
inline |
The constructor.
[in] | n | the number of parameters in total. |
[in] | num_dims | the number of dimensions in the process. |
[in] | graph | the graph object representing the dependency structure among the dimensions. |
[in] | tau | the location of each basis function. |
[in] | sigma | the bandwidth of each basis function. |
Definition at line 201 of file HawkesLearningTriggeringKernel.h.
void HawkesLearningTriggeringKernel::fit | ( | const std::vector< Sequence > & | data, |
const OPTION & | options | ||
) |
Maximum likelihood estimation for the model parameters.
[in] | data | vectors of observed sequences. |
[in] | options | data structure sotring different configuration for the optimization algorithm and the respective regularizations. |
Definition at line 230 of file HawkesLearningTriggeringKernel.cc.
|
protected |
Compute the negative loglikelihood.
[out] | objvalue | negative loglikelihood. |
[out] | gradient | gradient of the parameters. |
Definition at line 268 of file HawkesLearningTriggeringKernel.cc.
|
virtual |
Returns the gradient w.r.t. the model parameters on the k-th sequence.
[in] | k | sequence index. |
[out] | gradient | the gradient vector w.r.t. the model parameters. |
Implements IProcess.
Definition at line 392 of file HawkesLearningTriggeringKernel.cc.
|
protected |
Initialize the temporal features arrayK and arrayG from the input sequences.
[in] | data | input collection of sequences |
Definition at line 26 of file HawkesLearningTriggeringKernel.cc.
|
protected |
Initialize all helper variables.
Definition at line 13 of file HawkesLearningTriggeringKernel.cc.
|
protected |
Initialize the temporal features arrayK and arrayG from the input sequences where the dependency structure among the dimensions is given.
[in] | data | input collection of sequences |
Definition at line 97 of file HawkesLearningTriggeringKernel.cc.
|
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.
[in] | t | the given time. |
[in] | data | the given sequence of the past events until time t. |
[out] | intensity_dim | a 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. |
Implements IProcess.
Definition at line 399 of file HawkesLearningTriggeringKernel.cc.
|
virtual |
Returns the integral of the intensity function \(\int_{a}^b\lambda^*(\tau)d\tau\) where \(a = lower\) and \(b = upper\).
[in] | lower | starting point of the integral. |
[in] | upper | ending point of the integral. |
[in] | data | sequence of past events. |
Implements IProcess.
Definition at line 413 of file HawkesLearningTriggeringKernel.cc.
|
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.
t | the starting time. |
L | the duration. |
data | the given sequence of the past events until time t. |
intensity_upper_dim | a column vector of size num_dims_ storing the upper bound of the intensity function on each dimension from time t to t + L. |
Implements IProcess.
Definition at line 406 of file HawkesLearningTriggeringKernel.cc.
|
virtual |
Negative loglikelihood of general Hawkes process.
\begin{align} - \sum_{n=1}^D\bigg\{\frac{1}{C}\sum_{c=1}^C\bigg(\sum_{i = 1}^{n_c}\bigg(\log(\lambda^n_0 + \sum_{m=1}^D\sum_{t^m_{j,c}<t^n_{i,c}}\underbrace{\boldsymbol{k}_{mn}(t^n_{i,c} - t^m_{j,c})^\top}_{\text{arrayK[n][c][i](m,:)}}\boldsymbol{\alpha}_{mn})\bigg) - T_c\lambda_0^n - \sum_{m=1}^D\sum_{t^m_{j,c} < T_c}\underbrace{\int_{t^m_j}^{T_c}\boldsymbol{k}_{mn}(t - t^m_{j,c})^\top dt}_{\text{arrayG[n][c](m,:)}}\cdot\boldsymbol{\alpha}_{mn}\bigg)\bigg\}. \end{align}
[out] | objvalue | negative loglikelihood. |
[out] | gradient | gradient of the parameters. |
Implements IProcess.
Definition at line 333 of file HawkesLearningTriggeringKernel.cc.
void HawkesLearningTriggeringKernel::PlotTriggeringKernel | ( | const unsigned & | dim_m, |
const unsigned & | dim_n, | ||
const double & | T, | ||
const double & | delta | ||
) |
Visualize the learned triggering kernel between dimension dim_m and dimension dim_n.
dim_m | the influential dimension |
dim_n | the influenced dimension |
T | observation window. |
delta | unit step size to partition the time line. |
Definition at line 424 of file HawkesLearningTriggeringKernel.cc.
|
protected |
Post process the learned dependency structure.
Definition at line 171 of file HawkesLearningTriggeringKernel.cc.
|
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.
[in] | data | the sequence of past events. |
[in] | num_simulations | number of simulations we use to calculate the sample average. |
Implements IProcess.
Definition at line 419 of file HawkesLearningTriggeringKernel.cc.
|
protected |
Temporal features derived from the integral of the intensity.
The log-likelihood of observing a collection of C sequences can be derived as the following:
\begin{align} \sum_{n=1}^D\bigg\{\frac{1}{C}\sum_{c=1}^C\bigg(\sum_{i = 1}^{n_c}\bigg(\log(\lambda^n_0 + \sum_{m=1}^D\sum_{t^m_{j,c}<t^n_{i,c}}\underbrace{\boldsymbol{k}_{mn}(t^n_{i,c} - t^m_{j,c})^\top}_{\text{arrayK[n][c][i](m,:)}}\boldsymbol{\alpha}_{mn})\bigg) - T_c\lambda_0^n - \sum_{m=1}^D\sum_{t^m_{j,c} < T_c}\underbrace{\int_{t^m_j}^{T_c}\boldsymbol{k}_{mn}(t - t^m_{j,c})^\top dt}_{\text{arrayG[n][c](m,:)}}\cdot\boldsymbol{\alpha}_{mn}\bigg)\bigg\}. \end{align}
arrayG[n][c] is a \(D\) by num_rbfs_ matrix which stores the integral of the influence from dimension \(m\) to dimension \(n\) in the sequence c.
Definition at line 137 of file HawkesLearningTriggeringKernel.h.
|
protected |
Temporal features associated with the intensity function.
The log-likelihood of observing a collection of C sequences can be derived as the following:
\begin{align} \sum_{n=1}^D\bigg\{\frac{1}{C}\sum_{c=1}^C\bigg(\sum_{i = 1}^{n_c}\bigg(\log(\lambda^n_0 + \sum_{m=1}^D\sum_{t^m_{j,c}<t^n_{i,c}}\underbrace{\boldsymbol{k}_{mn}(t^n_{i,c} - t^m_{j,c})^\top}_{\text{arrayK[n][c][i](m,:)}}\boldsymbol{\alpha}_{mn})\bigg) - T_c\lambda_0^n - \sum_{m=1}^D\sum_{t^m_{j,c} < T_c}\underbrace{\int_{t^m_j}^{T_c}\boldsymbol{k}_{mn}(t - t^m_{j,c})^\top dt}_{\text{arrayG[n][c](m,:)}}\cdot\boldsymbol{\alpha}_{mn}\bigg)\bigg\}. \end{align}
arrayK[n][c][i] is a \(D\) by num_rbfs_ matrix which stores the influence from the event at the time \(t^m_{j,c}\) on dimension \(m\) to the event at the time \(t^n_{i,c}\) on dimension \(n\) in the sequence c.
Definition at line 126 of file HawkesLearningTriggeringKernel.h.
|
protected |
Helper variable defined as \(\text{erfc}\bigg(\frac{tau\_}{\sqrt{2}\text{sigma_}}\bigg)\).
Definition at line 115 of file HawkesLearningTriggeringKernel.h.
|
protected |
A graph object represents the dependency structure among the dimensions.
Definition at line 177 of file HawkesLearningTriggeringKernel.h.
|
protected |
Total number of RBF basis functions.
Definition at line 99 of file HawkesLearningTriggeringKernel.h.
|
protected |
Total number of observed sequences.
Definition at line 145 of file HawkesLearningTriggeringKernel.h.
|
protected |
A column vector of length \(C\) which is the total number of sequences. Each component records the observation window in the respective sequence.
Definition at line 141 of file HawkesLearningTriggeringKernel.h.
|
protected |
A configuration object which saves the optimization options.
Definition at line 173 of file HawkesLearningTriggeringKernel.h.
|
protected |
Definition at line 89 of file HawkesLearningTriggeringKernel.h.
|
protected |
A column vector of size num_rbfs_. Each component of sigma_ is the bandwidth of the respective RBF basis function.
Definition at line 103 of file HawkesLearningTriggeringKernel.h.
|
protected |
Helper variable defined as \(0.5\sqrt{2\pi}\)sigma_.
Definition at line 111 of file HawkesLearningTriggeringKernel.h.
|
protected |
Helper variable defined as \(\sqrt{2}\)sigma_.
Definition at line 107 of file HawkesLearningTriggeringKernel.h.
|
protected |
A column vector of size num_rbfs_. Each component of tau_ is the location of the respective RBF basis function.
Definition at line 94 of file HawkesLearningTriggeringKernel.h.