Multivariate Point Process Package
0.1
|
HawkesGeneralKernel implements the multivariate Hawkes process with customized triggering kernels. More...
#include "include/HawkesGeneralKernel.h"
Classes | |
struct | OPTION |
Options used to configure the fitting of the general Hawkes Process with customized triggering kernels. More... | |
Public Types | |
enum | Regularizer { L1, L22, NUCLEAR, NONE } |
Supported regularizations used to fit the standard Hawkes Process. More... | |
enum | RegCoef { LAMBDA0, LAMBDA } |
Regularization coefficients. More... | |
enum | OptMethod { SGD, PLBFGS } |
Optimization algorithms used to fit the standard Hawkes Process. More... | |
Public Member Functions | |
HawkesGeneralKernel (const unsigned &n, const unsigned &num_dims, const std::vector< std::vector< TriggeringKernel * > > &triggeringkernels) | |
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... | |
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 | RestoreOptionToDefault () |
restore to the default optimization configuration More... | |
Protected Member Functions inherited from IProcess | |
void | InitializeDimension (const std::vector< Sequence > &data) |
Protected Attributes | |
std::vector< std::vector< Eigen::MatrixXd > > | arrayK |
the temporal features associated with the intensity More... | |
std::vector< Eigen::MatrixXd > | arrayG |
summation of the integral of the triggering kernels More... | |
std::vector< std::vector< TriggeringKernel * > > | triggeringkernels_ |
a \(D\) by \(D\) grid of triggering-kernels 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... | |
SimpleRNG | RNG_ |
internal implementation for random number generator More... | |
unsigned | num_sequences_ |
total number of observed sequences More... | |
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... | |
HawkesGeneralKernel implements the multivariate Hawkes process with customized triggering kernels.
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 generally defined as the following:
\begin{align} \lambda^n(t) = \lambda_0^n + \sum_{m=1}^D\alpha_{mn}\sum_{t^m_j < t}\gamma(t - t^m_j), \end{align}
where \(\lambda_0^n\geq 0\) is the base intensity, \(D\) is the number of dimensions, \(\alpha_{mn}\geq 0\), and the triggering kernel \(\gamma(t - t^m_j)\) captures the extent to which an event on dimension m at the time \(t^m_j\) can trigger an event on dimension n in the near future. Normally, in the standard Hawkes process, we have \(\gamma(t - t^m_j) = \exp(-\beta_{mn}(t - t^m_j))\). However, in more general cases, the form of the triggering kernel can be formulated to catpure the phenomena of interest. The collection of \(\{\alpha_{mn}\}\) can be represented as a matrix \(\mathbf{A}(m,n) = \alpha_{mn}\), and the collection of \(\{\lambda_0^n\}\) can be represented as a column vector \(\boldsymbol{\lambda}_0\).
Definition at line 25 of file HawkesGeneralKernel.h.
Optimization algorithms used to fit the standard Hawkes Process.
Enumerator | |
---|---|
SGD |
stochastic gradient descend. |
PLBFGS |
Definition at line 125 of file HawkesGeneralKernel.h.
Regularization coefficients.
Enumerator | |
---|---|
LAMBDA0 |
Regularization coefficient for \(\|\boldsymbol{\lambda}_0\|\) |
LAMBDA |
Regularization coefficient for \(\|\mathbf{A}\|\) |
Definition at line 111 of file HawkesGeneralKernel.h.
Supported regularizations used to fit the standard Hawkes Process.
Definition at line 90 of file HawkesGeneralKernel.h.
|
inline |
The constructor.
[in] | n | the number of parameters in total. |
[in] | num_dims | the number of dimensions in the process. |
[in] | triggeringkernels | a D-by-D grid of triggering-kernels. |
Definition at line 173 of file HawkesGeneralKernel.h.
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 70 of file HawkesGeneralKernel.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 189 of file HawkesGeneralKernel.cc.
|
protected |
initialize the temporal features arrayK and arrayG from the input sequences
[in] | data | input collection of sequences |
Definition at line 17 of file HawkesGeneralKernel.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 196 of file HawkesGeneralKernel.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 266 of file HawkesGeneralKernel.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 231 of file HawkesGeneralKernel.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\alpha_{mn}\underbrace{\sum_{t^m_{i,c}<t^n_{i,c}}\gamma_{mn}(t^n_{i,c} - t^m_{i,c})}_{\text{arrayK[n][c](i,m)}})\bigg) - T_c\lambda_0^n - \sum_{m=1}^D\alpha_{mn}\underbrace{\sum_{t^m_{j,c} < T_c}\int_{t^m_j}^{T_c}\gamma_{mn}(t - t^m_{j,c})dt)}_{\text{arrayG[n](c,m)}}\bigg)\bigg\}. \end{align}
[out] | objvalue | negative loglikelihood. |
[out] | gradient | gradient of the parameters. |
Implements IProcess.
Definition at line 83 of file HawkesGeneralKernel.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 309 of file HawkesGeneralKernel.cc.
|
protected |
restore to the default optimization configuration
Definition at line 321 of file HawkesGeneralKernel.cc.
|
protected |
summation of the integral of the triggering kernels
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\alpha_{mn}\underbrace{\sum_{t^m_{i,c}<t^n_{i,c}}\gamma_{mn}(t^n_{i,c} - t^m_{i,c})}_{\text{arrayK[n][c](i,m)}})\bigg) - T_c\lambda_0^n - \sum_{m=1}^D\alpha_{mn}\underbrace{\sum_{t^m_{j,c} < T_c}\int_{t^m_j}^{T_c}\gamma_{mn}(t - t^m_{j,c})dt)}_{\text{arrayG[n](c,m)}}\bigg)\bigg\}. \end{align}
arrayG[n] is an \(C\) by \(D\) matrix where \(n_c\) is the number of events on the nth dimension in the sequence c. arrayK[n][c](i,m) stores the cumulative influence of the past events on dimension \(m\) in the sequence \(c\) to the occurence of the \(i\)th event on dimension \(n\).
Definition at line 52 of file HawkesGeneralKernel.h.
|
protected |
the temporal features associated with 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\alpha_{mn}\underbrace{\sum_{t^m_{j,c}<t^n_{i,c}}\gamma_{mn}(t^n_{i,c} - t^m_{j,c})}_{\text{arrayK[n][c](i,m)}})\bigg) - T_c\lambda_0^n - \sum_{m=1}^D\alpha_{mn}\underbrace{\sum_{t^m_{j,c} < T_c}\int_{t^m_j}^{T_c}\gamma_{mn}(t - t^m_{j,c})dt)}_{\text{arrayG[n](c,m)}}\bigg)\bigg\}. \end{align}
arrayK[n][c] is an \(n_c\) by \(D\) matrix where \(n_c\) is the number of events on the nth dimension in the sequence c. arrayK[n][c](i,m) stores the cumulative influence of the past events on dimension \(m\) in the sequence \(c\) to the occurence of the \(i\)th event on dimension \(n\).
Definition at line 40 of file HawkesGeneralKernel.h.
|
protected |
total number of observed sequences
Definition at line 73 of file HawkesGeneralKernel.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 63 of file HawkesGeneralKernel.h.
|
protected |
A configuration object which saves the optimization options.
Definition at line 163 of file HawkesGeneralKernel.h.
|
protected |
internal implementation for random number generator
Definition at line 68 of file HawkesGeneralKernel.h.
|
protected |
a \(D\) by \(D\) grid of triggering-kernels
Definition at line 58 of file HawkesGeneralKernel.h.