Multivariate Point Process Package  0.1
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
IProcess Interface Referenceabstract

IProcess defines a general interface for each specific point process. More...

#include "include/Process.h"

Inheritance diagram for IProcess:
HawkesGeneralKernel HawkesLearningTriggeringKernel LowRankHawkesProcess PlainHawkes PlainTerminating Poisson SelfInhibitingProcess TerminatingProcessLearningTriggeringKernel

Public Member Functions

 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...
 
virtual void NegLoglikelihood (double &objvalue, Eigen::VectorXd &Gradient)=0
 Returns the negative loglikelihood value and the gradient vectors. Given a collection of sequences \(\mathcal{C} = \{\mathcal{S}^c\}\) where \(\mathcal{S}^c = \{t^{c,d}_i\}^{d=1\dotso D}_{i=1\dotso n^c_d}\), and \(n^c_d\) is the number of events on dimension d in sequence \(\mathcal{S}^c\). The negative loglihood of observing \(\mathcal{C}\) is defined as

\begin{align} \ell(\mathcal{C}) = - \sum_{d=1}^D\frac{1}{|\mathcal{C}|}\sum_{c=1}^{|\mathcal{C}|}\Bigg(\sum_{i=1}^{n^c_d}\lambda^*_d(t_i^{c,d}) - \int_{0}^{T_c}\lambda^*_d(\tau)d\tau\Bigg), \end{align}

where \(T_c\) is the observation window of the sequence \(\mathcal{S}^c\). More...

 
virtual double Intensity (const double &t, const Sequence &data, Eigen::VectorXd &intensity_dim)=0
 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)=0
 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)=0
 Returns the integral of the intensity function \(\int_{a}^b\lambda^*(\tau)d\tau\) where \(a = lower\) and \(b = upper\). More...
 
virtual void Gradient (const unsigned &k, Eigen::VectorXd &gradient)=0
 Returns the gradient w.r.t. the model parameters on the k-th sequence. 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...
 
virtual double PredictNextEventTime (const Sequence &data, const unsigned &num_simulations)=0
 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...
 

Protected Member Functions

void InitializeDimension (const std::vector< Sequence > &data)
 

Protected Attributes

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

IProcess defines a general interface for each specific point process.

The IProcess class provides the basic interface for a general point process, which can be used for general fitting, simulation, and model checking. Customized point processes are supported by simply inheriting from this class.

Definition at line 20 of file Process.h.

Constructor & Destructor Documentation

IProcess::IProcess ( const unsigned &  n,
const unsigned &  num_dims 
)
inline

The constructor.

Parameters
[in]nthe number of parameters in total.
[in]num_dimsthe number of dimensions in the process.

Definition at line 54 of file Process.h.

Member Function Documentation

unsigned IProcess::GetNumDims ( )
inline

Return the number of dimensions in the process.

Returns
the number of dimensions in the process.

Definition at line 71 of file Process.h.

const Eigen::VectorXd& IProcess::GetParameters ( )
inline

Return the column vector of model parameters.

Returns
the column vector of model parameters.

Definition at line 65 of file Process.h.

virtual void IProcess::Gradient ( const unsigned &  k,
Eigen::VectorXd &  gradient 
)
pure 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.

Implemented in PlainHawkes, TerminatingProcessLearningTriggeringKernel, HawkesLearningTriggeringKernel, HawkesGeneralKernel, PlainTerminating, LowRankHawkesProcess, SelfInhibitingProcess, and Poisson.

void IProcess::InitializeDimension ( const std::vector< Sequence > &  data)
protected

Assign each event of each sequence in the variable data to the respective dimesion by initializing the varaible all_timestamp_per_dimension_.

Parameters
[in]datacollection of sequences.

Definition at line 8 of file Process.cc.

virtual double IProcess::Intensity ( const double &  t,
const Sequence data,
Eigen::VectorXd &  intensity_dim 
)
pure 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.

