particles.state_space_models

State-space models as Python objects.

Overview

This module defines:

  1. the StateSpaceModel class, which lets you define a state-space model as a Python object;

  2. FeynmanKac classes that automatically define the Bootstrap, guided or auxiliary Feynman-Kac models associated to a given state-space model;

  3. several standard state-space models (stochastic volatility, bearings-only tracking, and so on).

The recommended import is:

from particles import state_space_models as ssms

For more details on state-space models and their properties, see Chapters 2 and 4 of the book.

Defining a state-space model

Consider the following (simplified) stochastic volatility model:

\[\begin{split}Y_t|X_t=x_t &\sim N(0, e^{x_t}) \\ X_t|X_{t-1}=x_{t-1} &\sim N(0, \rho x_{t-1}) \\ X_0 &\sim N(0, \sigma^2 / (1 - \rho^2))\end{split}\]

To define this particular model, we sub-class StateSpaceModel as follows:

import numpy as np
from particles import distributions as dists

class SimplifiedStochVol(ssms.StateSpaceModel):
    default_params = {'sigma': 1., 'rho': 0.8}  # optional
    def PY(self, t, xp, x):  # dist of Y_t at time t, given X_t and X_{t-1}
        return dists.Normal(scale=np.exp(x))
    def PX(self, t, xp):  # dist of X_t at time t, given X_{t-1}
        return dists.Normal(loc=self.mu + self.rho * (xp - self.mu),
                            scale=self.sigma)
    def PX0(self):  # dist of X_0
        return dists.Normal(scale=self.sigma / np.sqrt(1. - self.rho**2))

Then we define a particular object (model) by instantiating this class:

my_stoch_vol_model = SimplifiedStochVol(sigma=0.3, rho=0.9)

Hopefully, the code above is fairly transparent, but here are some noteworthy details:

  • probability distributions are defined through ProbDist objects, which are defined in module distributions. Most basic probability distributions are defined there; see module distributions for more details.

  • The class above actually defines a parametric class of models; in particular, self.sigma and self.rho are attributes of this class that are set when we define object my_stoch_vol_model. Default values for these parameters may be defined in a dictionary called default_params. When this dictionary is defined, any un-defined parameter will be replaced by its default value:

    default_stoch_vol_model = SimplifiedStochVol()  # sigma=1., rho=0.8
    
  • There is no need to define a __init__() method, as it is already defined by the parent class. (This parent __init__() simply takes care of the default parameters, and may be overridden if needed.)

Now that our state-space model is properly defined, what can we do with it? First, we may simulate states and data from it:

x, y = my_stoch_vol_model.simulate(20)

This generates two lists of length 20: a list of states, X_0, …, X_{19} and a list of observations (data-points), Y_0, …, Y_{19}.

Associated Feynman-Kac models

Now that our state-space model is defined, we obtain the associated Bootstrap Feynman-Kac model as follows:

my_fk_model = ssms.Bootstrap(ssm=my_stoch_vol_model, data=y)

That’s it! You are now able to run a bootstrap filter for this model:

my_alg = particles.SMC(fk=my_fk_model, N=200)
my_alg.run()

In case you are not clear about what are Feynman-Kac models, and how one may associate a Feynman-Kac model to a given state-space model, see Chapter 5 of the book.

To generate a guided Feynman-Kac model, we must provide proposal kernels (that is, Markov kernels that define how we simulate particles X_t at time t, given an ancestor X_{t-1}):

class StochVol_with_prop(StochVol):
    def proposal0(self, data):
        return dists.Normal(scale = self.sigma)
    def proposal(t, xp, data):  # a silly proposal
        return dists.Normal(loc=rho * xp + data[t], scale=self.sigma)

my_second_ssm = StochVol_with_prop(sigma=0.3)
my_better_fk_model = ssms.GuidedPF(ssm=my_second_ssm, data=y)
# then run a SMC as above

Voilà! You have now implemented a guided filter.

Of course, the proposal distribution above does not make much sense; we use it to illustrate how proposals may be defined. Note in particular that it depends on data, an object that represents the complete dataset. Hence the proposal kernel at time t may depend on y_t but also y_{t-1}, or any other datapoint.

For auxiliary particle filters (APF), one must in addition specify auxiliary functions, that is the (log of) functions \(\eta_t\) that modify the resampling probabilities (see Section 10.3.3 in the book):

class StochVol_with_prop_and_aux_func(StochVol_with_prop):
    def logeta(self, t, x, data):
        "Log of auxiliary function eta_t at time t"
        return -(x-data[t])**2

my_third_ssm = StochVol_with_prop_and_aux_func()
apf_fk_model = ssms.AuxiliaryPF(ssm=my_third_ssm, data=y)

Again, this particular choice does not make much sense, and is just given to show how to define an auxiliary function.

Already implemented state-space models

This module implements a few basic state-space models that are often used as numerical examples:

Class

Comments

StochVol

Basic, univariate stochastic volatility model

StochVolLeverage

Univariate stochastic volatility model with leverage

MVStochVol

Multivariate stochastic volatility model

BearingsOnly

Bearings-only tracking

Gordon_etal

Popular toy model often used as a benchmark

DiscreteCox

A discrete Cox model (Y_t|X_t is Poisson)

ThetaLogistic

Theta-logistic model from Population Ecology

Note

Linear Gaussian state-space models are implemented in module kalman; similarly hidden Markov models (state-space models with a finite state-space) are implemented in module hmm.

Classes

APFMixin()

AuxiliaryBootstrap([ssm, data])

Base class for auxiliary bootstrap particle filters

AuxiliaryPF([ssm, data])

Auxiliary particle filter for a given state-space model.

BearingsOnly(**kwargs)

Bearings-only tracking SSM.

Bootstrap([ssm, data])

Bootstrap Feynman-Kac formalism of a given state-space model.

DiscreteCox(**kwargs)

A discrete Cox model.

Gordon_etal(**kwargs)

Popular toy example that appeared initially in Gordon et al (1993).

GuidedPF([ssm, data])

Guided filter for a given state-space model.

MVStochVol(**kwargs)

Multivariate stochastic volatility model.

StateSpaceModel(**kwargs)

Base class for state-space models.

StochVol(**kwargs)

Univariate stochastic volatility model.

StochVolLeverage(**kwargs)

Univariate stochastic volatility model with leverage effect.

ThetaLogistic(**kwargs)

Theta-Logistic state-space model (used in Ecology).