pyselfi.lotkavolterra.selfi module

SELFI for inference of the Lotka-Volterra model: main class.

class pyselfi.lotkavolterra.selfi.lotkavolterra_selfi(fname, pool_prefix, pool_suffix, prior, blackbox, theta_0, Ne, Ns, Delta_theta, phi_obs, LVsimulator)[source]

Bases: selfi

Main class for Lotka-Volterra inference with SELFI. Child of the main class of the SELFI code.

Variables
  • fname (str) – name and address of the selfi output file on disk

  • pool_prefix (str) – address and prefix of the simulation pools

  • pool_suffix (str) – suffix of the simulation pools

  • prior (prior) – the prior, defined as explained in the SELFI documentation

  • blackbox (blackbox) – the blackbox simulator, defined as explained in the SELFI documentation

  • theta_0 (array, double, dimension=S) – expansion point in parameter space

  • Ne (int) – number of simulations at the expansion point, to compute the covariance matrix

  • Ns (int) – number of simulations per expansion direction in parameter space. (Ne may be larger than Ns, in which case only the first Ns simulations at the expansion point are used to compute the gradient)

  • Delta_theta (double) – step size for finite differencing to compute the blackbox gradient in parameter space

  • phi_obs (array, double, dimension=P) – the vector of observed data

  • LVsimulator (LVsimulator) – a Lotka-Volterra simulator (instance of the class LVsimulator)

FisherRao_distance(omega_tilde1, omega_tilde2)[source]

Computes the Fisher-Rao distance between two compressed summaries \(\tilde{\boldsymbol{\omega}}_1\) and \(\tilde{\boldsymbol{\omega}}_2\).

Parameters
  • omega_tilde1 (array, double, dimension=4) – first vector of compressed summaries

  • omega_tilde2 (array, double, dimension=4) – second vector of compressed summaries

Returns

dist – Fisher-Rao distance between omega_tilde1 and omega_tilde2, using the Fisher matrix of the problem

Return type

double

Mahalanobis_ms_check(theta=None)[source]

Computes a quantitative check for model misspecification: the Mahalanobis distance between the expansion point (prior mean) \(\theta_\mathrm{0}\) and a point \(\theta\) (usually the posterior mean \(\gamma\)), using the prior inverse covariance matrix \(\textbf{S}\).

Parameters

theta (array, double, dimension=S, optional, default=posterior.mean) – a point in parameter space

Returns

dist – Mahalanobis distance between prior and theta/posterior

Return type

double

compute_Fisher_matrix()[source]

Computes the Fisher matrix of the problem, once all necessary quantities are here.

compute_grad_f_omega(t, t_s, tres, Delta_omega)[source]

Computes the gradient of the observations with respect to the input parameters \(\boldsymbol{\omega}\), \(\nabla_\boldsymbol{\omega} \boldsymbol{f}\). First, compute \(\nabla_\boldsymbol{\omega} \boldsymbol{\theta}\) using compute_dtheta_domega() of the LVsimulator, then combine with the usual SELFI blackbox gradient \(\nabla \boldsymbol{f}\) so that \(\nabla_\boldsymbol{\omega} \boldsymbol{f} = \nabla \boldsymbol{f} \cdot \nabla_\boldsymbol{\omega} \boldsymbol{\theta}\).

Parameters
  • t (array, double, dimension=n) – array of time values where we approximate X and Y values timestep at each iteration is given by t[n+1] - t[n].

  • t_s (array, double, dimension=S) – array of the support time values where theta is defined

  • tres (double) – time resolution, defined such that: t = np.linspace(0,tmax-1,tmax*tres+1) t_s = np.arange(tmax)

  • Delta_omega (double) – step size for finite differences

compute_grad_f_omega_direct(t, t_s, tres, Delta_omega, pool_prefix, N)[source]

Computes the gradient of the observations with respect to the input parameters \(\boldsymbol{\omega}\), \(\nabla_\boldsymbol{\omega} \boldsymbol{f}\), using direct simulations (instead of the SELFI framework).

Parameters
  • t (array, double, dimension=n) – array of time values where we approximate X and Y values timestep at each iteration is given by t[n+1] - t[n].

  • t_s (array, double, dimension=S) – array of the support time values where theta is defined

  • tres (double) – time resolution, defined such that: t = np.linspace(0,tmax-1,tmax*tres+1) t_s = np.arange(tmax)

  • Delta_omega (double) – step size for finite differences

  • pool_prefix (str) – prefix for the filename of simulation pools

  • N (int) – number of realizations required