Implemented in PlainHawkes, TerminatingProcessLearningTriggeringKernel, HawkesLearningTriggeringKernel, PlainTerminating, HawkesGeneralKernel, LowRankHawkesProcess, SelfInhibitingProcess, and Poisson.

virtual double IProcess::IntensityIntegral ( const double &  lower,
const double &  upper,
const Sequence data 
)
pure 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\).

Implemented in PlainHawkes, TerminatingProcessLearningTriggeringKernel, HawkesLearningTriggeringKernel, PlainTerminating, HawkesGeneralKernel, LowRankHawkesProcess, SelfInhibitingProcess, and Poisson.

virtual double IProcess::IntensityUpperBound ( const double &  t,
const double &  L,
const Sequence data,
Eigen::VectorXd &  intensity_upper_dim 
)
pure 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].

Implemented in PlainHawkes, TerminatingProcessLearningTriggeringKernel, HawkesLearningTriggeringKernel, PlainTerminating, HawkesGeneralKernel, LowRankHawkesProcess, SelfInhibitingProcess, and Poisson.

virtual void IProcess::NegLoglikelihood ( double &  objvalue,
Eigen::VectorXd &  Gradient 
)
pure virtual

Returns the negative loglikelihood value and the gradient vectors. Given a collection of sequences \(\mathcal{C} = \{\mathcal{S}^c\}\) where \(\mathcal{S}^c = \{t^{c,d}_i\}^{d=1\dotso D}_{i=1\dotso n^c_d}\), and \(n^c_d\) is the number of events on dimension d in sequence \(\mathcal{S}^c\). The negative loglihood of observing \(\mathcal{C}\) is defined as

\begin{align} \ell(\mathcal{C}) = - \sum_{d=1}^D\frac{1}{|\mathcal{C}|}\sum_{c=1}^{|\mathcal{C}|}\Bigg(\sum_{i=1}^{n^c_d}\lambda^*_d(t_i^{c,d}) - \int_{0}^{T_c}\lambda^*_d(\tau)d\tau\Bigg), \end{align}

where \(T_c\) is the observation window of the sequence \(\mathcal{S}^c\).

Parameters
[out]objvaluethe objective function value, which is the negative log-likelihood value evaluated on the given set of sequences using the current model parameters.
[out]Gradientthe gradient vector w.r.t. the model parameters.

Implemented in PlainHawkes, TerminatingProcessLearningTriggeringKernel, HawkesLearningTriggeringKernel, HawkesGeneralKernel, PlainTerminating, LowRankHawkesProcess, SelfInhibitingProcess, and Poisson.

void IProcess::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.

Warning
Currently, it supports the ploting of up to four dimensions in total.
Parameters
[in]dataan input sequence of events.

Definition at line 26 of file Process.cc.

void IProcess::PlotIntensityFunction ( const Sequence data,
const unsigned &  dim_id 
)

Plots the intensity function and the associated event points of the dimension dim_id.

Parameters
[in]dataan input sequence of events.
[in]dim_idthe index of the dimension we want to plot.

Definition at line 87 of file Process.cc.

virtual double IProcess::PredictNextEventTime ( const Sequence data,
const unsigned &  num_simulations 
)
pure 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.

Implemented in PlainHawkes, TerminatingProcessLearningTriggeringKernel, HawkesLearningTriggeringKernel, PlainTerminating, HawkesGeneralKernel, LowRankHawkesProcess, SelfInhibitingProcess, and Poisson.

void IProcess::SetParameters ( const Eigen::VectorXd &  v)
inline

Set the model parameters.

Parameters
[in]vA column vector storing the new values for the model parameters.

Definition at line 77 of file Process.h.

Member Data Documentation

std::vector<std::vector<std::vector<double> > > IProcess::all_timestamp_per_dimension_
protected

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.

Definition at line 38 of file Process.h.

unsigned IProcess::num_dims_
protected

The total number of dimensions of the process.

Definition at line 33 of file Process.h.

Eigen::VectorXd IProcess::parameters_
protected

A column vector represents all model parameters of the process.

Definition at line 28 of file Process.h.


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