birdepy.probability()

birdepy.probability(z0, zt, t, param, model='Verhulst 1', method='expm', **options)

Transition probabilities for continuous-time birth-and-death processes.

Parameters:
  • z0 (array_like) – States of birth-and-death process at time 0.

  • zt (array_like) – States of birth-and-death process at time(s) t.

  • t (array_like) – Elapsed time(s) (if has size greater than 1, then must be increasing).

  • param (array_like) – The parameters governing the evolution of the birth-and-death process. Array of real elements of size (n,), where ‘n’ is the number of parameters. These must be in the order given here).

  • model (string, optional) –

    Model specifying birth and death rates of process (see here). Should be one of:

    • ’Verhulst’ (default)

    • ’Ricker’

    • ’Beverton-Holt’

    • ’Hassell’

    • ’MS-S’

    • ’Moran’

    • ’pure-birth’

    • ’pure-death’

    • ’Poisson’

    • ’linear’

    • ’linear-migration’

    • ’MS-S’

    • ’M/M/inf’

    • ’loss-system’

    • ’custom’

    If set to ‘custom’, then arguments b_rate and d_rate must also be specified. See here for more information.

  • method (string, optional) –

    Transition probability approximation method. Should be one of (alphabetical ordering):

    • ’da’ (see here)

    • ’Erlang’(see here)

    • ’expm’ (default) (see here)

    • ’gwa’ (see here)

    • ’gwasa’ (see here)

    • ’ilt’ (see here)

    • ’oua’ (see here)

    • ’sim’ (see here)

    • ’uniform’ (see here)

  • options (dict, optional) – A dictionary of method specific options.

Returns:

transition_probability – An array of transition probabilities. If t has size bigger than 1, then the first coordinate corresponds to t, the second coordinate corresponds to z0 and the third coordinate corresponds to zt; for example if z0=[1,3,5,10], zt=[5,8] and t=[1,2,3], then transition_probability[2,0,1] corresponds to P(Z(3)=8 | Z(0)=1). If t has size 1 the first coordinate corresponds to z0 and the second coordinate corresponds to zt.

Return type:

ndarray

Examples

Approximate transition probability for a Verhulst model using various methods:

for method in ['da', 'Erlang', 'expm', 'gwa', 'gwasa', 'ilt', 'oua', 'uniform']:
    print(method, 'approximation:',
    bd.probability(19, 27, 1.0, [0.5, 0.3, 0.01, 0.01], model='Verhulst', method=method)[0][0])

Outputs:

da approximation: 0.016040426614336103
Erlang approximation: 0.0161337966847677
expm approximation: 0.016191045449304713
gwa approximation: 0.014646030484734228
gwasa approximation: 0.014622270048744283
ilt approximation: 0.01618465415009876
oua approximation: 0.021627234315268227
uniform approximation: 0.016191045442910168

Notes

Methods for computing transition probabilities and models are also described in [1]. If you use this function for published work, then please cite this paper.

For a text book treatment on the theory of birth-and-death processes see [2].

References

[View birdepy.probability() source code]

birdepy.probability(method=’da’)

birdepy.probability_da.probability_da(z0, zt, t, param, b_rate, d_rate, h_fun, k)

Transition probabilities for continuous-time birth-and-death processes using the diffusion approximation method.

To use this function call birdepy.probability() with method set to ‘da’:

bd.probability(z0, zt, t, param, method='da', k=10)

The parameters associated with this method (listed below) can be accessed using kwargs in birdepy.probability(). See documentation of birdepy.probability() for the main arguments.

Parameters:

k (int, optional) – 2**k + 1 equally spaced samples between observation times are used to compute integrals.

Examples

Approximate transition probability for a Verhulst model using diffusion approximation:

import birdepy as bd
bd.probability(19, 27, 1.0, [0.5, 0.3, 0.01, 0.01], model='Verhulst', method='da')[0][0])

