nifty8.re.gauss_markov module#

class GaussMarkovProcess(process: Callable, x0: float | jax.Array | nifty8.re.model.LazyModel, dt: float | jax.Array, name='xi', N_steps: int = None, **kwargs)[source]#

Bases: Model

__init__(process: Callable, x0: float | Array | LazyModel, dt: float | Array, name='xi', N_steps: int = None, **kwargs)[source]#

Wrap a callable and associate it with a domain.

Parameters:
  • call (callable) – Method acting on objects of type domain.

  • domain (tree-like structure of ShapeWithDtype, optional) – PyTree of objects with a shape and dtype attribute. Inferred from init if not specified.

  • target (tree-like structure of ShapeWithDtype, optional) – PyTree of objects with a shape and dtype attribute akin to the output of call. Inferrred from call and domain if not set.

  • init (callable, optional) – Initialization method taking a PRNG key as first argument and creating an object of type domain. Inferred from domain assuming a white normal prior if not set.

  • white_init (bool, optional) – If True, the domain is set to a white normal prior. Defaults to False.

IntegratedWienerProcess(x0: tuple | Array | LazyModel, sigma: tuple | float | Array | LazyModel, dt: float | Array, name: str = 'iwp', asperity: tuple | float | Array | LazyModel = None, N_steps: int = None)[source]#

Implements the Integrated Wiener Process.

The generalized IWP in continuous time takes the form:

d/dt x_t = y_t + sigma * asperity xi^1_t , \
d/dt y_t = sigma * xi^2_t

where xi^i_t are continuous time white noise processes. This is a standard IWP in x in the case that asperity is zero (None). If it is non-zero, a WP with relative strength set by asperity is added to the dynamics of the IWP.

For information on the parameters, please refer to WienerProcess. The asperity parameter may be defined analogously to sigma.

This is also the process used in the CorrelatedField to describe the deviations of the power spectra from power-laws on a double logarithmic scale. See also re.correlated_field.CorrelatedFieldMaker and references there for further information.

Notes:#

sigma and asperity may also be sequences. See notes on WienerProcess for further information.

OrnsteinUhlenbeckProcess(sigma: tuple | float | Array | LazyModel, gamma: tuple | float | Array | LazyModel, dt: float | Array, name: str = 'oup', x0: tuple | float | LazyModel = None, N_steps: int = None)[source]#

Implements the Ornstein Uhlenbeck process (OUP).

The stochastic differential equation of the OUP takes the form:

d/dt x_t + gamma x_t = sigma xi_t

where xi_t is continuous time white noise.

For information on the parameters, please refer to WienerProcess. The gamma parameter may be defined analogously to sigma. Unlike the WP and IWP, the OUP is a proper process, i.E. allows for a proper steady state solution. Therefore, x0 may not be provided at all, in which case it defaults to a random variable drawn from the steady state distribution of the OUP.

Notes:#

sigma and gamma may also be sequences. See notes on WienerProcess for further information.

WienerProcess(x0: tuple | float | LazyModel, sigma: tuple | float | Array | LazyModel, dt: float | Array, name: str = 'wp', N_steps: int = None)[source]#

Implements the Wiener process (WP).

The WP in continuous time takes the form:

d/dt x_t = sigma xi_t ,

where xi_t is continuous time white noise.

Parameters:#

x0: tuple, float, or LazyModel

Initial position of the WP. Can be passed as a fixed value, or a generative Model. Passing a tuple is a shortcut to set a normal prior with mean and std equal to the first and second entry of the tuple respectively on x0.

sigma: tuple, float, Array, LazyModel

Standard deviation of the WP. Analogously to x0 may also be passed on as a model. May also be passed as a sequence of length equal to dt in which case a different sigma is used for each time interval.

dt: float or Array of float

Stepsizes of the process. In case it is a single float, N_steps must be provided to indicate the number of steps taken.

name: str

Name of the key corresponding to the parameters of the WP. Default wp.

N_steps: int (optional)

Option to set the number of steps in case dt is a scalar.

Notes:#

In case sigma is time dependent, i.E. passed on as a sequence of length equal to xi, it is assumed to be constant within each timebin, i.E. sigma_t = sigma_i for t_i <= t < t_{i+1}.

discrete_gauss_markov_process(xi: Array, x0: Array, drift: Array, diffamp: Array)[source]#

Generator for a Gauss-Markov process (GMP).

Given the discrete transition probabilities via the drift and diffamp matrices, this function produces a series res of the form:

res_{i+1} = drift_i @ res_i + diffamp_i @ xi_i

which corresponds to a random realization of a GMP in case xi is standard normal distributed.

Parameters:#

xi: Array

Input random variables for the GMP.

x0: Array

Initial state of the process. The first entry of the result res_0 is always equal to x0.

drift: Array

Matrix or sequence of matrices that describe the mean of the transition probabilities of the GMP.

diffamp: Array

Matrix or sequence of matrices that correspond to the amplitude of the transition probabilities. For the GMP they imply that their matrix product diffamp_i @ diffamp_i.T is equal to the covariance of the transition probability iff xi is standard normal.

returns:
  • Sequence of vectors of the Gauss-Markov series. The sequence is returned

  • as an array where the first axis corresponds to the steps of the

  • sequence and the second axis are the entries of the vector.

Notes:#

In case the sequence xi has length N (i.E. xi.shape[0] = N) the result sequence is of length N+1 as x0 is returned as the first entry of the final sequence, followed by N transitions.

integrated_wiener_process(xi: Array, x0: Array, sigma: Array, dt: Array, asperity: float | Array = None)[source]#

Implements the (generalized) Integrated Wiener process (IWP).

ornstein_uhlenbeck_process(xi: Array, x0: float, sigma: float | Array, gamma: float | Array, dt: float | Array)[source]#

Implements the Ornstein Uhlenbeck process (OUP).

scalar_gauss_markov_process(xi, x0, drift, diffamp)[source]#

Simple wrapper of discrete_gm_general for 1D scalar processes.

wiener_process(xi: Array, x0: float | Array, sigma: float | Array, dt: float | Array)[source]#

Implements the Wiener process (WP).