load(fname)[source]

Loads the data compression information from the SELFI input file.

Parameters

fname (str) – input filename

loglikelihood_hyperparams(theta_fid, alpha_norm, t_smooth, t_chaos)[source]

Returns the log pdf of the likelihood for hyperparameters \(\alpha_\mathrm{norm}\), \(t_\mathrm{smooth}\), \(t_\mathrm{chaos}\), defined as the sum of the log posterior pdf at \(\boldsymbol{\theta}_\mathrm{fid,i}\) for \(0 \leq i \leq N_\mathrm{fid}-1\):

\(-2 \log p(\alpha_\mathrm{norm}, t_\mathrm{smooth}, t_\mathrm{chaos}) \equiv \sum_{i=1}^{N_\mathrm{fid}} [ \log \left| 2\pi\boldsymbol{\Gamma} \right| + (\boldsymbol{\theta}_\mathrm{fid,i}-\boldsymbol{\gamma})^\intercal \boldsymbol{\Gamma}^{-1} (\boldsymbol{\theta}_\mathrm{fid,i}-\boldsymbol{\gamma}) ]\).

Parameters
  • theta_fiducial (array, double, dimension=(N_fid,S)) – set of fiducial points in parameter space to optimize the prior

  • alpha_norm (double) – alpha_norm value where to evaluate the likelihood

  • t_smooth (double) – t_smooth value where to evaluate the likelihood

  • t_chaos (double) – t_chaos value where to evaluate the likelihood

Returns

logpdf – log pdf of the likelihood at the evaluation point

Return type

double

logposterior_hyperparameters(theta_fid, alpha_norm, t_smooth, t_chaos)[source]

Returns the log pdf of the posterior for hyperparameters \(\alpha_\mathrm{norm}\), \(t_\mathrm{smooth}\), \(t_\mathrm{chaos}\). In this case, it is the same as loglikelihood_hyperparams(), but it would be possible add a log prior here.

Parameters
  • theta_fiducial (array, double, dimension=(N_fid,S)) – set of fiducial points in parameter space to optimize the prior

  • alpha_norm (double) – alpha_norm value where to evaluate the likelihood

  • t_smooth (double) – t_smooth value where to evaluate the likelihood

  • t_chaos (double) – t_chaos value where to evaluate the likelihood

Returns

logpdf – log pdf of the posterior at the evaluation point

Return type

double

optimize_prior(theta_fid, x0=None, alpha_norm_min=1e-10, alpha_norm_max=20, t_smooth_min=1e-10, t_smooth_max=20, t_chaos_min=1e-10, t_chaos_max=20, method='L-BFGS-B', options={'disp': False, 'eps': 1e-06, 'ftol': 1e-10, 'gtol': 1e-10, 'maxiter': 50})[source]

Optimizes the SELFI prior. See section II.E. in Leclercq et al. (2019).

Parameters
  • theta_fid (array, double, dimension=S) – fiducial point in parameter space to optimize the prior

  • x0 (array, double, dimension=3, optional, default=(prior.get_alpha_norm(),prior.get_t_smooth(),prior.get_t_chaos())) – starting point in parameter space

  • alpha_norm_min (double, optional, default=1e-10) – boundary in parameter space

  • alpha_norm_max (double, optional, default=20) – boundary in parameter space

  • t_smooth_min (double, optional, default=1e-10) – boundary in parameter space

  • t_smooth_max (double, optional, default=20) – boundary in parameter space

  • t_chaos_min (double, optional, default=1e-10) – boundary in parameter space

  • t_chaos_max (double, optional, default=20) – boundary in parameter space

  • method (str, optional, default=’L-BFGS-B’) – optimization method. See documentation of scipy.optimize.minimize

  • options (dictionary, optional, default={‘maxiter’:50, ‘ftol’:1e-10, ‘gtol’:1e-10, ‘eps’:1e-6, ‘disp’:False}) – optimization options. See documentation of scipy.optimize.minimize

save()[source]

Saves the data compression information to the SELFI output file.

score_compression(phi)[source]

Compress a data vector \(\boldsymbol{\Phi}\) using score compression (the optimal quasi-maximum likelihood estimator).

Parameters

phi (array, double, dimension=P) – input data vector to be compressed

Returns

omega_tilde – output compressed summaries

Return type

array, double, dimension=4