Outputs:

0.002266101391343583

Notes

Estimation algorithms and models are also described in [1]. If you use this function for published work, then please cite this paper.

For a text book treatment on the theory of birth-and-death processes see [2].

For a method related to this one see [3].

References

[View birdepy.probability(method=”da”) source code]

birdepy.probability(method=’Erlang’)

birdepy.probability_Erlang.probability_Erlang(z0, zt, t, param, b_rate, d_rate, z_trunc, k)

Transition probabilities for continuous-time birth-and-death processes using the Erlangization method.

To use this function call birdepy.probability with method set to ‘Erlang’:

bd.probability(z0, zt, t, param, method='Erlang', z_trunc=(), k=1502)

The parameters associated with this method (listed below) can be accessed using kwargs in birdepy.probability(). See documentation of birdepy.probability() for the main arguments.

Parameters:
  • z_trunc (array_like, optional) – Truncation thresholds, i.e., minimum and maximum states of process considered. Array of real elements of size (2,) by default z_trunc=[z_min, z_max] where z_min=max(0, min(z0, zt) - 100) and z_max=max(z0, zt) + 100

  • k (int, optional) – Number of terms to include in approximation to probability.

Examples

Approximate transition probability for a Verhulst model using Erlangization:

import birdepy as bd
bd.probability(19, 27, 1.0, [0.5, 0.3, 0.02, 0.01], model='Verhulst', method='Erlang')[0][0]

Outputs:

0.002731464736623327

Notes

Estimation algorithms and models are also described in [1]. If you use this function for published work, then please cite this paper.

For a text book treatment on the theory of birth-and-death processes see [2].

For a method related to this one see [3].

References

[View birdepy.probability(method=”Erlang”) source code]

birdepy.probability(method=’expm’)

birdepy.probability_expm.probability_expm(z0, zt, t, param, b_rate, d_rate, z_trunc)

Transition probabilities for continuous-time birth-and-death processes using the matrix exponential method.

To use this function call birdepy.probability() with method set to ‘expm’:

birdepy.probability(z0, zt, t, param, method='expm', z_trunc=())

The parameters associated with this method (listed below) can be accessed using kwargs in birdepy.probability(). See documentation of birdepy.probability() for the main arguments.

Parameters:

z_trunc (array_like, optional) – Truncation thresholds, i.e., minimum and maximum states of process considered. Array of real elements of size (2,) by default z_trunc=[z_min, z_max] where z_min=max(0, min(z0, zt) - 100) and z_max=max(z0, zt) + 100

Examples

Approximate transition probability for a Verhulst model using a matrix exponential:

import birdepy as bd
bd.probability(19, 27, 1.0, [0.5, 0.3, 0.02, 0.01], model='Verhulst', method='expm')[0][0]

Outputs:

0.0027414224836612463

References

[View birdepy.probability(method=”expm”) source code]

birdepy.probability(method=’gwa’)

birdepy.probability_gwa.probability_gwa(z0, zt, t, param, b_rate, d_rate, anchor)

Transition probabilities for continuous-time birth-and-death processes using the Galton-Watson approximation method.

To use this function call birdepy.probability() with method set to ‘gwa’:

birdepy.probability(z0, zt, t, param, method='gwa', anchor='midpoint')

The parameters associated with this method (listed below) can be accessed using kwargs in birdepy.probability(). See documentation of birdepy.probability() for the main arguments.

Parameters:

anchor (string, optional) – Determines which state z is used to determine the linear approximation. Should be one of: ‘initial’ (z0 is used), ‘midpoint’ (default, 0.5*(z0+zt) is used), ‘terminal’ (zt is used), ‘max’ (max(z0, zt) is used), or ‘min’ (min(z0, zt) is used).

Examples

Approximate transition probability for a Verhulst model using a linear approximation:

import birdepy as bd
bd.probability(19, 27, 1.0, [0.5, 0.3, 0.02, 0.01], model='Verhulst', method='gwa')[0][0]

Outputs:

0.00227651766770292

Notes

Estimation algorithms and models are also described in [1]. If you use this function for published work, then please cite this paper.

For a text book treatment on the theory of birth-and-death processes see [2].

References

[View birdepy.probability(method=”gwa”) source code]

birdepy.probability(method=’gwasa’)

birdepy.probability_gwasa.probability_gwasa(z0, zt, t, param, b_rate, d_rate, anchor)

Transition probabilities for continuous-time birth-and-death processes using the Galton-Watson saddlepoint approximation method.

To use this function call birdepy.probability() with method set to ‘gwasa’:

birdepy.probability(z0, zt, t, param, method='gwasa', anchor='midpoint')

The parameters associated with this method (listed below) can be accessed using kwargs in birdepy.probability(). See documentation of birdepy.probability() for the main arguments.

Parameters:

anchor (string, optional) – Determines which state z is used to determine the linear approximation. Should be one of: ‘initial’ (z0 is used), ‘midpoint’ (default, 0.5*(z0+zt) is used), ‘terminal’ (zt is used), ‘max’ (max(z0, zt) is used), or ‘min’ (min(z0, zt) is used).

Examples

Approximate transition probability for a Verhulst model using the saddlepoint approximation to a linear approximation:

import birdepy as bd
bd.probability(19, 27, 1.0, [0.5, 0.3, 0.02, 0.01], model='Verhulst', method='gwasa')[0][0]

Outputs:

0.002271944691896704

Notes

Estimation algorithms and models are also described in [1]. If you use this function for published work, then please cite this paper.

For a text book treatment on the theory of birth-and-death processes see [2].

References

[View birdepy.probability(method=”gwasa”) source code]

birdepy.probability(method=’ilt’)

birdepy.probability_ilt.probability_ilt(z0, zt, t, param, b_rate, d_rate, lentz_eps, laplace_method, k)

Transition probabilities for continuous-time birth-and-death processes using the inverse Laplace transform method.

To use this function call birdepy.probability() with method set to ‘ilt’:

birdepy.probability(z0, zt, t, param, method='ilt', eps=1e-6,
                    laplace_method='cme-talbot', precision=100, k=25)

The parameters associated with this method (listed below) can be accessed using kwargs in birdepy.probability(). See documentation of birdepy.probability() for the main arguments.

Parameters:
  • lentz_eps (scalar, optional) – Termination threshold for Lentz algorithm computation of Laplace domain functions.

  • laplace_method (string, optional) – Numerical inverse Laplace transform algorithm to use. Should be one of: ‘cme-talbot’ (default), ‘cme’, ‘euler’, ‘gaver’, ‘talbot’, ‘stehfest’, ‘dehoog’, ‘cme-mp’ or ‘gwr’.

  • precision (int, optional) – Numerical precision (only used for methods that invoke mpmath).

  • k (int, optional) – Maximum number of terms used for Laplace transform numerical inversion. Only applicable if argument ‘laplace_method’ is set to ‘cme-talbot’, ‘cme’, ‘euler’, ‘gaver’ or ‘cme-mp’. See https://www.inverselaplace.org/ for more information.

Examples

Approximate transition probability for a Verhulst model using numerical inverse Laplace transform:

import birdepy as bd
bd.probability(19, 27, 1.0, [0.5, 0.3, 0.02, 0.01], model='Verhulst', method='ilt')[0][0]

Outputs:

0.0027403264310572615

Notes

Estimation algorithms and models are also described in [1]. If you use this function for published work, then please cite this paper.

For a text book treatment on the theory of birth-and-death processes see [2].

For more details on this method see [3].

References

[View birdepy.probability(method=”ilt”) source code]

birdepy.probability(method=’oua’)

birdepy.probability_oua.probability_oua(z0, zt, t, param, b_rate, d_rate, h_fun, zf_bld)

Transition probabilities for continuous-time birth-and-death processes using the Ornstein-Uhlenbeck approximation method.

To use this function call birdepy.probability() with method set to ‘oua’:

birdepy.probability(z0, zt, t, param, method='oua')

This function does not have any arguments which are not already described by the documentation of birdepy.probability()

Examples

Approximate transition probability for a Verhulst model using Ornstein–Uhlenbeck approximation:

import birdepy as bd
bd.probability(19, 27, 1.0, [0.5, 0.3, 0.02, 0.01], model='Verhulst', method='oua')[0][0]

Outputs:

0.0018882966813798246

Notes

Estimation algorithms and models are also described in [1]. If you use this function for published work, then please cite this paper.

For a text book treatment on the theory of birth-and-death processes see [2].

For a method related to this one see [3].

References

[View birdepy.probability(method=”oua”) source code]

birdepy.probability(method=’sim’)

birdepy.probability_sim.probability_sim(z0, zt, t, param, b_rate, d_rate, k, sim_method, tau, seed)

Transition probabilities for continuous-time birth-and-death processes using crude Monte Carlo simulation.

To use this function call birdepy.probability() with method set to ‘sim’:

birdepy.probability(z0, zt, t, param, method='sim', k=10**5, seed=None)

The parameters associated with this method (listed below) can be accessed using kwargs in birdepy.probability(). See documentation of birdepy.probability() for the main arguments.

Parameters:
  • k (int, optional) – Number of samples used to generate each probability estimate. The total number of samples used will be z0.size * k.

  • seed (int, Generator, optional) – If seed is not specified the random numbers are generated according to np.random.default_rng(). If seed is an int, random numbers are generated according to np.random.default_rng(seed). If seed is a Generator, then that object is used. See here for more information.

  • sim_method (string, optional) –

    Simulation algorithm used to generate samples (see here). Should be one of:

    • ’exact’ (default)

    • ’ea’

    • ’ma’

    • ’gwa’

Examples

>>> import birdepy as bd
>>> bd.probability(19, 27, 1.0, [0.5, 0.3, 0.02, 0.01], model='Verhulst', method='sim', k=10**5,
...                seed=2021)[0][0]
0.00294

Notes

Estimation algorithms and models are also described in [1]. If you use this function for published work, then please cite this paper.

For a text book treatment on the theory of birth-and-death processes see [2].

For more details on this method see [3].

References

[View birdepy.probability(method=”sim”) source code]

birdepy.probability(method=’uniform’)

birdepy.probability_uniform.probability_uniform(z0, zt, t, param, b_rate, d_rate, z_trunc, k)

Transition probabilities for continuous-time birth-and-death processes using the uniformization method.

To use this function call birdepy.probability() with method set to ‘uniform’:

birdepy.probability(z0, zt, t, param, method='uniform', k=1000, z_trunc=())

The parameters associated with this method (listed below) can be accessed using kwargs in birdepy.probability(). See documentation of birdepy.probability() for the main arguments.

Parameters:
  • z_trunc (array_like, optional) – Truncation thresholds, i.e., minimum and maximum states of process considered. Array of real elements of size (2,) by default z_trunc=[z_min, z_max] where z_min=max(0, min(z0, zt) - 100) and z_max=max(z0, zt) + 100.

  • k (int, optional) – Number of terms to include in approximation to probability.

Examples

Approximate transition probability for a Verhulst model using uniformization:

import birdepy as bd
bd.probability(19, 27, 1.0, [0.5, 0.3, 0.02, 0.01], model='Verhulst', method='uniform')[0][0]

Outputs:

0.002741422482539626

Notes

Estimation algorithms and models are also described in [1]. If you use this function for published work, then please cite this paper.

For a text book treatment on the theory of birth-and-death processes see [2].

For more information on this method see [3] and [4].

References

[View birdepy.probability(method=”uniform”) source code]